Public Member Functions | |
ReducedIntegral (int spatialDim, const CellType &maxCellType, int dim, const CellType &cellType, const QuadratureFamily &quad, bool isInternalBdry, const ParametrizedCurve &globalCurve, const Mesh &mesh, int verb) | |
ReducedIntegral (int spatialDim, const CellType &maxCellType, int dim, const CellType &cellType, const BasisFamily &testBasis, int alpha, int testDerivOrder, const QuadratureFamily &quad, bool isInternalBdry, const ParametrizedCurve &globalCurve, const Mesh &mesh, int verb) | |
ReducedIntegral (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, bool isInternalBdry, const ParametrizedCurve &globalCurve, const Mesh &mesh, int verb) | |
virtual | ~ReducedIntegral () |
void | transform (const CellJacobianBatch &JTrans, const CellJacobianBatch &JVol, const Array< int > &isLocalFlag, const Array< int > &facetNum, const RCP< Array< int > > &cellLIDs, const double *const coeffs, RCP< Array< double > > &A) const |
virtual void | transformZeroForm (const CellJacobianBatch &JTrans, const CellJacobianBatch &JVol, const Array< int > &isLocalFlag, const Array< int > &facetIndex, const RCP< Array< int > > &cellLIDs, const double *const coeffs, RCP< Array< double > > &A) const |
virtual void | transformTwoForm (const CellJacobianBatch &JTrans, const CellJacobianBatch &JVol, const Array< int > &facetIndex, const RCP< Array< int > > &cellLIDs, const double *const coeffs, RCP< Array< double > > &A) const |
void | transformOneForm (const CellJacobianBatch &JTrans, const CellJacobianBatch &JVol, const Array< int > &facetIndex, const RCP< Array< int > > &cellLIDs, const double *const coeffs, RCP< Array< double > > &A) const |
Static Protected Member Functions | |
static void | addFlops (const double &flops) |
Private Member Functions | |
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 Private Member Functions | |
static double & | totalFlops () |
Private Attributes | |
Array< Array< double > > | W_ |
Definition at line 45 of file SundanceReducedIntegral.hpp.
ReducedIntegral::ReducedIntegral | ( | int | spatialDim, | |
const CellType & | maxCellType, | |||
int | dim, | |||
const CellType & | cellType, | |||
const QuadratureFamily & | quad, | |||
bool | isInternalBdry, | |||
const ParametrizedCurve & | globalCurve, | |||
const Mesh & | mesh, | |||
int | verb | |||
) |
Construct a zero-form to be computed by reference integration with coefficients that are piecewise constant
Definition at line 83 of file SundanceReducedIntegral.cpp.
References Sundance::ElementIntegral::describe(), Sundance::QuadratureFamily::getPoints(), Playa::Out::os(), Sundance::ElementIntegral::setupVerb(), SUNDANCE_MSG1, and W_.
ReducedIntegral::ReducedIntegral | ( | int | spatialDim, | |
const CellType & | maxCellType, | |||
int | dim, | |||
const CellType & | cellType, | |||
const BasisFamily & | testBasis, | |||
int | alpha, | |||
int | testDerivOrder, | |||
const QuadratureFamily & | quad, | |||
bool | isInternalBdry, | |||
const ParametrizedCurve & | globalCurve, | |||
const Mesh & | mesh, | |||
int | verb | |||
) |
Construct a one form to be computed by reference integration with coefficients that are piecewise constant
Definition at line 121 of file SundanceReducedIntegral.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::BasisFamily::order(), Playa::Out::os(), Sundance::BasisFamily::refEval(), Sundance::ElementIntegral::setupVerb(), SUNDANCE_MSG1, SUNDANCE_MSG2, SUNDANCE_MSG4, value(), and W_.
ReducedIntegral::ReducedIntegral | ( | 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, | |||
bool | isInternalBdry, | |||
const ParametrizedCurve & | globalCurve, | |||
const Mesh & | mesh, | |||
int | verb | |||
) |
Construct a two-form to be computed by reference integration with coefficients that are piecewise constant
Definition at line 257 of file SundanceReducedIntegral.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::BasisFamily::order(), Playa::Out::os(), Sundance::BasisFamily::refEval(), Sundance::ElementIntegral::setupVerb(), SUNDANCE_MSG1, SUNDANCE_MSG2, SUNDANCE_MSG4, value(), and W_.
virtual Sundance::ReducedIntegral::~ReducedIntegral | ( | ) | [inline, virtual] |
virtual dtor
Definition at line 94 of file SundanceReducedIntegral.hpp.
static void Sundance::ReducedIntegral::addFlops | ( | const double & | flops | ) | [inline, static, protected] |
Reimplemented from Sundance::ElementIntegral.
Definition at line 170 of file SundanceReducedIntegral.hpp.
References totalFlops().
Referenced by ReducedIntegral(), transformOneForm(), transformTwoForm(), and transformZeroForm().
static double& Sundance::ReducedIntegral::totalFlops | ( | ) | [inline, static, private] |
Reimplemented from Sundance::ElementIntegral.
Definition at line 166 of file SundanceReducedIntegral.hpp.
Referenced by addFlops().
void Sundance::ReducedIntegral::transform | ( | const CellJacobianBatch & | JTrans, | |
const CellJacobianBatch & | JVol, | |||
const Array< int > & | isLocalFlag, | |||
const Array< int > & | facetNum, | |||
const RCP< Array< int > > & | cellLIDs, | |||
const double *const | coeffs, | |||
RCP< Array< double > > & | A | |||
) | const [inline] |
Definition at line 97 of file SundanceReducedIntegral.hpp.
References Sundance::ElementIntegral::order(), transformOneForm(), transformTwoForm(), and transformZeroForm().
Referenced by Sundance::IntegralGroup::evaluate().
void ReducedIntegral::transformOneForm | ( | const CellJacobianBatch & | JTrans, | |
const CellJacobianBatch & | JVol, | |||
const Array< int > & | facetIndex, | |||
const RCP< Array< int > > & | cellLIDs, | |||
const double *const | coeffs, | |||
RCP< Array< double > > & | A | |||
) | const |
Definition at line 491 of file SundanceReducedIntegral.cpp.
References addFlops(), Sundance::ElementIntegral::alpha(), Sundance::ElementIntegral::createOneFormTransformationMatrix(), Sundance::CellJacobianBatch::detJ(), dgemm_(), Sundance::ElementIntegral::G(), Sundance::ElementIntegral::globalCurve(), Sundance::ElementIntegral::integrationVerb(), Sundance::ElementIntegral::nFacetCases(), Sundance::ElementIntegral::nNodes(), Sundance::ElementIntegral::nRefDerivTest(), Sundance::CellJacobianBatch::numCells(), Sundance::ElementIntegral::order(), reduced1IntegrationTimer(), SUNDANCE_MSG1, SUNDANCE_MSG3, Sundance::ElementIntegral::testDerivOrder(), Sundance::ElementIntegral::transformVerb(), and W_.
Referenced by transform().
void ReducedIntegral::transformTwoForm | ( | const CellJacobianBatch & | JTrans, | |
const CellJacobianBatch & | JVol, | |||
const Array< int > & | facetIndex, | |||
const RCP< Array< int > > & | cellLIDs, | |||
const double *const | coeffs, | |||
RCP< Array< double > > & | A | |||
) | const [virtual] |
Definition at line 599 of file SundanceReducedIntegral.cpp.
References addFlops(), Sundance::ElementIntegral::alpha(), Sundance::ElementIntegral::beta(), Sundance::ElementIntegral::createTwoFormTransformationMatrix(), Sundance::CellJacobianBatch::detJ(), dgemm_(), Sundance::ElementIntegral::G(), Sundance::ElementIntegral::globalCurve(), Sundance::ElementIntegral::integrationVerb(), Sundance::ElementIntegral::nFacetCases(), Sundance::ElementIntegral::nNodes(), Sundance::ElementIntegral::nRefDerivTest(), Sundance::ElementIntegral::nRefDerivUnk(), Sundance::CellJacobianBatch::numCells(), Sundance::ElementIntegral::order(), reduced2IntegrationTimer(), SUNDANCE_MSG1, SUNDANCE_MSG2, Sundance::ElementIntegral::testDerivOrder(), Sundance::ElementIntegral::transformVerb(), Sundance::ElementIntegral::unkDerivOrder(), and W_.
Referenced by transform().
void ReducedIntegral::transformZeroForm | ( | const CellJacobianBatch & | JTrans, | |
const CellJacobianBatch & | JVol, | |||
const Array< int > & | isLocalFlag, | |||
const Array< int > & | facetIndex, | |||
const RCP< Array< int > > & | cellLIDs, | |||
const double *const | coeffs, | |||
RCP< Array< double > > & | A | |||
) | const [virtual] |
Definition at line 423 of file SundanceReducedIntegral.cpp.
References addFlops(), Sundance::CellJacobianBatch::detJ(), Sundance::ElementIntegral::globalCurve(), Sundance::ElementIntegral::integrationVerb(), Sundance::CellJacobianBatch::numCells(), Sundance::ElementIntegral::order(), reduced0IntegrationTimer(), SUNDANCE_MSG1, and W_.
Referenced by transform().
const double& Sundance::ReducedIntegral::value | ( | int | facetCase, | |
int | testDerivDir, | |||
int | testNode | |||
) | const [inline, private] |
Definition at line 162 of file SundanceReducedIntegral.hpp.
References Sundance::ElementIntegral::nNodesTest(), and W_.
double& Sundance::ReducedIntegral::value | ( | int | facetCase, | |
int | testDerivDir, | |||
int | testNode | |||
) | [inline, private] |
Definition at line 158 of file SundanceReducedIntegral.hpp.
References Sundance::ElementIntegral::nNodesTest(), and W_.
const double& Sundance::ReducedIntegral::value | ( | int | facetCase, | |
int | testDerivDir, | |||
int | testNode, | |||
int | unkDerivDir, | |||
int | unkNode | |||
) | const [inline, private] |
Definition at line 148 of file SundanceReducedIntegral.hpp.
References Sundance::ElementIntegral::nNodes(), Sundance::ElementIntegral::nNodesUnk(), Sundance::ElementIntegral::nRefDerivUnk(), and W_.
double& Sundance::ReducedIntegral::value | ( | int | facetCase, | |
int | testDerivDir, | |||
int | testNode, | |||
int | unkDerivDir, | |||
int | unkNode | |||
) | [inline, private] |
Definition at line 141 of file SundanceReducedIntegral.hpp.
References Sundance::ElementIntegral::nNodes(), Sundance::ElementIntegral::nNodesUnk(), Sundance::ElementIntegral::nRefDerivUnk(), and W_.
Referenced by ReducedIntegral().
Array<Array<double> > Sundance::ReducedIntegral::W_ [private] |
Definition at line 174 of file SundanceReducedIntegral.hpp.
Referenced by ReducedIntegral(), transformOneForm(), transformTwoForm(), transformZeroForm(), and value().