MFEM  v4.4.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-2022, 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  virtual 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 /** @brief A piecewise coefficient with the pieces keyed off the element
129  attribute numbers.
130 
131  A value of zero will be returned for any missing attribute numbers.
132 
133  This object will not assume ownership of any Coefficient objects
134  passed to it. Consequently, the caller must ensure that the
135  individual Coefficient objects are not deleted while this
136  PWCoefficient is still in use.
137 
138  \note The keys may either be domain attribute numbers or boundary
139  attribute numbers. If the PWCoefficient is used with a domain
140  integrator the keys are assumed to be domain attribute
141  numbers. Similarly, if the PWCoefficient is used with a boundary
142  integrator the keys are assumed to be boundary attribute numbers.
143 */
145 {
146 private:
147  /** Internal data structure to store pointers to the appropriate
148  coefficients for different regions of the mesh. The keys used
149  in the map are the mesh attribute numbers (either element
150  attribute or boundary element attribute depending upon
151  context). The values returned for any missing attributes will
152  be zero. The coefficient pointers may be NULL in which case a
153  value of zero is returned.
154 
155  The Coefficient objects contained in this map are NOT owned by
156  this PWCoefficient object. This means that they will not be
157  deleted when this object is deleted also the caller must ensure
158  that the various Coefficient objects are not deleted while this
159  PWCoefficient is still needed.
160  */
161  std::map<int, Coefficient*> pieces;
162 
163  /** Convenience function to check for compatible array lengths,
164  loop over the arrays, and add their attribute/Coefficient pairs
165  to the internal data structure.
166  */
167  void InitMap(const Array<int> & attr,
168  const Array<Coefficient*> & coefs);
169 
170 public:
171 
172  /// Constructs a piecewise coefficient
173  explicit PWCoefficient() {}
174 
175  /// Construct the coefficient using arrays describing the pieces
176  /** \param attr - an array of attribute numbers for each piece
177  \param coefs - the corresponding array of Coefficient pointers
178  Any missing attributes or NULL coefficient pointers will result in a
179  value of zero being returned for that attribute.
180 
181  \note Ownership of the Coefficient objects will NOT be
182  transferred to this object.
183  */
184  PWCoefficient(const Array<int> & attr,
185  const Array<Coefficient*> & coefs)
186  { InitMap(attr, coefs); }
187 
188  /// Set the time for time dependent coefficients
189  virtual void SetTime(double t);
190 
191  /// Replace a set of coefficients
192  void UpdateCoefficients(const Array<int> & attr,
193  const Array<Coefficient*> & coefs)
194  { InitMap(attr, coefs); }
195 
196  /// Replace a single Coefficient for a particular attribute
197  void UpdateCoefficient(int attr, Coefficient & coef)
198  { pieces[attr] = &coef; }
199 
200  /// Remove a single Coefficient for a particular attribute
201  void ZeroCoefficient(int attr)
202  { pieces.erase(attr); }
203 
204  /// Evaluate the coefficient.
205  virtual double Eval(ElementTransformation &T,
206  const IntegrationPoint &ip);
207 };
208 
209 /// A general function coefficient
211 {
212 protected:
213  std::function<double(const Vector &)> Function;
214  std::function<double(const Vector &, double)> TDFunction;
215 
216 public:
217  /// Define a time-independent coefficient from a std function
218  /** \param F time-independent std::function */
219  FunctionCoefficient(std::function<double(const Vector &)> F)
220  : Function(std::move(F))
221  { }
222 
223  /// Define a time-dependent coefficient from a std function
224  /** \param TDF time-dependent function */
225  FunctionCoefficient(std::function<double(const Vector &, double)> TDF)
226  : TDFunction(std::move(TDF))
227  { }
228 
229  /// (DEPRECATED) Define a time-independent coefficient from a C-function
230  /** @deprecated Use the method where the C-function, @a f, uses a const
231  Vector argument instead of Vector. */
232  MFEM_DEPRECATED FunctionCoefficient(double (*f)(Vector &))
233  {
234  Function = reinterpret_cast<double(*)(const Vector&)>(f);
235  TDFunction = NULL;
236  }
237 
238  /// (DEPRECATED) Define a time-dependent coefficient from a C-function
239  /** @deprecated Use the method where the C-function, @a tdf, uses a const
240  Vector argument instead of Vector. */
241  MFEM_DEPRECATED FunctionCoefficient(double (*tdf)(Vector &, double))
242  {
243  Function = NULL;
244  TDFunction = reinterpret_cast<double(*)(const Vector&,double)>(tdf);
245  }
246 
247  /// Evaluate the coefficient at @a ip.
248  virtual double Eval(ElementTransformation &T,
249  const IntegrationPoint &ip);
250 };
251 
252 class GridFunction;
253 
254 /// Coefficient defined by a GridFunction. This coefficient is mesh dependent.
256 {
257 private:
258  const GridFunction *GridF;
259  int Component;
260 
261 public:
262  GridFunctionCoefficient() : GridF(NULL), Component(1) { }
263  /** Construct GridFunctionCoefficient from a given GridFunction, and
264  optionally specify a component to use if it is a vector GridFunction. */
265  GridFunctionCoefficient (const GridFunction *gf, int comp = 1)
266  { GridF = gf; Component = comp; }
267 
268  /// Set the internal GridFunction
269  void SetGridFunction(const GridFunction *gf) { GridF = gf; }
270 
271  /// Get the internal GridFunction
272  const GridFunction * GetGridFunction() const { return GridF; }
273 
274  /// Evaluate the coefficient at @a ip.
275  virtual double Eval(ElementTransformation &T,
276  const IntegrationPoint &ip);
277 };
278 
279 
280 /** @brief A coefficient that depends on 1 or 2 parent coefficients and a
281  transformation rule represented by a C-function.
282 
283  \f$ C(x,t) = T(Q1(x,t)) \f$ or \f$ C(x,t) = T(Q1(x,t), Q2(x,t)) \f$
284 
285  where T is the transformation rule, and Q1/Q2 are the parent coefficients.*/
287 {
288 private:
289  Coefficient * Q1;
290  Coefficient * Q2;
291  double (*Transform1)(double);
292  double (*Transform2)(double,double);
293 
294 public:
295  TransformedCoefficient (Coefficient * q,double (*F)(double))
296  : Q1(q), Transform1(F) { Q2 = 0; Transform2 = 0; }
298  double (*F)(double,double))
299  : Q1(q1), Q2(q2), Transform2(F) { Transform1 = 0; }
300 
301  /// Set the time for internally stored coefficients
302  void SetTime(double t);
303 
304  /// Evaluate the coefficient at @a ip.
305  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
306 };
307 
308 /** @brief Delta function coefficient optionally multiplied by a weight
309  coefficient and a scaled time dependent C-function.
310 
311  \f$ F(x,t) = w(x,t) s T(t) d(x - xc) \f$
312 
313  where w is the optional weight coefficient, @a s is a scale factor
314  T is an optional time-dependent function and d is a delta function.
315 
316  WARNING this cannot be used as a normal coefficient. The usual Eval
317  method is disabled. */
319 {
320 protected:
321  double center[3], scale, tol;
323  int sdim;
324  double (*tdf)(double);
325 
326 public:
327 
328  /// Construct a unit delta function centered at (0.0,0.0,0.0)
330  {
331  center[0] = center[1] = center[2] = 0.; scale = 1.; tol = 1e-12;
332  weight = NULL; sdim = 0; tdf = NULL;
333  }
334 
335  /// Construct a delta function scaled by @a s and centered at (x,0.0,0.0)
336  DeltaCoefficient(double x, double s)
337  {
338  center[0] = x; center[1] = 0.; center[2] = 0.; scale = s; tol = 1e-12;
339  weight = NULL; sdim = 1; tdf = NULL;
340  }
341 
342  /// Construct a delta function scaled by @a s and centered at (x,y,0.0)
343  DeltaCoefficient(double x, double y, double s)
344  {
345  center[0] = x; center[1] = y; center[2] = 0.; scale = s; tol = 1e-12;
346  weight = NULL; sdim = 2; tdf = NULL;
347  }
348 
349  /// Construct a delta function scaled by @a s and centered at (x,y,z)
350  DeltaCoefficient(double x, double y, double z, double s)
351  {
352  center[0] = x; center[1] = y; center[2] = z; scale = s; tol = 1e-12;
353  weight = NULL; sdim = 3; tdf = NULL;
354  }
355 
356  /// Set the time for internally stored coefficients
357  void SetTime(double t);
358 
359  /// Set the center location of the delta function.
360  void SetDeltaCenter(const Vector& center);
361 
362  /// Set the scale value multiplying the delta function.
363  void SetScale(double s_) { scale = s_; }
364 
365  /// Set a time-dependent function that multiplies the Scale().
366  void SetFunction(double (*f)(double)) { tdf = f; }
367 
368  /** @brief Set the tolerance used during projection onto GridFunction to
369  identify the Mesh vertex where the Center() of the delta function
370  lies. (default 1e-12)*/
371  void SetTol(double tol_) { tol = tol_; }
372 
373  /// Set a weight Coefficient that multiplies the DeltaCoefficient.
374  /** The weight Coefficient multiplies the value returned by EvalDelta() but
375  not the value returned by Scale().
376  The weight Coefficient is also used as the L2-weight function when
377  projecting the DeltaCoefficient onto a GridFunction, so that the weighted
378  integral of the projection is exactly equal to the Scale(). */
379  void SetWeight(Coefficient *w) { weight = w; }
380 
381  /// Return a pointer to a c-array representing the center of the delta
382  /// function.
383  const double *Center() { return center; }
384 
385  /** @brief Return the scale factor times the optional time dependent
386  function. Returns \f$ s T(t) \f$ with \f$ T(t) = 1 \f$ when
387  not set by the user. */
388  double Scale() { return tdf ? (*tdf)(GetTime())*scale : scale; }
389 
390  /// Return the tolerance used to identify the mesh vertices
391  double Tol() { return tol; }
392 
393  /// See SetWeight() for description of the weight Coefficient.
394  Coefficient *Weight() { return weight; }
395 
396  /// Write the center of the delta function into @a center.
398 
399  /// The value of the function assuming we are evaluating at the delta center.
400  virtual double EvalDelta(ElementTransformation &T, const IntegrationPoint &ip);
401  /** @brief A DeltaFunction cannot be evaluated. Calling this method will
402  cause an MFEM error, terminating the application. */
403  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
404  { mfem_error("DeltaCoefficient::Eval"); return 0.; }
405  virtual ~DeltaCoefficient() { delete weight; }
406 };
407 
408 /** @brief Derived coefficient that takes the value of the parent coefficient
409  for the active attributes and is zero otherwise. */
411 {
412 private:
413  Coefficient *c;
414  Array<int> active_attr;
415 
416 public:
417  /** @brief Construct with a parent coefficient and an array with
418  ones marking the attributes on which this coefficient should be
419  active. */
421  { c = &c_; attr.Copy(active_attr); }
422 
423  /// Set the time for internally stored coefficients
424  void SetTime(double t);
425 
426  /// Evaluate the coefficient at @a ip.
427  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
428  { return active_attr[T.Attribute-1] ? c->Eval(T, ip, GetTime()) : 0.0; }
429 };
430 
431 /// Base class for vector Coefficients that optionally depend on time and space.
433 {
434 protected:
435  int vdim;
436  double time;
437 
438 public:
439  /// Initialize the VectorCoefficient with vector dimension @a vd.
440  VectorCoefficient(int vd) { vdim = vd; time = 0.; }
441 
442  /// Set the time for time dependent coefficients
443  virtual void SetTime(double t) { time = t; }
444 
445  /// Get the time for time dependent coefficients
446  double GetTime() { return time; }
447 
448  /// Returns dimension of the vector.
449  int GetVDim() { return vdim; }
450 
451  /** @brief Evaluate the vector coefficient in the element described by @a T
452  at the point @a ip, storing the result in @a V. */
453  /** @note When this method is called, the caller must make sure that the
454  IntegrationPoint associated with @a T is the same as @a ip. This can be
455  achieved by calling T.SetIntPoint(&ip). */
456  virtual void Eval(Vector &V, ElementTransformation &T,
457  const IntegrationPoint &ip) = 0;
458 
459  /** @brief Evaluate the vector coefficient in the element described by @a T
460  at all points of @a ir, storing the result in @a M. */
461  /** The dimensions of @a M are GetVDim() by ir.GetNPoints() and they must be
462  set by the implementation of this method.
463 
464  The general implementation provided by the base class (using the Eval
465  method for one IntegrationPoint at a time) can be overloaded for more
466  efficient implementation.
467 
468  @note The IntegrationPoint associated with @a T is not used, and this
469  method will generally modify this IntegrationPoint associated with @a T.
470  */
471  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
472  const IntegrationRule &ir);
473 
474  virtual ~VectorCoefficient() { }
475 };
476 
477 
478 /// Vector coefficient that is constant in space and time.
480 {
481 private:
482  Vector vec;
483 public:
484  /// Construct the coefficient with constant vector @a v.
486  : VectorCoefficient(v.Size()), vec(v) { }
488 
489  /// Evaluate the vector coefficient at @a ip.
490  virtual void Eval(Vector &V, ElementTransformation &T,
491  const IntegrationPoint &ip) { V = vec; }
492 
493  /// Return a reference to the constant vector in this class.
494  const Vector& GetVec() { return vec; }
495 };
496 
497 /** @brief A piecewise vector-valued coefficient with the pieces keyed off the
498  element attribute numbers.
499 
500  A value of zero will be returned for any missing attribute numbers.
501 
502  This object will not assume ownership of any VectorCoefficient
503  objects passed to it. Consequently, the caller must ensure that
504  the individual VectorCoefficient objects are not deleted while
505  this PWVectorCoefficient is still in use.
506 
507  \note The keys may either be domain attribute numbers or boundary
508  attribute numbers. If the PWVectorCoefficient is used with a
509  domain integrator the keys are assumed to be domain attribute
510  numbers. Similarly, if the PWVectorCoefficient is used with a
511  boundary integrator the keys are assumed to be boundary attribute
512  numbers.
513 */
515 {
516 private:
517  /** Internal data structure to store pointers to the appropriate
518  coefficients for different regions of the mesh. The keys used
519  in the map are the mesh attribute numbers (either element
520  attribute or boundary element attribute depending upon
521  context). The values returned for any missing attributes will
522  be zero. The coefficient pointers may be NULL in which case a
523  value of zero is returned.
524 
525  The VectorCoefficient objects contained in this map are NOT
526  owned by this PWVectorCoefficient object. This means that they
527  will not be deleted when this object is deleted also the caller
528  must ensure that the various VectorCoefficient objects are not
529  deleted while this PWVectorCoefficient is still needed.
530  */
531  std::map<int, VectorCoefficient*> pieces;
532 
533  /** Convenience function to check for compatible array lengths,
534  loop over the arrays, and add their attribute/VectorCoefficient
535  pairs to the internal data structure.
536  */
537  void InitMap(const Array<int> & attr,
538  const Array<VectorCoefficient*> & coefs);
539 
540 public:
541 
542  /// Constructs a piecewise vector coefficient of dimension vd
543  explicit PWVectorCoefficient(int vd): VectorCoefficient(vd) {}
544 
545  /// Construct the coefficient using arrays describing the pieces
546  /** \param vd - dimension of the vector-valued result
547  \param attr - an array of attribute numbers for each piece
548  \param coefs - the corresponding array of VectorCoefficient pointers
549  Any missing attributes or NULL coefficient pointers will result in a
550  zero vector being returned for that attribute.
551 
552  \note Ownership of the VectorCoefficient objects will NOT be
553  transferred to this object.
554  */
555  PWVectorCoefficient(int vd, const Array<int> & attr,
556  const Array<VectorCoefficient*> & coefs)
557  : VectorCoefficient(vd) { InitMap(attr, coefs); }
558 
559  /// Set the time for time dependent coefficients
560  virtual void SetTime(double t);
561 
562  /// Replace a set of coefficients
563  void UpdateCoefficients(const Array<int> & attr,
564  const Array<VectorCoefficient*> & coefs)
565  { InitMap(attr, coefs); }
566 
567  /// Replace a single Coefficient for a particular attribute
568  void UpdateCoefficient(int attr, VectorCoefficient & coef);
569 
570  /// Remove a single VectorCoefficient for a particular attribute
571  void ZeroCoefficient(int attr)
572  { pieces.erase(attr); }
573 
574  /// Evaluate the coefficient.
575  virtual void Eval(Vector &V, ElementTransformation &T,
576  const IntegrationPoint &ip);
578 };
579 
580 /// A general vector function coefficient
582 {
583 private:
584  std::function<void(const Vector &, Vector &)> Function;
585  std::function<void(const Vector &, double, Vector &)> TDFunction;
586  Coefficient *Q;
587 
588 public:
589  /// Define a time-independent vector coefficient from a std function
590  /** \param dim - the size of the vector
591  \param F - time-independent function
592  \param q - optional scalar Coefficient to scale the vector coefficient */
594  std::function<void(const Vector &, Vector &)> F,
595  Coefficient *q = nullptr)
596  : VectorCoefficient(dim), Function(std::move(F)), Q(q)
597  { }
598 
599  /// Define a time-dependent vector coefficient from a std function
600  /** \param dim - the size of the vector
601  \param TDF - time-dependent function
602  \param q - optional scalar Coefficient to scale the vector coefficient */
604  std::function<void(const Vector &, double, Vector &)> TDF,
605  Coefficient *q = nullptr)
606  : VectorCoefficient(dim), TDFunction(std::move(TDF)), Q(q)
607  { }
608 
610  /// Evaluate the vector coefficient at @a ip.
611  virtual void Eval(Vector &V, ElementTransformation &T,
612  const IntegrationPoint &ip);
613 
615 };
616 
617 /** @brief Vector coefficient defined by an array of scalar coefficients.
618  Coefficients that are not set will evaluate to zero in the vector. This
619  object takes ownership of the array of coefficients inside it and deletes
620  them at object destruction. */
622 {
623 private:
624  Array<Coefficient*> Coeff;
625  Array<bool> ownCoeff;
626 
627 public:
628  /** @brief Construct vector of dim coefficients. The actual coefficients
629  still need to be added with Set(). */
630  explicit VectorArrayCoefficient(int dim);
631 
632  /// Set the time for internally stored coefficients
633  void SetTime(double t);
634 
635  /// Returns i'th coefficient.
636  Coefficient* GetCoeff(int i) { return Coeff[i]; }
637 
638  /// Returns the entire array of coefficients.
639  Coefficient **GetCoeffs() { return Coeff; }
640 
641  /// Sets coefficient in the vector.
642  void Set(int i, Coefficient *c, bool own=true);
643 
644  /// Evaluates i'th component of the vector of coefficients and returns the
645  /// value.
646  double Eval(int i, ElementTransformation &T, const IntegrationPoint &ip)
647  { return Coeff[i] ? Coeff[i]->Eval(T, ip, GetTime()) : 0.0; }
648 
650  /** @brief Evaluate the coefficient. Each element of vector V comes from the
651  associated array of scalar coefficients. */
652  virtual void Eval(Vector &V, ElementTransformation &T,
653  const IntegrationPoint &ip);
654 
655  /// Destroys vector coefficient.
656  virtual ~VectorArrayCoefficient();
657 };
658 
659 /// Vector coefficient defined by a vector GridFunction
661 {
662 protected:
664 
665 public:
666  /** @brief Construct an empty coefficient. Calling Eval() before the grid
667  function is set will cause a segfault. */
669 
670  /** @brief Construct the coefficient with grid function @a gf. The
671  grid function is not owned by the coefficient. */
673 
674  /** @brief Set the grid function for this coefficient. Also sets the Vector
675  dimension to match that of the @a gf. */
676  void SetGridFunction(const GridFunction *gf);
677 
678  /// Returns a pointer to the grid function in this Coefficient
679  const GridFunction * GetGridFunction() const { return GridFunc; }
680 
681  /// Evaluate the vector coefficient at @a ip.
682  virtual void Eval(Vector &V, ElementTransformation &T,
683  const IntegrationPoint &ip);
684 
685  /** @brief Evaluate the vector coefficients at all of the locations in the
686  integration rule and write the vectors into the columns of matrix @a
687  M. */
688  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
689  const IntegrationRule &ir);
690 
692 };
693 
694 /// Vector coefficient defined as the Gradient of a scalar GridFunction
696 {
697 protected:
699 
700 public:
701 
702  /** @brief Construct the coefficient with a scalar grid function @a gf. The
703  grid function is not owned by the coefficient. */
705 
706  ///Set the scalar grid function.
707  void SetGridFunction(const GridFunction *gf);
708 
709  ///Get the scalar grid function.
710  const GridFunction * GetGridFunction() const { return GridFunc; }
711 
712  /// Evaluate the gradient vector coefficient at @a ip.
713  virtual void Eval(Vector &V, ElementTransformation &T,
714  const IntegrationPoint &ip);
715 
716  /** @brief Evaluate the gradient vector coefficient at all of the locations
717  in the integration rule and write the vectors into columns of matrix @a
718  M. */
719  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
720  const IntegrationRule &ir);
721 
723 };
724 
725 /// Vector coefficient defined as the Curl of a vector GridFunction
727 {
728 protected:
730 
731 public:
732  /** @brief Construct the coefficient with a vector grid function @a gf. The
733  grid function is not owned by the coefficient. */
735 
736  /// Set the vector grid function.
737  void SetGridFunction(const GridFunction *gf);
738 
739  /// Get the vector grid function.
740  const GridFunction * GetGridFunction() const { return GridFunc; }
741 
743  /// Evaluate the vector curl coefficient at @a ip.
744  virtual void Eval(Vector &V, ElementTransformation &T,
745  const IntegrationPoint &ip);
746 
748 };
749 
750 /// Scalar coefficient defined as the Divergence of a vector GridFunction
752 {
753 protected:
755 
756 public:
757  /** @brief Construct the coefficient with a vector grid function @a gf. The
758  grid function is not owned by the coefficient. */
760 
761  /// Set the vector grid function.
762  void SetGridFunction(const GridFunction *gf) { GridFunc = gf; }
763 
764  /// Get the vector grid function.
765  const GridFunction * GetGridFunction() const { return GridFunc; }
766 
767  /// Evaluate the scalar divergence coefficient at @a ip.
768  virtual double Eval(ElementTransformation &T,
769  const IntegrationPoint &ip);
770 
772 };
773 
774 /** @brief Vector coefficient defined by a scalar DeltaCoefficient and a
775  constant vector direction.
776 
777  WARNING this cannot be used as a normal coefficient. The usual Eval method
778  is disabled. */
780 {
781 protected:
784 
785 public:
786  /// Construct with a vector of dimension @a vdim_.
788  : VectorCoefficient(vdim_), dir(vdim_), d() { }
789 
790  /** @brief Construct with a Vector object representing the direction and a
791  unit delta function centered at (0.0,0.0,0.0) */
793  : VectorCoefficient(dir_.Size()), dir(dir_), d() { }
794 
795  /** @brief Construct with a Vector object representing the direction and a
796  delta function scaled by @a s and centered at (x,0.0,0.0) */
797  VectorDeltaCoefficient(const Vector& dir_, double x, double s)
798  : VectorCoefficient(dir_.Size()), dir(dir_), d(x,s) { }
799 
800  /** @brief Construct with a Vector object representing the direction and a
801  delta function scaled by @a s and centered at (x,y,0.0) */
802  VectorDeltaCoefficient(const Vector& dir_, double x, double y, double s)
803  : VectorCoefficient(dir_.Size()), dir(dir_), d(x,y,s) { }
804 
805  /** @brief Construct with a Vector object representing the direction and a
806  delta function scaled by @a s and centered at (x,y,z) */
807  VectorDeltaCoefficient(const Vector& dir_, double x, double y, double z,
808  double s)
809  : VectorCoefficient(dir_.Size()), dir(dir_), d(x,y,z,s) { }
810 
811  /// Set the time for internally stored coefficients
812  void SetTime(double t);
813 
814  /// Replace the associated DeltaCoefficient with a new DeltaCoefficient.
815  /** The new DeltaCoefficient cannot have a specified weight Coefficient, i.e.
816  DeltaCoefficient::Weight() should return NULL. */
817  void SetDeltaCoefficient(const DeltaCoefficient& d_) { d = d_; }
818 
819  /// Return the associated scalar DeltaCoefficient.
821 
822  void SetScale(double s) { d.SetScale(s); }
823  void SetDirection(const Vector& d_);
824 
825  void SetDeltaCenter(const Vector& center) { d.SetDeltaCenter(center); }
826  void GetDeltaCenter(Vector& center) { d.GetDeltaCenter(center); }
827 
828  /** @brief Return the specified direction vector multiplied by the value
829  returned by DeltaCoefficient::EvalDelta() of the associated scalar
830  DeltaCoefficient. */
831  virtual void EvalDelta(Vector &V, ElementTransformation &T,
832  const IntegrationPoint &ip);
833 
835  /** @brief A VectorDeltaFunction cannot be evaluated. Calling this method
836  will cause an MFEM error, terminating the application. */
837  virtual void Eval(Vector &V, ElementTransformation &T,
838  const IntegrationPoint &ip)
839  { mfem_error("VectorDeltaCoefficient::Eval"); }
841 };
842 
843 /** @brief Derived vector coefficient that has the value of the parent vector
844  where it is active and is zero otherwise. */
846 {
847 private:
849  Array<int> active_attr;
850 
851 public:
852  /** @brief Construct with a parent vector coefficient and an array of zeros
853  and ones representing the attributes for which this coefficient should be
854  active. */
856  : VectorCoefficient(vc.GetVDim())
857  { c = &vc; attr.Copy(active_attr); }
858 
859  /// Set the time for internally stored coefficients
860  void SetTime(double t);
861 
862  /// Evaluate the vector coefficient at @a ip.
863  virtual void Eval(Vector &V, ElementTransformation &T,
864  const IntegrationPoint &ip);
865 
866  /** @brief Evaluate the vector coefficient at all of the locations in the
867  integration rule and write the vectors into the columns of matrix @a
868  M. */
869  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
870  const IntegrationRule &ir);
871 };
872 
874 
875 /// Base class for Matrix Coefficients that optionally depend on time and space.
877 {
878 protected:
879  int height, width;
880  double time;
881  bool symmetric; // deprecated
882 
883 public:
884  /// Construct a dim x dim matrix coefficient.
885  explicit MatrixCoefficient(int dim, bool symm=false)
886  { height = width = dim; time = 0.; symmetric = symm; }
887 
888  /// Construct a h x w matrix coefficient.
889  MatrixCoefficient(int h, int w, bool symm=false) :
890  height(h), width(w), time(0.), symmetric(symm) { }
891 
892  /// Set the time for time dependent coefficients
893  virtual void SetTime(double t) { time = t; }
894 
895  /// Get the time for time dependent coefficients
896  double GetTime() { return time; }
897 
898  /// Get the height of the matrix.
899  int GetHeight() const { return height; }
900 
901  /// Get the width of the matrix.
902  int GetWidth() const { return width; }
903 
904  /// For backward compatibility get the width of the matrix.
905  int GetVDim() const { return width; }
906 
907  /** @deprecated Use SymmetricMatrixCoefficient instead */
908  bool IsSymmetric() const { return symmetric; }
909 
910  /** @brief Evaluate the matrix coefficient in the element described by @a T
911  at the point @a ip, storing the result in @a K. */
912  /** @note When this method is called, the caller must make sure that the
913  IntegrationPoint associated with @a T is the same as @a ip. This can be
914  achieved by calling T.SetIntPoint(&ip). */
915  virtual void Eval(DenseMatrix &K, ElementTransformation &T,
916  const IntegrationPoint &ip) = 0;
917 
918  /// (DEPRECATED) Evaluate a symmetric matrix coefficient.
919  /** @brief Evaluate the upper triangular entries of the matrix coefficient
920  in the symmetric case, similarly to Eval. Matrix entry (i,j) is stored
921  in K[j - i + os_i] for 0 <= i <= j < width, os_0 = 0,
922  os_{i+1} = os_i + width - i. That is, K = {M(0,0), ..., M(0,w-1),
923  M(1,1), ..., M(1,w-1), ..., M(w-1,w-1) with w = width.
924  @deprecated Use Eval() instead. */
926  const IntegrationPoint &ip)
927  { mfem_error("MatrixCoefficient::EvalSymmetric"); }
928 
929  virtual ~MatrixCoefficient() { }
930 };
931 
932 
933 /// A matrix coefficient that is constant in space and time.
935 {
936 private:
937  DenseMatrix mat;
938 public:
939  ///Construct using matrix @a m for the constant.
941  : MatrixCoefficient(m.Height(), m.Width()), mat(m) { }
943  /// Evaluate the matrix coefficient at @a ip.
945  const IntegrationPoint &ip) { M = mat; }
946 };
947 
948 
949 /** @brief A piecewise matrix-valued coefficient with the pieces keyed off the
950  element attribute numbers.
951 
952  A value of zero will be returned for any missing attribute numbers.
953 
954  This object will not assume ownership of any MatrixCoefficient
955  objects passed to it. Consequently, the caller must ensure that
956  the individual MatrixCoefficient objects are not deleted while
957  this PWMatrixCoefficient is still in use.
958 
959  \note The keys may either be domain attribute numbers or boundary
960  attribute numbers. If the PWMatrixCoefficient is used with a
961  domain integrator the keys are assumed to be domain attribute
962  numbers. Similarly, if the PWMatrixCoefficient is used with a
963  boundary integrator the keys are assumed to be boundary attribute
964  numbers.
965 */
967 {
968 private:
969  /** Internal data structure to store pointers to the appropriate
970  coefficients for different regions of the mesh. The keys used
971  in the map are the mesh attribute numbers (either element
972  attribute or boundary element attribute depending upon
973  context). The values returned for any missing attributes will
974  be zero. The coefficient pointers may be NULL in which case a
975  value of zero is returned.
976 
977  The MatrixCoefficient objects contained in this map are NOT
978  owned by this PWMatrixCoefficient object. This means that they
979  will not be deleted when this object is deleted also the caller
980  must ensure that the various MatrixCoefficient objects are not
981  deleted while this PWMatrixCoefficient is still needed.
982  */
983  std::map<int, MatrixCoefficient*> pieces;
984 
985  /** Convenience function to check for compatible array lengths,
986  loop over the arrays, and add their attribute/MatrixCoefficient
987  pairs to the internal data structure.
988  */
989  void InitMap(const Array<int> & attr,
990  const Array<MatrixCoefficient*> & coefs);
991 
992 public:
993 
994  /// Constructs a piecewise matrix coefficient of dimension dim by dim
995  explicit PWMatrixCoefficient(int dim, bool symm = false)
996  : MatrixCoefficient(dim, symm) {}
997 
998  /// Constructs a piecewise matrix coefficient of dimension h by w
999  explicit PWMatrixCoefficient(int h, int w, bool symm = false)
1000  : MatrixCoefficient(h, w, symm) {}
1001 
1002  /// Construct the coefficient using arrays describing the pieces
1003  /** \param dim - size of the square matrix-valued result
1004  \param attr - an array of attribute numbers for each piece
1005  \param coefs - the corresponding array of MatrixCoefficient pointers
1006  \param symm - true if the result will be symmetric, false otherwise
1007  Any missing attributes or NULL coefficient pointers will result in a
1008  zero matrix being returned.
1009 
1010  \note Ownership of the MatrixCoefficient objects will NOT be
1011  transferred to this object.
1012  */
1013  PWMatrixCoefficient(int dim, const Array<int> & attr,
1014  const Array<MatrixCoefficient*> & coefs,
1015  bool symm=false)
1016  : MatrixCoefficient(dim, symm) { InitMap(attr, coefs); }
1017 
1018  /// Construct the coefficient using arrays describing the pieces
1019  /** \param h - height of the matrix-valued result
1020  \param w - width of the matrix-valued result
1021  \param attr - an array of attribute numbers for each piece
1022  \param coefs - the corresponding array of MatrixCoefficient pointers
1023  \param symm - true if the result will be symmetric, false otherwise
1024  Any missing attributes or NULL coefficient pointers will result in a
1025  zero matrix being returned for that attribute.
1026 
1027  \note Ownership of the MatrixCoefficient objects will NOT be
1028  transferred to this object.
1029  */
1030  PWMatrixCoefficient(int h, int w, const Array<int> & attr,
1031  const Array<MatrixCoefficient*> & coefs,
1032  bool symm=false)
1033  : MatrixCoefficient(h, w, symm) { InitMap(attr, coefs); }
1034 
1035  /// Set the time for time dependent coefficients
1036  virtual void SetTime(double t);
1037 
1038  /// Replace a set of coefficients
1039  void UpdateCoefficients(const Array<int> & attr,
1040  const Array<MatrixCoefficient*> & coefs)
1041  { InitMap(attr, coefs); }
1042 
1043  /// Replace a single coefficient for a particular attribute
1044  void UpdateCoefficient(int attr, MatrixCoefficient & coef);
1045 
1046  /// Remove a single MatrixCoefficient for a particular attribute
1047  void ZeroCoefficient(int attr)
1048  { pieces.erase(attr); }
1049 
1050  /// Evaluate the coefficient.
1051  virtual void Eval(DenseMatrix &K, ElementTransformation &T,
1052  const IntegrationPoint &ip);
1053 };
1054 
1055 /** @brief A matrix coefficient with an optional scalar coefficient multiplier
1056  \a q. The matrix function can either be represented by a std function or
1057  a constant matrix provided when constructing this object. */
1059 {
1060 private:
1061  std::function<void(const Vector &, DenseMatrix &)> Function;
1062  std::function<void(const Vector &, Vector &)> SymmFunction; // deprecated
1063  std::function<void(const Vector &, double, DenseMatrix &)> TDFunction;
1064 
1065  Coefficient *Q;
1066  DenseMatrix mat;
1067 
1068 public:
1069  /// Define a time-independent square matrix coefficient from a std function
1070  /** \param dim - the size of the matrix
1071  \param F - time-independent function
1072  \param q - optional scalar Coefficient to scale the matrix coefficient */
1074  std::function<void(const Vector &, DenseMatrix &)> F,
1075  Coefficient *q = nullptr)
1076  : MatrixCoefficient(dim), Function(std::move(F)), Q(q), mat(0)
1077  { }
1078 
1079  /// Define a constant matrix coefficient times a scalar Coefficient
1080  /** \param m - constant matrix
1081  \param q - optional scalar Coefficient to scale the matrix coefficient */
1083  : MatrixCoefficient(m.Height(), m.Width()), Q(&q), mat(m)
1084  { }
1085 
1086  /** @brief Define a time-independent symmetric square matrix coefficient from
1087  a std function */
1088  /** \param dim - the size of the matrix
1089  \param SymmF - function used in EvalSymmetric
1090  \param q - optional scalar Coefficient to scale the matrix coefficient
1091  @deprecated Use another constructor without setting SymmFunction. */
1093  std::function<void(const Vector &, Vector &)> SymmF,
1094  Coefficient *q = NULL)
1095  : MatrixCoefficient(dim, true), SymmFunction(std::move(SymmF)), Q(q), mat(0)
1096  { }
1097 
1098  /// Define a time-dependent square matrix coefficient from a std function
1099  /** \param dim - the size of the matrix
1100  \param TDF - time-dependent function
1101  \param q - optional scalar Coefficient to scale the matrix coefficient */
1103  std::function<void(const Vector &, double, DenseMatrix &)> TDF,
1104  Coefficient *q = nullptr)
1105  : MatrixCoefficient(dim), TDFunction(std::move(TDF)), Q(q)
1106  { }
1107 
1108  /// Set the time for internally stored coefficients
1109  void SetTime(double t);
1110 
1111  /// Evaluate the matrix coefficient at @a ip.
1112  virtual void Eval(DenseMatrix &K, ElementTransformation &T,
1113  const IntegrationPoint &ip);
1114 
1115  /// (DEPRECATED) Evaluate the symmetric matrix coefficient at @a ip.
1116  /** @deprecated Use Eval() instead. */
1117  virtual void EvalSymmetric(Vector &K, ElementTransformation &T,
1118  const IntegrationPoint &ip);
1119 
1121 };
1122 
1123 
1124 /** @brief Matrix coefficient defined by a matrix of scalar coefficients.
1125  Coefficients that are not set will evaluate to zero in the vector. The
1126  coefficient is stored as a flat Array with indexing (i,j) -> i*width+j. */
1128 {
1129 private:
1130  Array<Coefficient *> Coeff;
1131  Array<bool> ownCoeff;
1132 
1133 public:
1134  /** @brief Construct a coefficient matrix of dimensions @a dim * @a dim. The
1135  actual coefficients still need to be added with Set(). */
1136  explicit MatrixArrayCoefficient (int dim);
1137 
1138  /// Set the time for internally stored coefficients
1139  void SetTime(double t);
1140 
1141  /// Get the coefficient located at (i,j) in the matrix.
1142  Coefficient* GetCoeff (int i, int j) { return Coeff[i*width+j]; }
1143 
1144  /** @brief Set the coefficient located at (i,j) in the matrix. By default by
1145  default this will take ownership of the Coefficient passed in, but this
1146  can be overridden with the @a own parameter. */
1147  void Set(int i, int j, Coefficient * c, bool own=true);
1148 
1149  /// Evaluate coefficient located at (i,j) in the matrix using integration
1150  /// point @a ip.
1151  double Eval(int i, int j, ElementTransformation &T, const IntegrationPoint &ip)
1152  { return Coeff[i*width+j] ? Coeff[i*width+j] -> Eval(T, ip, GetTime()) : 0.0; }
1153 
1154  /// Evaluate the matrix coefficient @a ip.
1155  virtual void Eval(DenseMatrix &K, ElementTransformation &T,
1156  const IntegrationPoint &ip);
1157 
1158  virtual ~MatrixArrayCoefficient();
1159 };
1160 
1161 
1162 /** @brief Derived matrix coefficient that has the value of the parent matrix
1163  coefficient where it is active and is zero otherwise. */
1165 {
1166 private:
1167  MatrixCoefficient *c;
1168  Array<int> active_attr;
1169 
1170 public:
1171  /** @brief Construct with a parent matrix coefficient and an array of zeros
1172  and ones representing the attributes for which this coefficient should be
1173  active. */
1175  : MatrixCoefficient(mc.GetHeight(), mc.GetWidth())
1176  { c = &mc; attr.Copy(active_attr); }
1177 
1178  /// Set the time for internally stored coefficients
1179  void SetTime(double t);
1180 
1181  /// Evaluate the matrix coefficient at @a ip.
1182  virtual void Eval(DenseMatrix &K, ElementTransformation &T,
1183  const IntegrationPoint &ip);
1184 };
1185 
1186 /// Coefficients based on sums, products, or other functions of coefficients.
1187 ///@{
1188 /** @brief Scalar coefficient defined as the linear combination of two scalar
1189  coefficients or a scalar and a scalar coefficient */
1191 {
1192 private:
1193  double aConst;
1194  Coefficient * a;
1195  Coefficient * b;
1196 
1197  double alpha;
1198  double beta;
1199 
1200 public:
1201  /// Constructor with one coefficient. Result is alpha_ * A + beta_ * B
1203  double alpha_ = 1.0, double beta_ = 1.0)
1204  : aConst(A), a(NULL), b(&B), alpha(alpha_), beta(beta_) { }
1205 
1206  /// Constructor with two coefficients. Result is alpha_ * A + beta_ * B.
1208  double alpha_ = 1.0, double beta_ = 1.0)
1209  : aConst(0.0), a(&A), b(&B), alpha(alpha_), beta(beta_) { }
1210 
1211  /// Set the time for internally stored coefficients
1212  void SetTime(double t);
1213 
1214  /// Reset the first term in the linear combination as a constant
1215  void SetAConst(double A) { a = NULL; aConst = A; }
1216  /// Return the first term in the linear combination
1217  double GetAConst() const { return aConst; }
1218 
1219  /// Reset the first term in the linear combination
1220  void SetACoef(Coefficient &A) { a = &A; }
1221  /// Return the first term in the linear combination
1222  Coefficient * GetACoef() const { return a; }
1223 
1224  /// Reset the second term in the linear combination
1225  void SetBCoef(Coefficient &B) { b = &B; }
1226  /// Return the second term in the linear combination
1227  Coefficient * GetBCoef() const { return b; }
1228 
1229  /// Reset the factor in front of the first term in the linear combination
1230  void SetAlpha(double alpha_) { alpha = alpha_; }
1231  /// Return the factor in front of the first term in the linear combination
1232  double GetAlpha() const { return alpha; }
1233 
1234  /// Reset the factor in front of the second term in the linear combination
1235  void SetBeta(double beta_) { beta = beta_; }
1236  /// Return the factor in front of the second term in the linear combination
1237  double GetBeta() const { return beta; }
1238 
1239  /// Evaluate the coefficient at @a ip.
1240  virtual double Eval(ElementTransformation &T,
1241  const IntegrationPoint &ip)
1242  {
1243  return alpha * ((a == NULL ) ? aConst : a->Eval(T, ip) )
1244  + beta * b->Eval(T, ip);
1245  }
1246 };
1247 
1248 
1249 /// Base class for symmetric matrix coefficients that optionally depend on time and space.
1251 {
1252 protected:
1253  int dim;
1254  double time;
1255 
1256 public:
1257  /// Construct a dim x dim matrix coefficient.
1258  explicit SymmetricMatrixCoefficient(int dimension)
1259  { dim = dimension; time = 0.; }
1260 
1261  /// Set the time for time dependent coefficients
1262  virtual void SetTime(double t) { time = t; }
1263 
1264  /// Get the time for time dependent coefficients
1265  double GetTime() { return time; }
1266 
1267  /// Get the size of the matrix.
1268  int GetSize() const { return dim; }
1269 
1270  /** @brief Evaluate the matrix coefficient in the element described by @a T
1271  at the point @a ip, storing the result in @a K. */
1272  /** @note When this method is called, the caller must make sure that the
1273  IntegrationPoint associated with @a T is the same as @a ip. This can be
1274  achieved by calling T.SetIntPoint(&ip). */
1275  virtual void Eval(DenseSymmetricMatrix &K, ElementTransformation &T,
1276  const IntegrationPoint &ip) = 0;
1277 
1279 };
1280 
1281 
1282 /// A matrix coefficient that is constant in space and time.
1284 {
1285 private:
1287 
1288 public:
1289  ///Construct using matrix @a m for the constant.
1291  : SymmetricMatrixCoefficient(m.Height()), mat(m) { }
1293  /// Evaluate the matrix coefficient at @a ip.
1295  const IntegrationPoint &ip) { M = mat; }
1296 };
1297 
1298 
1299 /** @brief A matrix coefficient with an optional scalar coefficient multiplier
1300  \a q. The matrix function can either be represented by a std function or
1301  a constant matrix provided when constructing this object. */
1303 {
1304 private:
1305  std::function<void(const Vector &, DenseSymmetricMatrix &)> Function;
1306  std::function<void(const Vector &, double, DenseSymmetricMatrix &)> TDFunction;
1307 
1308  Coefficient *Q;
1310 
1311 public:
1312  /// Define a time-independent symmetric matrix coefficient from a std function
1313  /** \param dim - the size of the matrix
1314  \param F - time-independent function
1315  \param q - optional scalar Coefficient to scale the matrix coefficient */
1317  std::function<void(const Vector &, DenseSymmetricMatrix &)> F,
1318  Coefficient *q = nullptr)
1319  : SymmetricMatrixCoefficient(dim), Function(std::move(F)), Q(q), mat(0)
1320  { }
1321 
1322  /// Define a constant matrix coefficient times a scalar Coefficient
1323  /** \param m - constant matrix
1324  \param q - optional scalar Coefficient to scale the matrix coefficient */
1326  Coefficient &q)
1327  : SymmetricMatrixCoefficient(m.Height()), Q(&q), mat(m)
1328  { }
1329 
1330  /// Define a time-dependent square matrix coefficient from a std function
1331  /** \param dim - the size of the matrix
1332  \param TDF - time-dependent function
1333  \param q - optional scalar Coefficient to scale the matrix coefficient */
1335  std::function<void(const Vector &, double, DenseSymmetricMatrix &)> TDF,
1336  Coefficient *q = nullptr)
1337  : SymmetricMatrixCoefficient(dim), TDFunction(std::move(TDF)), Q(q)
1338  { }
1339 
1340  /// Set the time for internally stored coefficients
1341  void SetTime(double t);
1342 
1343  /// Evaluate the matrix coefficient at @a ip.
1344  virtual void Eval(DenseSymmetricMatrix &K, ElementTransformation &T,
1345  const IntegrationPoint &ip);
1346 
1348 };
1349 
1350 
1351 /** @brief Scalar coefficient defined as the product of two scalar coefficients
1352  or a scalar and a scalar coefficient. */
1354 {
1355 private:
1356  double aConst;
1357  Coefficient * a;
1358  Coefficient * b;
1359 
1360 public:
1361  /// Constructor with one coefficient. Result is A * B.
1363  : aConst(A), a(NULL), b(&B) { }
1364 
1365  /// Constructor with two coefficients. Result is A * B.
1367  : aConst(0.0), a(&A), b(&B) { }
1368 
1369  /// Set the time for internally stored coefficients
1370  void SetTime(double t);
1371 
1372  /// Reset the first term in the product as a constant
1373  void SetAConst(double A) { a = NULL; aConst = A; }
1374  /// Return the first term in the product
1375  double GetAConst() const { return aConst; }
1376 
1377  /// Reset the first term in the product
1378  void SetACoef(Coefficient &A) { a = &A; }
1379  /// Return the first term in the product
1380  Coefficient * GetACoef() const { return a; }
1381 
1382  /// Reset the second term in the product
1383  void SetBCoef(Coefficient &B) { b = &B; }
1384  /// Return the second term in the product
1385  Coefficient * GetBCoef() const { return b; }
1386 
1387  /// Evaluate the coefficient at @a ip.
1388  virtual double Eval(ElementTransformation &T,
1389  const IntegrationPoint &ip)
1390  { return ((a == NULL ) ? aConst : a->Eval(T, ip) ) * b->Eval(T, ip); }
1391 };
1392 
1393 /** @brief Scalar coefficient defined as the ratio of two scalars where one or
1394  both scalars are scalar coefficients. */
1396 {
1397 private:
1398  double aConst;
1399  double bConst;
1400  Coefficient * a;
1401  Coefficient * b;
1402 
1403 public:
1404  /** Initialize a coefficient which returns A / B where @a A is a
1405  constant and @a B is a scalar coefficient */
1407  : aConst(A), bConst(1.0), a(NULL), b(&B) { }
1408  /** Initialize a coefficient which returns A / B where @a A and @a B are both
1409  scalar coefficients */
1411  : aConst(0.0), bConst(1.0), a(&A), b(&B) { }
1412  /** Initialize a coefficient which returns A / B where @a A is a
1413  scalar coefficient and @a B is a constant */
1415  : aConst(0.0), bConst(B), a(&A), b(NULL) { }
1416 
1417  /// Set the time for internally stored coefficients
1418  void SetTime(double t);
1419 
1420  /// Reset the numerator in the ratio as a constant
1421  void SetAConst(double A) { a = NULL; aConst = A; }
1422  /// Return the numerator of the ratio
1423  double GetAConst() const { return aConst; }
1424 
1425  /// Reset the denominator in the ratio as a constant
1426  void SetBConst(double B) { b = NULL; bConst = B; }
1427  /// Return the denominator of the ratio
1428  double GetBConst() const { return bConst; }
1429 
1430  /// Reset the numerator in the ratio
1431  void SetACoef(Coefficient &A) { a = &A; }
1432  /// Return the numerator of the ratio
1433  Coefficient * GetACoef() const { return a; }
1434 
1435  /// Reset the denominator in the ratio
1436  void SetBCoef(Coefficient &B) { b = &B; }
1437  /// Return the denominator of the ratio
1438  Coefficient * GetBCoef() const { return b; }
1439 
1440  /// Evaluate the coefficient
1441  virtual double Eval(ElementTransformation &T,
1442  const IntegrationPoint &ip)
1443  {
1444  double den = (b == NULL ) ? bConst : b->Eval(T, ip);
1445  MFEM_ASSERT(den != 0.0, "Division by zero in RatioCoefficient");
1446  return ((a == NULL ) ? aConst : a->Eval(T, ip) ) / den;
1447  }
1448 };
1449 
1450 /// Scalar coefficient defined as a scalar raised to a power
1452 {
1453 private:
1454  Coefficient * a;
1455 
1456  double p;
1457 
1458 public:
1459  /// Construct with a coefficient and a constant power @a p_. Result is A^p.
1461  : a(&A), p(p_) { }
1462 
1463  /// Set the time for internally stored coefficients
1464  void SetTime(double t);
1465 
1466  /// Reset the base coefficient
1467  void SetACoef(Coefficient &A) { a = &A; }
1468  /// Return the base coefficient
1469  Coefficient * GetACoef() const { return a; }
1470 
1471  /// Reset the exponent
1472  void SetExponent(double p_) { p = p_; }
1473  /// Return the exponent
1474  double GetExponent() const { return p; }
1475 
1476  /// Evaluate the coefficient at @a ip.
1477  virtual double Eval(ElementTransformation &T,
1478  const IntegrationPoint &ip)
1479  { return pow(a->Eval(T, ip), p); }
1480 };
1481 
1482 
1483 /// Scalar coefficient defined as the inner product of two vector coefficients
1485 {
1486 private:
1487  VectorCoefficient * a;
1488  VectorCoefficient * b;
1489 
1490  mutable Vector va;
1491  mutable Vector vb;
1492 public:
1493  /// Construct with the two vector coefficients. Result is \f$ A \cdot B \f$.
1495 
1496  /// Set the time for internally stored coefficients
1497  void SetTime(double t);
1498 
1499  /// Reset the first vector in the inner product
1500  void SetACoef(VectorCoefficient &A) { a = &A; }
1501  /// Return the first vector coefficient in the inner product
1502  VectorCoefficient * GetACoef() const { return a; }
1503 
1504  /// Reset the second vector in the inner product
1505  void SetBCoef(VectorCoefficient &B) { b = &B; }
1506  /// Return the second vector coefficient in the inner product
1507  VectorCoefficient * GetBCoef() const { return b; }
1508 
1509  /// Evaluate the coefficient at @a ip.
1510  virtual double Eval(ElementTransformation &T,
1511  const IntegrationPoint &ip);
1512 };
1513 
1514 /// Scalar coefficient defined as a cross product of two vectors in the xy-plane.
1516 {
1517 private:
1518  VectorCoefficient * a;
1519  VectorCoefficient * b;
1520 
1521  mutable Vector va;
1522  mutable Vector vb;
1523 
1524 public:
1525  /// Constructor with two vector coefficients. Result is \f$ A_x B_y - A_y * B_x; \f$.
1527 
1528  /// Set the time for internally stored coefficients
1529  void SetTime(double t);
1530 
1531  /// Reset the first vector in the product
1532  void SetACoef(VectorCoefficient &A) { a = &A; }
1533  /// Return the first vector of the product
1534  VectorCoefficient * GetACoef() const { return a; }
1535 
1536  /// Reset the second vector in the product
1537  void SetBCoef(VectorCoefficient &B) { b = &B; }
1538  /// Return the second vector of the product
1539  VectorCoefficient * GetBCoef() const { return b; }
1540 
1541  /// Evaluate the coefficient at @a ip.
1542  virtual double Eval(ElementTransformation &T,
1543  const IntegrationPoint &ip);
1544 };
1545 
1546 /// Scalar coefficient defined as the determinant of a matrix coefficient
1548 {
1549 private:
1550  MatrixCoefficient * a;
1551 
1552  mutable DenseMatrix ma;
1553 
1554 public:
1555  /// Construct with the matrix.
1557 
1558  /// Set the time for internally stored coefficients
1559  void SetTime(double t);
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 determinant coefficient at @a ip.
1567  virtual double Eval(ElementTransformation &T,
1568  const IntegrationPoint &ip);
1569 };
1570 
1571 /// Vector coefficient defined as the linear combination of two vectors
1573 {
1574 private:
1575  VectorCoefficient * ACoef;
1576  VectorCoefficient * BCoef;
1577 
1578  Vector A;
1579  Vector B;
1580 
1581  Coefficient * alphaCoef;
1582  Coefficient * betaCoef;
1583 
1584  double alpha;
1585  double beta;
1586 
1587  mutable Vector va;
1588 
1589 public:
1590  /** Constructor with no coefficients.
1591  To be used with the various "Set" methods */
1592  VectorSumCoefficient(int dim);
1593 
1594  /** Constructor with two vector coefficients.
1595  Result is alpha_ * A + beta_ * B */
1597  double alpha_ = 1.0, double beta_ = 1.0);
1598 
1599  /** Constructor with scalar coefficients.
1600  Result is alpha_ * A_ + beta_ * B_ */
1602  Coefficient &alpha_, Coefficient &beta_);
1603 
1604  /// Set the time for internally stored coefficients
1605  void SetTime(double t);
1606 
1607  /// Reset the first vector coefficient
1608  void SetACoef(VectorCoefficient &A_) { ACoef = &A_; }
1609  /// Return the first vector coefficient
1610  VectorCoefficient * GetACoef() const { return ACoef; }
1611 
1612  /// Reset the second vector coefficient
1613  void SetBCoef(VectorCoefficient &B_) { BCoef = &B_; }
1614  /// Return the second vector coefficient
1615  VectorCoefficient * GetBCoef() const { return BCoef; }
1616 
1617  /// Reset the factor in front of the first vector coefficient
1618  void SetAlphaCoef(Coefficient &A_) { alphaCoef = &A_; }
1619  /// Return the factor in front of the first vector coefficient
1620  Coefficient * GetAlphaCoef() const { return alphaCoef; }
1621 
1622  /// Reset the factor in front of the second vector coefficient
1623  void SetBetaCoef(Coefficient &B_) { betaCoef = &B_; }
1624  /// Return the factor in front of the second vector coefficient
1625  Coefficient * GetBetaCoef() const { return betaCoef; }
1626 
1627  /// Reset the first vector as a constant
1628  void SetA(const Vector &A_) { A = A_; ACoef = NULL; }
1629  /// Return the first vector constant
1630  const Vector & GetA() const { return A; }
1631 
1632  /// Reset the second vector as a constant
1633  void SetB(const Vector &B_) { B = B_; BCoef = NULL; }
1634  /// Return the second vector constant
1635  const Vector & GetB() const { return B; }
1636 
1637  /// Reset the factor in front of the first vector coefficient as a constant
1638  void SetAlpha(double alpha_) { alpha = alpha_; alphaCoef = NULL; }
1639  /// Return the factor in front of the first vector coefficient
1640  double GetAlpha() const { return alpha; }
1641 
1642  /// Reset the factor in front of the second vector coefficient as a constant
1643  void SetBeta(double beta_) { beta = beta_; betaCoef = NULL; }
1644  /// Return the factor in front of the second vector coefficient
1645  double GetBeta() const { return beta; }
1646 
1647  /// Evaluate the coefficient at @a ip.
1648  virtual void Eval(Vector &V, ElementTransformation &T,
1649  const IntegrationPoint &ip);
1651 };
1652 
1653 /// Vector coefficient defined as a product of scalar and vector coefficients.
1655 {
1656 private:
1657  double aConst;
1658  Coefficient * a;
1659  VectorCoefficient * b;
1660 
1661 public:
1662  /// Constructor with constant and vector coefficient. Result is A * B.
1664 
1665  /// Constructor with two coefficients. Result is A * B.
1667 
1668  /// Set the time for internally stored coefficients
1669  void SetTime(double t);
1670 
1671  /// Reset the scalar factor as a constant
1672  void SetAConst(double A) { a = NULL; aConst = A; }
1673  /// Return the scalar factor
1674  double GetAConst() const { return aConst; }
1675 
1676  /// Reset the scalar factor
1677  void SetACoef(Coefficient &A) { a = &A; }
1678  /// Return the scalar factor
1679  Coefficient * GetACoef() const { return a; }
1680 
1681  /// Reset the vector factor
1682  void SetBCoef(VectorCoefficient &B) { b = &B; }
1683  /// Return the vector factor
1684  VectorCoefficient * GetBCoef() const { return b; }
1685 
1686  /// Evaluate the coefficient at @a ip.
1687  virtual void Eval(Vector &V, ElementTransformation &T,
1688  const IntegrationPoint &ip);
1690 };
1691 
1692 /// Vector coefficient defined as a normalized vector field (returns v/|v|)
1694 {
1695 private:
1696  VectorCoefficient * a;
1697 
1698  double tol;
1699 
1700 public:
1701  /** @brief Return a vector normalized to a length of one
1702 
1703  This class evaluates the vector coefficient @a A and, if |A| > @a tol,
1704  returns the normalized vector A / |A|. If |A| <= @a tol, the zero
1705  vector is returned.
1706  */
1707  NormalizedVectorCoefficient(VectorCoefficient &A, double tol = 1e-6);
1708 
1709  /// Set the time for internally stored coefficients
1710  void SetTime(double t);
1711 
1712  /// Reset the vector coefficient
1713  void SetACoef(VectorCoefficient &A) { a = &A; }
1714  /// Return the vector coefficient
1715  VectorCoefficient * GetACoef() const { return a; }
1716 
1717  /// Evaluate the coefficient at @a ip.
1718  virtual void Eval(Vector &V, ElementTransformation &T,
1719  const IntegrationPoint &ip);
1721 };
1722 
1723 /// Vector coefficient defined as a cross product of two vectors
1725 {
1726 private:
1727  VectorCoefficient * a;
1728  VectorCoefficient * b;
1729 
1730  mutable Vector va;
1731  mutable Vector vb;
1732 
1733 public:
1734  /// Construct with the two coefficients. Result is A x B.
1736 
1737  /// Set the time for internally stored coefficients
1738  void SetTime(double t);
1739 
1740  /// Reset the first term in the product
1741  void SetACoef(VectorCoefficient &A) { a = &A; }
1742  /// Return the first term in the product
1743  VectorCoefficient * GetACoef() const { return a; }
1744 
1745  /// Reset the second term in the product
1746  void SetBCoef(VectorCoefficient &B) { b = &B; }
1747  /// Return the second term in the product
1748  VectorCoefficient * GetBCoef() const { return b; }
1749 
1750  /// Evaluate the coefficient at @a ip.
1751  virtual void Eval(Vector &V, ElementTransformation &T,
1752  const IntegrationPoint &ip);
1754 };
1755 
1756 /** @brief Vector coefficient defined as a product of a matrix coefficient and
1757  a vector coefficient. */
1759 {
1760 private:
1761  MatrixCoefficient * a;
1762  VectorCoefficient * b;
1763 
1764  mutable DenseMatrix ma;
1765  mutable Vector vb;
1766 
1767 public:
1768  /// Constructor with two coefficients. Result is A*B.
1770 
1771  /// Set the time for internally stored coefficients
1772  void SetTime(double t);
1773 
1774  /// Reset the matrix coefficient
1775  void SetACoef(MatrixCoefficient &A) { a = &A; }
1776  /// Return the matrix coefficient
1777  MatrixCoefficient * GetACoef() const { return a; }
1778 
1779  /// Reset the vector coefficient
1780  void SetBCoef(VectorCoefficient &B) { b = &B; }
1781  /// Return the vector coefficient
1782  VectorCoefficient * GetBCoef() const { return b; }
1783 
1784  /// Evaluate the vector coefficient at @a ip.
1785  virtual void Eval(Vector &V, ElementTransformation &T,
1786  const IntegrationPoint &ip);
1788 };
1789 
1790 /// Convenient alias for the MatrixVectorProductCoefficient
1792 
1793 /// Constant matrix coefficient defined as the identity of dimension d
1795 {
1796 private:
1797  int dim;
1798 
1799 public:
1800  /// Construct with the dimension of the square identity matrix.
1802  : MatrixCoefficient(d, d), dim(d) { }
1803 
1804  /// Evaluate the matrix coefficient at @a ip.
1805  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1806  const IntegrationPoint &ip);
1807 };
1808 
1809 /// Matrix coefficient defined as the linear combination of two matrices
1811 {
1812 private:
1813  MatrixCoefficient * a;
1814  MatrixCoefficient * b;
1815 
1816  double alpha;
1817  double beta;
1818 
1819  mutable DenseMatrix ma;
1820 
1821 public:
1822  /// Construct with the two coefficients. Result is alpha_ * A + beta_ * B.
1824  double alpha_ = 1.0, double beta_ = 1.0);
1825 
1826  /// Set the time for internally stored coefficients
1827  void SetTime(double t);
1828 
1829  /// Reset the first matrix coefficient
1830  void SetACoef(MatrixCoefficient &A) { a = &A; }
1831  /// Return the first matrix coefficient
1832  MatrixCoefficient * GetACoef() const { return a; }
1833 
1834  /// Reset the second matrix coefficient
1835  void SetBCoef(MatrixCoefficient &B) { b = &B; }
1836  /// Return the second matrix coefficient
1837  MatrixCoefficient * GetBCoef() const { return b; }
1838 
1839  /// Reset the factor in front of the first matrix coefficient
1840  void SetAlpha(double alpha_) { alpha = alpha_; }
1841  /// Return the factor in front of the first matrix coefficient
1842  double GetAlpha() const { return alpha; }
1843 
1844  /// Reset the factor in front of the second matrix coefficient
1845  void SetBeta(double beta_) { beta = beta_; }
1846  /// Return the factor in front of the second matrix coefficient
1847  double GetBeta() const { return beta; }
1848 
1849  /// Evaluate the matrix coefficient at @a ip.
1850  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1851  const IntegrationPoint &ip);
1852 };
1853 
1854 /// Matrix coefficient defined as the product of two matrices
1856 {
1857 private:
1858  MatrixCoefficient * a;
1859  MatrixCoefficient * b;
1860 
1861  mutable DenseMatrix ma;
1862  mutable DenseMatrix mb;
1863 
1864 public:
1865  /// Construct with the two coefficients. Result is A * B.
1867 
1868  /// Reset the first matrix coefficient
1869  void SetACoef(MatrixCoefficient &A) { a = &A; }
1870  /// Return the first matrix coefficient
1871  MatrixCoefficient * GetACoef() const { return a; }
1872 
1873  /// Reset the second matrix coefficient
1874  void SetBCoef(MatrixCoefficient &B) { b = &B; }
1875  /// Return the second matrix coefficient
1876  MatrixCoefficient * GetBCoef() const { return b; }
1877 
1878  /// Evaluate the matrix coefficient at @a ip.
1879  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1880  const IntegrationPoint &ip);
1881 };
1882 
1883 /** @brief Matrix coefficient defined as a product of a scalar coefficient and a
1884  matrix coefficient.*/
1886 {
1887 private:
1888  double aConst;
1889  Coefficient * a;
1890  MatrixCoefficient * b;
1891 
1892 public:
1893  /// Constructor with one coefficient. Result is A*B.
1895 
1896  /// Constructor with two coefficients. Result is A*B.
1898 
1899  /// Set the time for internally stored coefficients
1900  void SetTime(double t);
1901 
1902  /// Reset the scalar factor as a constant
1903  void SetAConst(double A) { a = NULL; aConst = A; }
1904  /// Return the scalar factor
1905  double GetAConst() const { return aConst; }
1906 
1907  /// Reset the scalar factor
1908  void SetACoef(Coefficient &A) { a = &A; }
1909  /// Return the scalar factor
1910  Coefficient * GetACoef() const { return a; }
1911 
1912  /// Reset the matrix factor
1913  void SetBCoef(MatrixCoefficient &B) { b = &B; }
1914  /// Return the matrix factor
1915  MatrixCoefficient * GetBCoef() const { return b; }
1916 
1917  /// Evaluate the matrix coefficient at @a ip.
1918  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1919  const IntegrationPoint &ip);
1920 };
1921 
1922 /// Matrix coefficient defined as the transpose a matrix coefficient
1924 {
1925 private:
1926  MatrixCoefficient * a;
1927 
1928 public:
1929  /// Construct with the matrix coefficient. Result is \f$ A^T \f$.
1931 
1932  /// Set the time for internally stored coefficients
1933  void SetTime(double t);
1934 
1935  /// Reset the matrix coefficient
1936  void SetACoef(MatrixCoefficient &A) { a = &A; }
1937  /// Return the matrix coefficient
1938  MatrixCoefficient * GetACoef() const { return a; }
1939 
1940  /// Evaluate the matrix coefficient at @a ip.
1941  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1942  const IntegrationPoint &ip);
1943 };
1944 
1945 /// Matrix coefficient defined as the inverse a matrix coefficient.
1947 {
1948 private:
1949  MatrixCoefficient * a;
1950 
1951 public:
1952  /// Construct with the matrix coefficient. Result is \f$ A^{-1} \f$.
1954 
1955  /// Set the time for internally stored coefficients
1956  void SetTime(double t);
1957 
1958  /// Reset the matrix coefficient
1959  void SetACoef(MatrixCoefficient &A) { a = &A; }
1960  /// Return the matrix coefficient
1961  MatrixCoefficient * GetACoef() const { return a; }
1962 
1963  /// Evaluate the matrix coefficient at @a ip.
1964  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1965  const IntegrationPoint &ip);
1966 };
1967 
1968 /// Matrix coefficient defined as the outer product of two vector coefficients.
1970 {
1971 private:
1972  VectorCoefficient * a;
1973  VectorCoefficient * b;
1974 
1975  mutable Vector va;
1976  mutable Vector vb;
1977 
1978 public:
1979  /// Construct with two vector coefficients. Result is \f$ A B^T \f$.
1981 
1982  /// Set the time for internally stored coefficients
1983  void SetTime(double t);
1984 
1985  /// Reset the first vector in the outer product
1986  void SetACoef(VectorCoefficient &A) { a = &A; }
1987  /// Return the first vector coefficient in the outer product
1988  VectorCoefficient * GetACoef() const { return a; }
1989 
1990  /// Reset the second vector in the outer product
1991  void SetBCoef(VectorCoefficient &B) { b = &B; }
1992  /// Return the second vector coefficient in the outer product
1993  VectorCoefficient * GetBCoef() const { return b; }
1994 
1995  /// Evaluate the matrix coefficient at @a ip.
1996  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
1997  const IntegrationPoint &ip);
1998 };
1999 
2000 /** @brief Matrix coefficient defined as -a k x k x, for a vector k and scalar a
2001 
2002  This coefficient returns \f$a * (|k|^2 I - k \otimes k)\f$, where I is
2003  the identity matrix and \f$\otimes\f$ indicates the outer product. This
2004  can be evaluated for vectors of any dimension but in three
2005  dimensions it corresponds to computing the cross product with k twice.
2006 */
2008 {
2009 private:
2010  double aConst;
2011  Coefficient * a;
2012  VectorCoefficient * k;
2013 
2014  mutable Vector vk;
2015 
2016 public:
2019 
2020  /// Set the time for internally stored coefficients
2021  void SetTime(double t);
2022 
2023  /// Reset the scalar factor as a constant
2024  void SetAConst(double A) { a = NULL; aConst = A; }
2025  /// Return the scalar factor
2026  double GetAConst() const { return aConst; }
2027 
2028  /// Reset the scalar factor
2029  void SetACoef(Coefficient &A) { a = &A; }
2030  /// Return the scalar factor
2031  Coefficient * GetACoef() const { return a; }
2032 
2033  /// Reset the vector factor
2034  void SetKCoef(VectorCoefficient &K) { k = &K; }
2035  /// Return the vector factor
2036  VectorCoefficient * GetKCoef() const { return k; }
2037 
2038  /// Evaluate the matrix coefficient at @a ip.
2039  virtual void Eval(DenseMatrix &M, ElementTransformation &T,
2040  const IntegrationPoint &ip);
2041 };
2042 ///@}
2043 
2044 class QuadratureFunction;
2045 
2046 /** @brief Vector quadrature function coefficient which requires that the
2047  quadrature rules used for this vector coefficient be the same as those that
2048  live within the supplied QuadratureFunction. */
2050 {
2051 private:
2052  const QuadratureFunction &QuadF; //do not own
2053  int index;
2054 
2055 public:
2056  /// Constructor with a quadrature function as input
2058 
2059  /** Set the starting index within the QuadFunc that'll be used to project
2060  outwards as well as the corresponding length. The projected length should
2061  have the bounds of 1 <= length <= (length QuadFunc - index). */
2062  void SetComponent(int index_, int length_);
2063 
2064  const QuadratureFunction& GetQuadFunction() const { return QuadF; }
2065 
2067  virtual void Eval(Vector &V, ElementTransformation &T,
2068  const IntegrationPoint &ip);
2069 
2071 };
2072 
2073 /** @brief Quadrature function coefficient which requires that the quadrature
2074  rules used for this coefficient be the same as those that live within the
2075  supplied QuadratureFunction. */
2077 {
2078 private:
2079  const QuadratureFunction &QuadF;
2080 
2081 public:
2082  /// Constructor with a quadrature function as input
2084 
2085  const QuadratureFunction& GetQuadFunction() const { return QuadF; }
2086 
2087  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
2088 
2090 };
2091 
2092 /** @brief Compute the Lp norm of a function f.
2093  \f$ \| f \|_{Lp} = ( \int_\Omega | f |^p d\Omega)^{1/p} \f$ */
2094 double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh,
2095  const IntegrationRule *irs[]);
2096 
2097 /** @brief Compute the Lp norm of a vector function f = {f_i}_i=1...N.
2098  \f$ \| f \|_{Lp} = ( \sum_i \| f_i \|_{Lp}^p )^{1/p} \f$ */
2099 double ComputeLpNorm(double p, VectorCoefficient &coeff, Mesh &mesh,
2100  const IntegrationRule *irs[]);
2101 
2102 #ifdef MFEM_USE_MPI
2103 /** @brief Compute the global Lp norm of a function f.
2104  \f$ \| f \|_{Lp} = ( \int_\Omega | f |^p d\Omega)^{1/p} \f$ */
2105 double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh,
2106  const IntegrationRule *irs[]);
2107 
2108 /** @brief Compute the global Lp norm of a vector function f = {f_i}_i=1...N.
2109  \f$ \| f \|_{Lp} = ( \sum_i \| f_i \|_{Lp}^p )^{1/p} \f$ */
2110 double ComputeGlobalLpNorm(double p, VectorCoefficient &coeff, ParMesh &pmesh,
2111  const IntegrationRule *irs[]);
2112 #endif
2113 
2114 }
2115 
2116 #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.
Matrix coefficient defined as the product of two matrices.
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.
virtual void SetTime(double t)
Set the time for time dependent coefficients.
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.
A piecewise coefficient with the pieces keyed off the element attribute numbers.
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 UpdateCoefficients(const Array< int > &attr, const Array< VectorCoefficient * > &coefs)
Replace a set of coefficients.
void SetBCoef(VectorCoefficient &B_)
Reset the second vector coefficient.
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.
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.
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:521
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.
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.
void SetTime(double t)
Set the time for internally stored coefficients.
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.
void SetTime(double t)
Set the time for internally stored coefficients.
void UpdateCoefficients(const Array< int > &attr, const Array< Coefficient * > &coefs)
Replace a set of coefficients.
void SetTime(double t)
Set the time for internally stored coefficients.
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 ...
PWVectorCoefficient(int vd, const Array< int > &attr, const Array< VectorCoefficient * > &coefs)
Construct the coefficient using arrays describing the pieces.
Scalar coefficient defined as the product of two scalar coefficients or a scalar and a scalar coeffic...
void SetTime(double t)
Set the time for internally stored coefficients.
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:856
void SetTime(double t)
Set the time for internally stored coefficients.
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.
MatrixCoefficient * GetACoef() const
Return the first matrix coefficient.
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:199
void SetExponent(double p_)
Reset the exponent.
void SetAlphaCoef(Coefficient &A_)
Reset the factor in front of the first vector coefficient.
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...
void UpdateCoefficient(int attr, VectorCoefficient &coef)
Replace a single Coefficient for a particular attribute.
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.
void SetBCoef(MatrixCoefficient &B)
Reset the second matrix coefficient.
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.
MatrixCoefficient * GetBCoef() const
Return the second matrix 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.
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.
virtual void SetTime(double t)
Set the time for time dependent coefficients.
Definition: coefficient.hpp:48
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)
void SetTime(double t)
Set the time for internally stored coefficients.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
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.
void UpdateCoefficients(const Array< int > &attr, const Array< MatrixCoefficient * > &coefs)
Replace a set of coefficients.
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.
void ZeroCoefficient(int attr)
Remove a single VectorCoefficient for a particular attribute.
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.
virtual void SetTime(double t)
Set the time for time dependent coefficients.
Definition: coefficient.cpp:75
PWVectorCoefficient(int vd)
Constructs a piecewise vector coefficient of dimension vd.
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 .
void SetTime(double t)
Set the time for internally stored coefficients.
Coefficient * GetACoef() const
Return the first term in the linear combination.
void SetTime(double t)
Set the time for internally stored coefficients.
MatrixCoefficient(int dim, bool symm=false)
Construct a dim x dim matrix coefficient.
void SetTime(double t)
Set the time for internally stored coefficients.
void SetTime(double t)
Set the time for internally stored coefficients.
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.
void SetTime(double t)
Set the time for internally stored coefficients.
void UpdateCoefficient(int attr, Coefficient &coef)
Replace a single Coefficient for a particular attribute.
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
void SetTime(double t)
Set the time for internally stored coefficients.
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 SetTime(double t)
Set the time for internally stored coefficients.
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:154
void SetDeltaCenter(const Vector &center)
Coefficient * Weight()
See SetWeight() for description of the weight Coefficient.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the 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 ZeroCoefficient(int attr)
Remove a single MatrixCoefficient for a particular attribute.
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.
void UpdateCoefficient(int attr, MatrixCoefficient &coef)
Replace a single coefficient for a particular attribute.
virtual ~VectorArrayCoefficient()
Destroys vector coefficient.
void SetTime(double t)
Set the time for internally stored coefficients.
void SetACoef(VectorCoefficient &A_)
Reset the first vector coefficient.
void SetTime(double t)
Set the time for internally stored coefficients.
void SetBeta(double beta_)
Reset the factor in front of the second vector coefficient as a constant.
void SetTime(double t)
Set the time for internally stored coefficients.
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.
PWMatrixCoefficient(int dim, bool symm=false)
Constructs a piecewise matrix coefficient of dimension dim by dim.
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...
void SetTime(double t)
Set the time for internally stored coefficients.
const GridFunction * GetGridFunction() const
Get the vector grid function.
void SetComponent(int index_, int length_)
void SetBCoef(VectorCoefficient &B)
Reset the vector factor.
void operator=(double c)
Set the constants for all attributes to constant c.
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
Definition: fdual.hpp:543
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.
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.
virtual void SetTime(double t)
Set the time for time dependent coefficients.
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.
virtual void SetTime(double t)
Set the time for time dependent coefficients.
void SetBetaCoef(Coefficient &B_)
Reset the factor in front of the second vector coefficient.
Base class for Matrix Coefficients that optionally depend on time and space.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
Definition: coefficient.cpp:51
Coefficient * GetCoeff(int i)
Returns i&#39;th coefficient.
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...
MatrixProductCoefficient(MatrixCoefficient &A, MatrixCoefficient &B)
Construct with the two coefficients. Result is A * B.
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.
PowerCoefficient(Coefficient &A, double p_)
Construct with a coefficient and a constant power p_. Result is A^p.
PWMatrixCoefficient(int h, int w, bool symm=false)
Constructs a piecewise matrix coefficient of dimension h by w.
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.
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.
A piecewise matrix-valued coefficient with the pieces keyed off the element attribute numbers...
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.
Vector coefficient defined by a scalar DeltaCoefficient and a constant vector direction.
int GetVDim() const
For backward compatibility get the width of the matrix.
PWCoefficient()
Constructs a piecewise coefficient.
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
PWCoefficient(const Array< int > &attr, const Array< Coefficient * > &coefs)
Construct the coefficient using arrays describing the pieces.
VectorCoefficient * GetBCoef() const
Return the vector coefficient.
void SetTime(double t)
Set the time for internally stored coefficients.
virtual void SetTime(double t)
Set the time for time dependent coefficients.
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.
PWMatrixCoefficient(int h, int w, const Array< int > &attr, const Array< MatrixCoefficient * > &coefs, bool symm=false)
Construct the coefficient using arrays describing the pieces.
const Vector & GetB() const
Return the second vector constant.
Base class for symmetric matrix coefficients that optionally depend on time and space.
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.
void SetTime(double t)
Set the time for internally stored coefficients.
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
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 ZeroCoefficient(int attr)
Remove a single Coefficient for a particular attribute.
A piecewise vector-valued coefficient with the pieces keyed off the element attribute numbers...
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.
PWMatrixCoefficient(int dim, const Array< int > &attr, const Array< MatrixCoefficient * > &coefs, bool symm=false)
Construct the coefficient using arrays describing the pieces.
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)
void SetTime(double t)
Set the time for internally stored coefficients.
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))
virtual void SetTime(double t)
Set the time for time dependent coefficients.
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...
void SetTime(double t)
Set the time for internally stored coefficients.
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.
void SetTime(double t)
Set the time for internally stored coefficients.
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.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
Definition: coefficient.cpp:89
void SetBCoef(Coefficient &B)
Reset the denominator in the ratio.
void SetTime(double t)
Set the time for internally stored coefficients.
void SetACoef(MatrixCoefficient &A)
Reset the first matrix coefficient.
InverseMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
void SetTime(double t)
Set the time for internally stored coefficients.
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:757
void SetTime(double t)
Set the time for internally stored coefficients.
void SetTime(double t)
Set the time for internally stored coefficients.
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.
void SetTime(double t)
Set the time for internally stored coefficients.
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...
void SetTime(double t)
Set the time for internally stored coefficients.