37 int coarse_element = cf.
embeddings[fine_element].parent;
55 const int ne = qspace.
GetNE();
57 for (
int iel = 0; iel < ne; ++iel)
62 for (
int iq = 0; iq < ir.
Size(); ++iq)
67 values[iq_p] =
Eval(T, ip);
81 return (constants(att-1));
88 const bool compressed =
91 const int ne = qs.GetNE();
93 const int *attributes = [&]()
97 return qs.GetMesh()->GetElementAttributes().Read();
102 "Interior faces do not have attributes.");
103 return qs.GetMesh()->GetBdrFaceAttributes().Read();
107 MFEM_ABORT(
"Unsupported case.");
116 const int a = attributes[e];
117 const real_t elementConstant = d_c[
a - 1];
118 const int begin = compressed ? e*offsets[0] : offsets[e];
119 const int end = compressed ? (e+1)*offsets[0] : offsets[e+1];
120 for (
int i = begin; i < end; ++i)
122 d_qf[i] = elementConstant;
127void PWCoefficient::InitMap(
const Array<int> & attr,
130 MFEM_VERIFY(attr.
Size() == coefs.
Size(),
132 "Attribute and coefficient arrays have incompatible "
135 for (
int i=0; i<attr.
Size(); i++)
137 if (coefs[i] != NULL)
148 std::map<int, Coefficient*>::iterator
p = pieces.begin();
149 for (;
p != pieces.end();
p++)
151 if (
p->second != NULL)
153 p->second->SetTime(t);
162 std::map<int, Coefficient*>::const_iterator
p = pieces.find(att);
163 if (
p != pieces.end())
165 if (
p->second != NULL)
167 return p->second->Eval(T, ip);
202 return sqrt(transip[0] * transip[0] + transip[1] * transip[1]);
209 return atan2(transip[1], transip[0]);
216 return sqrt(transip * transip);
223 return atan2(transip[1], transip[0]);
230 return atan2(sqrt(transip[0] * transip[0] + transip[1] * transip[1]),
240 return GridF->
GetValue(T, ip, Component);
246 return GridF->
GetValue(*coarse_T, coarse_ip, Component);
284 MFEM_VERIFY(vcenter.
Size() <= 3,
285 "SetDeltaCenter::Maximum number of dim supported is 3")
286 for (
int i = 0; i < vcenter.
Size(); i++) {
center[i] = vcenter[i]; }
327 const int ne = qspace.
GetNE();
330 for (
int iel = 0; iel < ne; ++iel)
335 for (
int iq = 0; iq < ir.
Size(); ++iq)
346void PWVectorCoefficient::InitMap(
const Array<int> & attr,
349 MFEM_VERIFY(attr.
Size() == coefs.
Size(),
350 "PWVectorCoefficient: "
351 "Attribute and coefficient arrays have incompatible "
354 for (
int i=0; i<attr.
Size(); i++)
356 if (coefs[i] != NULL)
366 "PWVectorCoefficient::UpdateCoefficient: "
367 "VectorCoefficient has incompatible dimension.");
368 pieces[attr] = &coef;
375 std::map<int, VectorCoefficient*>::iterator
p = pieces.begin();
376 for (;
p != pieces.end();
p++)
378 if (
p->second != NULL)
380 p->second->SetTime(t);
389 std::map<int, VectorCoefficient*>::const_iterator
p = pieces.find(att);
390 if (
p != pieces.end())
392 if (
p->second != NULL)
394 p->second->Eval(V, T, ip);
421 Function(transip, V);
425 TDFunction(transip,
GetTime(), V);
436 for (
int i = 0; i <
dim; i++)
445 for (
int i = 0; i <
vdim; i++)
447 if (Coeff[i]) { Coeff[i]->SetTime(t); }
454 if (ownCoeff[i]) {
delete Coeff[i]; }
461 for (
int i = 0; i <
vdim; i++)
463 if (ownCoeff[i]) {
delete Coeff[i]; }
471 for (
int i = 0; i <
vdim; i++)
473 V(i) = this->
Eval(i, T, ip);
526 gf -> FESpace() ->
GetMesh() -> SpaceDimension() : 0)
534 gf -> FESpace() ->
GetMesh() -> SpaceDimension() : 0;
571 const int gf_vdim = fes.
GetVDim();
573 if (mesh.
GetNE() == 0) {
return; }
576 "All mesh elements must be the same type!");
580 "All mesh elements must use the same quadrature rule!");
696 const int ne = qspace.
GetNE();
698 for (
int iel = 0; iel < ne; ++iel)
703 for (
int iq = 0; iq < ir.
Size(); ++iq)
715void PWMatrixCoefficient::InitMap(
const Array<int> & attr,
718 MFEM_VERIFY(attr.
Size() == coefs.
Size(),
719 "PWMatrixCoefficient: "
720 "Attribute and coefficient arrays have incompatible "
723 for (
int i=0; i<attr.
Size(); i++)
725 if (coefs[i] != NULL)
735 "PWMatrixCoefficient::UpdateCoefficient: "
736 "MatrixCoefficient has incompatible height.");
738 "PWMatrixCoefficient::UpdateCoefficient: "
739 "MatrixCoefficient has incompatible width.");
743 "PWMatrixCoefficient::UpdateCoefficient: "
744 "MatrixCoefficient has incompatible symmetry.");
746 pieces[attr] = &coef;
753 std::map<int, MatrixCoefficient*>::iterator
p = pieces.begin();
754 for (;
p != pieces.end();
p++)
756 if (
p->second != NULL)
758 p->second->SetTime(t);
767 std::map<int, MatrixCoefficient*>::const_iterator
p = pieces.find(att);
768 if (
p != pieces.end())
770 if (
p->second != NULL)
772 p->second->Eval(K, T, ip);
800 "MatrixFunctionCoefficient is not symmetric");
804 SymmFunction(transip, Ksym);
808 for (
int i=0; i<
height; ++i)
810 for (
int j=i; j<
width; ++j)
812 const real_t Kij = Ksym[j - i + os];
814 if (j != i) { K(j,i) = Kij; }
824 Function(transip, K);
828 TDFunction(transip,
GetTime(), K);
847 "MatrixFunctionCoefficient is not symmetric");
858 SymmFunction(transip, K);
870 MFEM_VERIFY(vdim ==
height*(
height+1)/2,
"Wrong sizes.");
873 const int ne = qspace.
GetNE();
877 for (
int iel = 0; iel < ne; ++iel)
882 for (
int iq = 0; iq < ir.
Size(); ++iq)
897 for (
int j = 0; j <
width; ++j)
899 for (
int i = 0; i <
height; ++ i)
925 Function(transip, K);
929 TDFunction(transip,
GetTime(), K);
958 if (Coeff[i]) { Coeff[i]->SetTime(t); }
965 if (ownCoeff[i*
width+j]) {
delete Coeff[i*
width+j]; }
966 Coeff[i*
width+j] = c;
967 ownCoeff[i*
width+j] = own;
974 if (ownCoeff[i]) {
delete Coeff[i]; }
982 for (
int i = 0; i <
height; i++)
984 for (
int j = 0; j <
width; j++)
986 K(i,j) = this->
Eval(i, j, T, ip);
996 for (
int i = 0; i <
height; i++)
1005 for (
int i=0; i <
height; i++)
1007 if (Coeff[i]) { Coeff[i]->SetTime(t); }
1014 MFEM_ASSERT(i < height && i >= 0,
"Row "
1015 << i <<
" does not exist. " <<
1016 "Matrix height = " <<
height <<
".");
1017 if (ownCoeff[i]) {
delete Coeff[i]; }
1024 for (
int i=0; i <
height; i++)
1026 if (ownCoeff[i]) {
delete Coeff[i]; }
1034 MFEM_ASSERT(i < height && i >= 0,
"Row "
1035 << i <<
" does not exist. " <<
1036 "Matrix height = " <<
height <<
".");
1039 Coeff[i] ->
Eval(V, T, ip);
1053 for (
int i = 0; i <
height; i++)
1055 this->
Eval(i, V, T, ip);
1093 const real_t d_alpha_a = aConst*alpha;
1094 const real_t d_beta = beta;
1099 d_qf[i] = d_alpha_a + d_beta*d_qf[i];
1107 add(alpha, qf, beta, qf_b, qf);
1148 qf = aConst / bConst;
1182 :
a(&A),
b(&B), va(A.GetVDim()), vb(B.GetVDim())
1185 "InnerProductCoefficient: "
1186 "Arguments have incompatible dimensions.");
1207 "Incompatible vector coefficients: a->GetVDim(): "
1210 const int vdim = a->
GetVDim();
1211 MFEM_VERIFY(vdim >= 1,
"invalid vdim: " << vdim);
1215 auto dot_d = qf.
Write();
1223 auto a_d = qf_a.
Read();
1224 auto b_d = qf_b.
Read();
1228 const real_t *ai = a_d + i*vdim;
1229 const real_t *bi = b_d + i*vdim;
1230 real_t dot = ai[0]*bi[0];
1231 for (
int d = 1; d < vdim; d++)
1241 :
a(&A),
b(&B), va(A.GetVDim()), vb(B.GetVDim())
1244 "VectorRotProductCoefficient: "
1245 "Arguments must have dimension equal to two.");
1260 return va[0] * vb[1] - va[1] * vb[0];
1264 :
a(&A), ma(A.GetHeight(), A.GetWidth())
1267 "DeterminantCoefficient: "
1268 "Argument must be a square matrix.");
1285 :
a(&A), ma(A.GetHeight(), A.GetWidth())
1288 "TraceCoefficient: "
1289 "Argument must be a square matrix.");
1307 ACoef(NULL), BCoef(NULL),
1309 alphaCoef(NULL), betaCoef(NULL),
1310 alpha(1.0), beta(1.0)
1319 ACoef(&A_), BCoef(&B_),
1320 A(A_.GetVDim()), B(A_.GetVDim()),
1321 alphaCoef(NULL), betaCoef(NULL),
1325 "VectorSumCoefficient: "
1326 "Arguments must have the same dimension.");
1334 ACoef(&A_), BCoef(&B_),
1339 alpha(0.0), beta(0.0)
1342 "VectorSumCoefficient: "
1343 "Arguments must have the same dimension.");
1348 if (ACoef) { ACoef->
SetTime(t); }
1349 if (BCoef) { BCoef->
SetTime(t); }
1350 if (alphaCoef) { alphaCoef->
SetTime(t); }
1351 if (betaCoef) { betaCoef->
SetTime(t); }
1359 if ( ACoef) { ACoef->
Eval(A, T, ip); }
1360 if ( BCoef) { BCoef->
Eval(B, T, ip); }
1361 if (alphaCoef) { alpha = alphaCoef->
Eval(T, ip); }
1362 if ( betaCoef) { beta = betaCoef->
Eval(T, ip); }
1363 add(alpha, A, beta, B, V);
1388 real_t sa = (a == NULL) ? aConst : a->
Eval(T, ip);
1409 V *= (nv > tol) ? (1.0/nv) : 0.0;
1418 "VectorCrossProductCoefficient: "
1419 "Arguments must have dimension equal to three.");
1435 V[0] = va[1] * vb[2] - va[2] * vb[1];
1436 V[1] = va[2] * vb[0] - va[0] * vb[2];
1437 V[2] = va[0] * vb[1] - va[1] * vb[0];
1443 ma(A.GetHeight(), A.GetWidth()), vb(B.GetVDim())
1446 "MatrixVectorProductCoefficient: "
1447 "Arguments have incompatible dimensions.");
1471 for (
int d=0; d<dim; d++) { M(d,d) = 1.0; }
1479 ma(A.GetHeight(), A.GetWidth())
1482 "MatrixSumCoefficient: "
1483 "Arguments must have the same dimensions.");
1497 if ( beta != 1.0 ) { M *= beta; }
1506 ma(A.GetHeight(), A.GetWidth()),
1507 mb(B.GetHeight(), B.GetWidth())
1510 "MatrixProductCoefficient: "
1511 "Arguments must have compatible dimensions.");
1545 real_t sa = (a == NULL) ? aConst : a->
Eval(T, ip);
1572 "InverseMatrixCoefficient: "
1573 "Argument must be a square matrix.");
1594 "ExponentialMatrixCoefficient: "
1595 <<
"Argument must be a square 2x2 matrix."
1617 va(A.GetVDim()), vb(B.GetVDim())
1633 for (
int i=0; i<va.
Size(); i++)
1635 for (
int j=0; j<vb.
Size(); j++)
1637 M(i, j) = va[i] * vb[j];
1667 for (
int i=0; i<vk.Size(); i++)
1670 for (
int j=0; j<vk.Size(); j++)
1672 M(i, j) -= vk[i] * vk[j];
1675 M *= ((a == NULL ) ? aConst : a->
Eval(T, ip) );
1684 for (
int i = 0; i < mesh.
GetNE(); i++)
1691 tr->SetIntPoint(&ip);
1718 for (
int i = 0; i < mesh.
GetNE(); i++)
1725 tr->SetIntPoint(&ip);
1726 coeff.
Eval(vval, *tr, ip);
1729 for (
int idim(0); idim < vdim; ++idim)
1731 norm += ip.
weight * tr->Weight() * pow(fabs( vval(idim) ),
p);
1736 for (
int idim(0); idim < vdim; ++idim)
1738 val = fabs(vval(idim));
1800 MPI_Comm comm = pmesh.
GetComm();
1808 if (glob_norm < 0.0)
1810 glob_norm = -pow(-glob_norm, 1.0/
p);
1814 glob_norm = pow(glob_norm, 1.0/
p);
1832 MPI_Comm comm = pmesh.
GetComm();
1840 if (glob_norm < 0.0)
1842 glob_norm = -pow(-glob_norm, 1.0/
p);
1846 glob_norm = pow(glob_norm, 1.0/
p);
1865 MFEM_VERIFY(index_ >= 0,
"Index must be >= 0");
1866 MFEM_VERIFY(index_ < QuadF.
GetVDim(),
1867 "Index must be < QuadratureFunction length");
1870 MFEM_VERIFY(length_ > 0,
"Length must be > 0");
1871 MFEM_VERIFY(length_ <= QuadF.
GetVDim() - index,
1872 "Length must be <= (QuadratureFunction length - index)");
1886 if (el_idx < 0) { V = 0.0;
return; }
1899 for (
int i = 0; i <
vdim; i++)
1901 V(i) = temp(index + i);
1916 MFEM_VERIFY(qf.
GetVDim() == 1,
"QuadratureFunction's vdim must be 1");
1927 if (el_idx < 0) {
return 0.0; }
1941 :
Vector(), storage(storage_), vdim(0), qs(qs_), qf(NULL)
1994 MakeRef(qf_coeff->GetQuadFunction());
2012 else if (
auto *qf_coeff =
2015 MakeRef(qf_coeff->GetQuadFunction());
2032 else if (
auto *const_sym_coeff =
2042 const int width = coeff.
GetWidth();
2043 vdim = sym ? height*(height + 1)/2 : width*height;
2047 if (sym) { sym_coeff->ProjectSymmetric(*
qf); }
2062 MFEM_CONTRACT_VAR(qs2);
2063 MFEM_VERIFY(qs2 != NULL,
"Invalid QuadratureSpace.")
2082 for (
int iq = 0; iq < nq; ++iq)
2084 for (
int vd = 0; vd<
vdim; ++vd)
2086 (*this)[vd + iq*
vdim] = constant[vd];
2094 const int width = constant.
Width();
2095 const int height = constant.
Height();
2096 vdim = width*height;
2098 for (
int iq = 0; iq < nq; ++iq)
2100 for (
int j = 0; j < width; ++j)
2102 for (
int i = 0; i < height; ++i)
2104 const real_t val = transpose ? constant(j,i) : constant(i,j);
2105 (*this)[i + j*height + iq*
vdim] = val;
2114 const int height = constant.
Height();
2116 vdim = sym ? height*(height + 1)/2 : height*height;
2118 for (
int iq = 0; iq < nq; ++iq)
2120 for (
int vd = 0; vd <
vdim; ++vd)
2122 const real_t value = sym ? constant.
GetData()[vd] : constant(vd % height,
2124 (*this)[vd + iq*
vdim] = value;
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Size() const
Return the logical size of the array.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
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.
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) override
Fill the QuadratureFunction qf with the constant value.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
CrossCrossCoefficient(real_t A, VectorCoefficient &K)
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void SetGridFunction(const GridFunction *gf)
Set the vector grid function.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector curl coefficient at ip.
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...
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void GetDeltaCenter(Vector ¢er)
Write the center of the delta function into center.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
real_t Scale()
Return the scale factor times the optional time dependent function. Returns with when not set by th...
void SetDeltaCenter(const Vector ¢er)
Set the center location of the delta function.
virtual real_t EvalDelta(ElementTransformation &T, const IntegrationPoint &ip)
The value of the function assuming we are evaluating at the delta center.
Data type dense matrix using column-major storage.
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication.
void Transpose()
(*this) = (*this)^t
void SetRow(int r, const real_t *row)
void GetColumnReference(int c, Vector &col)
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
void Invert()
Replaces the current matrix with its inverse.
real_t Trace() const
Trace of a square matrix.
void UseExternalData(real_t *d, int h, int w)
Change the data array and the size of the DenseMatrix.
void Add(const real_t c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
void SetSize(int s)
Change the size of the DenseSymmetricMatrix to s x s.
void UseExternalData(real_t *d, int s)
Change the data array and the size of the DenseSymmetricMatrix.
real_t * GetData() const
Returns the matrix data array.
DeterminantCoefficient(MatrixCoefficient &A)
Construct with the matrix.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the determinant coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the scalar divergence coefficient at ip.
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
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
ExponentialMatrixCoefficient(MatrixCoefficient &A)
Construct the matrix coefficient. Result is .
Class representing the storage layout of a FaceQuadratureFunction.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Mesh * GetMesh() const
Returns the mesh.
int GetVDim() const
Returns the vector dimension of the finite element space.
std::function< real_t(const Vector &)> Function
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
std::function< real_t(const Vector &, real_t)> TDFunction
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the gradient vector coefficient at ip.
void SetGridFunction(const GridFunction *gf)
Set 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...
const GridFunction * GridFunc
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
Class for grid function - Vector with associated FE space.
virtual real_t GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
void GetGradients(ElementTransformation &tr, const IntegrationRule &ir, DenseMatrix &grad) const
Extension of GetGradient(...) for a collection of IntegrationPoints.
FiniteElementSpace * FESpace()
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
real_t GetDivergence(ElementTransformation &tr) const
void GetCurl(ElementTransformation &tr, Vector &curl) const
void GetGradient(ElementTransformation &tr, Vector &grad) const
Gradient of a scalar function at a quadrature point.
void GetVectorValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
InnerProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with the two vector coefficients. Result is .
Class for integration point with weight.
Class for an integration rule - an Array of IntegrationPoint.
int GetNPoints() const
Returns the number of the points in the integration rule.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
InverseMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
void SetTime(real_t t) override
Set the time for internally stored coefficients.
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) override
Set the time for internally stored coefficients.
void Eval(int i, Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
MatrixArrayVectorCoefficient(int dim)
Construct a coefficient matrix of dimensions dim * dim. The actual coefficients still need to be adde...
void SetTime(real_t t) override
Set the time for internally stored coefficients.
virtual ~MatrixArrayVectorCoefficient()
void Set(int i, VectorCoefficient *c, bool own=true)
Set the coefficient located at the i-th row of the matrix. By this will take ownership of the Coeffic...
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 void SetTime(real_t t)
Set the time for time dependent coefficients.
real_t 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 ...
int GetWidth() const
Get the width of the matrix.
int GetHeight() const
Get the height of the matrix.
A matrix coefficient that is constant in space and time.
void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void EvalSymmetric(Vector &K, ElementTransformation &T, const IntegrationPoint &ip) override
(DEPRECATED) Evaluate the symmetric matrix coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
MatrixProductCoefficient(MatrixCoefficient &A, MatrixCoefficient &B)
Construct with the two coefficients. Result is A * B.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
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 Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
MatrixVectorProductCoefficient(MatrixCoefficient &A, VectorCoefficient &B)
Constructor with two coefficients. Result is A*B.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
Element::Type GetElementType(int i) const
Returns the type of element i.
const CoarseFineTransformations & GetRefinementTransforms() const
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
int SpaceDimension() const
Dimension of the physical space containing the mesh.
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
NormalizedVectorCoefficient(VectorCoefficient &A, real_t tol=1e-6)
Return a vector normalized to a length of one.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
void SetTime(real_t t) override
Set the time for internally stored coefficients.
OuterProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with two vector coefficients. Result is .
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient.
void UpdateCoefficient(int attr, Coefficient &coef)
Replace a single Coefficient for a particular attribute.
void SetTime(real_t t) override
Set the time for time dependent coefficients.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient.
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf with the piecewise constant values.
void UpdateCoefficient(int attr, MatrixCoefficient &coef)
Replace a single coefficient for a particular attribute.
void SetTime(real_t t) override
Set the time for time dependent coefficients.
void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient.
void UpdateCoefficient(int attr, VectorCoefficient &coef)
Replace a single Coefficient for a particular attribute.
void SetTime(real_t t) override
Set the time for time dependent coefficients.
Class for parallel meshes.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
Quadrature function coefficient which requires that the quadrature rules used for this coefficient be...
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient in the element described by T at the point ip.
QuadratureFunctionCoefficient(const QuadratureFunction &qf)
Constructor with a quadrature function as input.
Represents values or vectors of values at quadrature points on a mesh.
void SetVDim(int vdim_)
Set the vector dimension, updating the size by calling Vector::SetSize().
QuadratureSpaceBase * GetSpace()
Get the associated QuadratureSpaceBase object.
void ProjectGridFunction(const GridFunction &gf)
Evaluate a grid function at each quadrature point.
void GetValues(int idx, Vector &values)
Return all values associated with mesh element idx in a Vector.
int GetVDim() const
Get the vector dimension.
const IntegrationRule & GetIntRule(int idx) const
Get the IntegrationRule associated with entity (element or face) idx.
Abstract base class for QuadratureSpace and FaceQuadratureSpace.
int GetNE() const
Return the number of entities.
virtual ElementTransformation * GetTransformation(int idx)=0
Get the (element or face) transformation of entity idx.
int GetSize() const
Return the total number of quadrature points.
const IntegrationRule & GetIntRule(int idx) const
Return the IntegrationRule associated with entity idx.
virtual int GetEntityIndex(const ElementTransformation &T) const =0
Returns the index in the quadrature space of the entity associated with the transformation T.
virtual int GetPermutedIndex(int idx, int iq) const =0
Returns the permuted index of the iq quadrature point in entity idx.
const Array< int > & Offsets(QSpaceOffsetStorage storage) const
Entity quadrature point offset array.
int GetOrder() const
Return the order of the quadrature rule(s) used by all elements.
Mesh * GetMesh() const
Returns the mesh.
Class representing the storage layout of a QuadratureFunction.
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
ScalarMatrixProductCoefficient(real_t A, MatrixCoefficient &B)
Constructor with one coefficient. Result is A*B.
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
ScalarVectorProductCoefficient(real_t A, VectorCoefficient &B)
Constructor with constant and vector coefficient. Result is A * B.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
Base class for symmetric matrix coefficients that optionally depend on time and space.
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 ...
DenseSymmetricMatrix mat_aux
Internal matrix used when evaluating this coefficient as a DenseMatrix.
A matrix coefficient that is constant in space and time.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void Eval(DenseSymmetricMatrix &K, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the trace coefficient at ip.
TraceCoefficient(MatrixCoefficient &A)
Construct with the matrix.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
TransposeMatrixCoefficient(MatrixCoefficient &A)
Construct with the matrix coefficient. Result is .
void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the matrix coefficient at ip.
virtual ~VectorArrayCoefficient()
Destroys vector coefficient.
VectorArrayCoefficient(int dim)
Construct vector of dim coefficients. The actual coefficients still need to be added with Set().
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void Set(int i, Coefficient *c, bool own=true)
Sets coefficient in the vector.
real_t Eval(int i, ElementTransformation &T, const IntegrationPoint &ip)
Base class for vector Coefficients that optionally depend on time and space.
int GetVDim()
Returns dimension of the vector.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
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.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
VectorCrossProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Construct with the two coefficients. Result is A x B.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
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...
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
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
VectorGridFunctionCoefficient()
Construct an empty coefficient. Calling Eval() before the grid function is set will cause a segfault.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
Vector quadrature function coefficient which requires that the quadrature rules used for this vector ...
void Project(QuadratureFunction &qf) override
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
VectorQuadratureFunctionCoefficient(const QuadratureFunction &qf)
Constructor with a quadrature function as input.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
void SetComponent(int index_, int length_)
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the vector coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
VectorRotProductCoefficient(VectorCoefficient &A, VectorCoefficient &B)
Constructor with two vector coefficients. Result is .
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient at ip.
void SetTime(real_t t) override
Set the time for internally stored coefficients.
VectorSumCoefficient(int dim)
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
real_t Norml2() const
Returns the l2 norm of the vector.
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
int Size() const
Returns the size of the vector.
void SetSize(int s)
Resize the vector to size s.
void Reciprocal()
(*this)(i) = 1.0 / (*this)(i)
virtual real_t * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Vector & operator=(const real_t *v)
Copy Size() entries from v.
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
int index(int i, int j, int nx, int ny)
real_t infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
real_t ComputeGlobalLpNorm(real_t p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
void add(const Vector &v1, const Vector &v2, Vector &v)
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.
real_t LpNormLoop(real_t p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
void forall(int N, lambda &&body)
ElementTransformation * RefinedToCoarse(Mesh &coarse_mesh, const ElementTransformation &T, const IntegrationPoint &ip, IntegrationPoint &coarse_ip)
real_t p(const Vector &x, real_t t)
MFEM_HOST_DEVICE real_t norm(const Complex &z)
Helper struct to convert a C++ type to an MPI type.