Intrepid2
Intrepid2_CellToolsDefRefToPhys.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), or
38 // Mauro Perego (mperego@sandia.gov)
39 //
40 // ************************************************************************
41 // @HEADER
42 
43 
49 #ifndef __INTREPID2_CELLTOOLS_DEF_REF_TO_PHYS_HPP__
50 #define __INTREPID2_CELLTOOLS_DEF_REF_TO_PHYS_HPP__
51 
52 // disable clang warnings
53 #if defined (__clang__) && !defined (__INTEL_COMPILER)
54 #pragma clang system_header
55 #endif
56 
57 namespace Intrepid2 {
58 
59  //============================================================================================//
60  // //
61  // Reference-to-physical frame mapping and its inverse //
62  // //
63  //============================================================================================//
64 
65  namespace FunctorCellTools {
69  template<typename physPointViewType,
70  typename worksetCellType,
71  typename basisValType>
73  physPointViewType _physPoints;
74  const worksetCellType _worksetCells;
75  const basisValType _basisVals;
76 
77  KOKKOS_INLINE_FUNCTION
78  F_mapToPhysicalFrame( physPointViewType physPoints_,
79  worksetCellType worksetCells_,
80  basisValType basisVals_ )
81  : _physPoints(physPoints_), _worksetCells(worksetCells_), _basisVals(basisVals_) {}
82 
83  KOKKOS_INLINE_FUNCTION
84  void operator()(const size_type iter) const {
85  size_type cell, pt;
86  unrollIndex( cell, pt,
87  _physPoints.extent(0),
88  _physPoints.extent(1),
89  iter );
90  auto phys = Kokkos::subdynrankview( _physPoints, cell, pt, Kokkos::ALL());
91  const auto dofs = Kokkos::subdynrankview( _worksetCells, cell, Kokkos::ALL(), Kokkos::ALL());
92 
93  const auto valRank = _basisVals.rank();
94  const auto val = ( valRank == 2 ? Kokkos::subdynrankview( _basisVals, Kokkos::ALL(), pt) :
95  Kokkos::subdynrankview( _basisVals, cell, Kokkos::ALL(), pt));
96 
97  const ordinal_type dim = phys.extent(0);
98  const ordinal_type cardinality = val.extent(0);
99 
100  for (ordinal_type i=0;i<dim;++i) {
101  phys(i) = 0;
102  for (ordinal_type bf=0;bf<cardinality;++bf)
103  phys(i) += dofs(bf, i)*val(bf);
104  }
105  }
106  };
107  }
108 
109 /*
110  template<typename SpT>
111  template<typename physPointValueType, class ...physPointProperties,
112  typename refPointValueType, class ...refPointProperties,
113  typename worksetCellValueType, class ...worksetCellProperties>
114  void
115  CellTools<SpT>::
116  mapToPhysicalFrame( Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
117  const Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
118  const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
119  const shards::CellTopology cellTopo ) {
120 
121  auto basis = createHGradBasis<refPointValueType,refPointValueType>(cellTopo);
122  mapToPhysicalFrame(physPoints,
123  refPoints,
124  worksetCell,
125  basis);
126  }
127 */
128 
129  template<typename SpT>
130  template<typename physPointValueType, class ...physPointProperties,
131  typename refPointValueType, class ...refPointProperties,
132  typename worksetCellValueType, class ...worksetCellProperties,
133  typename HGradBasisPtrType>
134  void
136  mapToPhysicalFrame( Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
137  const Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
138  const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
139  const HGradBasisPtrType basis ) {
140 #ifdef HAVE_INTREPID2_DEBUG
141  CellTools_mapToPhysicalFrameArgs( physPoints, refPoints, worksetCell, basis->getBaseCellTopology() );
142 #endif
143  const auto cellTopo = basis->getBaseCellTopology();
144  const auto numCells = worksetCell.extent(0);
145 
146  //points can be rank-2 (P,D), or rank-3 (C,P,D)
147  const auto refPointRank = refPoints.rank();
148  const auto numPoints = (refPointRank == 2 ? refPoints.extent(0) : refPoints.extent(1));
149  const auto basisCardinality = basis->getCardinality();
150  auto vcprop = Kokkos::common_view_alloc_prop(physPoints);
151 
152  typedef Kokkos::DynRankView<physPointValueType,physPointProperties...> physPointViewType;
153  typedef Kokkos::DynRankView<decltype(basis->getDummyOutputValue()),SpT> valViewType;
154 
155  valViewType vals;
156 
157  switch (refPointRank) {
158  case 2: {
159  // refPoints is (P,D): single set of ref. points is mapped to one or multiple physical cells
160  vals = valViewType(Kokkos::view_alloc("CellTools::mapToPhysicalFrame::vals", vcprop), basisCardinality, numPoints);
161  basis->getValues(vals,
162  refPoints,
163  OPERATOR_VALUE);
164  break;
165  }
166  case 3: {
167  // refPoints is (C,P,D): multiple sets of ref. points are mapped to matching number of physical cells.
168  //vals = valViewType("CellTools::mapToPhysicalFrame::vals", numCells, basisCardinality, numPoints);
169  vals = valViewType(Kokkos::view_alloc("CellTools::mapToPhysicalFrame::vals", vcprop), numCells, basisCardinality, numPoints);
170  for (size_type cell=0;cell<numCells;++cell)
171  basis->getValues(Kokkos::subdynrankview( vals, cell, Kokkos::ALL(), Kokkos::ALL() ),
172  Kokkos::subdynrankview( refPoints, cell, Kokkos::ALL(), Kokkos::ALL() ),
173  OPERATOR_VALUE);
174  break;
175  }
176  }
177 
178  typedef Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCellViewType;
180  typedef typename ExecSpace<typename worksetCellViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
181 
182  const auto loopSize = physPoints.extent(0)*physPoints.extent(1);
183  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
184  Kokkos::parallel_for( policy, FunctorType(physPoints, worksetCell, vals) );
185  }
186 
187  template<typename SpT>
188  template<typename refSubcellPointValueType, class ...refSubcellPointProperties,
189  typename paramPointValueType, class ...paramPointProperties>
190  void
192  mapToReferenceSubcell( Kokkos::DynRankView<refSubcellPointValueType,refSubcellPointProperties...> refSubcellPoints,
193  const Kokkos::DynRankView<paramPointValueType,paramPointProperties...> paramPoints,
194  const ordinal_type subcellDim,
195  const ordinal_type subcellOrd,
196  const shards::CellTopology parentCell ) {
197 #ifdef HAVE_INTREPID2_DEBUG
198  INTREPID2_TEST_FOR_EXCEPTION( !hasReferenceCell(parentCell), std::invalid_argument,
199  ">>> ERROR (Intrepid2::CellTools::mapToReferenceSubcell): the specified cell topology does not have a reference cell.");
200 
201  INTREPID2_TEST_FOR_EXCEPTION( subcellDim != 1 &&
202  subcellDim != 2, std::invalid_argument,
203  ">>> ERROR (Intrepid2::CellTools::mapToReferenceSubcell): method defined only for 1 and 2-dimensional subcells.");
204 
205  INTREPID2_TEST_FOR_EXCEPTION( subcellOrd < 0 ||
206  subcellOrd >= static_cast<ordinal_type>(parentCell.getSubcellCount(subcellDim)), std::invalid_argument,
207  ">>> ERROR (Intrepid2::CellTools::mapToReferenceSubcell): subcell ordinal out of range.");
208 
209  // refSubcellPoints is rank-2 (P,D1), D1 = cell dimension
210  INTREPID2_TEST_FOR_EXCEPTION( refSubcellPoints.rank() != 2, std::invalid_argument,
211  ">>> ERROR (Intrepid2::CellTools::mapToReferenceSubcell): refSubcellPoints must have rank 2.");
212  INTREPID2_TEST_FOR_EXCEPTION( refSubcellPoints.extent(1) != parentCell.getDimension(), std::invalid_argument,
213  ">>> ERROR (Intrepid2::CellTools::mapToReferenceSubcell): refSubcellPoints dimension (1) does not match to parent cell dimension.");
214 
215  // paramPoints is rank-2 (P,D2) with D2 = subcell dimension
216  INTREPID2_TEST_FOR_EXCEPTION( paramPoints.rank() != 2, std::invalid_argument,
217  ">>> ERROR (Intrepid2::CellTools::mapToReferenceSubcell): paramPoints must have rank 2.");
218  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(paramPoints.extent(1)) != subcellDim, std::invalid_argument,
219  ">>> ERROR (Intrepid2::CellTools::mapToReferenceSubcell): paramPoints dimension (1) does not match to subcell dimension.");
220 
221  // cross check: refSubcellPoints and paramPoints: dimension 0 must match
222  INTREPID2_TEST_FOR_EXCEPTION( refSubcellPoints.extent(0) < paramPoints.extent(0), std::invalid_argument,
223  ">>> ERROR (Intrepid2::CellTools::mapToReferenceSubcell): refSubcellPoints dimension (0) does not match to paramPoints dimension(0).");
224 #endif
225 
226  if(!isSubcellParametrizationSet_)
227  setSubcellParametrization();
228 
229  INTREPID2_TEST_FOR_EXCEPTION( subcellDim != 1 &&
230  subcellDim != 2, std::invalid_argument,
231  ">>> ERROR (Intrepid2::CellTools::mapToReferenceSubcell): method defined only for 1 and 2-subcells");
232 
233 
234  // Get the subcell map, i.e., the coefficients of the parametrization function for the subcell
235  subcellParamViewType subcellMap;
236  getSubcellParametrization( subcellMap,
237  subcellDim,
238  parentCell );
239 
240  // Apply the parametrization map to every point in parameter domain
241  mapToReferenceSubcell( refSubcellPoints, paramPoints, subcellMap, subcellDim, subcellOrd, parentCell.getDimension());
242  }
243 
244 
245  template<typename SpT>
246  template<typename refSubcellPointValueType, class ...refSubcellPointProperties,
247  typename paramPointValueType, class ...paramPointProperties>
248  void
249  KOKKOS_INLINE_FUNCTION
251  mapToReferenceSubcell( Kokkos::DynRankView<refSubcellPointValueType,refSubcellPointProperties...> refSubcellPoints,
252  const Kokkos::DynRankView<paramPointValueType,paramPointProperties...> paramPoints,
253  const subcellParamViewType subcellMap,
254  const ordinal_type subcellDim,
255  const ordinal_type subcellOrd,
256  const ordinal_type parentCellDim ) {
257 
258  const ordinal_type numPts = paramPoints.extent(0);
259 
260  // Apply the parametrization map to every point in parameter domain
261  switch (subcellDim) {
262  case 2: {
263  for (ordinal_type pt=0;pt<numPts;++pt) {
264  const auto u = paramPoints(pt, 0);
265  const auto v = paramPoints(pt, 1);
266 
267  // map_dim(u,v) = c_0(dim) + c_1(dim)*u + c_2(dim)*v because both Quad and Tri ref faces are affine!
268  for (ordinal_type i=0;i<parentCellDim;++i)
269  refSubcellPoints(pt, i) = subcellMap(subcellOrd, i, 0) + ( subcellMap(subcellOrd, i, 1)*u +
270  subcellMap(subcellOrd, i, 2)*v );
271  }
272  break;
273  }
274  case 1: {
275  for (ordinal_type pt=0;pt<numPts;++pt) {
276  const auto u = paramPoints(pt, 0);
277  for (ordinal_type i=0;i<parentCellDim;++i)
278  refSubcellPoints(pt, i) = subcellMap(subcellOrd, i, 0) + ( subcellMap(subcellOrd, i, 1)*u );
279  }
280  break;
281  }
282  default: {}
283  }
284  }
285 }
286 
287 #endif
Functor for mapping reference points to physical frame see Intrepid2::CellTools for more...
static void mapToPhysicalFrame(Kokkos::DynRankView< physPointValueType, physPointProperties...> physPoints, const Kokkos::DynRankView< refPointValueType, refPointProperties...> refPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties...> worksetCell, const HGradBasisPtrType basis)
Computes F, the reference-to-physical frame map.
static void mapToReferenceSubcell(Kokkos::DynRankView< refSubcellPointValueType, refSubcellPointProperties...> refSubcellPoints, const Kokkos::DynRankView< paramPointValueType, paramPointProperties...> paramPoints, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell)
Computes parameterization maps of 1- and 2-subcells of reference cells.