Intrepid2
Intrepid2_ProjectedGeometry.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Intrepid2 Package
4//
5// Copyright 2007 NTESS and the Intrepid2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
15#ifndef Intrepid2_ProjectedGeometry_h
16#define Intrepid2_ProjectedGeometry_h
17
18#include "Intrepid2_ScalarView.hpp"
19
20#include "Intrepid2_Basis.hpp"
24
26
27namespace Intrepid2
28{
32 template<int spaceDim, typename PointScalar, typename DeviceType>
34 {
35 public:
36 using ViewType = ScalarView< PointScalar, DeviceType>;
37 using ConstViewType = ScalarView<const PointScalar, DeviceType>;
38 using BasisPtr = Teuchos::RCP< Basis<DeviceType,PointScalar,PointScalar> >;
39
49 template<class ExactGeometry, class ExactGeometryGradient>
50 static void projectOntoHGRADBasis(ViewType projectedBasisNodes, BasisPtr targetHGradBasis, CellGeometry<PointScalar,spaceDim,DeviceType> flatCellGeometry,
51 const ExactGeometry &exactGeometry, const ExactGeometryGradient &exactGeometryGradient)
52 {
53 const ordinal_type numCells = flatCellGeometry.extent_int(0); // (C,N,D)
54
55 INTREPID2_TEST_FOR_EXCEPTION(spaceDim != targetHGradBasis->getBaseCellTopology().getDimension(), std::invalid_argument, "spaceDim must match the cell topology on which target basis is defined");
56 INTREPID2_TEST_FOR_EXCEPTION(projectedBasisNodes.rank() != 3, std::invalid_argument, "projectedBasisNodes must have shape (C,F,D)");
57 INTREPID2_TEST_FOR_EXCEPTION(projectedBasisNodes.extent_int(0) != numCells, std::invalid_argument, "cell counts must match in projectedBasisNodes and cellNodesToMap");
58 INTREPID2_TEST_FOR_EXCEPTION(projectedBasisNodes.extent_int(1) != targetHGradBasis->getCardinality(), std::invalid_argument, "projectedBasisNodes must have shape (C,F,D)");
59 INTREPID2_TEST_FOR_EXCEPTION(projectedBasisNodes.extent_int(2) != spaceDim, std::invalid_argument, "projectedBasisNodes must have shape (C,F,D)");
60
61 using ExecutionSpace = typename DeviceType::execution_space;
64
65 ProjectionStruct projectionStruct;
66 ordinal_type targetQuadratureDegree(targetHGradBasis->getDegree()), targetDerivativeQuadratureDegree(targetHGradBasis->getDegree());
67 projectionStruct.createHGradProjectionStruct(targetHGradBasis, targetQuadratureDegree, targetDerivativeQuadratureDegree);
68
69 auto evaluationPointsRefSpace = projectionStruct.getAllEvalPoints();
70 auto evaluationGradPointsRefSpace = projectionStruct.getAllDerivEvalPoints();
71 const ordinal_type numPoints = evaluationPointsRefSpace.extent(0);
72 const ordinal_type numGradPoints = evaluationGradPointsRefSpace.extent(0);
73
74 auto elementOrientations = flatCellGeometry.getOrientations();
75
76// printFunctor1(elementOrientations, std::cout);
77
78 // the evaluation points are all still in reference space; map to physical space:
79 ViewType evaluationPoints ("ProjectedGeometry evaluation points (value)", numCells, numPoints, spaceDim);
80 ViewType evaluationGradPoints("ProjectedGeometry evaluation points (gradient)", numCells, numGradPoints, spaceDim);
81
83 BasisPtr hgradLinearBasisForFlatGeometry = flatCellGeometry.basisForNodes();
84 if (numPoints > 0)
85 {
86 CellTools::mapToPhysicalFrame(evaluationPoints, evaluationPointsRefSpace, flatCellGeometry, hgradLinearBasisForFlatGeometry);
87 }
88 if (numGradPoints > 0)
89 {
90 CellTools::mapToPhysicalFrame(evaluationGradPoints, evaluationGradPointsRefSpace, flatCellGeometry, hgradLinearBasisForFlatGeometry);
91 }
92
93 auto refData = flatCellGeometry.getJacobianRefData(evaluationGradPoints);
94
95 // evaluate, transform, and project in each component
96 auto policy = Kokkos::MDRangePolicy<ExecutionSpace,Kokkos::Rank<2>>({0,0}, {numCells,numPoints});
97 auto gradPolicy = Kokkos::MDRangePolicy<ExecutionSpace,Kokkos::Rank<3>>({0,0,0},{numCells,numGradPoints,spaceDim});
98
99 ViewType evaluationValues ("exact geometry values", numCells, numPoints);
100 ViewType evaluationGradients ("exact geometry gradients", numCells, numGradPoints, spaceDim);
101
102// printView(evaluationPoints, std::cout, "evaluationPoints");
103
104 for (int comp=0; comp<spaceDim; comp++)
105 {
106 Kokkos::parallel_for("evaluate geometry function for projection", policy,
107 KOKKOS_LAMBDA (const int &cellOrdinal, const int &pointOrdinal) {
108 Kokkos::Array<PointScalar,spaceDim> point;
109 for (int d=0; d<spaceDim; d++)
110 {
111 point[d] = evaluationPoints(cellOrdinal,pointOrdinal,d);
112 }
113 evaluationValues(cellOrdinal,pointOrdinal) = exactGeometry(point,comp);
114 });
115
116// printView(evaluationValues, std::cout, "evaluationValues");
117
118 // projection occurs in ref space, so we need to apply inverse of the pullback
119 // HGRADtransformVALUE is identity, so evaluationValues above is correct
120 // HGRADtransformGRAD is multiplication by inverse of Jacobian, so here we want to multiply by Jacobian
121
122 auto gradPointsJacobians = flatCellGeometry.allocateJacobianData(evaluationGradPoints);
123 flatCellGeometry.setJacobian(gradPointsJacobians,evaluationGradPoints,refData);
124
125 Kokkos::parallel_for("evaluate geometry gradients for projection", gradPolicy,
126 KOKKOS_LAMBDA (const int &cellOrdinal, const int &pointOrdinal, const int &d2) {
127 Kokkos::Array<PointScalar,spaceDim> point;
128 for (int d=0; d<spaceDim; d++)
129 {
130 point[d] = evaluationGradPoints(cellOrdinal,pointOrdinal,d);
131 }
132 evaluationGradients(cellOrdinal,pointOrdinal,d2) = exactGeometryGradient(point,comp,d2);
133 });
134
135 // apply Jacobian
136 Data<PointScalar,DeviceType> gradientData(evaluationGradients);
137 auto transformedGradientData = Data<PointScalar,DeviceType>::allocateMatVecResult(gradPointsJacobians,gradientData);
138
139 transformedGradientData.storeMatVec(gradPointsJacobians,gradientData);
140
141 auto projectedBasisNodesForComp = Kokkos::subview(projectedBasisNodes,Kokkos::ALL(),Kokkos::ALL(),comp);
142
143 ProjectionTools::getHGradBasisCoeffs(projectedBasisNodesForComp,
144 evaluationValues,
145 transformedGradientData.getUnderlyingView(),
146 elementOrientations,
147 targetHGradBasis.get(),
148 &projectionStruct);
149 }
150 }
151 };
152}
153
154#endif /* Intrepid2_ProjectedGeometry_h */
Header file for the abstract base class Intrepid2::Basis.
Allows definition of cell geometry information, including uniform and curvilinear mesh definition,...
Header file for the Intrepid2::CellTools class.
Header file for the Intrepid2::ProjectionTools.
Utility methods for Intrepid2 unit tests.
CellGeometry provides the nodes for a set of cells; has options that support efficient definition of ...
BasisPtr basisForNodes() const
H^1 Basis used in the reference-to-physical transformation. Linear for straight-edged geometry; highe...
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 ScalarView< PointScalar, DeviceType > &points) const
Computes reference-space data for the specified points, to be used in setJacobian().
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 ...
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...
A stateless class for operations on cell data. Provides methods for:
static void mapToPhysicalFrame(PhysPointValueType physPoints, const RefPointValueType refPoints, const WorksetType worksetCell, const HGradBasisPtrType basis)
Computes F, the reference-to-physical frame map.
Wrapper around a Kokkos::View that allows data that is constant or repeating in various logical dimen...
static Data< DataScalar, DeviceType > allocateMatVecResult(const Data< DataScalar, DeviceType > &matData, const Data< DataScalar, DeviceType > &vecData, const bool transposeMatrix=false)
Allows generation of geometry degrees of freedom based on a provided map from straight-edged mesh dom...
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.
view_type getAllDerivEvalPoints(EvalPointsType type=TARGET) const
Returns all the evaluation points for basis/target derivatives in the cell.
view_type getAllEvalPoints(EvalPointsType type=TARGET) const
Returns the basis/target evaluation points in the cell.
A class providing static members to perform projection-based interpolations:
static void getHGradBasisCoeffs(BasisCoeffsViewType basisCoeffs, const TargetValueViewType targetAtEvalPoints, const TargetGradViewType targetGradAtGradEvalPoints, const OrientationViewType cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
Computes the basis coefficients of the HGrad projection of the target function.