MFEM v2.0
|
#include <fe.hpp>
Public Member Functions | |
VectorFiniteElement (int D, int G, int Do, int O, int F=FunctionSpace::Pk) | |
Protected Member Functions | |
void | CalcVShape_RT (ElementTransformation &Trans, DenseMatrix &shape) const |
void | CalcVShape_ND (ElementTransformation &Trans, DenseMatrix &shape) const |
void | Project_RT (const double *nk, const Array< int > &d2n, VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const |
void | Project_RT (const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const |
void | ProjectGrad_RT (const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const |
void | ProjectCurl_RT (const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const |
void | Project_ND (const double *tk, const Array< int > &d2t, VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const |
void | Project_ND (const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const |
void | ProjectGrad_ND (const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const |
void | LocalInterpolation_RT (const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &I) const |
void | LocalInterpolation_ND (const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &I) const |
Protected Attributes | |
DenseMatrix | Jinv |
DenseMatrix | vshape |
Private Member Functions | |
virtual void | CalcShape (const IntegrationPoint &ip, Vector &shape) const |
Overrides the scalar CalcShape function to print an error. | |
virtual void | CalcDShape (const IntegrationPoint &ip, DenseMatrix &dshape) const |
Overrides the scalar CalcDShape function to print an error. |
VectorFiniteElement::VectorFiniteElement | ( | int | D, |
int | G, | ||
int | Do, | ||
int | O, | ||
int | F = FunctionSpace::Pk |
||
) | [inline] |
Definition at line 270 of file fe.hpp.
References FiniteElement::RangeType, and FiniteElement::VECTOR.
void VectorFiniteElement::CalcDShape | ( | const IntegrationPoint & | ip, |
DenseMatrix & | dshape | ||
) | const [private, virtual] |
Overrides the scalar CalcDShape function to print an error.
Implements FiniteElement.
Definition at line 250 of file fe.cpp.
References mfem_error().
void VectorFiniteElement::CalcShape | ( | const IntegrationPoint & | ip, |
Vector & | shape | ||
) | const [private, virtual] |
Overrides the scalar CalcShape function to print an error.
Implements FiniteElement.
Definition at line 243 of file fe.cpp.
References mfem_error().
void VectorFiniteElement::CalcVShape_ND | ( | ElementTransformation & | Trans, |
DenseMatrix & | shape | ||
) | const [protected] |
Definition at line 270 of file fe.cpp.
References CalcInverse(), FiniteElement::CalcVShape(), FiniteElement::Dim, FiniteElement::Dof, ElementTransformation::GetIntPoint(), DenseMatrix::Height(), ElementTransformation::Jacobian(), Jinv, Mult(), DenseMatrix::SetSize(), vshape, and DenseMatrix::Width().
Referenced by ND_TriangleElement::CalcVShape(), ND_TetrahedronElement::CalcVShape(), ND_QuadrilateralElement::CalcVShape(), ND_HexahedronElement::CalcVShape(), Nedelec1TetFiniteElement::CalcVShape(), and Nedelec1HexFiniteElement::CalcVShape().
void VectorFiniteElement::CalcVShape_RT | ( | ElementTransformation & | Trans, |
DenseMatrix & | shape | ||
) | const [protected] |
Definition at line 257 of file fe.cpp.
References FiniteElement::CalcVShape(), FiniteElement::Dim, FiniteElement::Dof, ElementTransformation::GetIntPoint(), ElementTransformation::Jacobian(), MultABt(), vshape, and ElementTransformation::Weight().
Referenced by RT_TetrahedronElement::CalcVShape(), RT_TriangleElement::CalcVShape(), RT_HexahedronElement::CalcVShape(), RT_QuadrilateralElement::CalcVShape(), RT0TetFiniteElement::CalcVShape(), RT1HexFiniteElement::CalcVShape(), RT0HexFiniteElement::CalcVShape(), RT2QuadFiniteElement::CalcVShape(), RT2TriangleFiniteElement::CalcVShape(), RT1QuadFiniteElement::CalcVShape(), RT1TriangleFiniteElement::CalcVShape(), RT0QuadFiniteElement::CalcVShape(), and RT0TriangleFiniteElement::CalcVShape().
void VectorFiniteElement::LocalInterpolation_ND | ( | const double * | tk, |
const Array< int > & | d2t, | ||
ElementTransformation & | Trans, | ||
DenseMatrix & | I | ||
) | const [protected] |
Definition at line 491 of file fe.cpp.
References FiniteElement::CalcVShape(), FiniteElement::Dim, FiniteElement::Dof, Geometries, FiniteElement::GeomType, Geometry::GetCenter(), IntegrationRule::IntPoint(), ElementTransformation::Jacobian(), DenseMatrix::Mult, FiniteElement::Nodes, IntegrationPoint::Set3(), ElementTransformation::SetIntPoint(), ElementTransformation::Transform(), and vshape.
Referenced by ND_TriangleElement::GetLocalInterpolation(), ND_TetrahedronElement::GetLocalInterpolation(), ND_QuadrilateralElement::GetLocalInterpolation(), and ND_HexahedronElement::GetLocalInterpolation().
void VectorFiniteElement::LocalInterpolation_RT | ( | const double * | nk, |
const Array< int > & | d2n, | ||
ElementTransformation & | Trans, | ||
DenseMatrix & | I | ||
) | const [protected] |
Definition at line 457 of file fe.cpp.
References CalcAdjugateTranspose(), FiniteElement::CalcVShape(), FiniteElement::Dim, FiniteElement::Dof, Geometries, FiniteElement::GeomType, Geometry::GetCenter(), IntegrationRule::IntPoint(), ElementTransformation::Jacobian(), Jinv, DenseMatrix::Mult, FiniteElement::Nodes, IntegrationPoint::Set3(), ElementTransformation::SetIntPoint(), ElementTransformation::Transform(), and vshape.
Referenced by RT_TetrahedronElement::GetLocalInterpolation(), RT_TriangleElement::GetLocalInterpolation(), RT_HexahedronElement::GetLocalInterpolation(), and RT_QuadrilateralElement::GetLocalInterpolation().
void VectorFiniteElement::Project_ND | ( | const double * | tk, |
const Array< int > & | d2t, | ||
VectorCoefficient & | vc, | ||
ElementTransformation & | Trans, | ||
Vector & | dofs | ||
) | const [protected] |
Definition at line 389 of file fe.cpp.
References FiniteElement::Dim, FiniteElement::Dof, VectorCoefficient::Eval(), VectorCoefficient::GetVDim(), DenseMatrix::InnerProduct(), IntegrationRule::IntPoint(), ElementTransformation::Jacobian(), FiniteElement::Nodes, and ElementTransformation::SetIntPoint().
Referenced by ND_TriangleElement::Project(), ND_TetrahedronElement::Project(), ND_QuadrilateralElement::Project(), and ND_HexahedronElement::Project().
void VectorFiniteElement::Project_ND | ( | const double * | tk, |
const Array< int > & | d2t, | ||
const FiniteElement & | fe, | ||
ElementTransformation & | Trans, | ||
DenseMatrix & | I | ||
) | const [protected] |
Definition at line 406 of file fe.cpp.
References FiniteElement::CalcShape(), FiniteElement::Dim, FiniteElement::Dof, FiniteElement::GetDof(), FiniteElement::GetRangeType(), IntegrationRule::IntPoint(), ElementTransformation::Jacobian(), mfem_error(), DenseMatrix::Mult, FiniteElement::Nodes, FiniteElement::SCALAR, ElementTransformation::SetIntPoint(), and DenseMatrix::SetSize().
void VectorFiniteElement::Project_RT | ( | const double * | nk, |
const Array< int > & | d2n, | ||
const FiniteElement & | fe, | ||
ElementTransformation & | Trans, | ||
DenseMatrix & | I | ||
) | const [protected] |
Definition at line 311 of file fe.cpp.
References CalcAdjugateTranspose(), FiniteElement::CalcShape(), FiniteElement::Dim, FiniteElement::Dof, FiniteElement::GetDof(), FiniteElement::GetRangeType(), IntegrationRule::IntPoint(), ElementTransformation::Jacobian(), Jinv, mfem_error(), FiniteElement::Nodes, FiniteElement::SCALAR, ElementTransformation::SetIntPoint(), and DenseMatrix::SetSize().
void VectorFiniteElement::Project_RT | ( | const double * | nk, |
const Array< int > & | d2n, | ||
VectorCoefficient & | vc, | ||
ElementTransformation & | Trans, | ||
Vector & | dofs | ||
) | const [protected] |
Definition at line 289 of file fe.cpp.
References CalcAdjugate(), FiniteElement::Dim, FiniteElement::Dof, VectorCoefficient::Eval(), DenseMatrix::InnerProduct(), IntegrationRule::IntPoint(), ElementTransformation::Jacobian(), Jinv, FiniteElement::Nodes, and ElementTransformation::SetIntPoint().
Referenced by RT_TetrahedronElement::Project(), RT_TriangleElement::Project(), RT_HexahedronElement::Project(), and RT_QuadrilateralElement::Project().
void VectorFiniteElement::ProjectCurl_RT | ( | const double * | nk, |
const Array< int > & | d2n, | ||
const FiniteElement & | fe, | ||
ElementTransformation & | Trans, | ||
DenseMatrix & | curl | ||
) | const [protected] |
Definition at line 372 of file fe.cpp.
References FiniteElement::CalcCurlShape(), FiniteElement::Dim, FiniteElement::Dof, FiniteElement::GetDof(), IntegrationRule::IntPoint(), FiniteElement::Nodes, and DenseMatrix::SetSize().
Referenced by RT_TetrahedronElement::ProjectCurl(), and RT_HexahedronElement::ProjectCurl().
void VectorFiniteElement::ProjectGrad_ND | ( | const double * | tk, |
const Array< int > & | d2t, | ||
const FiniteElement & | fe, | ||
ElementTransformation & | Trans, | ||
DenseMatrix & | grad | ||
) | const [protected] |
Definition at line 440 of file fe.cpp.
References FiniteElement::CalcDShape(), FiniteElement::Dim, FiniteElement::Dof, FiniteElement::GetDim(), FiniteElement::GetDof(), IntegrationRule::IntPoint(), FiniteElement::Nodes, and DenseMatrix::SetSize().
Referenced by ND_TriangleElement::ProjectGrad(), ND_TetrahedronElement::ProjectGrad(), ND_QuadrilateralElement::ProjectGrad(), and ND_HexahedronElement::ProjectGrad().
void VectorFiniteElement::ProjectGrad_RT | ( | const double * | nk, |
const Array< int > & | d2n, | ||
const FiniteElement & | fe, | ||
ElementTransformation & | Trans, | ||
DenseMatrix & | grad | ||
) | const [protected] |
Definition at line 349 of file fe.cpp.
References FiniteElement::CalcDShape(), FiniteElement::Dim, FiniteElement::Dof, FiniteElement::GetDim(), FiniteElement::GetDof(), IntegrationRule::IntPoint(), mfem_error(), FiniteElement::Nodes, and DenseMatrix::SetSize().
Referenced by RT_TriangleElement::ProjectGrad(), and RT_QuadrilateralElement::ProjectGrad().
DenseMatrix VectorFiniteElement::Jinv [mutable, protected] |
Definition at line 222 of file fe.hpp.
Referenced by CalcVShape_ND(), RT0TetFiniteElement::GetLocalInterpolation(), RT1HexFiniteElement::GetLocalInterpolation(), RT0HexFiniteElement::GetLocalInterpolation(), RT2QuadFiniteElement::GetLocalInterpolation(), RT1QuadFiniteElement::GetLocalInterpolation(), RT1TriangleFiniteElement::GetLocalInterpolation(), RT0QuadFiniteElement::GetLocalInterpolation(), RT0TriangleFiniteElement::GetLocalInterpolation(), LocalInterpolation_RT(), RT0TetFiniteElement::Project(), RT1HexFiniteElement::Project(), RT0HexFiniteElement::Project(), RT2QuadFiniteElement::Project(), RT1QuadFiniteElement::Project(), RT1TriangleFiniteElement::Project(), RT0QuadFiniteElement::Project(), RT0TriangleFiniteElement::Project(), and Project_RT().
DenseMatrix VectorFiniteElement::vshape [mutable, protected] |
Definition at line 223 of file fe.hpp.
Referenced by CalcVShape_ND(), CalcVShape_RT(), RT0TetFiniteElement::GetLocalInterpolation(), RT1HexFiniteElement::GetLocalInterpolation(), RT0HexFiniteElement::GetLocalInterpolation(), Nedelec1TetFiniteElement::GetLocalInterpolation(), Nedelec1HexFiniteElement::GetLocalInterpolation(), RT2QuadFiniteElement::GetLocalInterpolation(), RT1QuadFiniteElement::GetLocalInterpolation(), RT1TriangleFiniteElement::GetLocalInterpolation(), RT0QuadFiniteElement::GetLocalInterpolation(), RT0TriangleFiniteElement::GetLocalInterpolation(), LocalInterpolation_ND(), and LocalInterpolation_RT().