Intrepid2
Intrepid2_ProjectedGeometry.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid2 Package
5 // Copyright (2007) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Kyungjoo Kim (kyukim@sandia.gov),
38 // Mauro Perego (mperego@sandia.gov), or
39 // Nate Roberts (nvrober@sandia.gov),
40 //
41 // ************************************************************************
42 // @HEADER
43 
49 #ifndef Intrepid2_ProjectedGeometry_h
50 #define Intrepid2_ProjectedGeometry_h
51 
52 #include "Intrepid2_ScalarView.hpp"
53 
54 #include "Intrepid2_Basis.hpp"
56 #include "Intrepid2_CellTools.hpp"
58 
59 #include "Intrepid2_TestUtils.hpp"
60 
61 namespace Intrepid2
62 {
66  template<int spaceDim, typename PointScalar, typename DeviceType>
68  {
69  public:
70  using ViewType = ScalarView< PointScalar, DeviceType>;
71  using ConstViewType = ScalarView<const PointScalar, DeviceType>;
72  using BasisPtr = Teuchos::RCP< Basis<DeviceType,PointScalar,PointScalar> >;
73 
83  template<class ExactGeometry, class ExactGeometryGradient>
84  static void projectOntoHGRADBasis(ViewType projectedBasisNodes, BasisPtr targetHGradBasis, CellGeometry<PointScalar,spaceDim,DeviceType> flatCellGeometry,
85  const ExactGeometry &exactGeometry, const ExactGeometryGradient &exactGeometryGradient)
86  {
87  const ordinal_type numCells = flatCellGeometry.extent_int(0); // (C,N,D)
88 
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)");
94 
95  using ExecutionSpace = typename DeviceType::execution_space;
96  using ProjectionTools = Experimental::ProjectionTools<ExecutionSpace>; // TODO: when ProjectionTools supports it, replace template argument with DeviceType
97  using ProjectionStruct = Experimental::ProjectionStruct<ExecutionSpace,PointScalar>; // TODO: when ProjectionTools supports it, replace template argument with DeviceType
98 
99  ProjectionStruct projectionStruct;
100  ordinal_type targetQuadratureDegree(targetHGradBasis->getDegree()), targetDerivativeQuadratureDegree(targetHGradBasis->getDegree());
101  projectionStruct.createHGradProjectionStruct(targetHGradBasis, targetQuadratureDegree, targetDerivativeQuadratureDegree);
102 
103  const ordinal_type numPoints = projectionStruct.getNumTargetEvalPoints();
104  const ordinal_type numGradPoints = projectionStruct.getNumTargetDerivEvalPoints();
105 
106  ViewType evaluationPointsRefSpace ("ProjectedGeometry evaluation points ref space (value)", numCells, numPoints, spaceDim);
107  ViewType evaluationGradPointsRefSpace("ProjectedGeometry evaluation points ref space (gradient)", numCells, numGradPoints, spaceDim);
108 
109  auto elementOrientations = flatCellGeometry.getOrientations();
110  ProjectionTools::getHGradEvaluationPoints(evaluationPointsRefSpace, evaluationGradPointsRefSpace, elementOrientations, targetHGradBasis.get(), &projectionStruct);
111 
112 // printFunctor1(elementOrientations, std::cout);
113 
114  // the evaluation points are all still in reference space; map to physical space:
115  ViewType evaluationPoints ("ProjectedGeometry evaluation points (value)", numCells, numPoints, spaceDim);
116  ViewType evaluationGradPoints("ProjectedGeometry evaluation points (gradient)", numCells, numGradPoints, spaceDim);
117 
119  BasisPtr hgradLinearBasisForFlatGeometry = flatCellGeometry.basisForNodes();
120  if (numPoints > 0)
121  {
122  CellTools::mapToPhysicalFrame(evaluationPoints, evaluationPointsRefSpace, flatCellGeometry, hgradLinearBasisForFlatGeometry);
123  }
124  if (numGradPoints > 0)
125  {
126  CellTools::mapToPhysicalFrame(evaluationGradPoints, evaluationGradPointsRefSpace, flatCellGeometry, hgradLinearBasisForFlatGeometry);
127  }
128 
129  auto refData = flatCellGeometry.getJacobianRefData(evaluationGradPoints);
130 
131  // evaluate, transform, and project in each component
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});
134 
135  ViewType evaluationValues ("exact geometry values", numCells, numPoints);
136  ViewType evaluationGradients ("exact geometry gradients", numCells, numGradPoints, spaceDim);
137 
138 // printView(evaluationPoints, std::cout, "evaluationPoints");
139 
140  for (int comp=0; comp<spaceDim; comp++)
141  {
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++)
146  {
147  point[d] = evaluationPoints(cellOrdinal,pointOrdinal,d);
148  }
149  evaluationValues(cellOrdinal,pointOrdinal) = exactGeometry(point,comp);
150  });
151 
152 // printView(evaluationValues, std::cout, "evaluationValues");
153 
154  // projection occurs in ref space, so we need to apply inverse of the pullback
155  // HGRADtransformVALUE is identity, so evaluationValues above is correct
156  // HGRADtransformGRAD is multiplication by inverse of Jacobian, so here we want to multiply by Jacobian
157 
158  auto gradPointsJacobians = flatCellGeometry.allocateJacobianData(evaluationGradPoints);
159  flatCellGeometry.setJacobian(gradPointsJacobians,evaluationGradPoints,refData);
160 
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++)
165  {
166  point[d] = evaluationGradPoints(cellOrdinal,pointOrdinal,d);
167  }
168  evaluationGradients(cellOrdinal,pointOrdinal,d2) = exactGeometryGradient(point,comp,d2);
169  });
170 
171  // apply Jacobian
172  Data<PointScalar,DeviceType> gradientData(evaluationGradients);
173  auto transformedGradientData = Data<PointScalar,DeviceType>::allocateMatVecResult(gradPointsJacobians,gradientData);
174 
175  transformedGradientData.storeMatVec(gradPointsJacobians,gradientData);
176 
177  auto projectedBasisNodesForComp = Kokkos::subview(projectedBasisNodes,Kokkos::ALL(),Kokkos::ALL(),comp);
178 
179  ProjectionTools::getHGradBasisCoeffs(projectedBasisNodesForComp,
180  evaluationValues,
181  transformedGradientData.getUnderlyingView(),
182  evaluationPointsRefSpace,
183  evaluationGradPointsRefSpace,
184  elementOrientations,
185  targetHGradBasis.get(),
186  &projectionStruct);
187  }
188  }
189  };
190 }
191 
192 #endif /* Intrepid2_ProjectedGeometry_h */
CellGeometry provides the nodes for a set of cells; has options that support efficient definition of ...
A stateless class for operations on cell data. Provides methods for:
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...
Header file for the Intrepid2::Experimental::ProjectionTools.
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 ...
static void mapToPhysicalFrame(Kokkos::DynRankView< physPointValueType, physPointProperties... > physPoints, const Kokkos::DynRankView< refPointValueType, refPointProperties... > refPoints, const WorksetType worksetCell, const HGradBasisPtrType basis)
Computes F, the reference-to-physical frame map.
A class providing static members to perform projection-based interpolations:
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...
Header file for the Intrepid2::CellTools class.
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.