MFEM
v4.1.0
Finite element discretization library
|
#include <eltrans.hpp>
Public Member Functions | |
void | SetFE (const FiniteElement *FE) |
const FiniteElement * | GetFE () const |
DenseMatrix & | GetPointMat () |
Read and write access to the underlying point matrix describing the transformation. More... | |
void | FinalizeTransformation () |
void | SetIdentityTransformation (Geometry::Type GeomType) |
virtual void | Transform (const IntegrationPoint &, Vector &) |
virtual void | Transform (const IntegrationRule &, DenseMatrix &) |
virtual void | Transform (const DenseMatrix &matrix, DenseMatrix &result) |
Transform columns of 'matrix', store result in 'result'. More... | |
virtual int | Order () |
virtual int | OrderJ () |
virtual int | OrderW () |
virtual int | OrderGrad (const FiniteElement *fe) |
Order of adj(J)^t.grad(fi) More... | |
virtual int | TransformBack (const Vector &v, IntegrationPoint &ip) |
Transform a point pt from physical space to a point ip in reference space. More... | |
virtual | ~IsoparametricTransformation () |
Public Member Functions inherited from mfem::ElementTransformation | |
ElementTransformation () | |
void | SetIntPoint (const IntegrationPoint *ip) |
const IntegrationPoint & | GetIntPoint () |
const DenseMatrix & | Jacobian () |
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint(). More... | |
const DenseMatrix & | Hessian () |
double | Weight () |
const DenseMatrix & | AdjugateJacobian () |
const DenseMatrix & | InverseJacobian () |
Geometry::Type | GetGeometryType () const |
Return the Geometry::Type of the reference element. More... | |
int | GetDimension () const |
Return the dimension of the reference element. More... | |
int | GetSpaceDim () const |
Get the dimension of the target (physical) space. More... | |
virtual | ~ElementTransformation () |
Additional Inherited Members | |
Public Attributes inherited from mfem::ElementTransformation | |
int | Attribute |
int | ElementNo |
Protected Types inherited from mfem::ElementTransformation | |
enum | StateMasks { JACOBIAN_MASK = 1, WEIGHT_MASK = 2, ADJUGATE_MASK = 4, INVERSE_MASK = 8, HESSIAN_MASK = 16 } |
Protected Member Functions inherited from mfem::ElementTransformation | |
double | EvalWeight () |
const DenseMatrix & | EvalAdjugateJ () |
const DenseMatrix & | EvalInverseJ () |
Protected Attributes inherited from mfem::ElementTransformation | |
const IntegrationPoint * | IntPoint |
DenseMatrix | dFdx |
DenseMatrix | adjJ |
DenseMatrix | invJ |
DenseMatrix | d2Fdx2 |
double | Wght |
int | EvalState |
Geometry::Type | geom |
int | space_dim |
Definition at line 291 of file eltrans.hpp.
|
inlinevirtual |
Definition at line 341 of file eltrans.hpp.
|
inline |
Definition at line 322 of file eltrans.hpp.
|
inline |
Definition at line 308 of file eltrans.hpp.
|
inline |
Read and write access to the underlying point matrix describing the transformation.
The dimensions of the matrix are space-dim x dof. The transformation is defined as
x=F(xh)=P.phi(xh),
where xh (x hat) is the reference point, x is the corresponding physical point, P is the point matrix, and phi(xh) is the column-vector of all basis functions evaluated at xh. The columns of P represent the control points in physical space defining the transformation.
Definition at line 321 of file eltrans.hpp.
|
inlinevirtual |
Implements mfem::ElementTransformation.
Definition at line 330 of file eltrans.hpp.
|
virtual |
Order of adj(J)^t.grad(fi)
Implements mfem::ElementTransformation.
Definition at line 464 of file eltrans.cpp.
|
virtual |
Implements mfem::ElementTransformation.
Definition at line 436 of file eltrans.cpp.
|
virtual |
Implements mfem::ElementTransformation.
Definition at line 450 of file eltrans.cpp.
|
inline |
Definition at line 307 of file eltrans.hpp.
void mfem::IsoparametricTransformation::SetIdentityTransformation | ( | Geometry::Type | GeomType | ) |
Definition at line 370 of file eltrans.cpp.
|
virtual |
Implements mfem::ElementTransformation.
Definition at line 483 of file eltrans.cpp.
|
virtual |
Implements mfem::ElementTransformation.
Definition at line 493 of file eltrans.cpp.
|
virtual |
Transform columns of 'matrix', store result in 'result'.
Implements mfem::ElementTransformation.
Definition at line 519 of file eltrans.cpp.
|
inlinevirtual |
Transform a point pt from physical space to a point ip in reference space.
Attempt to find the IntegrationPoint that is transformed into the given point in physical space. If the inversion fails a non-zero value is returned. This method is not 100 percent reliable for non-linear transformations.
Implements mfem::ElementTransformation.
Definition at line 335 of file eltrans.hpp.