12 #ifndef MFEM_COEFFICIENT
13 #define MFEM_COEFFICIENT
15 #include "../config/config.hpp"
16 #include "../linalg/linalg.hpp"
175 { GridF = gf; Component = comp; }
189 double (*Transform1)(double);
190 double (*Transform2)(double,double);
194 : Q1(q), Transform1(F) { Q2 = 0; Transform2 = 0; }
196 double (*F)(
double,
double))
197 : Q1(q1), Q2(q2), Transform2(F) { Transform1 = 0; }
261 {
mfem_error(
"DeltaCoefficient::Eval");
return 0.; }
274 { c = &_c; attr.
Copy(active_attr); }
388 {
return Coeff[i] ? Coeff[i]->Eval(T, ip,
GetTime()) : 0.0; }
521 {
mfem_error(
"VectorDeltaCoefficient::Eval"); }
535 { c = &vc; attr.
Copy(active_attr); }
665 { c = &mc; attr.
Copy(active_attr); }
686 double _alpha = 1.0,
double _beta = 1.0)
687 : a(&A), b(&B), alpha(_alpha), beta(_beta) { }
692 {
return alpha * a->
Eval(T, ip) + beta * b->
Eval(T, ip); }
709 {
return a->
Eval(T, ip) * b->
Eval(T, ip); }
728 {
return pow(a->
Eval(T, ip), p); }
796 double _alpha = 1.0,
double _beta = 1.0);
884 double _alpha = 1.0,
double _beta = 1.0);
Vector coefficient defined as a cross product of two vectors.
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 ...
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. ...
void GetDeltaCenter(Vector ¢er)
void SetTol(double _tol)
Set the tolerance used during projection onto GridFunction to identifying the Mesh vertex where the C...
Class for an integration rule - an Array of IntegrationPoint.
Vector coefficient defined by an array of scalar coefficients.
Class for grid function - Vector with associated FE space.
GridFunction * GetGridFunction() const
void Set(int i, int j, Coefficient *c, bool own=true)
IdentityMatrixCoefficient(int d)
PWConstCoefficient(Vector &c)
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
void SetWeight(Coefficient *w)
Set a weight Coefficient that multiplies the DeltaCoefficient.
Subclass constant coefficient.
Coefficient ** GetCoeffs()
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 ...
TransposeMatrixCoefficient(MatrixCoefficient &A)
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 ~MatrixCoefficient()
VectorCrossProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
FunctionCoefficient(double(*tdf)(const Vector &, double))
Define a time-dependent coefficient from a C-function.
Scalar coefficient defined as the inner product of two vector coefficients.
virtual ~VectorFunctionCoefficient()
void SetSize(int s)
Resize the vector to size s.
double Eval(int i, ElementTransformation &T, const IntegrationPoint &ip)
Evaluates i'th component of the vector.
Vector coefficient defined as a matrix vector product.
void SetDeltaCenter(const Vector ¢er)
double Scale()
Return the scale set by SetScale() multiplied by the time-dependent function specified by SetFunction...
FunctionCoefficient(double(*f)(Vector &))
(DEPRECATED) Define a time-independent coefficient from a C-function
Scalar coefficient defined as a cross product of two vectors in 2D.
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
Scalar coefficient defined as the product of two scalar coefficients.
SumCoefficient(Coefficient &A, Coefficient &B, double _alpha=1.0, double _beta=1.0)
void Copy(Array ©) const
Create a copy of the current array.
GradientGridFunctionCoefficient(GridFunction *gf)
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
double & operator()(int i)
Member function to access or modify the value of the i-th constant.
Data type dense matrix using column-major storage.
Delta function coefficient.
int Size() const
Returns the size of the vector.
virtual ~CurlGridFunctionCoefficient()
void SetFunction(double(*f)(double))
Set a time-dependent function that multiplies the Scale().
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 coefficient.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
double Tol()
See SetTol() for description of the tolerance parameter.
GridFunction * GetGridFunction() const
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 ...
Vector coefficient defined as the sum of two vector coefficients.
virtual ~MatrixArrayCoefficient()
Matrix coefficient defined as a product of a scalar and a matrix.
VectorFunctionCoefficient(int dim, void(*TDF)(const Vector &, double, Vector &), Coefficient *q=NULL)
Construct a time-dependent vector coefficient from a C-function.
virtual ~VectorGridFunctionCoefficient()
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
A DeltaFunction cannot be evaluated. Calling this method will cause an MFEM error, terminating the application.
VectorDeltaCoefficient(const Vector &_dir, double x, double y, double s)
DeltaCoefficient(double x, double y, double z, double s)
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate coefficient.
Matrix coefficient defined as the sum of two matrix coefficients.
DeltaCoefficient & GetDeltaCoefficient()
Return the associated scalar DeltaCoefficient.
DeltaCoefficient(double x, double s)
DeterminantCoefficient(MatrixCoefficient &A)
double(* Function)(const Vector &)
Matrix coefficient defined as the inverse a matrix.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
A VectorDeltaFunction cannot be evaluated. Calling this method will cause an MFEM error...
void SetGridFunction(GridFunction *gf)
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
Coefficient * GetCoeff(int i, int j)
Matrix coefficient defined as the outer product of two vectors.
double Eval(int i, int j, ElementTransformation &T, const IntegrationPoint &ip)
VectorCoefficient defined on a subset of domain or boundary attributes.
InnerProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
MatrixFunctionCoefficient(int dim, void(*F)(const Vector &, DenseMatrix &), Coefficient *q=NULL)
Construct a time-independent square matrix coefficient from a C-function.
double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
VectorConstantCoefficient(const Vector &v)
ConstantCoefficient(double c=1.0)
c is value of constant function
VectorDeltaCoefficient(int _vdim)
Coefficient defined on a subset of domain or boundary attributes.
virtual ~VectorCoefficient()
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
void SetDeltaCenter(const Vector ¢er)
Coefficient * Weight()
See SetWeight() for description of the weight Coefficient.
ScalarVectorProductCoefficient(Coefficient &A, VectorCoefficient &B)
VectorDeltaCoefficient(const Vector &_dir)
MatrixRestrictedCoefficient(MatrixCoefficient &mc, Array< int > &attr)
void SetDirection(const Vector &_d)
VectorDeltaCoefficient(const Vector &_dir, double x, double y, double z, double s)
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)
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
MatrixFunctionCoefficient(int dim, void(*TDF)(const Vector &, double, DenseMatrix &), Coefficient *q=NULL)
Construct a time-dependent square matrix coefficient from a C-function.
void Set(int i, Coefficient *c, bool own=true)
Sets coefficient in the vector.
void SetDeltaCoefficient(const DeltaCoefficient &_d)
Replace the associated DeltaCoeficient with a new DeltaCoeficient.
virtual ~DeltaCoefficient()
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 ~VectorArrayCoefficient()
Destroys vector coefficient.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
int GetVDim()
Returns dimension of the vector.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
GridFunction * GetGridFunction() const
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
VectorSumCoefficient(VectorCoefficient &A, VectorCoefficient &B, double _alpha=1.0, double _beta=1.0)
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
void operator=(double c)
Set domain constants equal to the same constant c.
PowerCoefficient(Coefficient &A, double _p)
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
void SetGridFunction(GridFunction *gf)
void GetDeltaCenter(Vector ¢er)
MatrixFunctionCoefficient(const DenseMatrix &m, Coefficient &q)
Construct a constant matrix coefficient times a scalar Coefficient.
Scalar coefficient defined as the Divergence of a vector GridFunction.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
MatVecCoefficient(MatrixCoefficient &A, VectorCoefficient &B)
Base class Coefficient that may optionally depend on time.
Vector coefficient defined as a product of a scalar and a vector.
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
FunctionCoefficient(double(*f)(const Vector &))
Define a time-independent coefficient from a C-function.
void UpdateConstants(Vector &c)
Update constants.
GridFunction * GetGridFunction() const
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient function.
Coefficient * GetCoeff(int i)
Returns i'th coefficient.
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 EvalDelta(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Return the specified direction vector multiplied by the value returned by DeltaCoefficient::EvalDelta...
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
void SetGridFunction(GridFunction *gf)
PWConstCoefficient(int NumOfSubD=0)
Constructs a piecewise constant coefficient in NumOfSubD subdomains.
virtual double EvalDelta(ElementTransformation &T, const IntegrationPoint &ip)
Return the Scale() multiplied by the weight Coefficient, if any.
virtual ~DivergenceGridFunctionCoefficient()
MatrixCoefficient defined on a subset of domain or boundary attributes.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
void SetGridFunction(GridFunction *gf)
VectorDeltaCoefficient: DeltaCoefficient with a direction.
MatrixArrayCoefficient(int dim)
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 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 ...
class for piecewise constant coefficient
Class for integration point with weight.
virtual ~GradientGridFunctionCoefficient()
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
void SetGridFunction(GridFunction *gf)
MatrixSumCoefficient(MatrixCoefficient &A, MatrixCoefficient &B, double _alpha=1.0, double _beta=1.0)
FunctionCoefficient(double(*tdf)(Vector &, double))
(DEPRECATED) Define a time-dependent coefficient from a C-function
VectorCoefficient(int vd)
int GetNConst()
Returns the number of constants.
Matrix coefficient defined as the transpose a matrix.
GridFunctionCoefficient()
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
VectorDeltaCoefficient(const Vector &_dir, double x, double s)
DeltaCoefficient(double x, double y, double s)
Scalar coefficient defined as a scalar raised to a power.
DivergenceGridFunctionCoefficient(GridFunction *gf)
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.
VectorFunctionCoefficient(int dim, void(*F)(const Vector &, Vector &), Coefficient *q=NULL)
Construct a time-independent vector coefficient from a C-function.
Scalar coefficient defined as the determinant of a matrix coefficient.
class for C-function coefficient
ProductCoefficient(Coefficient &A, Coefficient &B)
Vector coefficient defined by a vector GridFunction.
MatrixCoefficient(int dim)
VectorRestrictedCoefficient(VectorCoefficient &vc, Array< int > &attr)
VectorRotProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient.
MatrixConstantCoefficient(const DenseMatrix &m)
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
GridFunctionCoefficient(GridFunction *gf, int comp=1)
CurlGridFunctionCoefficient(GridFunction *gf)
Matrix coefficient defined as the identity of dimension d.
GridFunction * GetGridFunction() const
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
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 ...
double(* TDFunction)(const Vector &, double)
virtual ~VectorDeltaCoefficient()
Class for parallel meshes.
InverseMatrixCoefficient(MatrixCoefficient &A)
Coefficients based on sums and products of other coefficients.
MatrixCoefficient(int h, int w)
virtual ~MatrixFunctionCoefficient()
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
ScalarMatrixProductCoefficient(Coefficient &A, MatrixCoefficient &B)
VectorGridFunctionCoefficient()
RestrictedCoefficient(Coefficient &_c, Array< int > &attr)