22 MFEM_ABORT(
"Not supported.");
28 mfem_error(
"LinearFormIntegrator::AssembleRHSElementVect(...)");
35 mfem_error(
"LinearFormIntegrator::AssembleRHSElementVect(...)");
66 add(elvect, ip.
weight * val, shape, elvect);
73 MFEM_ASSERT(
delta != NULL,
"coefficient must be DeltaCoefficient");
104 Q.
Eval(Qvec, Tr, ip);
114 MFEM_ASSERT(
vec_delta != NULL,
"coefficient must be VectorDeltaCoefficient");
124 dshape.
Mult(Qvec, elvect);
139 int intorder = oa * el.
GetOrder() + ob;
152 add(elvect, ip.
weight * val, shape, elvect);
168 int intorder = oa * el.
GetOrder() + ob;
186 add(elvect, val, shape, elvect);
204 int intorder = oa * el.
GetOrder() + ob;
221 Q.
Eval(Qvec, Tr, ip);
242 mfem_error(
"These methods make sense only in 2D problems.");
248 int intorder = oa * el.
GetOrder() + ob;
258 tangent(0) = Jac(0,0);
259 tangent(1) = Jac(1,0);
261 Q.
Eval(Qvec, Tr, ip);
265 add(elvect, ip.
weight*(Qvec*tangent), shape, elvect);
297 Q.
Eval (Qvec, Tr, ip);
299 for (
int k = 0; k < vdim; k++)
303 for (
int s = 0; s < dof; s++)
305 elvect(dof*k+s) += ip.
weight * cf * shape(s);
314 MFEM_ASSERT(
vec_delta != NULL,
"coefficient must be VectorDeltaCoefficient");
325 MultVWt(shape, Qvec, elvec_as_mat);
332 const int dof = el.
GetDof();
338 elvect.
SetSize(dof*(vdim/sdim));
358 Q.
Eval(Qvec, Tr, ip);
361 for (
int k = 0; k < vdim/sdim; k++)
363 for (
int d=0; d < sdim; ++d) { part_x(d) = Qvec(k*sdim+d); }
364 dshape.
Mult(part_x, pelvect);
365 for (
int s = 0; s < dof; ++s) { elvect(s+k*dof) += pelvect(s); }
373 MFEM_ABORT(
"Not implemented!");
403 for (
int k = 0; k < vdim; k++)
404 for (
int s = 0; s < dof; s++)
406 elvect(dof*k+s) += vec(k) * shape(s);
444 for (
int k = 0; k < vdim; k++)
446 for (
int s = 0; s < dof; s++)
448 elvect(dof*k+s) += vec(k) * shape(s);
482 QF.
Eval (vec, Tr, ip);
491 MFEM_ASSERT(
vec_delta != NULL,
"coefficient must be VectorDeltaCoefficient");
501 vshape.
Mult(vec, elvect);
509 int n=(spaceDim == 3)? spaceDim : 1;
529 QF->
Eval(vec, Tr, ip);
532 curlshape.
AddMult (vec, elvect);
541 "coefficient must be VectorDeltaCoefficient");
543 int n=(spaceDim == 3)? spaceDim : 1;
550 curlshape.
Mult(vec, elvect);
577 add(elvect, ip.
weight * val, divshape, elvect);
584 MFEM_ASSERT(
delta != NULL,
"coefficient must be DeltaCoefficient");
613 nor *= Sign * ip.
weight * F -> Eval (Tr, ip);
614 for (
int j = 0; j < dof; j++)
615 for (
int k = 0; k <
dim; k++)
617 elvect(dof*k+j) += nor(k) * shape(j);
635 int intorder = oa * el.
GetOrder() + ob;
648 val *= F->
Eval(Tr, ip);
651 elvect.
Add(val, shape);
679 F.
Eval(Fvec, Tr, ip);
684 elvect.
Add(val, shape);
698 MFEM_VERIFY(vdim == 2,
"VectorFEBoundaryTangentLFIntegrator "
699 "must be called with vector basis functions of dimension 2.");
707 int intorder = oa * el.
GetOrder() + ob;
718 f.
Eval(f_loc, Tr, ip);
727 f_hat(0) = J(0,0) * f_loc(0) + J(1,0) * f_loc(1);
737 f_hat(0) = -f_hat(0);
747 mfem_error(
"BoundaryFlowIntegrator::AssembleRHSElementVect\n"
748 " is not implemented as boundary integrator!\n"
749 " Use LinearForm::AddBdrFaceIntegrator instead of\n"
750 " LinearForm::AddBoundaryIntegrator.");
756 int dim, ndof, order;
757 real_t un, w, vu_data[3], nor_data[3];
796 nor(0) = 2*eip.
x - 1.0;
804 w = 0.5*alpha*un - beta*fabs(un);
806 elvect.
Add(w, shape);
813 mfem_error(
"DGDirichletLFIntegrator::AssembleRHSElementVect");
820 bool kappa_is_nonzero = (
kappa != 0.);
862 nor(0) = 2*eip.
x - 1.0;
894 if (kappa_is_nonzero)
904 mfem_error(
"DGElasticityDirichletLFIntegrator::AssembleRHSElementVect");
910 MFEM_ASSERT(Tr.
Elem2No < 0,
"interior boundary is not supported");
912#ifdef MFEM_THREAD_SAFE
924 const int ndofs = el.
GetDof();
925 const int nvdofs =
dim*ndofs;
967 nor(0) = 2*eip.
x - 1.0;
1023 for (
int im = 0, i = 0; im <
dim; ++im)
1028 for (
int idof = 0; idof < ndofs; ++idof, ++i)
1046 for (
int i = 0; i < n; i++)
1048 elvect(i) = dist(generator);
1053 if (!save_factors || !L[iel])
1085 const int vdim = vqfc.
GetVDim();
1086 const int ndofs = fe.
GetDof();
1091 for (
int q = 0; q < nqp; q++)
1096 vqfc.
Eval(temp, Tr, ip);
1098 for (
int ind = 0; ind < vdim; ind++)
1100 for (
int nd = 0; nd < ndofs; nd++)
1102 elvect(nd + ind * ndofs) += w * shape(nd) * temp(ind);
1117 const int ndofs = fe.
GetDof();
1121 for (
int q = 0; q < nqp; q++)
1128 shape *= (w * temp);
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
bool Factor(int m, real_t TOL=0.0) override
Compute the Cholesky factorization of the current matrix.
void LMult(int m, int n, real_t *X) const
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
virtual real_t EvalDelta(ElementTransformation &T, const IntegrationPoint &ip)
The value of the function assuming we are evaluating at the delta center.
VectorDeltaCoefficient * vec_delta
Data type dense matrix using column-major storage.
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication.
void MultTranspose(const real_t *x, real_t *y) const
Multiply a vector with the transpose matrix.
real_t * Data() const
Returns the matrix data array.
void AddMult(const Vector &x, Vector &y, const real_t a=1.0) const override
y += a * A.x
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Abstract class for all finite elements.
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
int GetRangeDim() const
Returns the vector dimension for vector-valued finite elements, which is also the dimension of the in...
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
int GetDim() const
Returns the reference space dimension for the finite element.
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
int Space() const
Returns the type of FunctionSpace on the element.
void CalcPhysDivShape(ElementTransformation &Trans, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in physical space at the po...
void CalcPhysVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Equivalent to the CalcVShape() method with the same arguments.
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
virtual void CalcPhysCurlShape(ElementTransformation &Trans, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in physical space at the point de...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in physical space at the point ...
@ Pk
Polynomials of order k.
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.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
const IntegrationRule * GetIntegrationRule() const
Equivalent to GetIntRule, but retained for backward compatibility with applications.
const IntegrationRule * IntRule
void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat) override
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 Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
real_t Eval(ElementTransformation &T, const IntegrationPoint &ip) override
Evaluate the coefficient in the element described by T at the point ip.
const QuadratureFunction & GetQuadFunction() const
QuadratureSpaceBase * GetSpace()
Get the associated QuadratureSpaceBase object.
void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect) override
const IntegrationRule & GetIntRule(int idx) const
Return the IntegrationRule associated with entity idx.
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
int GetVDim()
Returns dimension of the vector.
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 EvalDelta(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Return the specified direction vector multiplied by the value returned by DeltaCoefficient::EvalDelta...
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
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 ...
const QuadratureFunction & GetQuadFunction() const
void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect) override
Vector & Set(const real_t a, const Vector &x)
(*this) = a * x
void SetSize(int s)
Resize the vector to size s.
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Vector & Add(const real_t a, const Vector &Va)
(*this) += a * Va
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
void CalcOrtho(const DenseMatrix &J, Vector &n)
void Swap(Array< T > &, Array< T > &)
void mfem_error(const char *msg)
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
void add(const Vector &v1, const Vector &v2, Vector &v)
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
real_t p(const Vector &x, real_t t)