Intrepid2
|
Defines expert-level interfaces for the evaluation of functions and operators in physical space (supported for FE, FV, and FD methods) and FE reference space; in addition, provides several function transformation utilities. More...
#include <Intrepid2_FunctionSpaceTools.hpp>
Static Public Member Functions | |
template<typename outputValueType , class... outputProperties, typename inputValueType , class... inputProperties> | |
static void | HGRADtransformVALUE (Kokkos::DynRankView< outputValueType, outputProperties...> output, const Kokkos::DynRankView< inputValueType, inputProperties...> input) |
Transformation of a (scalar) value field in the H-grad space, defined at points on a reference cell, stored in the user-provided container input and indexed by (F,P), into the output container output, defined on cells in physical space and indexed by (C,F,P). More... | |
template<typename outputValValueType , class... outputValProperties, typename jacobianInverseValueType , class... jacobianInverseProperties, typename inputValValueType , class... inputValProperties> | |
static void | HGRADtransformGRAD (Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianInverseValueType, jacobianInverseProperties...> jacobianInverse, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals) |
Transformation of a gradient field in the H-grad space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P,D), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P,D). More... | |
template<typename outputValValueType , class... outputValProperties, typename jacobianInverseValueType , class... jacobianInverseProperties, typename inputValValueType , class... inputValProperties> | |
static void | HCURLtransformVALUE (Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianInverseValueType, jacobianInverseProperties...> jacobianInverse, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals) |
Transformation of a (vector) value field in the H-curl space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P,D), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P,D). More... | |
template<typename outputValValueType , class... outputValProperties, typename jacobianValueType , class... jacobianProperties, typename jacobianDetValueType , class... jacobianDetProperties, typename inputValValueType , class... inputValProperties> | |
static void | HCURLtransformCURL (Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianValueType, jacobianProperties...> jacobian, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals) |
Transformation of a 3D curl field in the H-curl space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P,D), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P,D). More... | |
template<typename outputValValueType , class... outputValProperties, typename jacobianDetValueType , class... jacobianDetProperties, typename inputValValueType , class... inputValProperties> | |
static void | HCURLtransformCURL (Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals) |
Transformation of a 2D curl field in the H-curl space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P,D), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P). More... | |
template<typename outputValValueType , class... outputValProperties, typename jacobianValueType , class... jacobianProperties, typename jacobianDetValueType , class... jacobianDetProperties, typename inputValValueType , class... inputValProperties> | |
static void | HGRADtransformCURL (Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianValueType, jacobianProperties...> jacobian, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals) |
Transformation of a 2D curl field in the H-grad space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P,D), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P,D). More... | |
template<typename outputValValueType , class... outputValProperties, typename jacobianValueType , class... jacobianProperties, typename jacobianDetValueType , class... jacobianDetProperties, typename inputValValueType , class... inputValProperties> | |
static void | HDIVtransformVALUE (Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianValueType, jacobianProperties...> jacobian, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals) |
Transformation of a (vector) value field in the H-div space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P,D), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P,D). More... | |
template<typename outputValValueType , class... outputValProperties, typename jacobianDetValueType , class... jacobianDetProperties, typename inputValValueType , class... inputValProperties> | |
static void | HDIVtransformDIV (Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals) |
Transformation of a divergence field in the H-div space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P). More... | |
template<typename outputValValueType , class... outputValProperties, typename jacobianDetValueType , class... jacobianDetProperties, typename inputValValueType , class... inputValProperties> | |
static void | HVOLtransformVALUE (Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals) |
Transformation of a (scalar) value field in the H-vol space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P). More... | |
template<typename outputValueValueType , class... outputValueProperties, typename leftValueValueType , class... leftValueProperties, typename rightValueValueType , class... rightValueProperties> | |
static void | integrate (Kokkos::DynRankView< outputValueValueType, outputValueProperties...> outputValues, const Kokkos::DynRankView< leftValueValueType, leftValueProperties...> leftValues, const Kokkos::DynRankView< rightValueValueType, rightValueProperties...> rightValues, const bool sumInto=false) |
Contracts leftValues and rightValues arrays on the point and possibly space dimensions and stores the result in outputValues; this is a generic, high-level integration routine that calls either FunctionSpaceTools::operatorIntegral, or FunctionSpaceTools::functionalIntegral, or FunctionSpaceTools::dataIntegral methods, depending on the rank of the outputValues array. More... | |
template<typename outputValValueType , class... outputValProperties, typename inputDetValueType , class... inputDetPropertes, typename inputWeightValueType , class... inputWeightPropertes> | |
static bool | computeCellMeasure (Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< inputDetValueType, inputDetPropertes...> inputDet, const Kokkos::DynRankView< inputWeightValueType, inputWeightPropertes...> inputWeights) |
Returns the weighted integration measures outputVals with dimensions (C,P) used for the computation of cell integrals, by multiplying absolute values of the user-provided cell Jacobian determinants inputDet with dimensions (C,P) with the user-provided integration weights inputWeights with dimensions (P). More... | |
template<typename outputValValueType , class... outputValProperties, typename inputJacValueType , class... inputJacProperties, typename inputWeightValueType , class... inputWeightPropertes, typename scratchValueType , class... scratchProperties> | |
static void | computeFaceMeasure (Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< inputJacValueType, inputJacProperties...> inputJac, const Kokkos::DynRankView< inputWeightValueType, inputWeightPropertes...> inputWeights, const int whichFace, const shards::CellTopology parentCell, const Kokkos::DynRankView< scratchValueType, scratchProperties...> scratch) |
Returns the weighted integration measures outputVals with dimensions (C,P) used for the computation of face integrals, based on the provided cell Jacobian array inputJac with dimensions (C,P,D,D) and the provided integration weights inputWeights with dimensions (P). More... | |
template<typename outputValValueType , class... outputValProperties, typename inputJacValueType , class... inputJacProperties, typename inputWeightValueType , class... inputWeightProperties, typename scratchValueType , class... scratchProperties> | |
static void | computeEdgeMeasure (Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< inputJacValueType, inputJacProperties...> inputJac, const Kokkos::DynRankView< inputWeightValueType, inputWeightProperties...> inputWeights, const int whichEdge, const shards::CellTopology parentCell, const Kokkos::DynRankView< scratchValueType, scratchProperties...> scratch) |
Returns the weighted integration measures outVals with dimensions (C,P) used for the computation of edge integrals, based on the provided cell Jacobian array inputJac with dimensions (C,P,D,D) and the provided integration weights inWeights with dimensions (P). More... | |
template<typename outputValValueType , class... outputValProperties, typename inputMeasureValueType , class... inputMeasureProperties, typename inputValValueType , class... inputValProperteis> | |
static void | multiplyMeasure (Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< inputMeasureValueType, inputMeasureProperties...> inputMeasure, const Kokkos::DynRankView< inputValValueType, inputValProperteis...> inputVals) |
Multiplies fields inputVals by weighted measures inputMeasure and returns the field array outputVals; this is a simple redirection to the call FunctionSpaceTools::scalarMultiplyDataField. More... | |
template<typename outputFieldValueType , class... outputFieldProperties, typename inputDataValueType , class... inputDataPropertes, typename inputFieldValueType , class... inputFieldProperties> | |
static void | scalarMultiplyDataField (Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataPropertes...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields, const bool reciprocal=false) |
Scalar multiplication of data and fields; please read the description below. More... | |
template<typename outputDataValuetype , class... outputDataProperties, typename inputDataLeftValueType , class... inputDataLeftProperties, typename inputDataRightValueType , class... inputDataRightProperties> | |
static void | scalarMultiplyDataData (Kokkos::DynRankView< outputDataValuetype, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight, const bool reciprocal=false) |
Scalar multiplication of data and data; please read the description below. More... | |
template<typename outputFieldValueType , class... outputFieldProperties, typename inputDataValueType , class... inputDataProperties, typename inputFieldValueType , class... inputFieldProperties> | |
static void | dotMultiplyDataField (Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields) |
Dot product of data and fields; please read the description below. More... | |
template<typename outputDataValueType , class... outputDataProperties, typename inputDataLeftValueType , class... inputDataLeftProperties, typename inputDataRightValueType , class... inputDataRightProperties> | |
static void | dotMultiplyDataData (Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight) |
Dot product of data and data; please read the description below. More... | |
template<typename outputFieldValueType , class... outputFieldProperties, typename inputDataValueType , class... inputDataProperties, typename inputFieldValueType , class... inputFieldProperties> | |
static void | vectorMultiplyDataField (Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields) |
Cross or outer product of data and fields; please read the description below. More... | |
template<typename outputDataValueType , class... outputDataProperties, typename inputDataLeftValueType , class... inputDataLeftProperties, typename inputDataRightValueType , class... inputDataRightProperties> | |
static void | vectorMultiplyDataData (Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight) |
Cross or outer product of data and data; please read the description below. More... | |
template<typename outputFieldValueType , class... outputFieldProperties, typename inputDataValueType , class... inputDataProperties, typename inputFieldValueType , class... inputFieldProperties> | |
static void | tensorMultiplyDataField (Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields, const char transpose= 'N') |
Matrix-vector or matrix-matrix product of data and fields; please read the description below. More... | |
template<typename outputDataValueType , class... outputDataProperties, typename inputDataLeftValueType , class... inputDataLeftProperties, typename inputDataRightValueType , class... inputDataRightProperties> | |
static void | tensorMultiplyDataData (Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight, const char transpose= 'N') |
Matrix-vector or matrix-matrix product of data and data; please read the description below. More... | |
template<typename inoutOperatorValueType , class... inoutOperatorProperties, typename fieldSignValueType , class... fieldSignProperties> | |
static void | applyLeftFieldSigns (Kokkos::DynRankView< inoutOperatorValueType, inoutOperatorProperties...> inoutOperator, const Kokkos::DynRankView< fieldSignValueType, fieldSignProperties...> fieldSigns) |
Applies left (row) signs, stored in the user-provided container fieldSigns and indexed by (C,L), to the operator inoutOperator indexed by (C,L,R). More... | |
template<typename inoutOperatorValueType , class... inoutOperatorProperties, typename fieldSignValueType , class... fieldSignProperties> | |
static void | applyRightFieldSigns (Kokkos::DynRankView< inoutOperatorValueType, inoutOperatorProperties...> inoutOperator, const Kokkos::DynRankView< fieldSignValueType, fieldSignProperties...> fieldSigns) |
Applies right (column) signs, stored in the user-provided container fieldSigns and indexed by (C,R), to the operator inoutOperator indexed by (C,L,R). More... | |
template<typename inoutFunctionValueType , class... inoutFunctionProperties, typename fieldSignValueType , class... fieldSignProperties> | |
static void | applyFieldSigns (Kokkos::DynRankView< inoutFunctionValueType, inoutFunctionProperties...> inoutFunction, const Kokkos::DynRankView< fieldSignValueType, fieldSignProperties...> fieldSigns) |
Applies field signs, stored in the user-provided container fieldSigns and indexed by (C,F), to the function inoutFunction indexed by (C,F), (C,F,P), (C,F,P,D1) or (C,F,P,D1,D2). More... | |
template<typename outputPointValueType , class... outputPointProperties, typename inputCoeffValueType , class... inputCoeffProperties, typename inputFieldValueType , class... inputFieldProperties> | |
static void | evaluate (Kokkos::DynRankView< outputPointValueType, outputPointProperties...> outputPointVals, const Kokkos::DynRankView< inputCoeffValueType, inputCoeffProperties...> inputCoeffs, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields) |
Computes point values outPointVals of a discrete function specified by the basis inFields and coefficients inCoeffs. More... | |
Defines expert-level interfaces for the evaluation of functions and operators in physical space (supported for FE, FV, and FD methods) and FE reference space; in addition, provides several function transformation utilities.
Definition at line 78 of file Intrepid2_FunctionSpaceTools.hpp.
|
static |
Applies field signs, stored in the user-provided container fieldSigns and indexed by (C,F), to the function inoutFunction indexed by (C,F), (C,F,P), (C,F,P,D1) or (C,F,P,D1,D2).
Returns
See Section Pullbacks for discussion of field signs.
inoutFunction | [in/out] - Input / output function array. |
fieldSigns | [in] - Right field signs. |
Definition at line 933 of file Intrepid2_FunctionSpaceToolsDef.hpp.
|
static |
Applies left (row) signs, stored in the user-provided container fieldSigns and indexed by (C,L), to the operator inoutOperator indexed by (C,L,R).
Mathematically, this method computes the matrix-matrix product
where is array of matrices indexed by cell number c and stored in the rank-3 array inoutOperator, and is array of left field signs indexed by cell number c and stored in the rank-2 container fieldSigns; see Section Pullbacks for discussion of field signs. This operation is required for operators generated by HCURL and HDIV-conforming vector-valued finite element basis functions; see Sections Pullbacks and Section Evaluation of finite element operators and functionals for applications of this method.
inoutOperator | [in/out] - Input / output operator array. |
fieldSigns | [in] - Left field signs. |
Definition at line 811 of file Intrepid2_FunctionSpaceToolsDef.hpp.
|
static |
Applies right (column) signs, stored in the user-provided container fieldSigns and indexed by (C,R), to the operator inoutOperator indexed by (C,L,R).
Mathematically, this method computes the matrix-matrix product
where is array of matrices indexed by cell number c and stored in the rank-3 container inoutOperator, and is array of right field signs indexed by cell number c and stored in the rank-2 container fieldSigns; see Section Pullbacks for discussion of field signs. This operation is required for operators generated by HCURL and HDIV-conforming vector-valued finite element basis functions; see Sections Pullbacks and Section Evaluation of finite element operators and functionals for applications of this method.
inoutOperator | [in/out] - Input / output operator array. |
fieldSigns | [in] - Right field signs. |
Definition at line 870 of file Intrepid2_FunctionSpaceToolsDef.hpp.
|
static |
Returns the weighted integration measures outputVals with dimensions (C,P) used for the computation of cell integrals, by multiplying absolute values of the user-provided cell Jacobian determinants inputDet with dimensions (C,P) with the user-provided integration weights inputWeights with dimensions (P).
Returns a rank-2 array (C, P) array such that
where is a cubature rule defined on a reference cell (a set of integration points and their associated weights; see Intrepid2::Cubature::getCubature for getting cubature rules on reference cells).
outputVals | [out] - Output array with weighted cell measures. |
inputDet | [in] - Input array containing determinants of cell Jacobians. |
inputWeights | [in] - Input integration weights. |
Definition at line 438 of file Intrepid2_FunctionSpaceToolsDef.hpp.
|
static |
Returns the weighted integration measures outVals with dimensions (C,P) used for the computation of edge integrals, based on the provided cell Jacobian array inputJac with dimensions (C,P,D,D) and the provided integration weights inWeights with dimensions (P).
Returns a rank-2 array (C, P) array such that
where:
outputVals | [out] - Output array with weighted edge measures. |
inputJac | [in] - Input array containing cell Jacobians. |
inputWeights | [in] - Input integration weights. |
whichEdge | [in] - Index of the edge subcell relative to the parent cell; defines the domain of integration. |
parentCell | [in] - Parent cell topology. |
scratch | [in] - Scratch space, sized like inputJac |
Definition at line 525 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::CellTools< ExecSpaceType >::getPhysicalEdgeTangents(), Intrepid2::ArrayTools< ExecSpaceType >::scalarMultiplyDataData(), and Intrepid2::RealSpaceTools< ExecSpaceType >::vectorNorm().
Referenced by Intrepid2::CubatureControlVolumeBoundary< ExecSpaceType, pointValueType, weightValueType >::getCubature().
|
static |
Returns the weighted integration measures outputVals with dimensions (C,P) used for the computation of face integrals, based on the provided cell Jacobian array inputJac with dimensions (C,P,D,D) and the provided integration weights inputWeights with dimensions (P).
Returns a rank-2 array (C, P) array such that
where:
outputVals | [out] - Output array with weighted face measures. |
inputJac | [in] - Input array containing cell Jacobians. |
inputWeights | [in] - Input integration weights. |
whichFace | [in] - Index of the face subcell relative to the parent cell; defines the domain of integration. |
parentCell | [in] - Parent cell topology. |
scratch | [in] - Scratch space, sized like inputJac |
Definition at line 478 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::CellTools< ExecSpaceType >::getPhysicalFaceNormals(), Intrepid2::ArrayTools< ExecSpaceType >::scalarMultiplyDataData(), and Intrepid2::RealSpaceTools< ExecSpaceType >::vectorNorm().
Referenced by Intrepid2::CubatureControlVolumeBoundary< ExecSpaceType, pointValueType, weightValueType >::getCubature().
|
static |
Dot product of data and data; please read the description below.
There are two use cases:
For input fields containers without a dimension index, this operation reduces to scalar multiplication.
outputData | [out] - Output (dot product) data array. |
inputDataLeft | [in] - Left input data array. |
inputDataRight | [in] - Right input data array. |
Definition at line 639 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::dotMultiplyDataData().
|
static |
Dot product of data and fields; please read the description below.
There are two use cases:
For input fields containers without a dimension index, this operation reduces to scalar multiplication.
outputFields | [out] - Output (dot product) fields array. |
inputData | [in] - Data array. |
inputFields | [in] - Input fields array. |
Definition at line 623 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::dotMultiplyDataField().
|
static |
Computes point values outPointVals of a discrete function specified by the basis inFields and coefficients inCoeffs.
The array inFields with dimensions (C,F,P), (C,F,P,D1), or (C,F,P,D1,D2) represents the signed, transformed field (basis) values at points in REFERENCE frame; the outPointVals array with dimensions (C,P), (C,P,D1), or (C,P,D1,D2), respectively, represents values of a discrete function at points in PHYSICAL frame. The array inCoeffs dimensioned (C,F) supplies the coefficients for the field (basis) array.
Returns rank-2,3 or 4 array such that
where is scalar, vector or tensor valued finite element basis defined on physical cell and are the field signs of the basis functions; see Section Pullbacks. This method implements the last step in a four step process; please see Section Evaluation of finite element fields for details about the first three steps that prepare the necessary data for this method.
outputPointVals | [out] - Output point values of a discrete function. |
inputCoeffs | [in] - Coefficients associated with the fields (basis) array. |
inputFields | [in] - Field (basis) values. |
Definition at line 1003 of file Intrepid2_FunctionSpaceToolsDef.hpp.
|
static |
Transformation of a 3D curl field in the H-curl space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P,D), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P,D).
Computes pullback of curls of HCURL functions for points in one or more physical cells that are images of a given set of points in the reference cell:
In this case and the user-provided container should contain the curls of the vector function set at the reference points:
The method returns
See Section Pullbacks for more details about pullbacks.
outputVals | [out] - Output array with transformed values |
jacobian | [in] - Input array containing cell Jacobians. |
jacobianDet | [in] - Input array containing cell Jacobian determinants. |
inputVals | [in] - Input array of reference HCURL curls. |
Definition at line 189 of file Intrepid2_FunctionSpaceToolsDef.hpp.
|
static |
Transformation of a 2D curl field in the H-curl space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P,D), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P).
Computes pullback of curls of HCURL functions for points in one or more physical cells that are images of a given set of points in the reference cell:
In this case and the user-provided container should contain the 2d curls of the vector function set at the reference points:
The method returns
See Section Pullbacks for more details about pullbacks.
outputVals | [out] - Output array with transformed values |
jacobianDet | [in] - Input array containing cell Jacobian determinants. |
inputVals | [in] - Input array of reference HCURL curls. |
Definition at line 207 of file Intrepid2_FunctionSpaceToolsDef.hpp.
|
static |
Transformation of a (vector) value field in the H-curl space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P,D), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P,D).
Computes pullback of HCURL functions for points in one or more physical cells that are images of a given set of points in the reference cell:
In this case and the user-provided container should contain the values of the vector function set at the reference points:
The method returns
See Section Pullbacks for more details about pullbacks.
outputVals | [out] - Output array with transformed values |
jacobianInverse | [in] - Input array containing cell Jacobian inverses. |
inputVals | [in] - Input array of reference HCURL values. |
Definition at line 174 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::matvecProductDataField().
|
static |
Transformation of a divergence field in the H-div space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P).
Computes pullback of the divergence of HDIV functions for points in one or more physical cells that are images of a given set of points in the reference cell:
In this case and the user-provided container should contain the divergencies of the vector function set at the reference points:
The method returns
See Section Pullbacks for more details about pullbacks.
outputVals | [out] - Output array with transformed values |
jacobianDet | [in] - Input array containing cell Jacobian determinants. |
inputVals | [in] - Input array of reference HDIV divergences. |
Definition at line 266 of file Intrepid2_FunctionSpaceToolsDef.hpp.
|
static |
Transformation of a (vector) value field in the H-div space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P,D), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P,D).
Computes pullback of HDIV functions for points in one or more physical cells that are images of a given set of points in the reference cell:
In this case and the user-provided container should contain the values of the vector function set at the reference points:
The method returns
See Section Pullbacks for more details about pullbacks.
outputVals | [out] - Output array with transformed values |
jacobian | [in] - Input array containing cell Jacobians. |
jacobianDet | [in] - Input array containing cell Jacobian determinants. |
inputVals | [in] - Input array of reference HDIV values. |
Definition at line 250 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::matvecProductDataField(), and Intrepid2::ArrayTools< ExecSpaceType >::scalarMultiplyDataField().
|
static |
Transformation of a 2D curl field in the H-grad space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P,D), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P,D).
Computes pullback of curls of 2D HGRAD functions for points in one or more physical cells that are images of a given set of points in the reference cell:
In this case and the user-provided container should contain the curls of the vector function set at the reference points:
The method returns
See Section Pullbacks for more details about pullbacks.
outputVals | [out] - Output array with transformed values |
jacobian | [in] - Input array containing cell Jacobians. |
jacobianDet | [in] - Input array containing cell Jacobian determinants. |
inputVals | [in] - Input array of reference HDIV values. |
Definition at line 228 of file Intrepid2_FunctionSpaceToolsDef.hpp.
|
static |
Transformation of a gradient field in the H-grad space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P,D), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P,D).
Computes pullback of gradients of HGRAD functions for points in one or more physical cells that are images of a given set of points in the reference cell:
In this case and the user-provided container should contain the gradients of the function set at the reference points:
The method returns
See Section Pullbacks for more details about pullbacks.
outputVals | [out] - Output array with transformed values |
jacobianInverse | [in] - Input array containing cell Jacobian inverses. |
inputVals | [in] - Input array of reference HGRAD gradients. |
Definition at line 124 of file Intrepid2_FunctionSpaceToolsDef.hpp.
|
static |
Transformation of a (scalar) value field in the H-grad space, defined at points on a reference cell, stored in the user-provided container input and indexed by (F,P), into the output container output, defined on cells in physical space and indexed by (C,F,P).
Computes pullback of HGRAD functions for points in one or more physical cells that are images of a given set of points in the reference cell:
In this case and the user-provided container should contain the values of the function set at the reference points:
The method returns
i.e., it simply replicates the values in the user-provided container to every cell. See Section Pullbacks for more details about pullbacks.
output | [out] - Output array with transformed values |
input | [in] - Input array of reference HGRAD values. |
Definition at line 60 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::RealSpaceTools< ExecSpaceType >::clone(), and Intrepid2::ArrayTools< ExecSpaceType >::cloneFields().
|
static |
Transformation of a (scalar) value field in the H-vol space, defined at points on a reference cell, stored in the user-provided container inputVals and indexed by (F,P), into the output container outputVals, defined on cells in physical space and indexed by (C,F,P).
Computes pullback of HVOL functions for points in one or more physical cells that are images of a given set of points in the reference cell:
In this case and the user-provided container should contain the values of the functions in the set at the reference points:
The method returns
See Section Pullbacks for more details about pullbacks.
outputVals | [out] - Output array with transformed values |
jacobianDet | [in] - Input array containing cell Jacobian determinants. |
inputVals | [in] - Input array of reference HVOL values. |
Definition at line 280 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::scalarMultiplyDataField().
|
static |
Contracts leftValues and rightValues arrays on the point and possibly space dimensions and stores the result in outputValues; this is a generic, high-level integration routine that calls either FunctionSpaceTools::operatorIntegral, or FunctionSpaceTools::functionalIntegral, or FunctionSpaceTools::dataIntegral methods, depending on the rank of the outputValues array.
outputValues | [out] - Output array. |
leftValues | [in] - Left input array. |
rightValues | [in] - Right input array. |
sumInto | [in] - If TRUE, sum into given output array, otherwise overwrite it. Default: FALSE. |
Definition at line 294 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::contractDataDataScalar(), Intrepid2::ArrayTools< ExecSpaceType >::contractDataDataTensor(), Intrepid2::ArrayTools< ExecSpaceType >::contractDataDataVector(), Intrepid2::ArrayTools< ExecSpaceType >::contractDataFieldScalar(), Intrepid2::ArrayTools< ExecSpaceType >::contractDataFieldTensor(), Intrepid2::ArrayTools< ExecSpaceType >::contractDataFieldVector(), Intrepid2::ArrayTools< ExecSpaceType >::contractFieldFieldScalar(), Intrepid2::ArrayTools< ExecSpaceType >::contractFieldFieldTensor(), and Intrepid2::ArrayTools< ExecSpaceType >::contractFieldFieldVector().
|
static |
Multiplies fields inputVals by weighted measures inputMeasure and returns the field array outputVals; this is a simple redirection to the call FunctionSpaceTools::scalarMultiplyDataField.
outputVals | [out] - Output array with scaled field values. |
inputMeasure | [in] - Input array containing weighted measures. |
inputVals | [in] - Input fields. |
Definition at line 571 of file Intrepid2_FunctionSpaceToolsDef.hpp.
|
static |
Scalar multiplication of data and data; please read the description below.
There are two use cases:
outputData | [out] - Output data array. |
inputDataLeft | [in] - Left (multiplying) data array. |
inputDataRight | [in] - Right (being multiplied) data array. |
reciprocal | [in] - If TRUE, divides input fields by the data (instead of multiplying). Default: FALSE. |
Definition at line 605 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::scalarMultiplyDataData().
|
static |
Scalar multiplication of data and fields; please read the description below.
There are two use cases:
outputFields | [out] - Output (product) fields array. |
inputData | [in] - Data (multiplying) array. |
inputFields | [in] - Input (being multiplied) fields array. |
reciprocal | [in] - If TRUE, divides input fields by the data (instead of multiplying). Default: FALSE. |
Definition at line 587 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::scalarMultiplyDataField().
|
static |
Matrix-vector or matrix-matrix product of data and data; please read the description below.
There are four use cases:
outputData | [out] - Output (matrix-vector product) data array. |
inputDataLeft | [in] - Left input data array. |
inputDataRight | [in] - Right input data array. |
transpose | [in] - If 'T', use transposed tensor; if 'N', no transpose. Default: 'N'. |
Definition at line 751 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::matmatProductDataData(), and Intrepid2::ArrayTools< ExecSpaceType >::matvecProductDataData().
|
static |
Matrix-vector or matrix-matrix product of data and fields; please read the description below.
There are four use cases:
outputFields | [out] - Output (matrix-vector or matrix-matrix product) fields array. |
inputData | [in] - Data array. |
inputFields | [in] - Input fields array. |
transpose | [in] - If 'T', use transposed left data tensor; if 'N', no transpose. Default: 'N'. |
Definition at line 717 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::matmatProductDataField(), and Intrepid2::ArrayTools< ExecSpaceType >::matvecProductDataField().
|
static |
Cross or outer product of data and data; please read the description below.
There are four use cases:
outputData | [out] - Output (cross or outer product) data array. |
inputDataLeft | [in] - Left input data array. |
inputDataRight | [in] - Right input data array. |
Definition at line 686 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::crossProductDataData(), and Intrepid2::ArrayTools< ExecSpaceType >::outerProductDataData().
|
static |
Cross or outer product of data and fields; please read the description below.
There are four use cases:
outputFields | [out] - Output (cross or outer product) fields array. |
inputData | [in] - Data array. |
inputFields | [in] - Input fields array. |
Definition at line 655 of file Intrepid2_FunctionSpaceToolsDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::crossProductDataField(), and Intrepid2::ArrayTools< ExecSpaceType >::outerProductDataField().