49 #ifndef Intrepid2_ProjectedGeometry_h 50 #define Intrepid2_ProjectedGeometry_h 52 #include "Intrepid2_ScalarView.hpp" 66 template<
int spaceDim,
typename Po
intScalar,
typename DeviceType>
70 using ViewType = ScalarView< PointScalar, DeviceType>;
71 using ConstViewType = ScalarView<const PointScalar, DeviceType>;
72 using BasisPtr = Teuchos::RCP< Basis<DeviceType,PointScalar,PointScalar> >;
83 template<
class ExactGeometry,
class ExactGeometryGradient>
85 const ExactGeometry &exactGeometry,
const ExactGeometryGradient &exactGeometryGradient)
87 const ordinal_type numCells = flatCellGeometry.
extent_int(0);
89 INTREPID2_TEST_FOR_EXCEPTION(spaceDim != targetHGradBasis->getBaseCellTopology().getDimension(), std::invalid_argument,
"spaceDim must match the cell topology on which target basis is defined");
90 INTREPID2_TEST_FOR_EXCEPTION(projectedBasisNodes.rank() != 3, std::invalid_argument,
"projectedBasisNodes must have shape (C,F,D)");
91 INTREPID2_TEST_FOR_EXCEPTION(projectedBasisNodes.extent_int(0) != numCells, std::invalid_argument,
"cell counts must match in projectedBasisNodes and cellNodesToMap");
92 INTREPID2_TEST_FOR_EXCEPTION(projectedBasisNodes.extent_int(1) != targetHGradBasis->getCardinality(), std::invalid_argument,
"projectedBasisNodes must have shape (C,F,D)");
93 INTREPID2_TEST_FOR_EXCEPTION(projectedBasisNodes.extent_int(2) != spaceDim, std::invalid_argument,
"projectedBasisNodes must have shape (C,F,D)");
95 using ExecutionSpace =
typename DeviceType::execution_space;
99 ProjectionStruct projectionStruct;
100 ordinal_type targetQuadratureDegree(targetHGradBasis->getDegree()), targetDerivativeQuadratureDegree(targetHGradBasis->getDegree());
103 const ordinal_type numPoints = projectionStruct.getNumTargetEvalPoints();
104 const ordinal_type numGradPoints = projectionStruct.getNumTargetDerivEvalPoints();
106 ViewType evaluationPointsRefSpace (
"ProjectedGeometry evaluation points ref space (value)", numCells, numPoints, spaceDim);
107 ViewType evaluationGradPointsRefSpace(
"ProjectedGeometry evaluation points ref space (gradient)", numCells, numGradPoints, spaceDim);
110 ProjectionTools::getHGradEvaluationPoints(evaluationPointsRefSpace, evaluationGradPointsRefSpace, elementOrientations, targetHGradBasis.get(), &projectionStruct);
115 ViewType evaluationPoints (
"ProjectedGeometry evaluation points (value)", numCells, numPoints, spaceDim);
116 ViewType evaluationGradPoints(
"ProjectedGeometry evaluation points (gradient)", numCells, numGradPoints, spaceDim);
119 BasisPtr hgradLinearBasisForFlatGeometry = flatCellGeometry.
basisForNodes();
124 if (numGradPoints > 0)
132 auto policy = Kokkos::MDRangePolicy<ExecutionSpace,Kokkos::Rank<2>>({0,0}, {numCells,numPoints});
133 auto gradPolicy = Kokkos::MDRangePolicy<ExecutionSpace,Kokkos::Rank<3>>({0,0,0},{numCells,numGradPoints,spaceDim});
135 ViewType evaluationValues (
"exact geometry values", numCells, numPoints);
136 ViewType evaluationGradients (
"exact geometry gradients", numCells, numGradPoints, spaceDim);
140 for (
int comp=0; comp<spaceDim; comp++)
142 Kokkos::parallel_for(
"evaluate geometry function for projection", policy,
143 KOKKOS_LAMBDA (
const int &cellOrdinal,
const int &pointOrdinal) {
144 Kokkos::Array<PointScalar,spaceDim> point;
145 for (
int d=0; d<spaceDim; d++)
147 point[d] = evaluationPoints(cellOrdinal,pointOrdinal,d);
149 evaluationValues(cellOrdinal,pointOrdinal) = exactGeometry(point,comp);
159 flatCellGeometry.
setJacobian(gradPointsJacobians,evaluationGradPoints,refData);
161 Kokkos::parallel_for(
"evaluate geometry gradients for projection", gradPolicy,
162 KOKKOS_LAMBDA (
const int &cellOrdinal,
const int &pointOrdinal,
const int &d2) {
163 Kokkos::Array<PointScalar,spaceDim> point;
164 for (
int d=0; d<spaceDim; d++)
166 point[d] = evaluationGradPoints(cellOrdinal,pointOrdinal,d);
168 evaluationGradients(cellOrdinal,pointOrdinal,d2) = exactGeometryGradient(point,comp,d2);
175 transformedGradientData.storeMatVec(gradPointsJacobians,gradientData);
177 auto projectedBasisNodesForComp = Kokkos::subview(projectedBasisNodes,Kokkos::ALL(),Kokkos::ALL(),comp);
179 ProjectionTools::getHGradBasisCoeffs(projectedBasisNodesForComp,
181 transformedGradientData.getUnderlyingView(),
182 evaluationPointsRefSpace,
183 evaluationGradPointsRefSpace,
185 targetHGradBasis.get(),
CellGeometry provides the nodes for a set of cells; has options that support efficient definition of ...
KOKKOS_INLINE_FUNCTION std::enable_if< std::is_integral< iType >::value, int >::type extent_int(const iType &r) const
Returns the logical extent of the container in the specified dimension as an int; the shape of CellGe...
Wrapper around a Kokkos::View that allows data that is constant or repeating in various logical dimen...
void setJacobian(Data< PointScalar, DeviceType > &jacobianData, const TensorPoints< PointScalar, DeviceType > &points, const Data< PointScalar, DeviceType > &refData, const int startCell=0, const int endCell=-1) const
Compute Jacobian values for the reference-to-physical transformation, and place them in the provided ...
Data< PointScalar, DeviceType > getJacobianRefData(const TensorPoints< PointScalar, DeviceType > &points) const
Computes reference-space data for the specified points, to be used in setJacobian().
static Data< DataScalar, DeviceType > allocateMatVecResult(const Data< DataScalar, DeviceType > &matData, const Data< DataScalar, DeviceType > &vecData)
Data< PointScalar, DeviceType > allocateJacobianData(const TensorPoints< PointScalar, DeviceType > &points, const int startCell=0, const int endCell=-1) const
Allocate a container into which Jacobians of the reference-to-physical mapping can be placed...
Data< Orientation, DeviceType > getOrientations()
Returns the orientations for all cells. Calls initializeOrientations() if it has not previously been ...
Allows generation of geometry degrees of freedom based on a provided map from straight-edged mesh dom...
Allows definition of cell geometry information, including uniform and curvilinear mesh definition...
BasisPtr basisForNodes() const
H^1 Basis used in the reference-to-physical transformation. Linear for straight-edged geometry; highe...
Utility methods for Intrepid2 unit tests.
static void projectOntoHGRADBasis(ViewType projectedBasisNodes, BasisPtr targetHGradBasis, CellGeometry< PointScalar, spaceDim, DeviceType > flatCellGeometry, const ExactGeometry &exactGeometry, const ExactGeometryGradient &exactGeometryGradient)
Generate geometry degrees of freedom based on a provided map from straight-edged mesh domain to curvi...
An helper class to compute the evaluation points and weights needed for performing projections...
void createHGradProjectionStruct(const BasisPtrType cellBasis, const ordinal_type targetCubDegree, const ordinal_type targetGradCubDegre)
Initialize the ProjectionStruct for HGRAD projections.
Header file for the abstract base class Intrepid2::Basis.