12#ifndef MFEM_COEFFICIENT
13#define MFEM_COEFFICIENT
26class QuadratureSpaceBase;
27class QuadratureFunction;
170 std::map<int, Coefficient*> pieces;
195 { InitMap(attr, coefs); }
203 { InitMap(attr, coefs); }
207 { pieces[attr] = &coef; }
211 { pieces.erase(attr); }
393 { GridF = gf; Component = comp; }
430 : Q1(q), Transform1(std::move(F)) { Q2 = 0; Transform2 = 0; }
433 : Q1(q1), Q2(q2), Transform2(std::move(F)) { Transform1 = 0; }
538 {
mfem_error(
"DeltaCoefficient::Eval");
return 0.; }
555 { c = &c_; attr.
Copy(active_attr); }
672 std::map<int, VectorCoefficient*> pieces;
706 { InitMap(attr, coefs); }
713 { pieces.erase(attr); }
804 {
return Coeff[i] ? Coeff[i]->Eval(T, ip,
GetTime()) : 0.0; }
1003 {
mfem_error(
"VectorDeltaCoefficient::Eval"); }
1021 { c = &vc; attr.
Copy(active_attr); }
1099 {
mfem_error(
"MatrixCoefficient::EvalSymmetric"); }
1157 std::map<int, MatrixCoefficient*> pieces;
1215 { InitMap(attr, coefs); }
1222 { pieces.erase(attr); }
1236 std::function<void(
const Vector &,
Vector &)> SymmFunction;
1352 { c = &mc; attr.
Copy(active_attr); }
1380 : aConst(A),
a(NULL),
b(&B),
alpha(alpha_),
beta(beta_) { }
1385 : aConst(0.0),
a(&A),
b(&B),
alpha(alpha_),
beta(beta_) { }
1419 return alpha * ((
a == NULL ) ? aConst :
a->Eval(T, ip) )
1420 +
beta *
b->Eval(T, ip);
1556 : aConst(A),
a(NULL),
b(&B) { }
1560 : aConst(0.0),
a(&A),
b(&B) { }
1583 {
return ((
a == NULL ) ? aConst :
a->Eval(T, ip) ) *
b->Eval(T, ip); }
1600 : aConst(A), bConst(1.0),
a(NULL),
b(&B) { }
1604 : aConst(0.0), bConst(1.0),
a(&A),
b(&B) { }
1608 : aConst(0.0), bConst(B),
a(&A),
b(NULL) { }
1637 real_t den = (
b == NULL ) ? bConst :
b->Eval(T, ip);
1638 MFEM_ASSERT(den != 0.0,
"Division by zero in RatioCoefficient");
1639 return ((
a == NULL ) ? aConst :
a->Eval(T, ip) ) / den;
1672 {
return pow(
a->Eval(T, ip),
p); }
2303 return int(
a) & int(
b);
void Copy(Array ©) const
Create a copy of the internal array to the provided copy.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
CartesianCoefficient(int comp_)
comp_ index of the desired component (0 -> x, 1 -> y, 2 -> z)
Scalar coefficient which returns the x-component of the evaluation point.
Scalar coefficient which returns the y-component of the evaluation point.
Scalar coefficient which returns the z-component of the evaluation point.
Class to represent a coefficient evaluated at quadrature points.
void SetConstant(real_t constant)
Set this vector to the given constant.
int vdim
Number of values per quadrature point.
void Project(Coefficient &coeff)
Evaluate the given Coefficient at the quadrature points defined by qs.
int GetVDim() const
Return the number of values per quadrature point.
CoefficientVector(QuadratureSpaceBase &qs_, CoefficientStorage storage_=CoefficientStorage::FULL)
Create an empty CoefficientVector.
QuadratureFunction * qf
Internal QuadratureFunction (owned, may be NULL).
CoefficientStorage storage
Storage optimizations (see CoefficientStorage).
QuadratureSpaceBase & qs
Associated QuadratureSpaceBase.
void ProjectTranspose(MatrixCoefficient &coeff)
Project the transpose of coeff.
void MakeRef(const QuadratureFunction &qf_)
Make this vector a reference to the given QuadratureFunction.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
real_t GetTime()
Get the time for time dependent coefficients.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip, real_t t)
Evaluate the coefficient in the element described by T at the point ip at time t.
virtual void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
A coefficient that is constant across space and time.
void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf with the constant value.
ConstantCoefficient(real_t c=1.0)
c is value of constant function
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
Matrix coefficient defined as -a k x k x, for a vector k and scalar a.
void SetKCoef(VectorCoefficient &K)
Reset the vector factor.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
CrossCrossCoefficient(real_t A, VectorCoefficient &K)
void SetACoef(Coefficient &A)
Reset the scalar factor.
void SetAConst(real_t A)
Reset the scalar factor as a constant.
real_t GetAConst() const
Return the scalar factor.
Coefficient * GetACoef() const
Return the scalar factor.
VectorCoefficient * GetKCoef() const
Return the vector factor.
Vector coefficient defined as the Curl of a vector GridFunction.
void SetGridFunction(const GridFunction *gf)
Set the vector grid function.
const GridFunction * GetGridFunction() const
Get the vector grid function.
const GridFunction * GridFunc
CurlGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a vector grid function gf. The grid function is not owned by the coeff...
virtual ~CurlGridFunctionCoefficient()
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector curl coefficient at ip.
CylindricalAzimuthalCoefficient()
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
CylindricalRadialCoefficient()
Delta function coefficient optionally multiplied by a weight coefficient and a scaled time dependent ...
void SetTime(real_t t)
Set the time for internally stored coefficients.
DeltaCoefficient(real_t x, real_t s)
Construct a delta function scaled by s and centered at (x,0.0,0.0)
DeltaCoefficient(real_t x, real_t y, real_t s)
Construct a delta function scaled by s and centered at (x,y,0.0)
void GetDeltaCenter(Vector ¢er)
Write the center of the delta function into center.
Coefficient * Weight()
See SetWeight() for description of the weight Coefficient.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
A DeltaFunction cannot be evaluated. Calling this method will cause an MFEM error,...
void SetTol(real_t tol_)
Set the tolerance used during projection onto GridFunction to identify the Mesh vertex where the Cent...
real_t Scale()
Return the scale factor times the optional time dependent function. Returns with when not set by th...
DeltaCoefficient()
Construct a unit delta function centered at (0.0,0.0,0.0)
void SetDeltaCenter(const Vector ¢er)
Set the center location of the delta function.
DeltaCoefficient(real_t x, real_t y, real_t z, real_t s)
Construct a delta function scaled by s and centered at (x,y,z)
void SetFunction(real_t(*f)(real_t))
Set a time-dependent function that multiplies the Scale().
virtual ~DeltaCoefficient()
void SetWeight(Coefficient *w)
Set a weight Coefficient that multiplies the DeltaCoefficient.
virtual real_t EvalDelta(ElementTransformation &T, const IntegrationPoint &ip)
The value of the function assuming we are evaluating at the delta center.
void SetScale(real_t s_)
Set the scale value multiplying the delta function.
real_t Tol()
Return the tolerance used to identify the mesh vertices.
Data type dense matrix using column-major storage.
Scalar coefficient defined as the determinant of a matrix coefficient.
DeterminantCoefficient(MatrixCoefficient &A)
Construct with the matrix.
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the determinant coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
Scalar coefficient defined as the Divergence of a vector GridFunction.
void SetGridFunction(const GridFunction *gf)
Set the vector grid function.
virtual ~DivergenceGridFunctionCoefficient()
const GridFunction * GetGridFunction() const
Get the vector grid function.
DivergenceGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a vector grid function gf. The grid function is not owned by the coeff...
const GridFunction * GridFunc
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the scalar divergence coefficient at ip.
A general function coefficient.
FunctionCoefficient(std::function< real_t(const Vector &, real_t)> TDF)
Define a time-dependent coefficient from a std function.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
std::function< real_t(const Vector &)> Function
std::function< real_t(const Vector &, real_t)> TDFunction
MFEM_DEPRECATED FunctionCoefficient(real_t(*f)(Vector &))
(DEPRECATED) Define a time-independent coefficient from a C-function
FunctionCoefficient(std::function< real_t(const Vector &)> F)
Define a time-independent coefficient from a std function.
MFEM_DEPRECATED FunctionCoefficient(real_t(*tdf)(Vector &, real_t))
(DEPRECATED) Define a time-dependent coefficient from a C-function
Vector coefficient defined as the Gradient of a scalar GridFunction.
void SetGridFunction(const GridFunction *gf)
Set the scalar grid function.
const GridFunction * GetGridFunction() const
Get the scalar grid function.
GradientGridFunctionCoefficient(const GridFunction *gf)
Construct the coefficient with a scalar grid function gf. The grid function is not owned by the coeff...
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the gradient vector coefficient at ip.
virtual ~GradientGridFunctionCoefficient()
const GridFunction * GridFunc
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
GridFunctionCoefficient()
const GridFunction * GetGridFunction() const
Get the internal GridFunction.
void SetGridFunction(const GridFunction *gf)
Set the internal GridFunction.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
virtual void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
GridFunctionCoefficient(const GridFunction *gf, int comp=1)
Class for grid function - Vector with associated FE space.
Constant matrix coefficient defined as the identity of dimension d.
IdentityMatrixCoefficient(int d)
Construct with the dimension of the square identity matrix.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
Scalar coefficient defined as the inner product of two vector coefficients.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
VectorCoefficient * GetACoef() const
Return the first vector coefficient in the inner product.
VectorCoefficient * GetBCoef() const
Return the second vector coefficient in the inner product.
void SetBCoef(VectorCoefficient &B)
Reset the second vector in the inner product.
void SetACoef(VectorCoefficient &A)
Reset the first vector in the inner product.
InnerProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with the two vector coefficients. Result is .
void SetTime(real_t t)
Set the time for internally stored coefficients.
Class for integration point with weight.
Class for an integration rule - an Array of IntegrationPoint.
Matrix coefficient defined as the inverse a matrix coefficient.
InverseMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
Matrix coefficient defined by a matrix of scalar coefficients. Coefficients that are not set will eva...
Coefficient * GetCoeff(int i, int j)
Get the coefficient located at (i,j) in the matrix.
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...
MatrixArrayCoefficient(int dim)
Construct a coefficient matrix of dimensions dim * dim. The actual coefficients still need to be adde...
virtual ~MatrixArrayCoefficient()
real_t Eval(int i, int j, ElementTransformation &T, const IntegrationPoint &ip)
void SetTime(real_t t)
Set the time for internally stored coefficients.
Base class for Matrix Coefficients that optionally depend on time and space.
virtual void Project(QuadratureFunction &qf, bool transpose=false)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points....
virtual ~MatrixCoefficient()
MatrixCoefficient(int dim, bool symm=false)
Construct a dim x dim matrix coefficient.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
real_t GetTime()
Get the time for time dependent coefficients.
int GetVDim() const
For backward compatibility get the width of the matrix.
MatrixCoefficient(int h, int w, bool symm=false)
Construct a h x w matrix coefficient.
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 ...
int GetWidth() const
Get the width of the matrix.
virtual void EvalSymmetric(Vector &K, ElementTransformation &T, const IntegrationPoint &ip)
(DEPRECATED) Evaluate a symmetric matrix coefficient.
int GetHeight() const
Get the height of the matrix.
A matrix coefficient that is constant in space and time.
MatrixConstantCoefficient(const DenseMatrix &m)
Construct using matrix m for the constant.
const DenseMatrix & GetMatrix()
Return a reference to the constant matrix.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
A matrix coefficient with an optional scalar coefficient multiplier q. The matrix function can either...
MatrixFunctionCoefficient(int dim, std::function< void(const Vector &, DenseMatrix &)> F, Coefficient *q=nullptr)
Define a time-independent square matrix coefficient from a std function.
virtual void EvalSymmetric(Vector &K, ElementTransformation &T, const IntegrationPoint &ip)
(DEPRECATED) Evaluate the symmetric matrix coefficient at ip.
MatrixFunctionCoefficient(int dim, std::function< void(const Vector &, real_t, DenseMatrix &)> TDF, Coefficient *q=nullptr)
Define a time-dependent square matrix coefficient from a std function.
MatrixFunctionCoefficient(const DenseMatrix &m, Coefficient &q)
Define a constant matrix coefficient times a scalar Coefficient.
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
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.
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual ~MatrixFunctionCoefficient()
Matrix coefficient defined as the product of two matrices.
MatrixProductCoefficient(MatrixCoefficient &A, MatrixCoefficient &B)
Construct with the two coefficients. Result is A * B.
MatrixCoefficient * GetACoef() const
Return the first matrix coefficient.
MatrixCoefficient * GetBCoef() const
Return the second matrix coefficient.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetBCoef(MatrixCoefficient &B)
Reset the second matrix coefficient.
void SetACoef(MatrixCoefficient &A)
Reset the first matrix coefficient.
Derived matrix coefficient that has the value of the parent matrix coefficient where it is active and...
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
MatrixRestrictedCoefficient(MatrixCoefficient &mc, Array< int > &attr)
Construct with a parent matrix coefficient and an array of zeros and ones representing the attributes...
void SetTime(real_t t)
Set the time for internally stored coefficients.
Matrix coefficient defined as the linear combination of two matrices.
MatrixCoefficient * GetBCoef() const
Return the second matrix coefficient.
void SetBeta(real_t beta_)
Reset the factor in front of the second matrix coefficient.
real_t GetBeta() const
Return the factor in front of the second matrix coefficient.
void SetAlpha(real_t alpha_)
Reset the factor in front of the first matrix coefficient.
MatrixSumCoefficient(MatrixCoefficient &A, MatrixCoefficient &B, real_t alpha_=1.0, real_t beta_=1.0)
Construct with the two coefficients. Result is alpha_ * A + beta_ * B.
void SetACoef(MatrixCoefficient &A)
Reset the first matrix coefficient.
void SetTime(real_t t)
Set the time for internally stored coefficients.
MatrixCoefficient * GetACoef() const
Return the first matrix coefficient.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
real_t GetAlpha() const
Return the factor in front of the first matrix coefficient.
void SetBCoef(MatrixCoefficient &B)
Reset the second matrix coefficient.
Vector coefficient defined as a product of a matrix coefficient and a vector coefficient.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
MatrixVectorProductCoefficient(MatrixCoefficient &A, VectorCoefficient &B)
Constructor with two coefficients. Result is A*B.
void SetTime(real_t t)
Set the time for internally stored coefficients.
void SetBCoef(VectorCoefficient &B)
Reset the vector coefficient.
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
VectorCoefficient * GetBCoef() const
Return the vector coefficient.
Vector coefficient defined as a normalized vector field (returns v/|v|)
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
VectorCoefficient * GetACoef() const
Return the vector coefficient.
void SetACoef(VectorCoefficient &A)
Reset the vector coefficient.
NormalizedVectorCoefficient(VectorCoefficient &A, real_t tol=1e-6)
Return a vector normalized to a length of one.
Matrix coefficient defined as the outer product of two vector coefficients.
VectorCoefficient * GetACoef() const
Return the first vector coefficient in the outer product.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetACoef(VectorCoefficient &A)
Reset the first vector in the outer product.
void SetTime(real_t t)
Set the time for internally stored coefficients.
void SetBCoef(VectorCoefficient &B)
Reset the second vector in the outer product.
OuterProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with two vector coefficients. Result is .
VectorCoefficient * GetBCoef() const
Return the second vector coefficient in the outer product.
A piecewise coefficient with the pieces keyed off the element attribute numbers.
PWCoefficient(const Array< int > &attr, const Array< Coefficient * > &coefs)
Construct the coefficient using arrays describing the pieces.
PWCoefficient()
Constructs a piecewise coefficient.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
void ZeroCoefficient(int attr)
Remove a single Coefficient for a particular attribute.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
void UpdateCoefficient(int attr, Coefficient &coef)
Replace a single Coefficient for a particular attribute.
void UpdateCoefficients(const Array< int > &attr, const Array< Coefficient * > &coefs)
Replace a set of coefficients.
A piecewise constant coefficient with the constants keyed off the element attribute numbers.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
int GetNConst()
Returns the number of constants representing different attributes.
real_t & operator()(int i)
Return a reference to the i-th constant.
void UpdateConstants(Vector &c)
Update the constants with vector c.
PWConstCoefficient(int NumOfSubD=0)
Constructs a piecewise constant coefficient in NumOfSubD subdomains.
PWConstCoefficient(Vector &c)
Construct the constant coefficient using a vector of constants.
void operator=(real_t c)
Set the constants for all attributes to constant c.
A piecewise matrix-valued coefficient with the pieces keyed off the element attribute numbers.
void UpdateCoefficient(int attr, MatrixCoefficient &coef)
Replace a single coefficient for a particular attribute.
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
PWMatrixCoefficient(int dim, bool symm=false)
Constructs a piecewise matrix coefficient of dimension dim by dim.
PWMatrixCoefficient(int dim, const Array< int > &attr, const Array< MatrixCoefficient * > &coefs, bool symm=false)
Construct the coefficient using arrays describing the pieces.
PWMatrixCoefficient(int h, int w, const Array< int > &attr, const Array< MatrixCoefficient * > &coefs, bool symm=false)
Construct the coefficient using arrays describing the pieces.
PWMatrixCoefficient(int h, int w, bool symm=false)
Constructs a piecewise matrix coefficient of dimension h by w.
void ZeroCoefficient(int attr)
Remove a single MatrixCoefficient for a particular attribute.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
void UpdateCoefficients(const Array< int > &attr, const Array< MatrixCoefficient * > &coefs)
Replace a set of coefficients.
A piecewise vector-valued coefficient with the pieces keyed off the element attribute numbers.
void UpdateCoefficients(const Array< int > &attr, const Array< VectorCoefficient * > &coefs)
Replace a set of coefficients.
PWVectorCoefficient(int vd, const Array< int > &attr, const Array< VectorCoefficient * > &coefs)
Construct the coefficient using arrays describing the pieces.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
void ZeroCoefficient(int attr)
Remove a single VectorCoefficient for a particular attribute.
void UpdateCoefficient(int attr, VectorCoefficient &coef)
Replace a single Coefficient for a particular attribute.
PWVectorCoefficient(int vd)
Constructs a piecewise vector coefficient of dimension vd.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
Class for parallel meshes.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
virtual ~PositionVectorCoefficient()
PositionVectorCoefficient(int dim)
Scalar coefficient defined as a scalar raised to a power.
real_t GetExponent() const
Return the exponent.
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
void SetExponent(real_t p_)
Reset the exponent.
Coefficient * GetACoef() const
Return the base coefficient.
PowerCoefficient(Coefficient &A, real_t p_)
Construct with a coefficient and a constant power p_. Result is A^p.
void SetACoef(Coefficient &A)
Reset the base coefficient.
Scalar coefficient defined as the product of two scalar coefficients or a scalar and a scalar coeffic...
void SetBCoef(Coefficient &B)
Reset the second term in the product.
Coefficient * GetACoef() const
Return the first term in the product.
ProductCoefficient(real_t A, Coefficient &B)
Constructor with one coefficient. Result is A * B.
ProductCoefficient(Coefficient &A, Coefficient &B)
Constructor with two coefficients. Result is A * B.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
real_t GetAConst() const
Return the first term in the product.
Coefficient * GetBCoef() const
Return the second term in the product.
void SetACoef(Coefficient &A)
Reset the first term in the product.
void SetAConst(real_t A)
Reset the first term in the product as a constant.
void SetTime(real_t t)
Set the time for internally stored coefficients.
Quadrature function coefficient which requires that the quadrature rules used for this coefficient be...
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
QuadratureFunctionCoefficient(const QuadratureFunction &qf)
Constructor with a quadrature function as input.
const QuadratureFunction & GetQuadFunction() const
virtual ~QuadratureFunctionCoefficient()
virtual void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
Represents values or vectors of values at quadrature points on a mesh.
Abstract base class for QuadratureSpace and FaceQuadratureSpace.
Scalar coefficient defined as the ratio of two scalars where one or both scalars are scalar coefficie...
Coefficient * GetBCoef() const
Return the denominator of the ratio.
void SetBConst(real_t B)
Reset the denominator in the ratio as a constant.
void SetAConst(real_t A)
Reset the numerator in the ratio as a constant.
void SetBCoef(Coefficient &B)
Reset the denominator in the ratio.
RatioCoefficient(real_t A, Coefficient &B)
real_t GetBConst() const
Return the denominator of the ratio.
Coefficient * GetACoef() const
Return the numerator of the ratio.
RatioCoefficient(Coefficient &A, Coefficient &B)
void SetTime(real_t t)
Set the time for internally stored coefficients.
RatioCoefficient(Coefficient &A, real_t B)
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
real_t GetAConst() const
Return the numerator of the ratio.
void SetACoef(Coefficient &A)
Reset the numerator in the ratio.
Derived coefficient that takes the value of the parent coefficient for the active attributes and is z...
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
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 a product of a scalar coefficient and a matrix coefficient.
real_t GetAConst() const
Return the scalar factor.
void SetACoef(Coefficient &A)
Reset the scalar factor.
void SetBCoef(MatrixCoefficient &B)
Reset the matrix factor.
MatrixCoefficient * GetBCoef() const
Return the matrix factor.
ScalarMatrixProductCoefficient(real_t A, MatrixCoefficient &B)
Constructor with one coefficient. Result is A*B.
Coefficient * GetACoef() const
Return the scalar factor.
void SetAConst(real_t A)
Reset the scalar factor as a constant.
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
Vector coefficient defined as a product of scalar and vector coefficients.
real_t GetAConst() const
Return the scalar factor.
Coefficient * GetACoef() const
Return the scalar factor.
void SetBCoef(VectorCoefficient &B)
Reset the vector factor.
VectorCoefficient * GetBCoef() const
Return the vector factor.
ScalarVectorProductCoefficient(real_t A, VectorCoefficient &B)
Constructor with constant and vector coefficient. Result is A * B.
void SetACoef(Coefficient &A)
Reset the scalar factor.
void SetAConst(real_t A)
Reset the scalar factor as a constant.
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
SphericalAzimuthalCoefficient()
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
SphericalPolarCoefficient()
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
SphericalRadialCoefficient()
Scalar coefficient defined as the linear combination of two scalar coefficients or a scalar and a sca...
void SetBCoef(Coefficient &B)
Reset the second term in the linear combination.
Coefficient * GetACoef() const
Return the first term in the linear combination.
SumCoefficient(real_t A, Coefficient &B, real_t alpha_=1.0, real_t beta_=1.0)
Constructor with one coefficient. Result is alpha_ * A + beta_ * B.
void SetAlpha(real_t alpha_)
Reset the factor in front of the first term in the linear combination.
SumCoefficient(Coefficient &A, Coefficient &B, real_t alpha_=1.0, real_t beta_=1.0)
Constructor with two coefficients. Result is alpha_ * A + beta_ * B.
void SetACoef(Coefficient &A)
Reset the first term in the linear combination.
void SetBeta(real_t beta_)
Reset the factor in front of the second term in the linear combination.
real_t GetBeta() const
Return the factor in front of the second term in the linear combination.
real_t GetAlpha() const
Return the factor in front of the first term in the linear combination.
void SetTime(real_t t)
Set the time for internally stored coefficients.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
Coefficient * GetBCoef() const
Return the second term in the linear combination.
void SetAConst(real_t A)
Reset the first term in the linear combination as a constant.
real_t GetAConst() const
Return the first term in the linear combination.
Base class for symmetric matrix coefficients that optionally depend on time and space.
const DenseSymmetricMatrix & GetMatrix()
Return a reference to the constant matrix.
DenseSymmetricMatrix mat
Internal matrix used when evaluating this coefficient as a DenseMatrix.
virtual void ProjectSymmetric(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
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 as ...
SymmetricMatrixCoefficient(int dimension)
Construct a dim x dim matrix coefficient.
virtual ~SymmetricMatrixCoefficient()
int GetSize() const
Get the size of the matrix.
A matrix coefficient that is constant in space and time.
virtual void Eval(DenseSymmetricMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
SymmetricMatrixConstantCoefficient(const DenseSymmetricMatrix &m)
Construct using matrix m for the constant.
A matrix coefficient with an optional scalar coefficient multiplier q. The matrix function can either...
virtual void Eval(DenseSymmetricMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
virtual ~SymmetricMatrixFunctionCoefficient()
SymmetricMatrixFunctionCoefficient(const DenseSymmetricMatrix &m, Coefficient &q)
Define a constant matrix coefficient times a scalar Coefficient.
void SetTime(real_t t)
Set the time for internally stored coefficients.
SymmetricMatrixFunctionCoefficient(int dim, std::function< void(const Vector &, DenseSymmetricMatrix &)> F, Coefficient *q=nullptr)
Define a time-independent symmetric matrix coefficient from a std function.
SymmetricMatrixFunctionCoefficient(int dim, std::function< void(const Vector &, real_t, DenseSymmetricMatrix &)> TDF, Coefficient *q=nullptr)
Define a time-dependent square matrix coefficient from a std function.
Matrix coefficient defined as the transpose a matrix coefficient.
MatrixCoefficient * GetACoef() const
Return the matrix coefficient.
TransposeMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
void SetACoef(MatrixCoefficient &A)
Reset the matrix coefficient.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient at ip.
void SetTime(real_t t)
Set the time for internally stored coefficients.
Vector coefficient defined by an array of scalar coefficients. Coefficients that are not set will eva...
virtual ~VectorArrayCoefficient()
Destroys vector coefficient.
Coefficient ** GetCoeffs()
Returns the entire array of coefficients.
void SetTime(real_t t)
Set the time for internally stored coefficients.
VectorArrayCoefficient(int dim)
Construct vector of dim coefficients. The actual coefficients still need to be added with Set().
void Set(int i, Coefficient *c, bool own=true)
Sets coefficient in the vector.
real_t Eval(int i, ElementTransformation &T, const IntegrationPoint &ip)
Coefficient * GetCoeff(int i)
Returns i'th coefficient.
Base class for vector Coefficients that optionally depend on time and space.
virtual ~VectorCoefficient()
int GetVDim()
Returns dimension of the vector.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
VectorCoefficient(int vd)
Initialize the VectorCoefficient with vector dimension vd.
real_t GetTime()
Get the time for time dependent coefficients.
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 ...
virtual void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
Vector coefficient that is constant in space and time.
const Vector & GetVec() const
Return a reference to the constant vector in this class.
VectorConstantCoefficient(const Vector &v)
Construct the coefficient with constant vector v.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
Vector coefficient defined as a cross product of two vectors.
void SetBCoef(VectorCoefficient &B)
Reset the second term in the product.
VectorCoefficient * GetACoef() const
Return the first term in the product.
VectorCoefficient * GetBCoef() const
Return the second term in the product.
void SetTime(real_t t)
Set the time for internally stored coefficients.
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.
void SetACoef(VectorCoefficient &A)
Reset the first term in the product.
Vector coefficient defined by a scalar DeltaCoefficient and a constant vector direction.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
A VectorDeltaFunction cannot be evaluated. Calling this method will cause an MFEM error,...
void GetDeltaCenter(Vector ¢er)
VectorDeltaCoefficient(const Vector &dir_, real_t x, real_t y, real_t z, real_t s)
Construct with a Vector object representing the direction and a delta function scaled by s and center...
void SetTime(real_t t)
Set the time for internally stored coefficients.
VectorDeltaCoefficient(const Vector &dir_)
Construct with a Vector object representing the direction and a unit delta function centered at (0....
VectorDeltaCoefficient(const Vector &dir_, real_t x, real_t s)
Construct with a Vector object representing the direction and a delta function scaled by s and center...
virtual ~VectorDeltaCoefficient()
void SetDeltaCoefficient(const DeltaCoefficient &d_)
Replace the associated DeltaCoefficient with a new DeltaCoefficient.
void SetDirection(const Vector &d_)
virtual void EvalDelta(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Return the specified direction vector multiplied by the value returned by DeltaCoefficient::EvalDelta...
VectorDeltaCoefficient(const Vector &dir_, real_t x, real_t y, real_t s)
Construct with a Vector object representing the direction and a delta function scaled by s and center...
DeltaCoefficient & GetDeltaCoefficient()
Return the associated scalar DeltaCoefficient.
VectorDeltaCoefficient(int vdim_)
Construct with a vector of dimension vdim_.
void SetDeltaCenter(const Vector ¢er)
A general vector function coefficient.
VectorFunctionCoefficient(int dim, std::function< void(const Vector &, real_t, Vector &)> TDF, Coefficient *q=nullptr)
Define a time-dependent vector coefficient from a std function.
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(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
virtual ~VectorFunctionCoefficient()
Vector coefficient defined by a vector GridFunction.
void SetGridFunction(const GridFunction *gf)
Set the grid function for this coefficient. Also sets the Vector dimension to match that of the gf.
const GridFunction * GridFunc
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
virtual ~VectorGridFunctionCoefficient()
VectorGridFunctionCoefficient()
Construct an empty coefficient. Calling Eval() before the grid function is set will cause a segfault.
virtual void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
const GridFunction * GetGridFunction() const
Returns a pointer to the grid function in this Coefficient.
Vector quadrature function coefficient which requires that the quadrature rules used for this vector ...
virtual ~VectorQuadratureFunctionCoefficient()
VectorQuadratureFunctionCoefficient(const QuadratureFunction &qf)
Constructor with a quadrature function as input.
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 void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
const QuadratureFunction & GetQuadFunction() const
void SetComponent(int index_, int length_)
Derived vector coefficient that has the value of the parent vector where it is active and is zero oth...
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient at ip.
VectorRestrictedCoefficient(VectorCoefficient &vc, Array< int > &attr)
Construct with a parent vector coefficient and an array of zeros and ones representing the attributes...
void SetTime(real_t t)
Set the time for internally stored coefficients.
Scalar coefficient defined as a cross product of two vectors in the xy-plane.
VectorRotProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Constructor with two vector coefficients. Result is .
void SetBCoef(VectorCoefficient &B)
Reset the second vector in the product.
void SetACoef(VectorCoefficient &A)
Reset the first vector in the product.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
VectorCoefficient * GetBCoef() const
Return the second vector of the product.
void SetTime(real_t t)
Set the time for internally stored coefficients.
VectorCoefficient * GetACoef() const
Return the first vector of the product.
Vector coefficient defined as the linear combination of two vectors.
void SetB(const Vector &B_)
Reset the second vector as a constant.
real_t GetAlpha() const
Return the factor in front of the first vector coefficient.
void SetTime(real_t t)
Set the time for internally stored coefficients.
const Vector & GetB() const
Return the second vector constant.
void SetAlpha(real_t alpha_)
Reset the factor in front of the first vector coefficient as a constant.
const Vector & GetA() const
Return the first vector constant.
real_t GetBeta() const
Return the factor in front of the second vector coefficient.
void SetA(const Vector &A_)
Reset the first vector as a constant.
void SetBeta(real_t beta_)
Reset the factor in front of the second vector coefficient as a constant.
Coefficient * GetAlphaCoef() const
Return the factor in front of the first vector coefficient.
void SetBCoef(VectorCoefficient &B_)
Reset the second vector coefficient.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
VectorSumCoefficient(int dim)
void SetAlphaCoef(Coefficient &A_)
Reset the factor in front of the first vector coefficient.
Coefficient * GetBetaCoef() const
Return the factor in front of the second vector coefficient.
VectorCoefficient * GetBCoef() const
Return the second vector coefficient.
void SetBetaCoef(Coefficient &B_)
Reset the factor in front of the second vector coefficient.
VectorCoefficient * GetACoef() const
Return the first vector coefficient.
void SetACoef(VectorCoefficient &A_)
Reset the first vector coefficient.
int Size() const
Returns the size of the vector.
void SetSize(int s)
Resize the vector to size s.
constexpr int dimension
This example only works in 3D. Kernels for 2D are not implemented.
int operator&(CoefficientStorage a, CoefficientStorage b)
VectorCoefficient DiagonalMatrixCoefficient
void mfem_error(const char *msg)
real_t ComputeGlobalLpNorm(real_t p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
CartesianZCoefficient CylindricalZCoefficient
real_t ComputeLpNorm(real_t p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
Compute the Lp norm of a function f. .
CoefficientStorage
Flags that determine what storage optimizations to use in CoefficientVector.
@ SYMMETRIC
Store the triangular part of symmetric matrices.
@ CONSTANTS
Store constants using only vdim entries.
@ COMPRESSED
Enable all above compressions.
@ FULL
Store the coefficient as a full QuadratureFunction.
CoefficientStorage operator|(CoefficientStorage a, CoefficientStorage b)
std::function< real_t(const Vector &)> f(real_t mass_coeff)
MatrixVectorProductCoefficient MatVecCoefficient
Convenient alias for the MatrixVectorProductCoefficient.
real_t p(const Vector &x, real_t t)