Intrepid2
Intrepid2_CellToolsDefPhysToRef.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_PHYS_TO_REF_HPP__
50 #define __INTREPID2_CELLTOOLS_DEF_PHYS_TO_REF_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  // //
62  // Reference-to-physical frame mapping and its inverse //
63  // //
64  //============================================================================================//
65 
66  template<typename SpT>
67  template<typename refPointValueType, class ...refPointProperties,
68  typename physPointValueType, class ...physPointProperties,
69  typename worksetCellValueType, class ...worksetCellProperties>
70  void
72  mapToReferenceFrame( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
73  const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
74  const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
75  const shards::CellTopology cellTopo ) {
76 #ifdef HAVE_INTREPID2_DEBUG
77  CellTools_mapToReferenceFrameArgs(refPoints, physPoints, worksetCell, cellTopo);
78 #endif
79  typedef RealSpaceTools<SpT> rst;
80  typedef Kokkos::DynRankView<typename ScalarTraits<refPointValueType>::scalar_type,SpT> refPointViewSpType;
81 
82  const auto spaceDim = cellTopo.getDimension();
83  refPointViewSpType
84  cellCenter("CellTools::mapToReferenceFrame::cellCenter", spaceDim),
85  cellVertex("CellTools::mapToReferenceFrame::cellCenter", spaceDim);
86  getReferenceCellCenter(cellCenter, cellVertex, cellTopo);
87 
88  // Default: map (C,P,D) array of physical pt. sets to (C,P,D) array. Requires (C,P,D) initial guess.
89  const auto numCells = worksetCell.extent(0);
90  const auto numPoints = physPoints.extent(1);
91 
92  // init guess is created locally and non fad whatever refpoints type is
93  using result_layout = typename DeduceLayout< decltype(refPoints) >::result_layout;
94  using device_type = typename decltype(refPoints)::device_type;
95  auto vcprop = Kokkos::common_view_alloc_prop(refPoints);
96  using common_value_type = typename decltype(vcprop)::value_type;
97  Kokkos::DynRankView< common_value_type, result_layout, device_type > initGuess ( Kokkos::view_alloc("CellTools::mapToReferenceFrame::initGuess", vcprop), numCells, numPoints, spaceDim );
98  //refPointViewSpType initGuess("CellTools::mapToReferenceFrame::initGuess", numCells, numPoints, spaceDim);
99  rst::clone(initGuess, cellCenter);
100 
101  mapToReferenceFrameInitGuess(refPoints, initGuess, physPoints, worksetCell, cellTopo);
102  }
103 
104 
105  template<typename SpT>
106  template<typename refPointValueType, class ...refPointProperties,
107  typename initGuessValueType, class ...initGuessProperties,
108  typename physPointValueType, class ...physPointProperties,
109  typename worksetCellValueType, class ...worksetCellProperties,
110  typename HGradBasisPtrType>
111  void
113  mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
114  const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
115  const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
116  const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
117  const HGradBasisPtrType basis ) {
118 #ifdef HAVE_INTREPID2_DEBUG
119  CellTools_mapToReferenceFrameInitGuessArgs(refPoints, initGuess, physPoints, worksetCell,
120  basis->getBaseCellTopology());
121 
122 #endif
123  const auto cellTopo = basis->getBaseCellTopology();
124  const auto spaceDim = cellTopo.getDimension();
125 
126  // Default: map (C,P,D) array of physical pt. sets to (C,P,D) array.
127  // Requires (C,P,D) temp arrays and (C,P,D,D) Jacobians.
128  const auto numCells = worksetCell.extent(0);
129  const auto numPoints = physPoints.extent(1);
130 
131  typedef RealSpaceTools<SpT> rst;
132  const auto tol = tolerence();
133 
134  using result_layout = typename DeduceLayout< decltype(refPoints) >::result_layout;
135  using device_type = typename decltype(refPoints)::device_type;
136  auto vcprop = Kokkos::common_view_alloc_prop(refPoints);
137  typedef Kokkos::DynRankView<typename decltype(vcprop)::value_type, result_layout, device_type > viewType;
138 
139  // Temp arrays for Newton iterates and Jacobians. Resize according to rank of ref. point array
140  viewType xOld(Kokkos::view_alloc("CellTools::mapToReferenceFrameInitGuess::xOld", vcprop), numCells, numPoints, spaceDim);
141  viewType xTmp(Kokkos::view_alloc("CellTools::mapToReferenceFrameInitGuess::xTmp", vcprop), numCells, numPoints, spaceDim);
142 
143  // deep copy may not work with FAD but this is right thing to do as it can move data between devices
144  Kokkos::deep_copy(xOld, initGuess);
145 
146  // jacobian should select fad dimension between xOld and worksetCell as they are input; no front interface yet
147  auto vcpropJ = Kokkos::common_view_alloc_prop(refPoints, worksetCell);
148  typedef Kokkos::DynRankView<typename decltype(vcpropJ)::value_type, result_layout, device_type > viewTypeJ;
149  viewTypeJ jacobian(Kokkos::view_alloc("CellTools::mapToReferenceFrameInitGuess::jacobian", vcpropJ), numCells, numPoints, spaceDim, spaceDim);
150  viewTypeJ jacobianInv(Kokkos::view_alloc("CellTools::mapToReferenceFrameInitGuess::jacobianInv", vcpropJ), numCells, numPoints, spaceDim, spaceDim);
151 
152  typedef Kokkos::DynRankView<typename ScalarTraits<refPointValueType>::scalar_type,SpT> errorViewType;
153  errorViewType
154  xScalarTmp ("CellTools::mapToReferenceFrameInitGuess::xScalarTmp", numCells, numPoints, spaceDim),
155  errorPointwise("CellTools::mapToReferenceFrameInitGuess::errorPointwise", numCells, numPoints),
156  errorCellwise ("CellTools::mapToReferenceFrameInitGuess::errorCellwise", numCells);
157 
158  //auto errorPointwise = Kokkos::createDynRankView(xTmp, "CellTools::mapToReferenceFrameInitGuess::errorPointwise", numCells, numPoints);
159  //auto errorCellwise = Kokkos::createDynRankView(xTmp, "CellTools::mapToReferenceFrameInitGuess::errorCellwise", numCells);
160 
161  // Newton method to solve the equation F(refPoints) - physPoints = 0:
162  // refPoints = xOld - DF^{-1}(xOld)*(F(xOld) - physPoints) = xOld + DF^{-1}(xOld)*(physPoints - F(xOld))
163  for (ordinal_type iter=0;iter<Parameters::MaxNewton;++iter) {
164 
165  // Jacobians at the old iterates and their inverses.
166  setJacobian(jacobian, xOld, worksetCell, basis);
167  setJacobianInv(jacobianInv, jacobian);
168 
169  // The Newton step.
170  mapToPhysicalFrame(xTmp, xOld, worksetCell, basis); // xTmp <- F(xOld)
171  rst::subtract(xTmp, physPoints, xTmp); // xTmp <- physPoints - F(xOld)
172  rst::matvec(refPoints, jacobianInv, xTmp); // refPoints <- DF^{-1}( physPoints - F(xOld) )
173  rst::add(refPoints, xOld); // refPoints <- DF^{-1}( physPoints - F(xOld) ) + xOld
174 
175  // l2 error (Euclidean distance) between old and new iterates: |xOld - xNew|
176  rst::subtract(xTmp, xOld, refPoints);
177 
178  // extract values
179  rst::extractScalarValues(xScalarTmp, xTmp);
180  rst::vectorNorm(errorPointwise, xScalarTmp, NORM_TWO);
181 
182  // Average L2 error for a multiple sets of physical points: error is rank-2 (C,P) array
183  rst::vectorNorm(errorCellwise, errorPointwise, NORM_ONE);
184  const auto errorTotal = rst::Serial::vectorNorm(errorCellwise, NORM_ONE);
185 
186  // Stopping criterion:
187  if (errorTotal < tol)
188  break;
189 
190  // initialize next Newton step ( this is not device friendly )
191  Kokkos::deep_copy(xOld, refPoints);
192  }
193  }
194 
195 }
196 
197 #endif
static constexpr ordinal_type MaxNewton
Maximum number of Newton iterations used internally in methods such as computing the action of the in...
static void mapToReferenceFrameInitGuess(Kokkos::DynRankView< refPointValueType, refPointProperties...> refPoints, const Kokkos::DynRankView< initGuessValueType, initGuessProperties...> initGuess, const Kokkos::DynRankView< physPointValueType, physPointProperties...> physPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties...> worksetCell, const HGradBasisPtrType basis)
Computation of , the inverse of the reference-to-physical frame map using user-supplied initial guess...
static void mapToReferenceFrame(Kokkos::DynRankView< refPointValueType, refPointProperties...> refPoints, const Kokkos::DynRankView< physPointValueType, physPointProperties...> physPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties...> worksetCell, const shards::CellTopology cellTopo)
Computes , the inverse of the reference-to-physical frame map using a default initial guess...
Implementation of basic linear algebra functionality in Euclidean space.