12 #ifndef MFEM_BILININTEG
13 #define MFEM_BILININTEG
15 #include "../config/config.hpp"
76 Vector &flux,
int with_coef = 1) { }
96 { bfi = _bfi; own_bfi = _own_bfi; }
124 { bfi = _bfi; own_bfi = _own_bfi; }
142 { integrator = integ; own_integrator = own_integ; }
163 { integrators.Append(integ); }
208 return "MixedScalarIntegrator: "
209 "Trial and test spaces must both be scalar fields.";
232 #ifndef MFEM_THREAD_SAFE
267 :
same_calc_shape(false), Q(NULL), VQ(diag?NULL:&dq), DQ(diag?&dq:NULL),
282 return "MixedVectorIntegrator: "
283 "Trial and test spaces must both be vector fields";
309 #ifndef MFEM_THREAD_SAFE
335 bool _cross_2d =
false)
336 : VQ(&vq), transpose(_transpose) , cross_2d(_cross_2d) {}
342 return ((transpose &&
355 return "MixedScalarVectorIntegrator: "
356 "Trial space must be a vector field "
357 "and the test space must be a scalar field";
361 return "MixedScalarVectorIntegrator: "
362 "Trial space must be a scalar field "
363 "and the test space must be a vector field";
389 #ifndef MFEM_THREAD_SAFE
433 return (trial_fe.
GetDim() == 1 && test_fe.
GetDim() == 1 &&
440 return "MixedScalarDerivativeIntegrator: "
441 "Trial and test spaces must both be scalar fields in 1D "
442 "and the trial space must implement CaldDShape.";
468 return (trial_fe.
GetDim() == 1 && test_fe.
GetDim() == 1 &&
475 return "MixedScalarWeakDerivativeIntegrator: "
476 "Trial and test spaces must both be scalar fields in 1D "
477 "and the test space must implement CalcDShape with "
478 "map type \"VALUE\".";
512 return "MixedScalarDivergenceIntegrator: "
513 "Trial must be H(Div) and the test space must be a "
548 return "MixedVectorDivergenceIntegrator: "
549 "Trial must be H(Div) and the test space must be a "
585 return "MixedScalarWeakGradientIntegrator: "
586 "Trial space must be a scalar field "
587 "and the test space must be H(Div)";
614 return (trial_fe.
GetDim() == 2 && test_fe.
GetDim() == 2 &&
621 return "MixedScalarCurlIntegrator: "
622 "Trial must be H(Curl) and the test space must be a "
655 return (trial_fe.
GetDim() == 2 && test_fe.
GetDim() == 2 &&
662 return "MixedScalarWeakCurlIntegrator: "
663 "Trial space must be a scalar field "
664 "and the test space must be H(Curl)";
719 return "MixedDotProductIntegrator: "
720 "Trial space must be a vector field "
721 "and the test space must be a scalar field";
745 return "MixedWeakGradDotIntegrator: "
746 "Trial space must be a vector field "
747 "and the test space must be a vector field with a divergence";
768 return (trial_fe.
GetDim() == 3 && test_fe.
GetDim() == 3 &&
776 return "MixedWeakDivCrossIntegrator: "
777 "Trial space must be a vector field in 3D "
778 "and the test space must be a scalar field with a gradient";
813 return "MixedGradGradIntegrator: "
814 "Trial and test spaces must both be scalar fields "
815 "with a gradient operator.";
859 return "MixedCrossGradGradIntegrator: "
860 "Trial and test spaces must both be scalar fields "
861 "with a gradient operator.";
893 return (trial_fe.
GetDim() == 3 && test_fe.
GetDim() == 3 &&
902 return "MixedCurlCurlIntegrator"
903 "Trial and test spaces must both be vector fields in 3D "
930 return (trial_fe.
GetDim() == 3 && test_fe.
GetDim() == 3 &&
939 return "MixedCrossCurlCurlIntegrator: "
940 "Trial and test spaces must both be vector fields in 3D "
967 return (trial_fe.
GetDim() == 3 && test_fe.
GetDim() == 3 &&
976 return "MixedCrossCurlGradIntegrator"
977 "Trial space must be a vector field in 3D with a curl"
978 "and the test space must be a scalar field with a gradient";
1004 return (trial_fe.
GetDim() == 3 && test_fe.
GetDim() == 3 &&
1013 return "MixedCrossGradCurlIntegrator"
1014 "Trial space must be a scalar field in 3D with a gradient"
1015 "and the test space must be a vector field with a curl";
1042 return (trial_fe.
GetDim() == 3 && test_fe.
GetDim() == 3 &&
1050 return "MixedWeakCurlCrossIntegrator: "
1051 "Trial space must be a vector field in 3D "
1052 "and the test space must be a vector field with a curl";
1074 return (trial_fe.
GetDim() == 2 && test_fe.
GetDim() == 2 &&
1082 return "MixedScalarWeakCurlCrossIntegrator: "
1083 "Trial space must be a vector field in 2D "
1084 "and the test space must be a vector field with a curl";
1109 return (trial_fe.
GetDim() == 3 && test_fe.
GetDim() == 3 &&
1117 return "MixedCrossGradIntegrator: "
1118 "Trial space must be a scalar field with a gradient operator"
1119 " and the test space must be a vector field both in 3D.";
1146 return (trial_fe.
GetDim() == 3 && test_fe.
GetDim() == 3 &&
1154 return "MixedCrossCurlIntegrator: "
1155 "Trial space must be a vector field in 3D with a curl "
1156 "and the test space must be a vector field";
1178 return (trial_fe.
GetDim() == 2 && test_fe.
GetDim() == 2 &&
1186 return "MixedCrossCurlIntegrator: "
1187 "Trial space must be a vector field in 2D with a curl "
1188 "and the test space must be a vector field";
1212 return (trial_fe.
GetDim() == 2 && test_fe.
GetDim() == 2 &&
1220 return "MixedScalarCrossGradIntegrator: "
1221 "Trial space must be a scalar field in 2D with a gradient "
1222 "and the test space must be a scalar field";
1243 return (trial_fe.
GetDim() == 2 && test_fe.
GetDim() == 2 &&
1250 return "MixedScalarCrossProductIntegrator: "
1251 "Trial space must be a vector field in 2D "
1252 "and the test space must be a scalar field";
1268 return (trial_fe.
GetDim() == 2 && test_fe.
GetDim() == 2 &&
1275 return "MixedScalarWeakCrossProductIntegrator: "
1276 "Trial space must be a scalar field in 2D "
1277 "and the test space must be a vector field";
1305 return "MixedDirectionalDerivativeIntegrator: "
1306 "Trial space must be a scalar field with a gradient "
1307 "and the test space must be a scalar field";
1336 return "MixedGradDivIntegrator: "
1337 "Trial space must be a scalar field with a gradient"
1338 "and the test space must be a vector field with a divergence";
1373 return "MixedDivGradIntegrator: "
1374 "Trial space must be a vector field with a divergence"
1375 "and the test space must be a scalar field with a gradient";
1408 return "MixedScalarWeakDivergenceIntegrator: "
1409 "Trial space must be a scalar field "
1410 "and the test space must be a scalar field with a gradient";
1444 return "MixedVectorGradientIntegrator: "
1445 "Trial spaces must be H1 and the test space must be a "
1446 "vector field in 2D or 3D";
1479 return (trial_fe.
GetDim() == 3 && test_fe.
GetDim() == 3 &&
1486 return "MixedVectorCurlIntegrator: "
1487 "Trial space must be H(Curl) and the test space must be a "
1488 "vector field in 3D";
1518 return (trial_fe.
GetDim() == 3 && test_fe.
GetDim() == 3 &&
1525 return "MixedVectorWeakCurlIntegrator: "
1526 "Trial space must be vector field in 3D and the "
1527 "test space must be H(Curl)";
1563 return "MixedVectorWeakDivergenceIntegrator: "
1564 "Trial space must be vector field and the "
1565 "test space must be H1";
1582 Vector vec, pointflux, shape;
1583 #ifndef MFEM_THREAD_SAFE
1620 Vector &flux,
int with_coef = 1);
1631 #ifndef MFEM_THREAD_SAFE
1671 #ifndef MFEM_THREAD_SAFE
1673 Vector shape, vec2, BdFidxT;
1680 : Q(q) { alpha = a; }
1697 : Q(q) { alpha = a; }
1709 Vector shape, te_shape, vec;
1721 { Q = NULL; VQ = NULL; MQ = NULL; Q_order = 0; }
1727 : Q(&q) { VQ = NULL; MQ = NULL; Q_order = qo; }
1730 { VQ = NULL; MQ = NULL; Q_order = 0; }
1733 : VQ(&q) { Q = NULL; MQ = NULL; Q_order = qo; }
1736 : MQ(&q) { Q = NULL; VQ = NULL; Q_order = qo; }
1759 #ifndef MFEM_THREAD_SAFE
1781 #ifndef MFEM_THREAD_SAFE
1805 #ifndef MFEM_THREAD_SAFE
1848 #ifndef MFEM_THREAD_SAFE
1870 Vector &flux,
int with_coef);
1882 #ifndef MFEM_THREAD_SAFE
1883 DenseMatrix dshape_hat, dshape, curlshape, Jadj, grad_hat, grad;
1910 { Q = q; VQ = vq; MQ = mq; }
1912 #ifndef MFEM_THREAD_SAFE
1968 #ifndef MFEM_THREAD_SAFE
2017 double q_lambda, q_mu;
2020 #ifndef MFEM_THREAD_SAFE
2027 { lambda = &l; mu = &m; }
2031 { lambda = NULL; mu = &m; q_lambda = q_l; q_mu = q_m; }
2058 { rho = NULL; u = &_u; alpha = a; beta = b; }
2062 { rho = &_rho; u = &_u; alpha = a; beta = b; }
2175 double alpha_,
double kappa_)
2188 #ifndef MFEM_THREAD_SAFE
2210 const int dim,
const int row_ndofs,
const int col_ndofs,
2211 const int row_offset,
const int col_offset,
2212 const double jmatcoef,
const Vector &col_nL,
const Vector &col_nM,
2224 Vector face_shape, shape1, shape2;
2242 Vector face_shape, normal, shape1_n, shape2_n;
2284 { ran_fe.
Project(dom_fe, Trans, elmat); }
Abstract class for Finite Elements.
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
VectorCurlCurlIntegrator()
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &elmat)
Support for use in BilinearForm. Can be used only when appropriate.
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleElementMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &elmat)
Support for use in BilinearForm. Can be used only when appropriate.
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
VectorFEMassIntegrator(Coefficient &q)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
int GetDim() const
Returns the reference space dimension for the finite element.
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual const char * FiniteElementTypeFailureMessage() const
Class for an integration rule - an Array of IntegrationPoint.
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
MixedCrossCurlGradIntegrator(VectorCoefficient &vq)
DGTraceIntegrator(VectorCoefficient &_u, double a, double b)
Construct integrator with rho = 1.
MixedGradDivIntegrator(VectorCoefficient &vq)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
DGDiffusionIntegrator(Coefficient &q, const double s, const double k)
Integrator defining a sum of multiple Integrators.
SumIntegrator(int own_integs=1)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
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...
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
DiffusionIntegrator(Coefficient &q)
Construct a diffusion integrator with a scalar coefficient q.
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedScalarCrossProductIntegrator(VectorCoefficient &vq)
CurlCurlIntegrator(MatrixCoefficient &m)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
MixedVectorWeakDivergenceIntegrator(VectorCoefficient &dq)
virtual double ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
MixedVectorCurlIntegrator()
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...
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
ElasticityIntegrator(Coefficient &l, Coefficient &m)
VectorFEMassIntegrator(VectorCoefficient *_vq)
NormalTraceJumpIntegrator()
VectorDiffusionIntegrator(Coefficient &q)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Integrator for (curl u, curl v) for Nedelec elements.
MixedVectorMassIntegrator()
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
VectorCurlCurlIntegrator(Coefficient &q)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
MixedScalarMassIntegrator(Coefficient &q)
VectorMassIntegrator(Coefficient &q, const IntegrationRule *ir)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
alpha (q . grad u, v) using the "group" FE discretization
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
VectorDiffusionIntegrator()
int GetOrder() const
Returns the order of the finite element.
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Data type dense matrix using column-major storage.
DivDivIntegrator(Coefficient &q)
int Size() const
Returns the size of the vector.
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
MixedVectorWeakCurlIntegrator(MatrixCoefficient &mq)
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
virtual const char * FiniteElementTypeFailureMessage() const
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual ~TransposeIntegrator()
VectorFEDivergenceIntegrator()
int Space() const
Returns the type of space on each element.
void AddIntegrator(BilinearFormIntegrator *integ)
virtual const char * FiniteElementTypeFailureMessage() const
MixedVectorCurlIntegrator(VectorCoefficient &dq)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
MixedCrossProductIntegrator(VectorCoefficient &vq)
MixedScalarVectorIntegrator(VectorCoefficient &vq, bool _transpose=false, bool _cross_2d=false)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
MixedVectorGradientIntegrator(Coefficient &q)
(Q div u, div v) for RT elements
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual const char * FiniteElementTypeFailureMessage() const
MixedScalarCrossGradIntegrator(VectorCoefficient &vq)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef=1)
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
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...
virtual ~LumpedIntegrator()
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Implements CalcDivShape methods.
MixedCurlCurlIntegrator(MatrixCoefficient &mq)
VectorFEMassIntegrator(VectorCoefficient &vq)
VectorMassIntegrator()
Construct an integrator with coefficient 1.0.
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
VectorMassIntegrator(VectorCoefficient &q, int qo=0)
Construct an integrator with diagonal coefficient q.
VectorCrossProductInterpolator(VectorCoefficient &vc)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
MixedGradGradIntegrator(VectorCoefficient &dq)
Integrator that inverts the matrix assembled by another integrator.
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
ElasticityIntegrator(Coefficient &m, double q_l, double q_m)
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual const char * FiniteElementTypeFailureMessage() const
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
VectorInnerProductInterpolator(VectorCoefficient &vc)
MixedWeakDivCrossIntegrator(VectorCoefficient &vq)
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
MixedScalarWeakGradientIntegrator()
MixedVectorGradientIntegrator()
VectorDivergenceIntegrator(Coefficient &q)
MixedCrossCurlCurlIntegrator(VectorCoefficient &vq)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedVectorWeakCurlIntegrator(VectorCoefficient &dq)
MixedScalarWeakCurlIntegrator(Coefficient &q)
MixedCurlCurlIntegrator(Coefficient &q)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Implements CalcDShape methods.
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual const char * FiniteElementTypeFailureMessage() const
VectorFEDivergenceIntegrator(Coefficient &q)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
VectorFECurlIntegrator(Coefficient &q)
virtual void AssembleElementMatrix2(const FiniteElement &h1_fe, const FiniteElement &nd_fe, ElementTransformation &Trans, DenseMatrix &elmat)
MixedVectorWeakCurlIntegrator()
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual const char * FiniteElementTypeFailureMessage() const
DiffusionIntegrator()
Construct a diffusion integrator with coefficient Q = 1.
MixedVectorMassIntegrator(MatrixCoefficient &mq)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
MixedScalarDivergenceIntegrator(Coefficient &q)
MixedVectorIntegrator(MatrixCoefficient &mq)
Implements CalcCurlShape methods.
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
MixedVectorMassIntegrator(Coefficient &q)
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...
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedDirectionalDerivativeIntegrator(VectorCoefficient &vq)
MixedScalarMassIntegrator()
MixedScalarDerivativeIntegrator()
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
MixedCrossGradIntegrator(VectorCoefficient &vq)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
ScalarProductInterpolator(Coefficient &sc)
MixedGradGradIntegrator()
virtual const char * FiniteElementTypeFailureMessage() const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
MixedVectorCurlIntegrator(MatrixCoefficient &mq)
MixedCurlCurlIntegrator(VectorCoefficient &dq)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
virtual const char * FiniteElementTypeFailureMessage() const
BoundaryMassIntegrator(Coefficient &q)
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
MixedScalarWeakDerivativeIntegrator()
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the BilinearFormIntegrator.
MixedWeakGradDotIntegrator(VectorCoefficient &vq)
ScalarVectorProductInterpolator(Coefficient &sc)
virtual const char * FiniteElementTypeFailureMessage() const
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
MixedScalarWeakGradientIntegrator(Coefficient &q)
VectorDivergenceIntegrator(Coefficient *_q)
MixedVectorIntegrator(Coefficient &q)
DiffusionIntegrator(MatrixCoefficient &q)
Construct a diffusion integrator with a matrix coefficient q.
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual const char * FiniteElementTypeFailureMessage() const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedGradGradIntegrator(MatrixCoefficient &mq)
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
VectorFEMassIntegrator(MatrixCoefficient &mq)
MixedVectorMassIntegrator(VectorCoefficient &dq)
MixedScalarWeakCurlCrossIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Base class Coefficient that may optionally depend on time.
MixedScalarWeakCurlIntegrator()
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
MixedScalarCurlIntegrator()
virtual const char * FiniteElementTypeFailureMessage() const
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 ...
virtual const char * FiniteElementTypeFailureMessage() const
MixedCurlCurlIntegrator()
MixedVectorWeakDivergenceIntegrator(Coefficient &q)
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
MixedScalarWeakDivergenceIntegrator(VectorCoefficient &vq)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
DGDiffusionIntegrator(const double s, const double k)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
LumpedIntegrator(BilinearFormIntegrator *_bfi, int _own_bfi=1)
MixedCrossGradCurlIntegrator(VectorCoefficient &vq)
MixedVectorDivergenceIntegrator(VectorCoefficient &vq)
MassIntegrator(Coefficient &q, const IntegrationRule *ir=NULL)
Construct a mass integrator with coefficient q.
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual const char * FiniteElementTypeFailureMessage() const
MixedCrossGradGradIntegrator(VectorCoefficient &vq)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedVectorIntegrator(VectorCoefficient &dq, bool diag=true)
MixedScalarCrossCurlIntegrator(VectorCoefficient &vq)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Assemble an element matrix.
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
VectorFEWeakDivergenceIntegrator(Coefficient &q)
VectorMassIntegrator(Coefficient &q, int qo=0)
VectorFEMassIntegrator(Coefficient *_q)
DGElasticityIntegrator(double alpha_, double kappa_)
MixedScalarDerivativeIntegrator(Coefficient &q)
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
VectorFEMassIntegrator(MatrixCoefficient *_mq)
MixedDotProductIntegrator(VectorCoefficient &vq)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
MixedVectorCurlIntegrator(Coefficient &q)
MixedWeakCurlCrossIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual const char * FiniteElementTypeFailureMessage() const
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix2(const FiniteElement &nd_fe, const FiniteElement &rt_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
VectorMassIntegrator(MatrixCoefficient &q, int qo=0)
Construct an integrator with matrix coefficient q.
MixedScalarWeakCrossProductIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
GroupConvectionIntegrator(VectorCoefficient &q, double a=1.0)
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
MixedVectorGradientIntegrator(VectorCoefficient &dq)
static void AssembleBlock(const int dim, const int row_ndofs, const int col_ndofs, const int row_offset, const int col_offset, const double jmatcoef, const Vector &col_nL, const Vector &col_nM, const Vector &row_shape, const Vector &col_shape, const Vector &col_dshape_dnM, const DenseMatrix &col_dshape, DenseMatrix &elmat, DenseMatrix &jmat)
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
virtual ~InverseIntegrator()
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the BilinearFormIntegrator.
MixedVectorWeakDivergenceIntegrator(MatrixCoefficient &mq)
MixedScalarDivergenceIntegrator()
MixedVectorWeakDivergenceIntegrator()
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix2(const FiniteElement &rt_fe, const FiniteElement &l2_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
DGDiffusionIntegrator(MatrixCoefficient &q, const double s, const double k)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
Compute element energy: (1/2) (curl u, curl u)_E.
MixedCrossCurlIntegrator(VectorCoefficient &vq)
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual const char * FiniteElementTypeFailureMessage() const
MixedVectorProductIntegrator(VectorCoefficient &vq)
VectorScalarProductInterpolator(VectorCoefficient &vc)
virtual const char * FiniteElementTypeFailureMessage() const
ConvectionIntegrator(VectorCoefficient &q, double a=1.0)
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
MixedScalarCurlIntegrator(Coefficient &q)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
DerivativeIntegrator(Coefficient &q, int i)
DGTraceIntegrator(Coefficient &_rho, VectorCoefficient &_u, double a, double b)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
CurlCurlIntegrator(Coefficient &q)
Construct a bilinear form integrator for Nedelec elements.
virtual const char * FiniteElementTypeFailureMessage() const
MixedGradGradIntegrator(Coefficient &q)
MixedVectorGradientIntegrator(MatrixCoefficient &mq)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
MassIntegrator(const IntegrationRule *ir=NULL)
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
VectorFEWeakDivergenceIntegrator()
virtual double ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
Integrator for (Q u, v) for VectorFiniteElements.
MixedVectorWeakCurlIntegrator(Coefficient &q)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
VectorDivergenceIntegrator()
TransposeIntegrator(BilinearFormIntegrator *_bfi, int _own_bfi=1)
MixedDivGradIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
MixedScalarIntegrator(Coefficient &q)
DGElasticityIntegrator(Coefficient &lambda_, Coefficient &mu_, double alpha_, double kappa_)
Class for integrating (Q D_i(u), v); u and v are scalars.
InverseIntegrator(BilinearFormIntegrator *integ, int own_integ=1)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedScalarWeakDerivativeIntegrator(Coefficient &q)