MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
coefficient.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_COEFFICIENT
13 #define MFEM_COEFFICIENT
14 
15 #include <functional>
16 
17 #include "../config/config.hpp"
18 #include "../linalg/linalg.hpp"
19 #include "intrules.hpp"
20 #include "eltrans.hpp"
21 
22 namespace mfem
23 {
24 
25 class Mesh;
26 
27 #ifdef MFEM_USE_MPI
28 class ParMesh;
29 #endif
30 
31 
32 /** @brief Base class Coefficients that optionally depend on space and time.
33  These are used by the BilinearFormIntegrator, LinearFormIntegrator, and
34  NonlinearFormIntegrator classes to represent the physical coefficients in
35  the PDEs that are being discretized. This class can also be used in a more
36  general way to represent functions that don't necessarily belong to a FE
37  space, e.g., to project onto GridFunctions to use as initial conditions,
38  exact solutions, etc. See, e.g., ex4 or ex22 for these uses. */
40 {
41 protected:
42  double time;
43 
44 public:
45  Coefficient() { time = 0.; }
46 
47  /// Set the time for time dependent coefficients
48  void SetTime(double t) { time = t; }
49 
50  /// Get the time for time dependent coefficients
51  double GetTime() { return time; }
52 
53  /** @brief Evaluate the coefficient in the element described by @a T at the
54  point @a ip. */
55  /** @note When this method is called, the caller must make sure that the
56  IntegrationPoint associated with @a T is the same as @a ip. This can be
57  achieved by calling T.SetIntPoint(&ip). */
58  virtual double Eval(ElementTransformation &T,
59  const IntegrationPoint &ip) = 0;
60 
61  /** @brief Evaluate the coefficient in the element described by @a T at the
62  point @a ip at time @a t. */
63  /** @note When this method is called, the caller must make sure that the
64  IntegrationPoint associated with @a T is the same as @a ip. This can be
65  achieved by calling T.SetIntPoint(&ip). */
67  const IntegrationPoint &ip, double t)
68  {
69  SetTime(t);
70  return Eval(T, ip);
71  }
72 
73  virtual ~Coefficient() { }
74 };
75 
76 
77 /// A coefficient that is constant across space and time
79 {
80 public:
81  double constant;
82 
83  /// c is value of constant function
84  explicit ConstantCoefficient(double c = 1.0) { constant=c; }
85 
86  /// Evaluate the coefficient at @a ip.
87  virtual double Eval(ElementTransformation &T,
88  const IntegrationPoint &ip)
89  { return (constant); }
90 };
91 
92 /** @brief A piecewise constant coefficient with the constants keyed
93  off the element attribute numbers. */
95 {
96 private:
97  Vector constants;
98 
99 public:
100 
101  /// Constructs a piecewise constant coefficient in NumOfSubD subdomains
102  explicit PWConstCoefficient(int NumOfSubD = 0) : constants(NumOfSubD)
103  { constants = 0.0; }
104 
105  /// Construct the constant coefficient using a vector of constants.
106  /** @a c should be a vector defined by attributes, so for region with
107  attribute @a i @a c[i-1] is the coefficient in that region */
109  { constants.SetSize(c.Size()); constants=c; }
110 
111  /// Update the constants with vector @a c.
112  void UpdateConstants(Vector &c) { constants.SetSize(c.Size()); constants=c; }
113 
114  /// Return a reference to the i-th constant
115  double &operator()(int i) { return constants(i-1); }
116 
117  /// Set the constants for all attributes to constant @a c.
118  void operator=(double c) { constants = c; }
119 
120  /// Returns the number of constants representing different attributes.
121  int GetNConst() { return constants.Size(); }
122 
123  /// Evaluate the coefficient.
124  virtual double Eval(ElementTransformation &T,
125  const IntegrationPoint &ip);
126 };
127 
128 /// A general function coefficient
130 {
131 protected:
132  std::function<double(const Vector &)> Function;
133  std::function<double(const Vector &, double)> TDFunction;
134 
135 public:
136  /// Define a time-independent coefficient from a std function
137  /** \param F time-independent std::function */
138  FunctionCoefficient(std::function<double(const Vector &)> F)
139  : Function(std::move(F))
140  { }
141 
142  /// Define a time-dependent coefficient from a std function
143  /** \param TDF time-dependent function */
144  FunctionCoefficient(std::function<double(const Vector &, double)> TDF)
145  : TDFunction(std::move(TDF))
146  { }
147 
148  /// (DEPRECATED) Define a time-independent coefficient from a C-function
149  /** @deprecated Use the method where the C-function, @a f, uses a const
150  Vector argument instead of Vector. */
151  MFEM_DEPRECATED FunctionCoefficient(double (*f)(Vector &))
152  {
153  Function = reinterpret_cast<double(*)(const Vector&)>(f);
154  TDFunction = NULL;
155  }
156 
157  /// (DEPRECATED) Define a time-dependent coefficient from a C-function
158  /** @deprecated Use the method where the C-function, @a tdf, uses a const
159  Vector argument instead of Vector. */
160  MFEM_DEPRECATED FunctionCoefficient(double (*tdf)(Vector &, double))
161  {
162  Function = NULL;
163  TDFunction = reinterpret_cast<double(*)(const Vector&,double)>(tdf);
164  }
165 
166  /// Evaluate the coefficient at @a ip.
167  virtual double Eval(ElementTransformation &T,
168  const IntegrationPoint &ip);
169 };
170 
171 class GridFunction;
172 
173 /// Coefficient defined by a GridFunction. This coefficient is mesh dependent.
175 {
176 private:
177  const GridFunction *GridF;
178  int Component;
179 
180 public:
181  GridFunctionCoefficient() : GridF(NULL), Component(1) { }
182  /** Construct GridFunctionCoefficient from a given GridFunction, and
183  optionally specify a component to use if it is a vector GridFunction. */
184  GridFunctionCoefficient (const GridFunction *gf, int comp = 1)
185  { GridF = gf; Component = comp; }
186 
187  /// Set the internal GridFunction
188  void SetGridFunction(const GridFunction *gf) { GridF = gf; }
189 
190  /// Get the internal GridFunction
191  const GridFunction * GetGridFunction() const { return GridF; }
192 
193  /// Evaluate the coefficient at @a ip.
194  virtual double Eval(ElementTransformation &T,
195  const IntegrationPoint &ip);
196 };
197 
198 
199 /** @brief A coefficient that depends on 1 or 2 parent coefficients and a
200  transformation rule represented by a C-function.
201 
202  \f$ C(x,t) = T(Q1(x,t)) \f$ or \f$ C(x,t) = T(Q1(x,t), Q2(x,t)) \f$
203 
204  where T is the transformation rule, and Q1/Q2 are the parent coefficients.*/
206 {
207 private:
208  Coefficient * Q1;
209  Coefficient * Q2;
210  double (*Transform1)(double);
211  double (*Transform2)(double,double);
212 
213 public:
214  TransformedCoefficient (Coefficient * q,double (*F)(double))
215  : Q1(q), Transform1(F) { Q2 = 0; Transform2 = 0; }
217  double (*F)(double,double))
218  : Q1(q1), Q2(q2), Transform2(F) { Transform1 = 0; }
219 
220  /// Evaluate the coefficient at @a ip.
221  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
222 };
223 
224 /** @brief Delta function coefficient optionally multiplied by a weight
225  coefficient and a scaled time dependent C-function.
226 
227  \f$ F(x,t) = w(x,t) s T(t) d(x - xc) \f$
228 
229  where w is the optional weight coefficient, @a s is a scale factor
230  T is an optional time-dependent function and d is a delta function.
231 
232  WARNING this cannot be used as a normal coefficient. The usual Eval
233  method is disabled. */
235 {
236 protected:
237  double center[3], scale, tol;
239  int sdim;
240  double (*tdf)(double);
241 
242 public:
243 
244  /// Construct a unit delta function centered at (0.0,0.0,0.0)
246  {
247  center[0] = center[1] = center[2] = 0.; scale = 1.; tol = 1e-12;
248  weight = NULL; sdim = 0; tdf = NULL;
249  }
250 
251  /// Construct a delta function scaled by @a s and centered at (x,0.0,0.0)
252  DeltaCoefficient(double x, double s)
253  {
254  center[0] = x; center[1] = 0.; center[2] = 0.; scale = s; tol = 1e-12;
255  weight = NULL; sdim = 1; tdf = NULL;
256  }
257 
258  /// Construct a delta function scaled by @a s and centered at (x,y,0.0)
259  DeltaCoefficient(double x, double y, double s)
260  {
261  center[0] = x; center[1] = y; center[2] = 0.; scale = s; tol = 1e-12;
262  weight = NULL; sdim = 2; tdf = NULL;
263  }
264 
265  /// Construct a delta function scaled by @a s and centered at (x,y,z)
266  DeltaCoefficient(double x, double y, double z, double s)
267  {
268  center[0] = x; center[1] = y; center[2] = z; scale = s; tol = 1e-12;
269  weight = NULL; sdim = 3; tdf = NULL;
270  }
271 
272  /// Set the center location of the delta function.
273  void SetDeltaCenter(const Vector& center);
274 
275  /// Set the scale value multiplying the delta function.
276  void SetScale(double s_) { scale = s_; }
277 
278  /// Set a time-dependent function that multiplies the Scale().
279  void SetFunction(double (*f)(double)) { tdf = f; }
280 
281  /** @brief Set the tolerance used during projection onto GridFunction to
282  identify the Mesh vertex where the Center() of the delta function
283  lies. (default 1e-12)*/
284  void SetTol(double tol_) { tol = tol_; }
285 
286  /// Set a weight Coefficient that multiplies the DeltaCoefficient.
287  /** The weight Coefficient multiplies the value returned by EvalDelta() but
288  not the value returned by Scale().
289  The weight Coefficient is also used as the L2-weight function when
290  projecting the DeltaCoefficient onto a GridFunction, so that the weighted
291  integral of the projection is exactly equal to the Scale(). */
292  void SetWeight(Coefficient *w) { weight = w; }
293 
294  /// Return a pointer to a c-array representing the center of the delta
295  /// function.
296  const double *Center() { return center; }
297 
298  /** @brief Return the scale factor times the optional time dependent
299  function. Returns \f$ s T(t) \f$ with \f$ T(t) = 1 \f$ when
300  not set by the user. */
301  double Scale() { return tdf ? (*tdf)(GetTime())*scale : scale; }
302 
303  /// Return the tolerance used to identify the mesh vertices
304  double Tol() { return tol; }
305 
306  /// See SetWeight() for description of the weight Coefficient.
307  Coefficient *Weight() { return weight; }
308 
309  /// Write the center of the delta function into @a center.
311 
312  /// The value of the function assuming we are evaluating at the delta center.
313  virtual double EvalDelta(ElementTransformation &T, const IntegrationPoint &ip);
314  /** @brief A DeltaFunction cannot be evaluated. Calling this method will
315  cause an MFEM error, terminating the application. */
316  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
317  { mfem_error("DeltaCoefficient::Eval"); return 0.; }
318  virtual ~DeltaCoefficient() { delete weight; }
319 };
320 
321 /** @brief Derived coefficient that takes the value of the parent coefficient
322  for the active attributes and is zero otherwise. */
324 {
325 private:
326  Coefficient *c;
327  Array<int> active_attr;
328 
329 public:
330  /** @brief Construct with a parent coefficient and an array with
331  ones marking the attributes on which this coefficient should be
332  active. */
334  { c = &c_; attr.Copy(active_attr); }
335 
336  /// Evaluate the coefficient at @a ip.
337  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
338  { return active_attr[T.Attribute-1] ? c->Eval(T, ip, GetTime()) : 0.0; }
339 };
340 
341 /// Base class for vector Coefficients that optionally depend on time and space.
343 {
344 protected:
345  int vdim;
346  double time;
347 
348 public:
349  /// Initialize the VectorCoefficient with vector dimension @a vd.
350  VectorCoefficient(int vd) { vdim = vd; time = 0.; }
351 
352  /// Set the time for time dependent coefficients
353  void SetTime(double t) { time = t; }
354 
355  /// Get the time for time dependent coefficients
356  double GetTime() { return time; }
357 
358  /// Returns dimension of the vector.
359  int GetVDim() { return vdim; }
360 
361  /** @brief Evaluate the vector coefficient in the element described by @a T
362  at the point @a ip, storing the result in @a V. */
363  /** @note When this method is called, the caller must make sure that the
364  IntegrationPoint associated with @a T is the same as @a ip. This can be
365  achieved by calling T.SetIntPoint(&ip). */
366  virtual void Eval(Vector &V, ElementTransformation &T,
367  const IntegrationPoint &ip) = 0;
368 
369  /** @brief Evaluate the vector coefficient in the element described by @a T
370  at all points of @a ir, storing the result in @a M. */
371  /** The dimensions of @a M are GetVDim() by ir.GetNPoints() and they must be
372  set by the implementation of this method.
373 
374  The general implementation provided by the base class (using the Eval
375  method for one IntegrationPoint at a time) can be overloaded for more
376  efficient implementation.
377 
378  @note The IntegrationPoint associated with @a T is not used, and this
379  method will generally modify this IntegrationPoint associated with @a T.
380  */
381  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
382  const IntegrationRule &ir);
383 
384  virtual ~VectorCoefficient() { }
385 };
386 
387 
388 /// Vector coefficient that is constant in space and time.
390 {
391 private:
392  Vector vec;
393 public:
394  /// Construct the coefficient with constant vector @a v.
396  : VectorCoefficient(v.Size()), vec(v) { }
398 
399  /// Evaluate the vector coefficient at @a ip.
400  virtual void Eval(Vector &V, ElementTransformation &T,
401  const IntegrationPoint &ip) { V = vec; }
402 
403  /// Return a reference to the constant vector in this class.
404  const Vector& GetVec() { return vec; }
405 };
406 
407 /// A general vector function coefficient
409 {
410 private:
411  std::function<void(const Vector &, Vector &)> Function;
412  std::function<void(const Vector &, double, Vector &)> TDFunction;
413  Coefficient *Q;
414 
415 public:
416  /// Define a time-independent vector coefficient from a std function
417  /** \param dim - the size of the vector
418  \param F - time-independent function
419  \param q - optional scalar Coefficient to scale the vector coefficient */
421  std::function<void(const Vector &, Vector &)> F,
422  Coefficient *q = nullptr)
423  : VectorCoefficient(dim), Function(std::move(F)), Q(q)
424  { }
425 
426  /// Define a time-dependent vector coefficient from a std function
427  /** \param dim - the size of the vector
428  \param TDF - time-dependent function
429  \param q - optional scalar Coefficient to scale the vector coefficient */
431  std::function<void(const Vector &, double, Vector &)> TDF,
432  Coefficient *q = nullptr)
433  : VectorCoefficient(dim), TDFunction(std::move(TDF)), Q(q)
434  { }
435 
437  /// Evaluate the vector coefficient at @a ip.
438  virtual void Eval(Vector &V, ElementTransformation &T,
439  const IntegrationPoint &ip);
440 
442 };
443 
444 /** @brief Vector coefficient defined by an array of scalar coefficients.
445  Coefficients that are not set will evaluate to zero in the vector. This
446  object takes ownership of the array of coefficients inside it and deletes
447  them at object destruction. */
449 {
450 private:
451  Array<Coefficient*> Coeff;
452  Array<bool> ownCoeff;
453 
454 public:
455  /** @brief Construct vector of dim coefficients. The actual coefficients
456  still need to be added with Set(). */
457  explicit VectorArrayCoefficient(int dim);
458 
459  /// Returns i'th coefficient.
460  Coefficient* GetCoeff(int i) { return Coeff[i]; }
461 
462  /// Returns the entire array of coefficients.
463  Coefficient **GetCoeffs() { return Coeff; }
464 
465  /// Sets coefficient in the vector.
466  void Set(int i, Coefficient *c, bool own=true);
467 
468  /// Evaluates i'th component of the vector of coefficients and returns the
469  /// value.
470  double Eval(int i, ElementTransformation &T, const IntegrationPoint &ip)
471  { return Coeff[i] ? Coeff[i]->Eval(T, ip, GetTime()) : 0.0; }
472 
474  /** @brief Evaluate the coefficient. Each element of vector V comes from the
475  associated array of scalar coefficients. */
476  virtual void Eval(Vector &V, ElementTransformation &T,
477  const IntegrationPoint &ip);
478 
479  /// Destroys vector coefficient.
480  virtual ~VectorArrayCoefficient();
481 };
482 
483 /// Vector coefficient defined by a vector GridFunction
485 {
486 protected:
488 
489 public:
490  /** @brief Construct an empty coefficient. Calling Eval() before the grid
491  function is set will cause a segfault. */
493 
494  /** @brief Construct the coefficient with grid function @a gf. The
495  grid function is not owned by the coefficient. */
497 
498  /** @brief Set the grid function for this coefficient. Also sets the Vector
499  dimension to match that of the @a gf. */
500  void SetGridFunction(const GridFunction *gf);
501 
502  /// Returns a pointer to the grid function in this Coefficient
503  const GridFunction * GetGridFunction() const { return GridFunc; }
504 
505  /// Evaluate the vector coefficient at @a ip.
506  virtual void Eval(Vector &V, ElementTransformation &T,
507  const IntegrationPoint &ip);
508 
509  /** @brief Evaluate the vector coefficients at all of the locations in the
510  integration rule and write the vectors into the columns of matrix @a
511  M. */
512  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
513  const IntegrationRule &ir);
514 
516 };
517 
518 /// Vector coefficient defined as the Gradient of a scalar GridFunction
520 {
521 protected:
523 
524 public:
525 
526  /** @brief Construct the coefficient with a scalar grid function @a gf. The
527  grid function is not owned by the coefficient. */
529 
530  ///Set the scalar grid function.
531  void SetGridFunction(const GridFunction *gf);
532 
533  ///Get the scalar grid function.
534  const GridFunction * GetGridFunction() const { return GridFunc; }
535 
536  /// Evaluate the gradient vector coefficient at @a ip.
537  virtual void Eval(Vector &V, ElementTransformation &T,
538  const IntegrationPoint &ip);
539 
540  /** @brief Evaluate the gradient vector coefficient at all of the locations
541  in the integration rule and write the vectors into columns of matrix @a
542  M. */
543  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
544  const IntegrationRule &ir);
545 
547 };
548 
549 /// Vector coefficient defined as the Curl of a vector GridFunction
551 {
552 protected:
554 
555 public:
556  /** @brief Construct the coefficient with a vector grid function @a gf. The
557  grid function is not owned by the coefficient. */
559 
560  /// Set the vector grid function.
561  void SetGridFunction(const GridFunction *gf);
562 
563  /// Get the vector grid function.
564  const GridFunction * GetGridFunction() const { return GridFunc; }
565 
567  /// Evaluate the vector curl coefficient at @a ip.
568  virtual void Eval(Vector &V, ElementTransformation &T,
569  const IntegrationPoint &ip);
570 
572 };
573 
574 /// Scalar coefficient defined as the Divergence of a vector GridFunction
576 {
577 protected:
579 
580 public:
581  /** @brief Construct the coefficient with a vector grid function @a gf. The
582  grid function is not owned by the coefficient. */
584 
585  /// Set the vector grid function.
586  void SetGridFunction(const GridFunction *gf) { GridFunc = gf; }
587 
588  /// Get the vector grid function.
589  const GridFunction * GetGridFunction() const { return GridFunc; }
590 
591  /// Evaluate the scalar divergence coefficient at @a ip.
592  virtual double Eval(ElementTransformation &T,
593  const IntegrationPoint &ip);
594 
596 };
597 
598 /** @brief Vector coefficient defined by a scalar DeltaCoefficient and a
599  constant vector direction.
600 
601  WARNING this cannot be used as a normal coefficient. The usual Eval method
602  is disabled. */
604 {
605 protected:
608 
609 public:
610  /// Construct with a vector of dimension @a vdim_.
612  : VectorCoefficient(vdim_), dir(vdim_), d() { }
613 
614  /** @brief Construct with a Vector object representing the direction and a
615  unit delta function centered at (0.0,0.0,0.0) */
617  : VectorCoefficient(dir_.Size()), dir(dir_), d() { }
618 
619  /** @brief Construct with a Vector object representing the direction and a
620  delta function scaled by @a s and centered at (x,0.0,0.0) */
621  VectorDeltaCoefficient(const Vector& dir_, double x, double s)
622  : VectorCoefficient(dir_.Size()), dir(dir_), d(x,s) { }
623 
624  /** @brief Construct with a Vector object representing the direction and a
625  delta function scaled by @a s and centered at (x,y,0.0) */
626  VectorDeltaCoefficient(const Vector& dir_, double x, double y, double s)
627  : VectorCoefficient(dir_.Size()), dir(dir_), d(x,y,s) { }
628 
629  /** @brief Construct with a Vector object representing the direction and a
630  delta function scaled by @a s and centered at (x,y,z) */
631  VectorDeltaCoefficient(const Vector& dir_, double x, double y, double z,
632  double s)
633  : VectorCoefficient(dir_.Size()), dir(dir_), d(x,y,z,s) { }
634 
635  /// Replace the associated DeltaCoefficient with a new DeltaCoefficient.
636  /** The new DeltaCoefficient cannot have a specified weight Coefficient, i.e.
637  DeltaCoefficient::Weight() should return NULL. */
638  void SetDeltaCoefficient(const DeltaCoefficient& d_) { d = d_; }
639 
640  /// Return the associated scalar DeltaCoefficient.
642 
643  void SetScale(double s) { d.SetScale(s); }
644  void SetDirection(const Vector& d_);
645 
646  void SetDeltaCenter(const Vector& center) { d.SetDeltaCenter(center); }
647  void GetDeltaCenter(Vector& center) { d.GetDeltaCenter(center); }
648 
649  /** @brief Return the specified direction vector multiplied by the value
650  returned by DeltaCoefficient::EvalDelta() of the associated scalar
651  DeltaCoefficient. */
652  virtual void EvalDelta(Vector &V, ElementTransformation &T,
653  const IntegrationPoint &ip);
654 
656  /** @brief A VectorDeltaFunction cannot be evaluated. Calling this method
657  will cause an MFEM error, terminating the application. */
658  virtual void Eval(Vector &V, ElementTransformation &T,
659  const IntegrationPoint &ip)
660  { mfem_error("VectorDeltaCoefficient::Eval"); }
662 };
663 
664 /** @brief Derived vector coefficient that has the value of the parent vector
665  where it is active and is zero otherwise. */
667 {
668 private:
670  Array<int> active_attr;
671 
672 public:
673  /** @brief Construct with a parent vector coefficient and an array of zeros
674  and ones representing the attributes for which this coefficient should be
675  active. */
677  : VectorCoefficient(vc.GetVDim())
678  { c = &vc; attr.Copy(active_attr); }
679 
680  /// Evaluate the vector coefficient at @a ip.
681  virtual void Eval(Vector &V, ElementTransformation &T,
682  const IntegrationPoint &ip);
683 
684  /** @brief Evaluate the vector coefficient at all of the locations in the
685  integration rule and write the vectors into the columns of matrix @a
686  M. */
687  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
688  const IntegrationRule &ir);
689 };
690 
692 
693 /// Base class for Matrix Coefficients that optionally depend on time and space.
695 {
696 protected:
697  int height, width;
698  double time;
699  bool symmetric; // deprecated
700 
701 public:
702  /// Construct a dim x dim matrix coefficient.
703  explicit MatrixCoefficient(int dim, bool symm=false)
704  { height = width = dim; time = 0.; symmetric = symm; }
705 
706  /// Construct a h x w matrix coefficient.
707  MatrixCoefficient(int h, int w, bool symm=false) :
708  height(h), width(w), time(0.), symmetric(symm) { }
709 
710  /// Set the time for time dependent coefficients
711  void SetTime(double t) { time = t; }
712 
713  /// Get the time for time dependent coefficients
714  double GetTime() { return time; }
715 
716  /// Get the height of the matrix.
717  int GetHeight() const { return height; }
718 
719  /// Get the width of the matrix.
720  int GetWidth() const { return width; }
721 
722  /// For backward compatibility get the width of the matrix.
723  int GetVDim() const { return width; }
724 
725  /** @deprecated Use SymmetricMatrixCoefficient instead */
726  bool IsSymmetric() const { return symmetric; }
727 
728  /** @brief Evaluate the matrix coefficient in the element described by @a T
729  at the point @a ip, storing the result in @a K. */
730  /** @note When this method is called, the caller must make sure that the
731  IntegrationPoint associated with @a T is the same as @a ip. This can be
732  achieved by calling T.SetIntPoint(&ip). */
733  virtual void Eval(DenseMatrix &K, ElementTransformation &T,
734  const IntegrationPoint &ip) = 0;
735 
736  /// (DEPRECATED) Evaluate a symmetric matrix coefficient.
737  /** @brief Evaluate the upper triangular entries of the matrix coefficient
738  in the symmetric case, similarly to Eval. Matrix entry (i,j) is stored
739  in K[j - i + os_i] for 0 <= i <= j < width, os_0 = 0,
740  os_{i+1} = os_i + width - i. That is, K = {M(0,0), ..., M(0,w-1),
741  M(1,1), ..., M(1,w-1), ..., M(w-1,w-1) with w = width.
742  @deprecated Use Eval() instead. */
744  const IntegrationPoint &ip)
745  { mfem_error("MatrixCoefficient::EvalSymmetric"); }
746 
747  virtual ~MatrixCoefficient() { }
748 };
749 
750 
751 /// A matrix coefficient that is constant in space and time.
753 {
754 private:
755  DenseMatrix mat;
756 public:
757  ///Construct using matrix @a m for the constant.
759  : MatrixCoefficient(m.Height(), m.Width()), mat(m) { }
761  /// Evaluate the matrix coefficient at @a ip.
763  const IntegrationPoint &ip) { M = mat; }
764 };
765 
766 
767 /** @brief A matrix coefficient with an optional scalar coefficient multiplier
768  \a q. The matrix function can either be represented by a std function or
769  a constant matrix provided when constructing this object. */
771 {
772 private:
773  std::function<void(const Vector &, DenseMatrix &)> Function;
774  std::function<void(const Vector &, Vector &)> SymmFunction; // deprecated
775  std::function<void(const Vector &, double, DenseMatrix &)> TDFunction;
776 
777  Coefficient *Q;
778  DenseMatrix mat;
779 
780 public:
781  /// Define a time-independent square matrix coefficient from a std function
782  /** \param dim - the size of the matrix
783  \param F - time-independent function
784  \param q - optional scalar Coefficient to scale the matrix coefficient */
786  std::function<void(const Vector &, DenseMatrix &)> F,
787  Coefficient *q = nullptr)
788  : MatrixCoefficient(dim), Function(std::move(F)), Q(q), mat(0)
789  { }
790 
791  /// Define a constant matrix coefficient times a scalar Coefficient
792  /** \param m - constant matrix
793  \param q - optional scalar Coefficient to scale the matrix coefficient */
795  : MatrixCoefficient(m.Height(), m.Width()), Q(&q), mat(m)
796  { }
797 
798  /** @brief Define a time-independent symmetric square matrix coefficient from
799  a std function */
800  /** \param dim - the size of the matrix
801  \param SymmF - function used in EvalSymmetric
802  \param q - optional scalar Coefficient to scale the matrix coefficient
803  @deprecated Use another constructor without setting SymmFunction. */
805  std::function<void(const Vector &, Vector &)> SymmF,
806  Coefficient *q = NULL)
807  : MatrixCoefficient(dim, true), SymmFunction(std::move(SymmF)), Q(q), mat(0)
808  { }
809 
810  /// Define a time-dependent square matrix coefficient from a std function
811  /** \param dim - the size of the matrix
812  \param TDF - time-dependent function
813  \param q - optional scalar Coefficient to scale the matrix coefficient */
815  std::function<void(const Vector &, double, DenseMatrix &)> TDF,
816  Coefficient *q = nullptr)
817  : MatrixCoefficient(dim), TDFunction(std::move(TDF)), Q(q)
818  { }
819 
820  /// Evaluate the matrix coefficient at @a ip.
821  virtual void Eval(DenseMatrix &K, ElementTransformation &T,
822  const IntegrationPoint &ip);
823 
824  /// (DEPRECATED) Evaluate the symmetric matrix coefficient at @a ip.
825  /** @deprecated Use Eval() instead. */
826  virtual void EvalSymmetric(Vector &K, ElementTransformation &T,
827  const IntegrationPoint &ip);
828 
830 };
831 
832 
833 /** @brief Matrix coefficient defined by a matrix of scalar coefficients.
834  Coefficients that are not set will evaluate to zero in the vector. The
835  coefficient is stored as a flat Array with indexing (i,j) -> i*width+j. */
837 {
838 private:
839  Array<Coefficient *> Coeff;
840  Array<bool> ownCoeff;
841 
842 public:
843  /** @brief Construct a coefficient matrix of dimensions @a dim * @a dim. The
844  actual coefficients still need to be added with Set(). */
845  explicit MatrixArrayCoefficient (int dim);
846 
847  /// Get the coefficient located at (i,j) in the matrix.
848  Coefficient* GetCoeff (int i, int j) { return Coeff[i*width+j]; }
849 
850  /** @brief Set the coefficient located at (i,j) in the matrix. By default by
851  default this will take ownership of the Coefficient passed in, but this
852  can be overridden with the @a own parameter. */
853  void Set(int i, int j, Coefficient * c, bool own=true);
854 
855  /// Evaluate coefficient located at (i,j) in the matrix using integration
856  /// point @a ip.
857  double Eval(int i, int j, ElementTransformation &T, const IntegrationPoint &ip)
858  { return Coeff[i*width+j] ? Coeff[i*width+j] -> Eval(T, ip, GetTime()) : 0.0; }
859 
860  /// Evaluate the matrix coefficient @a ip.
861  virtual void Eval(DenseMatrix &K, ElementTransformation &T,
862  const IntegrationPoint &ip);
863 
864  virtual ~MatrixArrayCoefficient();
865 };
866 
867 
868 /** @brief Derived matrix coefficient that has the value of the parent matrix
869  coefficient where it is active and is zero otherwise. */
871 {
872 private:
874  Array<int> active_attr;
875 
876 public:
877  /** @brief Construct with a parent matrix coefficient and an array of zeros
878  and ones representing the attributes for which this coefficient should be
879  active. */
881  : MatrixCoefficient(mc.GetHeight(), mc.GetWidth())
882  { c = &mc; attr.Copy(active_attr); }
883 
884  /// Evaluate the matrix coefficient at @a ip.
885  virtual void Eval(DenseMatrix &K, ElementTransformation &T,
886  const IntegrationPoint &ip);
887 };
888 
889 /// Coefficients based on sums, products, or other functions of coefficients.
890 ///@{
891 /** @brief Scalar coefficient defined as the linear combination of two scalar
892  coefficients or a scalar and a scalar coefficient */
894 {
895 private:
896  double aConst;
897  Coefficient * a;
898  Coefficient * b;
899 
900  double alpha;
901  double beta;
902 
903 public:
904  /// Constructor with one coefficient. Result is alpha_ * A + beta_ * B
906  double alpha_ = 1.0, double beta_ = 1.0)
907  : aConst(A), a(NULL), b(&B), alpha(alpha_), beta(beta_) { }
908 
909  /// Constructor with two coefficients. Result is alpha_ * A + beta_ * B.
911  double alpha_ = 1.0, double beta_ = 1.0)
912  : aConst(0.0), a(&A), b(&B), alpha(alpha_), beta(beta_) { }
913 
914  /// Reset the first term in the linear combination as a constant
915  void SetAConst(double A) { a = NULL; aConst = A; }
916  /// Return the first term in the linear combination
917  double GetAConst() const { return aConst; }
918 
919  /// Reset the first term in the linear combination
920  void SetACoef(Coefficient &A) { a = &A; }
921  /// Return the first term in the linear combination
922  Coefficient * GetACoef() const { return a; }
923 
924  /// Reset the second term in the linear combination
925  void SetBCoef(Coefficient &B) { b = &B; }
926  /// Return the second term in the linear combination
927  Coefficient * GetBCoef() const { return b; }
928 
929  /// Reset the factor in front of the first term in the linear combination
930  void SetAlpha(double alpha_) { alpha = alpha_; }
931  /// Return the factor in front of the first term in the linear combination
932  double GetAlpha() const { return alpha; }
933 
934  /// Reset the factor in front of the second term in the linear combination
935  void SetBeta(double beta_) { beta = beta_; }
936  /// Return the factor in front of the second term in the linear combination
937  double GetBeta() const { return beta; }
938 
939  /// Evaluate the coefficient at @a ip.
940  virtual double Eval(ElementTransformation &T,
941  const IntegrationPoint &ip)
942  {
943  return alpha * ((a == NULL ) ? aConst : a->Eval(T, ip) )
944  + beta * b->Eval(T, ip);
945  }
946 };
947 
948 
949 /// Base class for symmetric matrix coefficients that optionally depend on time and space.
951 {
952 protected:
953  int dim;
954  double time;
955 
956 public:
957  /// Construct a dim x dim matrix coefficient.
958  explicit SymmetricMatrixCoefficient(int dimension)
959  { dim = dimension; time = 0.; }
960 
961  /// Set the time for time dependent coefficients
962  void SetTime(double t) { time = t; }
963 
964  /// Get the time for time dependent coefficients
965  double GetTime() { return time; }
966 
967  /// Get the size of the matrix.
968  int GetSize() const { return dim; }
969 
970  /** @brief Evaluate the matrix coefficient in the element described by @a T
971  at the point @a ip, storing the result in @a K. */
972  /** @note When this method is called, the caller must make sure that the
973  IntegrationPoint associated with @a T is the same as @a ip. This can be
974  achieved by calling T.SetIntPoint(&ip). */
975  virtual void Eval(DenseSymmetricMatrix &K, ElementTransformation &T,
976  const IntegrationPoint &ip) = 0;
977 
979 };
980 
981 
982 /// A matrix coefficient that is constant in space and time.
984 {
985 private:
987 
988 public:
989  ///Construct using matrix @a m for the constant.
991  : SymmetricMatrixCoefficient(m.Height()), mat(m) { }
993  /// Evaluate the matrix coefficient at @a ip.
995  const IntegrationPoint &ip) { M = mat; }
996 };
997 
998 
999 /** @brief A matrix coefficient with an optional scalar coefficient multiplier
1000  \a q. The matrix function can either be represented by a std function or
1001  a constant matrix provided when constructing this object. */
1003 {
1004 private:
1005  std::function<void(const Vector &, DenseSymmetricMatrix &)> Function;
1006  std::function<void(const Vector &, double, DenseSymmetricMatrix &)> TDFunction;
1007 
1008  Coefficient *Q;
1010 
1011 public:
1012  /// Define a time-independent symmetric matrix coefficient from a std function
1013  /** \param dim - the size of the matrix
1014  \param F - time-independent function
1015  \param q - optional scalar Coefficient to scale the matrix coefficient */
1017  std::function<void(const Vector &, DenseSymmetricMatrix &)> F,
1018  Coefficient *q = nullptr)
1019  : SymmetricMatrixCoefficient(dim), Function(std::move(F)), Q(q), mat(0)
1020  { }
1021 
1022  /// Define a constant matrix coefficient times a scalar Coefficient
1023  /** \param m - constant matrix
1024  \param q - optional scalar Coefficient to scale the matrix coefficient */
1026  Coefficient &q)
1027  : SymmetricMatrixCoefficient(m.Height()), Q(&q), mat(m)
1028  { }
1029 
1030  /// Define a time-dependent square matrix coefficient from a std function
1031  /** \param dim - the size of the matrix
1032  \param TDF - time-dependent function
1033  \param q - optional scalar Coefficient to scale the matrix coefficient */
1035  std::function<void(const Vector &, double, DenseSymmetricMatrix &)> TDF,
1036  Coefficient *q = nullptr)
1037  : SymmetricMatrixCoefficient(dim), TDFunction(std::move(TDF)), Q(q)
1038  { }
1039 
1040  /// Evaluate the matrix coefficient at @a ip.
1041  virtual void Eval(DenseSymmetricMatrix &K, ElementTransformation &T,
1042  const IntegrationPoint &ip);
1043 
1045 };
1046 
1047 
1048 /** @brief Scalar coefficient defined as the product of two scalar coefficients
1049  or a scalar and a scalar coefficient. */
1051 {
1052 private:
1053  double aConst;
1054  Coefficient * a;
1055  Coefficient * b;
1056 
1057 public:
1058  /// Constructor with one coefficient. Result is A * B.
1060  : aConst(A), a(NULL), b(&B) { }
1061 
1062  /// Constructor with two coefficients. Result is A * B.
1064  : aConst(0.0), a(&A), b(&B) { }
1065 
1066  /// Reset the first term in the product as a constant
1067  void SetAConst(double A) { a = NULL; aConst = A; }
1068  /// Return the first term in the product
1069  double GetAConst() const { return aConst; }
1070 
1071  /// Reset the first term in the product
1072  void SetACoef(Coefficient &A) { a = &A; }
1073  /// Return the first term in the product
1074  Coefficient * GetACoef() const { return a; }
1075 
1076  /// Reset the second term in the product
1077  void SetBCoef(Coefficient &B) { b = &B; }
1078  /// Return the second term in the product
1079  Coefficient * GetBCoef() const { return b; }
1080 
1081  /// Evaluate the coefficient at @a ip.
1082  virtual double Eval(ElementTransformation &T,
1083  const IntegrationPoint &ip)
1084  { return ((a == NULL ) ? aConst : a->Eval(T, ip) ) * b->Eval(T, ip); }
1085 };
1086 
1087 /** @brief Scalar coefficient defined as the ratio of two scalars where one or
1088  both scalars are scalar coefficients. */
1090 {
1091 private:
1092  double aConst;
1093  double bConst;
1094  Coefficient * a;
1095  Coefficient * b;
1096 
1097 public:
1098  /** Initialize a coefficient which returns A / B where @a A is a
1099  constant and @a B is a scalar coefficient */
1101  : aConst(A), bConst(1.0), a(NULL), b(&B) { }
1102  /** Initialize a coefficient which returns A / B where @a A and @a B are both
1103  scalar coefficients */
1105  : aConst(0.0), bConst(1.0), a(&A), b(&B) { }
1106  /** Initialize a coefficient which returns A / B where @a A is a
1107  scalar coefficient and @a B is a constant */
1109  : aConst(0.0), bConst(B), a(&A), b(NULL) { }
1110 
1111  /// Reset the numerator in the ratio as a constant
1112  void SetAConst(double A) { a = NULL; aConst = A; }
1113  /// Return the numerator of the ratio
1114  double GetAConst() const { return aConst; }
1115 
1116  /// Reset the denominator in the ratio as a constant
1117  void SetBConst(double B) { b = NULL; bConst = B; }
1118  /// Return the denominator of the ratio
1119  double GetBConst() const { return bConst; }
1120 
1121  /// Reset the numerator in the ratio
1122  void SetACoef(Coefficient &A) { a = &A; }
1123  /// Return the numerator of the ratio
1124  Coefficient * GetACoef() const { return a; }
1125 
1126  /// Reset the denominator in the ratio
1127  void SetBCoef(Coefficient &B) { b = &B; }
1128  /// Return the denominator of the ratio
1129  Coefficient * GetBCoef() const { return b; }
1130 
1131  /// Evaluate the coefficient
1132  virtual double Eval(ElementTransformation &T,
1133  const IntegrationPoint &ip)
1134  {
1135  double den = (b == NULL ) ? bConst : b->Eval(T, ip);
1136  MFEM_ASSERT(den != 0.0, "Division by zero in RatioCoefficient");
1137  return ((a == NULL ) ? aConst : a->Eval(T, ip) ) / den;
1138  }
1139 };
1140 
1141 /// Scalar coefficient defined as a scalar raised to a power
1143 {
1144 private:
1145  Coefficient * a;
1146 
1147  double p;
1148 
1149 public:
1150  /// Construct with a coefficient and a constant power @a p_. Result is A^p.
1152  : a(&A), p(p_) { }
1153 
1154  /// Reset the base coefficient
1155  void SetACoef(Coefficient &A) { a = &A; }
1156  /// Return the base coefficient
1157  Coefficient * GetACoef() const { return a; }
1158 
1159  /// Reset the exponent
1160  void SetExponent(double p_) { p = p_; }
1161  /// Return the exponent
1162  double GetExponent() const { return p; }
1163 
1164  /// Evaluate the coefficient at @a ip.
1165  virtual double Eval(ElementTransformation &T,
1166  const IntegrationPoint &ip)
1167  { return pow(a->Eval(T, ip), p); }
1168 };
1169 
1170 
1171 /// Scalar coefficient defined as the inner product of two vector coefficients
1173 {
1174 private:
1175  VectorCoefficient * a;
1176  VectorCoefficient * b;
1177 
1178  mutable Vector va;
1179  mutable Vector vb;
1180 public:
1181  /// Construct with the two vector coefficients. Result is \f$ A \cdot B \f$.
1183 
1184  /// Reset the first vector in the inner product
1185  void SetACoef(VectorCoefficient &A) { a = &A; }
1186  /// Return the first vector coefficient in the inner product
1187  VectorCoefficient * GetACoef() const { return a; }
1188 
1189  /// Reset the second vector in the inner product
1190  void SetBCoef(VectorCoefficient &B) { b = &B; }
1191  /// Return the second vector coefficient in the inner product
1192  VectorCoefficient * GetBCoef() const { return b; }
1193 
1194  /// Evaluate the coefficient at @a ip.
1195  virtual double Eval(ElementTransformation &T,
1196  const IntegrationPoint &ip);
1197 };
1198 
1199 /// Scalar coefficient defined as a cross product of two vectors in the xy-plane.
1201 {
1202 private:
1203  VectorCoefficient * a;
1204  VectorCoefficient * b;
1205 
1206  mutable Vector va;
1207  mutable Vector vb;
1208 
1209 public:
1210  /// Constructor with two vector coefficients. Result is \f$ A_x B_y - A_y * B_x; \f$.
1212 
1213  /// Reset the first vector in the product
1214  void SetACoef(VectorCoefficient &A) { a = &A; }
1215  /// Return the first vector of the product
1216  VectorCoefficient * GetACoef() const { return a; }
1217 
1218  /// Reset the second vector in the product
1219  void SetBCoef(VectorCoefficient &B) { b = &B; }
1220  /// Return the second vector of the product
1221  VectorCoefficient * GetBCoef() const { return b; }
1222 
1223  /// Evaluate the coefficient at @a ip.
1224  virtual double Eval(ElementTransformation &T,
1225  const IntegrationPoint &ip);
1226 };
1227 
1228 /// Scalar coefficient defined as the determinant of a matrix coefficient
1230 {
1231 private:
1232  MatrixCoefficient * a;
1233 
1234  mutable DenseMatrix ma;
1235 
1236 public:
1237  /// Construct with the matrix.
1239 
1240  /// Reset the matrix coefficient
1241  void SetACoef(MatrixCoefficient &A) { a = &A; }
1242  /// Return the matrix coefficient
1243  MatrixCoefficient * GetACoef() const { return a; }
1244 
1245  /// Evaluate the determinant coefficient at @a ip.
1246  virtual double Eval(ElementTransformation &T,
1247  const IntegrationPoint &ip);
1248 };
1249 
1250 /// Vector coefficient defined as the linear combination of two vectors
1252 {
1253 private:
1254  VectorCoefficient * ACoef;
1255  VectorCoefficient * BCoef;
1256 
1257  Vector A;
1258  Vector B;
1259 
1260  Coefficient * alphaCoef;
1261  Coefficient * betaCoef;
1262 
1263  double alpha;
1264  double beta;
1265 
1266  mutable Vector va;
1267 
1268 public:
1269  /** Constructor with no coefficients.
1270  To be used with the various "Set" methods */
1271  VectorSumCoefficient(int dim);
1272 
1273  /** Constructor with two vector coefficients.
1274  Result is alpha_ * A + beta_ * B */
1276  double alpha_ = 1.0, double beta_ = 1.0);
1277 
1278  /** Constructor with scalar coefficients.
1279  Result is alpha_ * A_ + beta_ * B_ */
1281  Coefficient &alpha_, Coefficient &beta_);
1282 
1283  /// Reset the first vector coefficient
1284  void SetACoef(VectorCoefficient &A) { ACoef = &A; }
1285  /// Return the first vector coefficient
1286  VectorCoefficient * GetACoef() const { return ACoef; }
1287 
1288  /// Reset the second vector coefficient
1289  void SetBCoef(VectorCoefficient &B) { BCoef = &B; }
1290  /// Return the second vector coefficient
1291  VectorCoefficient * GetBCoef() const { return BCoef; }
1292 
1293  /// Reset the factor in front of the first vector coefficient
1294  void SetAlphaCoef(Coefficient &A) { alphaCoef = &A; }
1295  /// Return the factor in front of the first vector coefficient
1296  Coefficient * GetAlphaCoef() const { return alphaCoef; }
1297 
1298  /// Reset the factor in front of the second vector coefficient
1299  void SetBetaCoef(Coefficient &B) { betaCoef = &B; }
1300  /// Return the factor in front of the second vector coefficient
1301  Coefficient * GetBetaCoef() const { return betaCoef; }
1302 
1303  /// Reset the first vector as a constant
1304  void SetA(const Vector &A_) { A = A_; ACoef = NULL; }
1305  /// Return the first vector constant
1306  const Vector & GetA() const { return A; }
1307 
1308  /// Reset the second vector as a constant
1309  void SetB(const Vector &B_) { B = B_; BCoef = NULL; }
1310  /// Return the second vector constant
1311  const Vector & GetB() const { return B; }
1312 
1313  /// Reset the factor in front of the first vector coefficient as a constant
1314  void SetAlpha(double alpha_) { alpha = alpha_; alphaCoef = NULL; }
1315  /// Return the factor in front of the first vector coefficient
1316  double GetAlpha() const { return alpha; }
1317 
1318  /// Reset the factor in front of the second vector coefficient as a constant
1319  void SetBeta(double beta_) { beta = beta_; betaCoef = NULL; }
1320  /// Return the factor in front of the second vector coefficient
1321  double GetBeta() const { return beta; }
1322 
1323  /// Evaluate the coefficient at @a ip.
1324  virtual void Eval(Vector &V, ElementTransformation &T,
1325  const IntegrationPoint &ip);
1327 };
1328 
1329 /// Vector coefficient defined as a product of scalar and vector coefficients.
1331 {
1332 private:
1333  double aConst;
1334  Coefficient * a;
1335  VectorCoefficient * b;
1336 
1337 public:
1338  /// Constructor with constant and vector coefficient. Result is A * B.
1340 
1341  /// Constructor with two coefficients. Result is A * B.
1343 
1344  /// Reset the scalar factor as a constant
1345  void SetAConst(double A) { a = NULL; aConst = A; }
1346  /// Return the scalar factor
1347  double GetAConst() const { return aConst; }
1348 
1349  /// Reset the scalar factor
1350  void SetACoef(Coefficient &A) { a = &A; }
1351  /// Return the scalar factor
1352  Coefficient * GetACoef() const { return a; }
1353 
1354  /// Reset the vector factor
1355  void SetBCoef(VectorCoefficient &B) { b = &B; }
1356  /// Return the vector factor
1357  VectorCoefficient * GetBCoef() const { return b; }
1358 
1359  /// Evaluate the coefficient at @a ip.
1360  virtual void Eval(Vector &V, ElementTransformation &T,
1361  const IntegrationPoint &ip);
1363 };
1364 
1365 /// Vector coefficient defined as a normalized vector field (returns v/|v|)
1367 {
1368 private:
1369  VectorCoefficient * a;
1370 
1371  double tol;
1372 
1373 public:
1374  /** @brief Return a vector normalized to a length of one
1375 
1376  This class evaluates the vector coefficient @a A and, if |A| > @a tol,
1377  returns the normalized vector A / |A|. If |A| <= @a tol, the zero
1378  vector is returned.
1379  */
1380  NormalizedVectorCoefficient(VectorCoefficient &A, double tol = 1e-6);
1381 
1382  /// Reset the vector coefficient
1383  void SetACoef(VectorCoefficient &A) { a = &A; }
1384  /// Return the vector coefficient
1385  VectorCoefficient * GetACoef() const { return a; }
1386 
1387  /// Evaluate the coefficient at @a ip.
1388  virtual void Eval(Vector &V, ElementTransformation &T,
1389  const IntegrationPoint &ip);
1391 };
1392 
1393 /// Vector coefficient defined as a cross product of two vectors
1395 {
1396 private:
1397  VectorCoefficient * a;
1398  VectorCoefficient * b;
1399 
1400  mutable Vector va;
1401  mutable Vector vb;
1402 
1403 public:
1404  /// Construct with the two coefficients. Result is A x B.
1406 
1407  /// Reset the first term in the product
1408  void SetACoef(VectorCoefficient &A) { a = &A; }
1409  /// Return the first term in the product
1410  VectorCoefficient * GetACoef() const { return a; }
1411 
1412  /// Reset the second term in the product
1413  void SetBCoef(VectorCoefficient &B) { b = &B; }
1414  /// Return the second term in the product
1415  VectorCoefficient * GetBCoef() const { return b; }
1416 
1417  /// Evaluate the coefficient at @a ip.
1418  virtual void Eval(Vector &V, ElementTransformation &T,
1419  const IntegrationPoint &ip);
1421 };
1422 
1423 /** @brief Vector coefficient defined as a product of a matrix coefficient and
1424  a vector coefficient. */
1426 {
1427 private:
1428  MatrixCoefficient * a;
1429  VectorCoefficient * b;
1430 
1431  mutable DenseMatrix ma;
1432  mutable Vector vb;
1433 
1434 public:
1435  /// Constructor with two coefficients. Result is A*B.
1437 
1438  /// Reset the matrix coefficient
1439  void SetACoef(MatrixCoefficient &A) { a = &A; }
1440  /// Return the matrix coefficient
1441  MatrixCoefficient * GetACoef() const { return a; }
1442 
1443  /// Reset the vector coefficient
1444  void SetBCoef(VectorCoefficient &B) { b = &B; }
1445  /// Return the vector coefficient
1446  VectorCoefficient * GetBCoef() const { return b; }
1447 
1448  /// Evaluate the vector coefficient at @a ip.
1449  virtual void Eval(Vector &V, ElementTransformation &T,
1450  const IntegrationPoint &ip);
1452 };
1453 
1454 /// Convenient alias for the MatrixVectorProductCoefficient
1456 
1457 /// Constant matrix coefficient defined as the identity of dimension d
1459 {
1460 private:
1461  int dim;
1462 
1463 public:
1464  /// Construct with the dimension of the square identity matrix.
1466  : MatrixCoefficient(d, d), dim(d) { }
1467 
1468  /// Evaluate the matrix coefficient at @a ip.
1469  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1470  const IntegrationPoint &ip);
1471 };
1472 
1473 /// Matrix coefficient defined as the linear combination of two matrices
1475 {
1476 private:
1477  MatrixCoefficient * a;
1478  MatrixCoefficient * b;
1479 
1480  double alpha;
1481  double beta;
1482 
1483  mutable DenseMatrix ma;
1484 
1485 public:
1486  /// Construct with the two coefficients. Result is alpha_ * A + beta_ * B.
1488  double alpha_ = 1.0, double beta_ = 1.0);
1489 
1490  /// Reset the first matrix coefficient
1491  void SetACoef(MatrixCoefficient &A) { a = &A; }
1492  /// Return the first matrix coefficient
1493  MatrixCoefficient * GetACoef() const { return a; }
1494 
1495  /// Reset the second matrix coefficient
1496  void SetBCoef(MatrixCoefficient &B) { b = &B; }
1497  /// Return the second matrix coefficient
1498  MatrixCoefficient * GetBCoef() const { return b; }
1499 
1500  /// Reset the factor in front of the first matrix coefficient
1501  void SetAlpha(double alpha_) { alpha = alpha_; }
1502  /// Return the factor in front of the first matrix coefficient
1503  double GetAlpha() const { return alpha; }
1504 
1505  /// Reset the factor in front of the second matrix coefficient
1506  void SetBeta(double beta_) { beta = beta_; }
1507  /// Return the factor in front of the second matrix coefficient
1508  double GetBeta() const { return beta; }
1509 
1510  /// Evaluate the matrix coefficient at @a ip.
1511  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1512  const IntegrationPoint &ip);
1513 };
1514 
1515 /** @brief Matrix coefficient defined as a product of a scalar coefficient and a
1516  matrix coefficient.*/
1518 {
1519 private:
1520  double aConst;
1521  Coefficient * a;
1522  MatrixCoefficient * b;
1523 
1524 public:
1525  /// Constructor with one coefficient. Result is A*B.
1527 
1528  /// Constructor with two coefficients. Result is A*B.
1530 
1531  /// Reset the scalar factor as a constant
1532  void SetAConst(double A) { a = NULL; aConst = A; }
1533  /// Return the scalar factor
1534  double GetAConst() const { return aConst; }
1535 
1536  /// Reset the scalar factor
1537  void SetACoef(Coefficient &A) { a = &A; }
1538  /// Return the scalar factor
1539  Coefficient * GetACoef() const { return a; }
1540 
1541  /// Reset the matrix factor
1542  void SetBCoef(MatrixCoefficient &B) { b = &B; }
1543  /// Return the matrix factor
1544  MatrixCoefficient * GetBCoef() const { return b; }
1545 
1546  /// Evaluate the matrix coefficient at @a ip.
1547  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1548  const IntegrationPoint &ip);
1549 };
1550 
1551 /// Matrix coefficient defined as the transpose a matrix coefficient
1553 {
1554 private:
1555  MatrixCoefficient * a;
1556 
1557 public:
1558  /// Construct with the matrix coefficient. Result is \f$ A^T \f$.
1560 
1561  /// Reset the matrix coefficient
1562  void SetACoef(MatrixCoefficient &A) { a = &A; }
1563  /// Return the matrix coefficient
1564  MatrixCoefficient * GetACoef() const { return a; }
1565 
1566  /// Evaluate the matrix coefficient at @a ip.
1567  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1568  const IntegrationPoint &ip);
1569 };
1570 
1571 /// Matrix coefficient defined as the inverse a matrix coefficient.
1573 {
1574 private:
1575  MatrixCoefficient * a;
1576 
1577 public:
1578  /// Construct with the matrix coefficient. Result is \f$ A^{-1} \f$.
1580 
1581  /// Reset the matrix coefficient
1582  void SetACoef(MatrixCoefficient &A) { a = &A; }
1583  /// Return the matrix coefficient
1584  MatrixCoefficient * GetACoef() const { return a; }
1585 
1586  /// Evaluate the matrix coefficient at @a ip.
1587  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1588  const IntegrationPoint &ip);
1589 };
1590 
1591 /// Matrix coefficient defined as the outer product of two vector coefficients.
1593 {
1594 private:
1595  VectorCoefficient * a;
1596  VectorCoefficient * b;
1597 
1598  mutable Vector va;
1599  mutable Vector vb;
1600 
1601 public:
1602  /// Construct with two vector coefficients. Result is \f$ A B^T \f$.
1604 
1605  /// Reset the first vector in the outer product
1606  void SetACoef(VectorCoefficient &A) { a = &A; }
1607  /// Return the first vector coefficient in the outer product
1608  VectorCoefficient * GetACoef() const { return a; }
1609 
1610  /// Reset the second vector in the outer product
1611  void SetBCoef(VectorCoefficient &B) { b = &B; }
1612  /// Return the second vector coefficient in the outer product
1613  VectorCoefficient * GetBCoef() const { return b; }
1614 
1615  /// Evaluate the matrix coefficient at @a ip.
1616  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1617  const IntegrationPoint &ip);
1618 };
1619 
1620 /** @brief Matrix coefficient defined as -a k x k x, for a vector k and scalar a
1621 
1622  This coefficient returns \f$a * (|k|^2 I - k \otimes k)\f$, where I is
1623  the identity matrix and \f$\otimes\f$ indicates the outer product. This
1624  can be evaluated for vectors of any dimension but in three
1625  dimensions it corresponds to computing the cross product with k twice.
1626 */
1628 {
1629 private:
1630  double aConst;
1631  Coefficient * a;
1632  VectorCoefficient * k;
1633 
1634  mutable Vector vk;
1635 
1636 public:
1639 
1640  /// Reset the scalar factor as a constant
1641  void SetAConst(double A) { a = NULL; aConst = A; }
1642  /// Return the scalar factor
1643  double GetAConst() const { return aConst; }
1644 
1645  /// Reset the scalar factor
1646  void SetACoef(Coefficient &A) { a = &A; }
1647  /// Return the scalar factor
1648  Coefficient * GetACoef() const { return a; }
1649 
1650  /// Reset the vector factor
1651  void SetKCoef(VectorCoefficient &K) { k = &K; }
1652  /// Return the vector factor
1653  VectorCoefficient * GetKCoef() const { return k; }
1654 
1655  /// Evaluate the matrix coefficient at @a ip.
1656  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1657  const IntegrationPoint &ip);
1658 };
1659 ///@}
1660 
1661 class QuadratureFunction;
1662 
1663 /** @brief Vector quadrature function coefficient which requires that the
1664  quadrature rules used for this vector coefficient be the same as those that
1665  live within the supplied QuadratureFunction. */
1667 {
1668 private:
1669  const QuadratureFunction &QuadF; //do not own
1670  int index;
1671 
1672 public:
1673  /// Constructor with a quadrature function as input
1675 
1676  /** Set the starting index within the QuadFunc that'll be used to project
1677  outwards as well as the corresponding length. The projected length should
1678  have the bounds of 1 <= length <= (length QuadFunc - index). */
1679  void SetComponent(int index_, int length_);
1680 
1681  const QuadratureFunction& GetQuadFunction() const { return QuadF; }
1682 
1684  virtual void Eval(Vector &V, ElementTransformation &T,
1685  const IntegrationPoint &ip);
1686 
1688 };
1689 
1690 /** @brief Quadrature function coefficient which requires that the quadrature
1691  rules used for this coefficient be the same as those that live within the
1692  supplied QuadratureFunction. */
1694 {
1695 private:
1696  const QuadratureFunction &QuadF;
1697 
1698 public:
1699  /// Constructor with a quadrature function as input
1701 
1702  const QuadratureFunction& GetQuadFunction() const { return QuadF; }
1703 
1704  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
1705 
1707 };
1708 
1709 /** @brief Compute the Lp norm of a function f.
1710  \f$ \| f \|_{Lp} = ( \int_\Omega | f |^p d\Omega)^{1/p} \f$ */
1711 double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh,
1712  const IntegrationRule *irs[]);
1713 
1714 /** @brief Compute the Lp norm of a vector function f = {f_i}_i=1...N.
1715  \f$ \| f \|_{Lp} = ( \sum_i \| f_i \|_{Lp}^p )^{1/p} \f$ */
1716 double ComputeLpNorm(double p, VectorCoefficient &coeff, Mesh &mesh,
1717  const IntegrationRule *irs[]);
1718 
1719 #ifdef MFEM_USE_MPI
1720 /** @brief Compute the global Lp norm of a function f.
1721  \f$ \| f \|_{Lp} = ( \int_\Omega | f |^p d\Omega)^{1/p} \f$ */
1722 double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh,
1723  const IntegrationRule *irs[]);
1724 
1725 /** @brief Compute the global Lp norm of a vector function f = {f_i}_i=1...N.
1726  \f$ \| f \|_{Lp} = ( \sum_i \| f_i \|_{Lp}^p )^{1/p} \f$ */
1727 double ComputeGlobalLpNorm(double p, VectorCoefficient &coeff, ParMesh &pmesh,
1728  const IntegrationRule *irs[]);
1729 #endif
1730 
1731 }
1732 
1733 #endif
FunctionCoefficient(std::function< double(const Vector &)> F)
Define a time-independent coefficient from a std function.
Vector coefficient defined as a cross product of two vectors.
double GetAConst() const
Return the first term in the product.
A coefficient that depends on 1 or 2 parent coefficients and a transformation rule represented by a C...
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
RatioCoefficient(Coefficient &A, Coefficient &B)
A matrix coefficient that is constant in space and time.
void SetACoef(Coefficient &A)
Reset the scalar factor.
double Eval(ElementTransformation &T, const IntegrationPoint &ip, double t)
Evaluate the coefficient in the element described by T at the point ip at time t. ...
Definition: coefficient.hpp:66
void GetDeltaCenter(Vector &center)
VectorCoefficient * GetACoef() const
Return the first vector coefficient in the inner product.
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
Vector coefficient defined by an array of scalar coefficients. Coefficients that are not set will eva...
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
Coefficient * GetBCoef() const
Return the second term in the product.
MatrixFunctionCoefficient(int dim, std::function< void(const Vector &, double, DenseMatrix &)> TDF, Coefficient *q=nullptr)
Define a time-dependent square matrix coefficient from a std function.
SymmetricMatrixFunctionCoefficient(int dim, std::function< void(const Vector &, DenseSymmetricMatrix &)> F, Coefficient *q=nullptr)
Define a time-independent symmetric matrix coefficient from a std function.
Coefficient * GetACoef() const
Return the scalar factor.
void Set(int i, int j, Coefficient *c, bool own=true)
Set the coefficient located at (i,j) in the matrix. By default by default this will take ownership of...
IdentityMatrixCoefficient(int d)
Construct with the dimension of the square identity matrix.
VectorCoefficient * GetACoef() const
Return the vector coefficient.
PWConstCoefficient(Vector &c)
Construct the constant coefficient using a vector of constants.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetWeight(Coefficient *w)
Set a weight Coefficient that multiplies the DeltaCoefficient.
Base class for vector Coefficients that optionally depend on time and space.
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
Coefficient ** GetCoeffs()
Returns the entire array of coefficients.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector curl coefficient at ip.
void SetAConst(double A)
Reset the first term in the product as a constant.
TransposeMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
double GetBConst() const
Return the denominator of the ratio.
void SetBCoef(VectorCoefficient &B)
Reset the second vector coefficient.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
double(* tdf)(double)
VectorCrossProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with the two coefficients. Result is A x B.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
Scalar coefficient defined as the inner product of two vector coefficients.
void SetTime(double t)
Set the time for time dependent coefficients.
NormalizedVectorCoefficient(VectorCoefficient &A, double tol=1e-6)
Return a vector normalized to a length of one.
void SetBCoef(MatrixCoefficient &B)
Reset the matrix factor.
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
void SetKCoef(VectorCoefficient &K)
Reset the vector factor.
MatrixVectorProductCoefficient(MatrixCoefficient &A, VectorCoefficient &B)
Constructor with two coefficients. Result is A*B.
double Eval(int i, ElementTransformation &T, const IntegrationPoint &ip)
VectorDeltaCoefficient(int vdim_)
Construct with a vector of dimension vdim_.
void SetDeltaCenter(const Vector &center)
Set the center location of the delta function.
Definition: coefficient.cpp:69
MatrixCoefficient(int h, int w, bool symm=false)
Construct a h x w matrix coefficient.
double Scale()
Return the scale factor times the optional time dependent function. Returns with when not set by th...
MatrixCoefficient * GetBCoef() const
Return the matrix factor.
Vector quadrature function coefficient which requires that the quadrature rules used for this vector ...
void SetACoef(Coefficient &A)
Reset the numerator in the ratio.
Vector coefficient that is constant in space and time.
Scalar coefficient defined as a cross product of two vectors in the xy-plane.
double GetTime()
Get the time for time dependent coefficients.
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
Scalar coefficient defined as the product of two scalar coefficients or a scalar and a scalar coeffic...
void SetBCoef(VectorCoefficient &B)
Reset the second term in the product.
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition: array.hpp:851
void SetBeta(double beta_)
Reset the factor in front of the second matrix coefficient.
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
Coefficient * GetBCoef() const
Return the denominator of the ratio.
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
double & operator()(int i)
Return a reference to the i-th constant.
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Delta function coefficient optionally multiplied by a weight coefficient and a scaled time dependent ...
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
void SetExponent(double p_)
Reset the exponent.
A matrix coefficient with an optional scalar coefficient multiplier q. The matrix function can either...
void SetFunction(double(*f)(double))
Set a time-dependent function that multiplies the Scale().
Scalar coefficient defined as the ratio of two scalars where one or both scalars are scalar coefficie...
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
CurlGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a vector grid function gf. The grid function is not owned by the coeff...
Matrix coefficient defined by a matrix of scalar coefficients. Coefficients that are not set will eva...
GradientGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a scalar grid function gf. The grid function is not owned by the coeff...
CrossCrossCoefficient(double A, VectorCoefficient &K)
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetB(const Vector &B_)
Reset the second vector as a constant.
void SetDeltaCoefficient(const DeltaCoefficient &d_)
Replace the associated DeltaCoefficient with a new DeltaCoefficient.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
double Tol()
Return the tolerance used to identify the mesh vertices.
SymmetricMatrixFunctionCoefficient(const DenseSymmetricMatrix &m, Coefficient &q)
Define a constant matrix coefficient times a scalar Coefficient.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetGridFunction(const GridFunction *gf)
Set the vector grid function.
Vector coefficient defined as the linear combination of two vectors.
void SetACoef(VectorCoefficient &A)
Reset the first vector coefficient.
MFEM_DEPRECATED FunctionCoefficient(double(*f)(Vector &))
(DEPRECATED) Define a time-independent coefficient from a C-function
Matrix coefficient defined as a product of a scalar coefficient and a matrix coefficient.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
Definition: coefficient.cpp:55
virtual void Eval(DenseSymmetricMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
A DeltaFunction cannot be evaluated. Calling this method will cause an MFEM error, terminating the application.
Coefficient * GetACoef() const
Return the scalar factor.
VectorCoefficient * GetKCoef() const
Return the vector factor.
void SetTime(double t)
Set the time for time dependent coefficients.
void SetBCoef(VectorCoefficient &B)
Reset the second vector in the inner product.
virtual void EvalSymmetric(Vector &K, ElementTransformation &T, const IntegrationPoint &ip)
(DEPRECATED) Evaluate a symmetric matrix coefficient.
Matrix coefficient defined as -a k x k x, for a vector k and scalar a.
DeltaCoefficient(double x, double y, double z, double s)
Construct a delta function scaled by s and centered at (x,y,z)
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
Definition: coefficient.cpp:31
double GetBeta() const
Return the factor in front of the second matrix coefficient.
MatrixFunctionCoefficient(int dim, std::function< void(const Vector &, DenseMatrix &)> F, Coefficient *q=nullptr)
Define a time-independent square matrix coefficient from a std function.
RestrictedCoefficient(Coefficient &c_, Array< int > &attr)
Construct with a parent coefficient and an array with ones marking the attributes on which this coeff...
Matrix coefficient defined as the linear combination of two matrices.
DeltaCoefficient & GetDeltaCoefficient()
Return the associated scalar DeltaCoefficient.
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
DeltaCoefficient(double x, double s)
Construct a delta function scaled by s and centered at (x,0.0,0.0)
DeterminantCoefficient(MatrixCoefficient &A)
Construct with the matrix.
VectorCoefficient * GetBCoef() const
Return the vector factor.
SumCoefficient(double A, Coefficient &B, double alpha_=1.0, double beta_=1.0)
Constructor with one coefficient. Result is alpha_ * A + beta_ * B.
Matrix coefficient defined as the inverse a matrix coefficient.
void SetACoef(Coefficient &A)
Reset the first term in the product.
ScalarMatrixProductCoefficient(double A, MatrixCoefficient &B)
Constructor with one coefficient. Result is A*B.
FunctionCoefficient(std::function< double(const Vector &, double)> TDF)
Define a time-dependent coefficient from a std function.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
A VectorDeltaFunction cannot be evaluated. Calling this method will cause an MFEM error...
Coefficient * GetACoef() const
Return the first term in the product.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
int GetHeight() const
Get the height of the matrix.
RatioCoefficient(double A, Coefficient &B)
Coefficient * GetCoeff(int i, int j)
Get the coefficient located at (i,j) in the matrix.
Matrix coefficient defined as the outer product of two vector coefficients.
void SetBCoef(Coefficient &B)
Reset the second term in the product.
MFEM_DEPRECATED FunctionCoefficient(double(*tdf)(Vector &, double))
(DEPRECATED) Define a time-dependent coefficient from a C-function
void SetACoef(MatrixCoefficient &A)
Reset the first matrix coefficient.
int GetSize() const
Get the size of the matrix.
double Eval(int i, int j, ElementTransformation &T, const IntegrationPoint &ip)
void SetAConst(double A)
Reset the scalar factor as a constant.
Derived vector coefficient that has the value of the parent vector where it is active and is zero oth...
Coefficient * GetACoef() const
Return the scalar factor.
double GetAConst() const
Return the scalar factor.
InnerProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with the two vector coefficients. Result is .
Coefficient * GetACoef() const
Return the first term in the linear combination.
MatrixCoefficient(int dim, bool symm=false)
Construct a dim x dim matrix coefficient.
Quadrature function coefficient which requires that the quadrature rules used for this coefficient be...
void SetAConst(double A)
Reset the scalar factor as a constant.
VectorCoefficient * GetBCoef() const
Return the second vector coefficient in the outer product.
double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
Compute the Lp norm of a function f. .
std::function< double(const Vector &, double)> TDFunction
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
VectorConstantCoefficient(const Vector &v)
Construct the coefficient with constant vector v.
ConstantCoefficient(double c=1.0)
c is value of constant function
Definition: coefficient.hpp:84
VectorCoefficient * GetBCoef() const
Return the second vector coefficient.
void SetBCoef(VectorCoefficient &B)
Reset the vector coefficient.
void SetBCoef(VectorCoefficient &B)
Reset the second vector in the outer product.
const QuadratureFunction & GetQuadFunction() const
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
A matrix coefficient with an optional scalar coefficient multiplier q. The matrix function can either...
const double * Center()
double GetTime()
Get the time for time dependent coefficients.
VectorCoefficient * GetACoef() const
Return the first vector of the product.
void SetA(const Vector &A_)
Reset the first vector as a constant.
Derived coefficient that takes the value of the parent coefficient for the active attributes and is z...
QuadratureFunctionCoefficient(QuadratureFunction &qf)
Constructor with a quadrature function as input.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:153
void SetDeltaCenter(const Vector &center)
Coefficient * Weight()
See SetWeight() for description of the weight Coefficient.
RatioCoefficient(Coefficient &A, double B)
MatrixFunctionCoefficient(int dim, std::function< void(const Vector &, Vector &)> SymmF, Coefficient *q=NULL)
Define a time-independent symmetric square matrix coefficient from a std function.
const Vector & GetVec()
Return a reference to the constant vector in this class.
MatrixRestrictedCoefficient(MatrixCoefficient &mc, Array< int > &attr)
Construct with a parent matrix coefficient and an array of zeros and ones representing the attributes...
std::function< double(const Vector &)> Function
Coefficient * GetBetaCoef() const
Return the factor in front of the second vector coefficient.
Vector coefficient defined as the Gradient of a scalar GridFunction.
Vector coefficient defined as the Curl of a vector GridFunction.
OuterProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with two vector coefficients. Result is .
const GridFunction * GetGridFunction() const
Get the internal GridFunction.
void SetGridFunction(const GridFunction *gf)
Set the vector grid function.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void Set(int i, Coefficient *c, bool own=true)
Sets coefficient in the vector.
void SetBeta(double beta_)
Reset the factor in front of the second term in the linear combination.
DeltaCoefficient()
Construct a unit delta function centered at (0.0,0.0,0.0)
void SetTol(double tol_)
Set the tolerance used during projection onto GridFunction to identify the Mesh vertex where the Cent...
void SetACoef(VectorCoefficient &A)
Reset the first vector in the inner product.
void SetAlpha(double alpha_)
Reset the factor in front of the first matrix coefficient.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
virtual ~VectorArrayCoefficient()
Destroys vector coefficient.
void SetBeta(double beta_)
Reset the factor in front of the second vector coefficient as a constant.
VectorFunctionCoefficient(int dim, std::function< void(const Vector &, Vector &)> F, Coefficient *q=nullptr)
Define a time-independent vector coefficient from a std function.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
int GetVDim()
Returns dimension of the vector.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetACoef(Coefficient &A)
Reset the scalar factor.
VectorCoefficient * GetACoef() const
Return the first vector coefficient in the outer product.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetACoef(VectorCoefficient &A)
Reset the first term in the product.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the determinant coefficient at ip.
VectorDeltaCoefficient(const Vector &dir_)
Construct with a Vector object representing the direction and a unit delta function centered at (0...
const GridFunction * GetGridFunction() const
Get the vector grid function.
void SetComponent(int index_, int length_)
void SetBCoef(VectorCoefficient &B)
Reset the vector factor.
void SetTime(double t)
Set the time for time dependent coefficients.
Definition: coefficient.hpp:48
void operator=(double c)
Set the constants for all attributes to constant c.
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
VectorCoefficient DiagonalMatrixCoefficient
A general vector function coefficient.
const GridFunction * GetGridFunction() const
Returns a pointer to the grid function in this Coefficient.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
Definition: coefficient.hpp:87
Coefficient * GetAlphaCoef() const
Return the factor in front of the first vector coefficient.
VectorCoefficient * GetBCoef() const
Return the second term in the product.
void GetDeltaCenter(Vector &center)
Write the center of the delta function into center.
Definition: coefficient.cpp:77
MatrixFunctionCoefficient(const DenseMatrix &m, Coefficient &q)
Define a constant matrix coefficient times a scalar Coefficient.
Scalar coefficient defined as the Divergence of a vector GridFunction.
MatrixVectorProductCoefficient MatVecCoefficient
Convenient alias for the MatrixVectorProductCoefficient.
void SetACoef(Coefficient &A)
Reset the base coefficient.
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
bool IsSymmetric() const
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
VectorDeltaCoefficient(const Vector &dir_, double x, double y, double z, double s)
Construct with a Vector object representing the direction and a delta function scaled by s and center...
void SetGridFunction(const GridFunction *gf)
Set the scalar grid function.
VectorDeltaCoefficient(const Vector &dir_, double x, double s)
Construct with a Vector object representing the direction and a delta function scaled by s and center...
const GridFunction * GetGridFunction() const
Get the scalar grid function.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
void SetBCoef(VectorCoefficient &B)
Reset the second vector in the product.
virtual ~Coefficient()
Definition: coefficient.hpp:73
VectorQuadratureFunctionCoefficient(QuadratureFunction &qf)
Constructor with a quadrature function as input.
Vector coefficient defined as a product of scalar and vector coefficients.
MatrixCoefficient * GetACoef() const
Return the first matrix coefficient.
double GetAConst() const
Return the first term in the linear combination.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void UpdateConstants(Vector &c)
Update the constants with vector c.
Base class for Matrix Coefficients that optionally depend on time and space.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
Definition: coefficient.cpp:24
Coefficient * GetCoeff(int i)
Returns i&#39;th coefficient.
void SetTime(double t)
Set the time for time dependent coefficients.
virtual void Eval(DenseSymmetricMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
virtual void EvalDelta(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Return the specified direction vector multiplied by the value returned by DeltaCoefficient::EvalDelta...
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
void SetAConst(double A)
Reset the first term in the linear combination as a constant.
PWConstCoefficient(int NumOfSubD=0)
Constructs a piecewise constant coefficient in NumOfSubD subdomains.
void SetAlphaCoef(Coefficient &A)
Reset the factor in front of the first vector coefficient.
PowerCoefficient(Coefficient &A, double p_)
Construct with a coefficient and a constant power p_. Result is A^p.
Coefficient * GetACoef() const
Return the numerator of the ratio.
virtual double EvalDelta(ElementTransformation &T, const IntegrationPoint &ip)
The value of the function assuming we are evaluating at the delta center.
Definition: coefficient.cpp:83
SymmetricMatrixFunctionCoefficient(int dim, std::function< void(const Vector &, double, DenseSymmetricMatrix &)> TDF, Coefficient *q=nullptr)
Define a time-dependent square matrix coefficient from a std function.
SymmetricMatrixConstantCoefficient(const DenseSymmetricMatrix &m)
Construct using matrix m for the constant.
Vector coefficient defined as a product of a matrix coefficient and a vector coefficient.
VectorCoefficient * GetACoef() const
Return the first vector coefficient.
const QuadratureFunction & GetQuadFunction() const
Derived matrix coefficient that has the value of the parent matrix coefficient where it is active and...
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
Definition: coefficient.cpp:49
Vector coefficient defined by a scalar DeltaCoefficient and a constant vector direction.
int GetVDim() const
For backward compatibility get the width of the matrix.
MatrixCoefficient * GetBCoef() const
Return the second matrix coefficient.
MatrixArrayCoefficient(int dim)
Construct a coefficient matrix of dimensions dim * dim. The actual coefficients still need to be adde...
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the gradient vector coefficient at ip.
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
Definition: coefficient.hpp:94
void SetACoef(VectorCoefficient &A)
Reset the vector coefficient.
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
Class for integration point with weight.
Definition: intrules.hpp:25
VectorCoefficient * GetBCoef() const
Return the vector coefficient.
double GetAConst() const
Return the scalar factor.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the scalar divergence coefficient at ip.
double GetTime()
Get the time for time dependent coefficients.
void SetAlpha(double alpha_)
Reset the factor in front of the first vector coefficient as a constant.
double GetBeta() const
Return the factor in front of the second term in the linear combination.
double GetAlpha() const
Return the factor in front of the first vector coefficient.
void SetAConst(double A)
Reset the numerator in the ratio as a constant.
void SetScale(double s_)
Set the scale value multiplying the delta function.
const Vector & GetB() const
Return the second vector constant.
Base class for symmetric matrix coefficients that optionally depend on time and space.
void SetBetaCoef(Coefficient &B)
Reset the factor in front of the second vector coefficient.
DivergenceGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a vector grid function gf. The grid function is not owned by the coeff...
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
SymmetricMatrixCoefficient(int dimension)
Construct a dim x dim matrix coefficient.
Coefficient * GetBCoef() const
Return the second term in the linear combination.
TransformedCoefficient(Coefficient *q, double(*F)(double))
const GridFunction * GetGridFunction() const
Get the vector grid function.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
Vector coefficient defined as a normalized vector field (returns v/|v|)
VectorCoefficient(int vd)
Initialize the VectorCoefficient with vector dimension vd.
int GetNConst()
Returns the number of constants representing different attributes.
int dim
Definition: ex24.cpp:53
void SetDirection(const Vector &d_)
void SetAConst(double A)
Reset the scalar factor as a constant.
VectorFunctionCoefficient(int dim, std::function< void(const Vector &, double, Vector &)> TDF, Coefficient *q=nullptr)
Define a time-dependent vector coefficient from a std function.
void SetACoef(VectorCoefficient &A)
Reset the first vector in the product.
void SetGridFunction(const GridFunction *gf)
Set the grid function for this coefficient. Also sets the Vector dimension to match that of the gf...
Matrix coefficient defined as the transpose a matrix coefficient.
void SetBConst(double B)
Reset the denominator in the ratio as a constant.
void SetBCoef(Coefficient &B)
Reset the second term in the linear combination.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
A matrix coefficient that is constant in space and time.
VectorCoefficient * GetBCoef() const
Return the second vector coefficient in the inner product.
double GetAConst() const
Return the scalar factor.
DeltaCoefficient(double x, double y, double s)
Construct a delta function scaled by s and centered at (x,y,0.0)
MatrixSumCoefficient(MatrixCoefficient &A, MatrixCoefficient &B, double alpha_=1.0, double beta_=1.0)
Construct with the two coefficients. Result is alpha_ * A + beta_ * B.
Scalar coefficient defined as a scalar raised to a power.
RefCoord t[3]
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
VectorArrayCoefficient(int dim)
Construct vector of dim coefficients. The actual coefficients still need to be added with Set()...
virtual void Eval(DenseSymmetricMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
double GetTime()
Get the time for time dependent coefficients.
Definition: coefficient.hpp:51
Scalar coefficient defined as the determinant of a matrix coefficient.
A general function coefficient.
void SetACoef(VectorCoefficient &A)
Reset the first vector in the outer product.
ProductCoefficient(Coefficient &A, Coefficient &B)
Constructor with two coefficients. Result is A * B.
Vector data type.
Definition: vector.hpp:60
TransformedCoefficient(Coefficient *q1, Coefficient *q2, double(*F)(double, double))
GridFunctionCoefficient(const GridFunction *gf, int comp=1)
Vector coefficient defined by a vector GridFunction.
int GetWidth() const
Get the width of the matrix.
VectorRestrictedCoefficient(VectorCoefficient &vc, Array< int > &attr)
Construct with a parent vector coefficient and an array of zeros and ones representing the attributes...
VectorRotProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Constructor with two vector coefficients. Result is .
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetACoef(Coefficient &A)
Reset the scalar factor.
RefCoord s[3]
MatrixConstantCoefficient(const DenseMatrix &m)
Construct using matrix m for the constant.
double GetAlpha() const
Return the factor in front of the first term in the linear combination.
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
double GetAlpha() const
Return the factor in front of the first matrix coefficient.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
VectorCoefficient * GetACoef() const
Return the first term in the product.
void SetAlpha(double alpha_)
Reset the factor in front of the first term in the linear combination.
double GetAConst() const
Return the numerator of the ratio.
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
Constant matrix coefficient defined as the identity of dimension d.
double GetExponent() const
Return the exponent.
virtual void EvalSymmetric(Vector &K, ElementTransformation &T, const IntegrationPoint &ip)
(DEPRECATED) Evaluate the symmetric matrix coefficient at ip.
Coefficient * GetACoef() const
Return the base coefficient.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
void SetGridFunction(const GridFunction *gf)
Set the internal GridFunction.
void SetBCoef(Coefficient &B)
Reset the denominator in the ratio.
InverseMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
Scalar coefficient defined as the linear combination of two scalar coefficients or a scalar and a sca...
ProductCoefficient(double A, Coefficient &B)
Constructor with one coefficient. Result is A * B.
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:734
VectorDeltaCoefficient(const Vector &dir_, double x, double y, double s)
Construct with a Vector object representing the direction and a delta function scaled by s and center...
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
const Vector & GetA() const
Return the first vector constant.
VectorCoefficient * GetBCoef() const
Return the second vector of the product.
void SetACoef(Coefficient &A)
Reset the first term in the linear combination.
double GetBeta() const
Return the factor in front of the second vector coefficient.
void SetBCoef(MatrixCoefficient &B)
Reset the second matrix coefficient.
SumCoefficient(Coefficient &A, Coefficient &B, double alpha_=1.0, double beta_=1.0)
Constructor with two coefficients. Result is alpha_ * A + beta_ * B.
ScalarVectorProductCoefficient(double A, VectorCoefficient &B)
Constructor with constant and vector coefficient. Result is A * B.
VectorGridFunctionCoefficient()
Construct an empty coefficient. Calling Eval() before the grid function is set will cause a segfault...