Public Member Functions | |
| RefIntegral (int spatialDim, const CellType &maxCellType, int dim, const CellType &cellType, const QuadratureFamily &quad_in, bool isInternalBdry, const ParametrizedCurve &globalCurve, const Mesh &mesh, int verb) | |
| RefIntegral (int spatialDim, const CellType &maxCellType, int dim, const CellType &cellType, const BasisFamily &testBasis, int alpha, int testDerivOrder, const QuadratureFamily &quad_in, bool isInternalBdry, const ParametrizedCurve &globalCurve, const Mesh &mesh, int verb) | |
| RefIntegral (int spatialDim, const CellType &maxCellType, int dim, const CellType &cellType, const BasisFamily &testBasis, int alpha, int testDerivOrder, const BasisFamily &unkBasis, int beta, int unkDerivOrder, const QuadratureFamily &quad_in, bool isInternalBdry, const ParametrizedCurve &globalCurve, const Mesh &mesh, int verb) | |
| virtual | ~RefIntegral () |
| void | print (std::ostream &os) const |
| void | transform (const CellJacobianBatch &JTrans, const CellJacobianBatch &JVol, const Array< int > &isLocalFlag, const Array< int > &facetNum, const RCP< Array< int > > &cellLIDs, const double &coeff, RCP< Array< double > > &A) const |
| void | transformTwoForm (const CellJacobianBatch &JTrans, const CellJacobianBatch &JVol, const Array< int > &facetNum, const RCP< Array< int > > &cellLIDs, const double &coeff, RCP< Array< double > > &A) const |
| void | transformOneForm (const CellJacobianBatch &JTrans, const CellJacobianBatch &JVol, const Array< int > &facetNum, const RCP< Array< int > > &cellLIDs, const double &coeff, RCP< Array< double > > &A) const |
| void | transformZeroForm (const CellJacobianBatch &JVol, const Array< int > &isLocalFlag, const RCP< Array< int > > &cellLIDs, const double &coeff, RCP< Array< double > > &A) const |
| double & | value (int facetCase, int testDerivDir, int testNode, int unkDerivDir, int unkNode) |
| const double & | value (int facetCase, int testDerivDir, int testNode, int unkDerivDir, int unkNode) const |
| double & | value (int facetCase, int testDerivDir, int testNode) |
| const double & | value (int facetCase, int testDerivDir, int testNode) const |
Static Public Member Functions | |
| static double & | totalFlops () |
Static Protected Member Functions | |
| static void | addFlops (const double &flops) |
Private Attributes | |
| Array< Array< double > > | W_ |
| Array< Array< Array< Array < double > > > > | W_ACI_F1_ |
| Array< Array< Array< Array < Array< Array< double > > > > > > | W_ACI_F2_ |
| QuadratureFamily | quad_ |
| Array< Point > | quadPts_ |
| Array< double > | quadWeights_ |
RefIntegral represents the integrals of a product of basis functions (or their derivatives) over a reference cell. This object can be created once and for all up front, and then reused to produce the integrals of constant-coefficient weak forms over simplicial cells by means of a linear transformation.
An instance of this object can represent either a one form,
or a two form,
In the notation for the two form, we have grouped together the index pairs
and
to indicate that we will think of this 4-tensor as a 2D array, with
combinations running down rows and
running across columns. We will consider the node number combination
to be a single index
, and
to be a single index
. Thus, the storage and use of one-forms and two-forms is essentially identical.
This object's job in life is to be multiplied with a transformation matrix
to produce an array of local element matrices
,
The index
is over cells (cells are processed in batches).
Physical storage is as a 1D vector stored in colum-major order.
Definition at line 88 of file SundanceRefIntegral.hpp.
| RefIntegral::RefIntegral | ( | int | spatialDim, |
| const CellType & | maxCellType, | ||
| int | dim, | ||
| const CellType & | cellType, | ||
| const QuadratureFamily & | quad_in, | ||
| bool | isInternalBdry, | ||
| const ParametrizedCurve & | globalCurve, | ||
| const Mesh & | mesh, | ||
| int | verb | ||
| ) |
Construct a reference zero-form
Definition at line 92 of file SundanceRefIntegral.cpp.
References Sundance::ElementIntegral::describe(), Sundance::QuadratureFamily::getPoints(), Sundance::QuadratureFamily::order(), Playa::Out::os(), quad_, quadWeights_, Sundance::ElementIntegral::setupVerb(), SUNDANCE_MSG1, and W_.
| RefIntegral::RefIntegral | ( | int | spatialDim, |
| const CellType & | maxCellType, | ||
| int | dim, | ||
| const CellType & | cellType, | ||
| const BasisFamily & | testBasis, | ||
| int | alpha, | ||
| int | testDerivOrder, | ||
| const QuadratureFamily & | quad_in, | ||
| bool | isInternalBdry, | ||
| const ParametrizedCurve & | globalCurve, | ||
| const Mesh & | mesh, | ||
| int | verb | ||
| ) |
Construct a reference one-form
Definition at line 137 of file SundanceRefIntegral.cpp.
References addFlops(), Sundance::ElementIntegral::assertLinearForm(), Sundance::ElementIntegral::chop(), Sundance::QuadratureType::createQuadFamily(), Sundance::ElementIntegral::describe(), Sundance::BasisFamily::dim(), Sundance::ElementIntegral::evalCellType(), Sundance::QuadratureType::findValidOrder(), Sundance::ElementIntegral::getQuad(), Playa::max(), Sundance::ElementIntegral::nFacetCases(), Sundance::ElementIntegral::nNodesTest(), Sundance::ElementIntegral::nRefDerivTest(), Sundance::QuadratureFamily::order(), Sundance::BasisFamily::order(), Playa::Out::os(), quad_, quadPts_, quadWeights_, Sundance::BasisFamily::refEval(), Sundance::ElementIntegral::setupVerb(), SUNDANCE_MSG1, SUNDANCE_MSG2, SUNDANCE_MSG4, value(), W_, and W_ACI_F1_.
| RefIntegral::RefIntegral | ( | int | spatialDim, |
| const CellType & | maxCellType, | ||
| int | dim, | ||
| const CellType & | cellType, | ||
| const BasisFamily & | testBasis, | ||
| int | alpha, | ||
| int | testDerivOrder, | ||
| const BasisFamily & | unkBasis, | ||
| int | beta, | ||
| int | unkDerivOrder, | ||
| const QuadratureFamily & | quad_in, | ||
| bool | isInternalBdry, | ||
| const ParametrizedCurve & | globalCurve, | ||
| const Mesh & | mesh, | ||
| int | verb | ||
| ) |
Construct a reference two-form
Definition at line 282 of file SundanceRefIntegral.cpp.
References addFlops(), Sundance::ElementIntegral::assertBilinearForm(), Sundance::ElementIntegral::chop(), Sundance::QuadratureType::createQuadFamily(), Sundance::ElementIntegral::describe(), Sundance::BasisFamily::dim(), Sundance::ElementIntegral::evalCellType(), Sundance::QuadratureType::findValidOrder(), Sundance::ElementIntegral::getQuad(), Playa::max(), Sundance::ElementIntegral::nFacetCases(), Sundance::ElementIntegral::nNodesTest(), Sundance::ElementIntegral::nNodesUnk(), Sundance::ElementIntegral::nRefDerivTest(), Sundance::ElementIntegral::nRefDerivUnk(), Sundance::QuadratureFamily::order(), Sundance::BasisFamily::order(), Playa::Out::os(), quad_, quadPts_, quadWeights_, Sundance::BasisFamily::refEval(), Sundance::ElementIntegral::setupVerb(), SUNDANCE_MSG1, SUNDANCE_MSG2, SUNDANCE_MSG4, value(), W_, and W_ACI_F2_.
| virtual Sundance::RefIntegral::~RefIntegral | ( | ) | [inline, virtual] |
virtual dtor
Definition at line 134 of file SundanceRefIntegral.hpp.
| static void Sundance::RefIntegral::addFlops | ( | const double & | flops | ) | [inline, static, protected] |
Reimplemented from Sundance::ElementIntegral.
Definition at line 206 of file SundanceRefIntegral.hpp.
References totalFlops().
Referenced by RefIntegral(), transformOneForm(), transformTwoForm(), and transformZeroForm().
| void Sundance::RefIntegral::print | ( | std::ostream & | os | ) | const |
| static double& Sundance::RefIntegral::totalFlops | ( | ) | [inline, static] |
Reimplemented from Sundance::ElementIntegral.
Definition at line 202 of file SundanceRefIntegral.hpp.
Referenced by addFlops().
| void Sundance::RefIntegral::transform | ( | const CellJacobianBatch & | JTrans, |
| const CellJacobianBatch & | JVol, | ||
| const Array< int > & | isLocalFlag, | ||
| const Array< int > & | facetNum, | ||
| const RCP< Array< int > > & | cellLIDs, | ||
| const double & | coeff, | ||
| RCP< Array< double > > & | A | ||
| ) | const [inline] |
Definition at line 140 of file SundanceRefIntegral.hpp.
References Sundance::ElementIntegral::order(), transformOneForm(), transformTwoForm(), and transformZeroForm().
Referenced by Sundance::IntegralGroup::evaluate().
| void RefIntegral::transformOneForm | ( | const CellJacobianBatch & | JTrans, |
| const CellJacobianBatch & | JVol, | ||
| const Array< int > & | facetNum, | ||
| const RCP< Array< int > > & | cellLIDs, | ||
| const double & | coeff, | ||
| RCP< Array< double > > & | A | ||
| ) | const |
Definition at line 573 of file SundanceRefIntegral.cpp.
References addFlops(), Sundance::ElementIntegral::alpha(), Sundance::ElementIntegral::cellType(), Sundance::ElementIntegral::chop(), Sundance::ElementIntegral::createOneFormTransformationMatrix(), Sundance::CellJacobianBatch::detJ(), dgemm_(), Sundance::ElementIntegral::dim(), Sundance::ElementIntegral::G(), Sundance::QuadratureFamily::getAdaptedWeights(), Sundance::ElementIntegral::globalCurve(), Sundance::ElementIntegral::integrationVerb(), Sundance::ElementIntegral::mesh(), Sundance::ElementIntegral::nFacetCases(), Sundance::ElementIntegral::nNodes(), Sundance::ElementIntegral::nNodesTest(), Sundance::ElementIntegral::nRefDerivTest(), Sundance::CellJacobianBatch::numCells(), Sundance::ElementIntegral::order(), quad_, quadPts_, quadWeights_, ref1IntegrationTimer(), SUNDANCE_MSG1, SUNDANCE_MSG3, Sundance::ElementIntegral::testDerivOrder(), Sundance::ElementIntegral::transformVerb(), W_, and W_ACI_F1_.
Referenced by transform().
| void RefIntegral::transformTwoForm | ( | const CellJacobianBatch & | JTrans, |
| const CellJacobianBatch & | JVol, | ||
| const Array< int > & | facetNum, | ||
| const RCP< Array< int > > & | cellLIDs, | ||
| const double & | coeff, | ||
| RCP< Array< double > > & | A | ||
| ) | const |
Definition at line 784 of file SundanceRefIntegral.cpp.
References addFlops(), Sundance::ElementIntegral::alpha(), Sundance::ElementIntegral::beta(), Sundance::ElementIntegral::cellType(), Sundance::ElementIntegral::chop(), Sundance::ElementIntegral::createTwoFormTransformationMatrix(), Sundance::CellJacobianBatch::detJ(), dgemm_(), Sundance::ElementIntegral::dim(), Sundance::ElementIntegral::G(), Sundance::QuadratureFamily::getAdaptedWeights(), Sundance::ElementIntegral::globalCurve(), Sundance::ElementIntegral::integrationVerb(), Sundance::ElementIntegral::mesh(), Sundance::ElementIntegral::nFacetCases(), Sundance::ElementIntegral::nNodes(), Sundance::ElementIntegral::nNodesTest(), Sundance::ElementIntegral::nNodesUnk(), Sundance::ElementIntegral::nRefDerivTest(), Sundance::ElementIntegral::nRefDerivUnk(), Sundance::CellJacobianBatch::numCells(), Sundance::ElementIntegral::order(), quad_, quadPts_, quadWeights_, ref2IntegrationTimer(), SUNDANCE_MSG1, SUNDANCE_MSG2, Sundance::ElementIntegral::testDerivOrder(), Sundance::ElementIntegral::transformVerb(), Sundance::ElementIntegral::unkDerivOrder(), W_, and W_ACI_F2_.
Referenced by transform().
| void RefIntegral::transformZeroForm | ( | const CellJacobianBatch & | JVol, |
| const Array< int > & | isLocalFlag, | ||
| const RCP< Array< int > > & | cellLIDs, | ||
| const double & | coeff, | ||
| RCP< Array< double > > & | A | ||
| ) | const |
Definition at line 458 of file SundanceRefIntegral.cpp.
References addFlops(), Sundance::ElementIntegral::cellType(), Sundance::ElementIntegral::chop(), Sundance::CellJacobianBatch::detJ(), Sundance::ElementIntegral::dim(), Sundance::QuadratureFamily::getAdaptedWeights(), Sundance::ElementIntegral::globalCurve(), Sundance::ElementIntegral::integrationVerb(), Sundance::ElementIntegral::mesh(), Sundance::CellJacobianBatch::numCells(), Sundance::ElementIntegral::order(), quad_, quadPts_, quadWeights_, ref0IntegrationTimer(), SUNDANCE_MSG1, and W_.
Referenced by transform().
| double& Sundance::RefIntegral::value | ( | int | facetCase, |
| int | testDerivDir, | ||
| int | testNode, | ||
| int | unkDerivDir, | ||
| int | unkNode | ||
| ) | [inline] |
Definition at line 177 of file SundanceRefIntegral.hpp.
References Sundance::ElementIntegral::nNodes(), Sundance::ElementIntegral::nNodesUnk(), Sundance::ElementIntegral::nRefDerivUnk(), and W_.
Referenced by RefIntegral().
| const double& Sundance::RefIntegral::value | ( | int | facetCase, |
| int | testDerivDir, | ||
| int | testNode, | ||
| int | unkDerivDir, | ||
| int | unkNode | ||
| ) | const [inline] |
Definition at line 184 of file SundanceRefIntegral.hpp.
References Sundance::ElementIntegral::nNodes(), Sundance::ElementIntegral::nNodesUnk(), Sundance::ElementIntegral::nRefDerivUnk(), and W_.
| double& Sundance::RefIntegral::value | ( | int | facetCase, |
| int | testDerivDir, | ||
| int | testNode | ||
| ) | [inline] |
Definition at line 194 of file SundanceRefIntegral.hpp.
References Sundance::ElementIntegral::nNodesTest(), and W_.
| const double& Sundance::RefIntegral::value | ( | int | facetCase, |
| int | testDerivDir, | ||
| int | testNode | ||
| ) | const [inline] |
Definition at line 198 of file SundanceRefIntegral.hpp.
References Sundance::ElementIntegral::nNodesTest(), and W_.
QuadratureFamily Sundance::RefIntegral::quad_ [private] |
The quadrature family needed for special integration
Definition at line 221 of file SundanceRefIntegral.hpp.
Referenced by RefIntegral(), transformOneForm(), transformTwoForm(), and transformZeroForm().
Array<Point> Sundance::RefIntegral::quadPts_ [private] |
The quadrature points
Definition at line 224 of file SundanceRefIntegral.hpp.
Referenced by RefIntegral(), transformOneForm(), transformTwoForm(), and transformZeroForm().
Array<double> Sundance::RefIntegral::quadWeights_ [private] |
The standard weights (in case of ACI these might change)
Definition at line 227 of file SundanceRefIntegral.hpp.
Referenced by RefIntegral(), transformOneForm(), transformTwoForm(), and transformZeroForm().
Array<Array<double> > Sundance::RefIntegral::W_ [private] |
Definition at line 210 of file SundanceRefIntegral.hpp.
Referenced by RefIntegral(), transformOneForm(), transformTwoForm(), transformZeroForm(), and value().
Array<Array<Array<Array<double> > > > Sundance::RefIntegral::W_ACI_F1_ [private] |
For ACI (ACI = Adaptive Cell Integration), store the reference integral values for one form The indexes facet, quadPoints, nRefDerivTest , nNodesTest
Definition at line 214 of file SundanceRefIntegral.hpp.
Referenced by RefIntegral(), and transformOneForm().
Array<Array<Array<Array<Array<Array<double> > > > > > Sundance::RefIntegral::W_ACI_F2_ [private] |
For ACI (ACI = Adaptive Cell Integration), store the reference integral values for two form The indexes facet, quadPoints, nRefDerivTest , nNodesTest , nRefDerivUnk , nNodesUnk
Definition at line 218 of file SundanceRefIntegral.hpp.
Referenced by RefIntegral(), and transformTwoForm().