MFEM v2.0
|
Container class for integration rules. More...
#include <intrules.hpp>
Public Member Functions | |
IntegrationRules (int refined=0) | |
Defines all integration rules. | |
const IntegrationRule & | Get (int GeomType, int Order) |
Returns an integration rule for given GeomType and Order. | |
void | Set (int GeomType, int Order, IntegrationRule &IntRule) |
void | SetOwnRules (int o) |
~IntegrationRules () | |
Destroys an IntegrationRules object. | |
Private Member Functions | |
void | PointIntegrationRules () |
void | SegmentIntegrationRules (int refined) |
void | TriangleIntegrationRules (int refined) |
void | SquareIntegrationRules () |
void | TetrahedronIntegrationRules (int refined) |
void | CubeIntegrationRules () |
Integration rules for reference cube. | |
void | DeleteIntRuleArray (Array< IntegrationRule * > &ir_array) |
Private Attributes | |
int | own_rules |
Array< IntegrationRule * > | PointIntRules |
Array< IntegrationRule * > | SegmentIntRules |
Array< IntegrationRule * > | TriangleIntRules |
Array< IntegrationRule * > | SquareIntRules |
Array< IntegrationRule * > | TetrahedronIntRules |
Array< IntegrationRule * > | CubeIntRules |
Container class for integration rules.
Definition at line 200 of file intrules.hpp.
IntegrationRules::IntegrationRules | ( | int | refined = 0 | ) | [explicit] |
Defines all integration rules.
Definition at line 175 of file intrules.cpp.
References CubeIntegrationRules(), own_rules, PointIntegrationRules(), SegmentIntegrationRules(), SquareIntegrationRules(), TetrahedronIntegrationRules(), and TriangleIntegrationRules().
IntegrationRules::~IntegrationRules | ( | ) |
Destroys an IntegrationRules object.
Definition at line 257 of file intrules.cpp.
References CubeIntRules, DeleteIntRuleArray(), own_rules, PointIntRules, SegmentIntRules, SquareIntRules, TetrahedronIntRules, and TriangleIntRules.
void IntegrationRules::CubeIntegrationRules | ( | ) | [private] |
Integration rules for reference cube.
Definition at line 792 of file intrules.cpp.
References CubeIntRules, SegmentIntRules, Array< T >::SetSize(), and Array< T >::Size().
Referenced by IntegrationRules().
void IntegrationRules::DeleteIntRuleArray | ( | Array< IntegrationRule * > & | ir_array | ) | [private] |
Definition at line 247 of file intrules.cpp.
References Array< T >::Size().
Referenced by ~IntegrationRules().
const IntegrationRule & IntegrationRules::Get | ( | int | GeomType, |
int | Order | ||
) |
Returns an integration rule for given GeomType and Order.
Definition at line 188 of file intrules.cpp.
References Geometry::CUBE, CubeIntRules, mfem_error(), Geometry::Name, Geometry::POINT, PointIntRules, Geometry::SEGMENT, SegmentIntRules, Array< T >::Size(), Geometry::SQUARE, SquareIntRules, Geometry::TETRAHEDRON, TetrahedronIntRules, Geometry::TRIANGLE, and TriangleIntRules.
Referenced by ElasticityIntegrator::AssembleElementMatrix(), VectorDiffusionIntegrator::AssembleElementMatrix(), DivDivIntegrator::AssembleElementMatrix(), VectorFEMassIntegrator::AssembleElementMatrix(), CurlCurlIntegrator::AssembleElementMatrix(), VectorMassIntegrator::AssembleElementMatrix(), ConvectionIntegrator::AssembleElementMatrix(), MassIntegrator::AssembleElementMatrix(), DiffusionIntegrator::AssembleElementMatrix(), VectorDivergenceIntegrator::AssembleElementMatrix2(), VectorFEMassIntegrator::AssembleElementMatrix2(), DerivativeIntegrator::AssembleElementMatrix2(), VectorFECurlIntegrator::AssembleElementMatrix2(), VectorFEDivergenceIntegrator::AssembleElementMatrix2(), MassIntegrator::AssembleElementMatrix2(), VectorFEBoundaryFluxLFIntegrator::AssembleRHSElementVect(), VectorBoundaryFluxLFIntegrator::AssembleRHSElementVect(), VectorFEDomainLFIntegrator::AssembleRHSElementVect(), VectorBoundaryLFIntegrator::AssembleRHSElementVect(), VectorDomainLFIntegrator::AssembleRHSElementVect(), BoundaryLFIntegrator::AssembleRHSElementVect(), DomainLFIntegrator::AssembleRHSElementVect(), DiffusionIntegrator::ComputeFluxEnergy(), GridFunction::ComputeH1Error(), GridFunction::ComputeL1Error(), GridFunction::ComputeL2Error(), GridFunction::ComputeMaxError(), GridFunction::ComputeW11Error(), Mesh::GetElementVolume(), and GridFunction::ProjectBdrCoefficientNormal().
void IntegrationRules::PointIntegrationRules | ( | ) | [private] |
Definition at line 270 of file intrules.cpp.
References PointIntRules, and Array< T >::SetSize().
Referenced by IntegrationRules().
void IntegrationRules::SegmentIntegrationRules | ( | int | refined | ) | [private] |
Definition at line 282 of file intrules.cpp.
References IntegrationRule::GaussianRule(), IntegrationRule::IntPoint(), SegmentIntRules, Array< T >::SetSize(), Array< T >::Size(), IntegrationPoint::weight, and IntegrationPoint::x.
Referenced by IntegrationRules().
void IntegrationRules::Set | ( | int | GeomType, |
int | Order, | ||
IntegrationRule & | IntRule | ||
) |
Definition at line 215 of file intrules.cpp.
References Geometry::CUBE, CubeIntRules, mfem_error(), Geometry::POINT, PointIntRules, Geometry::SEGMENT, SegmentIntRules, Array< T >::SetSize(), Array< T >::Size(), Geometry::SQUARE, SquareIntRules, Geometry::TETRAHEDRON, TetrahedronIntRules, Geometry::TRIANGLE, and TriangleIntRules.
void IntegrationRules::SetOwnRules | ( | int | o | ) | [inline] |
Definition at line 230 of file intrules.hpp.
References own_rules.
void IntegrationRules::SquareIntegrationRules | ( | ) | [private] |
Definition at line 700 of file intrules.cpp.
References SegmentIntRules, Array< T >::SetSize(), Array< T >::Size(), and SquareIntRules.
Referenced by IntegrationRules().
void IntegrationRules::TetrahedronIntegrationRules | ( | int | refined | ) | [private] |
Integration rules for reference tetrahedron {[0,0,0],[1,0,0],[0,1,0],[0,0,1]}
Definition at line 713 of file intrules.cpp.
References IntegrationRule::AddTetMidPoint(), IntegrationRule::AddTetPoints12(), IntegrationRule::AddTetPoints4(), IntegrationRule::AddTetPoints4b(), IntegrationRule::AddTetPoints6(), mfem_error(), Array< T >::SetSize(), Array< T >::Size(), and TetrahedronIntRules.
Referenced by IntegrationRules().
void IntegrationRules::TriangleIntegrationRules | ( | int | refined | ) | [private] |
Definition at line 339 of file intrules.cpp.
References IntegrationRule::AddTriMidPoint(), IntegrationRule::AddTriPoints3(), IntegrationRule::AddTriPoints3b(), IntegrationRule::AddTriPoints3R(), IntegrationRule::AddTriPoints6(), mfem_error(), Array< T >::SetSize(), and TriangleIntRules.
Referenced by IntegrationRules().
Array<IntegrationRule *> IntegrationRules::CubeIntRules [private] |
Definition at line 210 of file intrules.hpp.
Referenced by CubeIntegrationRules(), Get(), Set(), and ~IntegrationRules().
int IntegrationRules::own_rules [private] |
Definition at line 203 of file intrules.hpp.
Referenced by IntegrationRules(), SetOwnRules(), and ~IntegrationRules().
Array<IntegrationRule *> IntegrationRules::PointIntRules [private] |
Definition at line 205 of file intrules.hpp.
Referenced by Get(), PointIntegrationRules(), Set(), and ~IntegrationRules().
Array<IntegrationRule *> IntegrationRules::SegmentIntRules [private] |
Definition at line 206 of file intrules.hpp.
Referenced by CubeIntegrationRules(), Get(), SegmentIntegrationRules(), Set(), SquareIntegrationRules(), and ~IntegrationRules().
Array<IntegrationRule *> IntegrationRules::SquareIntRules [private] |
Definition at line 208 of file intrules.hpp.
Referenced by Get(), Set(), SquareIntegrationRules(), and ~IntegrationRules().
Array<IntegrationRule *> IntegrationRules::TetrahedronIntRules [private] |
Definition at line 209 of file intrules.hpp.
Referenced by Get(), Set(), TetrahedronIntegrationRules(), and ~IntegrationRules().
Array<IntegrationRule *> IntegrationRules::TriangleIntRules [private] |
Definition at line 207 of file intrules.hpp.
Referenced by Get(), Set(), TriangleIntegrationRules(), and ~IntegrationRules().