Intrepid2
Intrepid2_OrientationToolsDefCoeffMatrix_HDIV.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 
58 #ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HDIV_HPP__
59 #define __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HDIV_HPP__
60 
61 // disable clang warnings
62 #if defined (__clang__) && !defined (__INTEL_COMPILER)
63 #pragma clang system_header
64 #endif
65 
66 namespace Intrepid2 {
67 
68  template<typename ExecSpaceType,
69  typename outputValueType,
70  typename pointValueType>
71  class Basis_HDIV_TET_In_FEM;
72 
73  namespace Impl {
74 
75  template<typename OutputViewType,
76  typename subcellBasisType,
77  typename cellBasisType>
78  inline
79  void
81  getCoeffMatrix_HDIV(OutputViewType &output,
82  const subcellBasisType subcellBasis,
83  const cellBasisType cellBasis,
84  const ordinal_type subcellId,
85  const ordinal_type subcellOrt) {
86  typedef typename OutputViewType::execution_space space_type;
87  typedef typename OutputViewType::value_type value_type;
88 
89  // with shards, everything should be computed on host space
90  typedef typename
91  Kokkos::Impl::is_space<space_type>::host_mirror_space::execution_space host_space_type;
92 
93  typedef Kokkos::DynRankView<value_type,host_space_type> DynRankViewHostType;
94 
95  //
96  // Topology
97  //
98 
99  // populate points on a subcell and map to subcell
100  const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
101  const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
102 
103  const ordinal_type cellDim = cellTopo.getDimension();
104  const ordinal_type subcellDim = subcellTopo.getDimension();
105 
106  INTREPID2_TEST_FOR_EXCEPTION( subcellDim >= cellDim,
107  std::logic_error,
108  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
109  "cellDim must be greater than subcellDim.");
110 
111  const auto subcellBaseKey = subcellTopo.getBaseKey();
112  const auto cellBaseKey = cellTopo.getBaseKey();
113 
114  INTREPID2_TEST_FOR_EXCEPTION( subcellBaseKey != shards::Line<>::key &&
115  subcellBaseKey != shards::Quadrilateral<>::key &&
116  subcellBaseKey != shards::Triangle<>::key,
117  std::logic_error,
118  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
119  "subcellBasis must have line, quad, or triangle topology.");
120 
121  INTREPID2_TEST_FOR_EXCEPTION( cellBaseKey != shards::Quadrilateral<>::key &&
122  cellBaseKey != shards::Triangle<>::key &&
123  cellBaseKey != shards::Hexahedron<>::key &&
124  cellBaseKey != shards::Tetrahedron<>::key,
125  std::logic_error,
126  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
127  "cellBasis must have quad, triangle, hexhedron or tetrahedron topology.");
128 
129  //
130  // Function space
131  //
132  {
133  const bool isHDIV = cellBasis.getFunctionSpace() == FUNCTION_SPACE_HDIV;
134  INTREPID2_TEST_FOR_EXCEPTION( !isHDIV,
135  std::logic_error,
136  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
137  "cellBasis is not HDIV.");
138  {
139  const bool subcellBasisIsHGRAD = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
140  const bool subcellBasisIsHVOL = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HVOL;
141  const bool cellIsTri = cellBaseKey == shards::Triangle<>::key;
142  const bool cellIsTet = cellBaseKey == shards::Tetrahedron<>::key;
143  const bool cellIsHex = cellBaseKey == shards::Hexahedron<>::key;
144  const bool cellIsQuad = cellBaseKey == shards::Quadrilateral<>::key;
145 
146  const std::string subcellBasisName(subcellBasis.getName());
147  switch (subcellDim) {
148  case 1: {
149  //TODO: Hex, QUAD, TET and TRI element should have the same 1d basis
150  if (cellIsHex || cellIsQuad) {
151  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHGRAD,
152  std::logic_error,
153  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_DIV): "
154  "subcellBasis function space (1d) is not consistent to cellBasis, which should be open line hgrad, order -1.");
155  } else if (cellIsTet || cellIsTri) {
156  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHVOL,
157  std::logic_error,
158  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_DIV): "
159  "subcellBasis function space (1d) is not consistent to cellBasis, which should be HVOL line, order -1.");
160  }
161  break;
162  }
163  case 2: {
164  if (subcellBaseKey == shards::Quadrilateral<>::key) {
165  // quad face basis is tensor product of open line basis functions
166  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHGRAD,
167  std::logic_error,
168  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
169  "subcellBasis function space is not compatible, which should be open line hgrad, order -1.");
170  } else if (subcellBaseKey == shards::Triangle<>::key) {
171  // triangle face basis comes from HVOL basis
172  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHVOL,
173  std::logic_error,
174  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
175  "subcellBasis function space is not compatible, which should HVOL, order-1.");
176  }
177  break;
178  }
179  }
180  }
181  }
182 
183  //
184  // Collocation points
185  //
186  const ordinal_type numCellBasis = cellBasis.getCardinality();
187  const ordinal_type numSubcellBasis = subcellBasis.getCardinality();
188 
189  const ordinal_type ordSubcell = cellBasis.getDofOrdinal(subcellDim, subcellId, 0);
190  INTREPID2_TEST_FOR_EXCEPTION( ordSubcell == -1,
191  std::logic_error,
192  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
193  "Invalid subcellId returns -1 ordSubcell.");
194 
195  const ordinal_type ndofSubcell = cellBasis.getDofTag(ordSubcell)(3);
196 
197  // reference points on a subcell
198  DynRankViewHostType refPtsSubcell("refPtsSubcell", ndofSubcell, subcellDim);
199  DynRankViewHostType subcellDofCoords("subcellDofCoords", numSubcellBasis, subcellDim);
200  subcellBasis.getDofCoords(subcellDofCoords);
201  for(ordinal_type i=0; i<ndofSubcell; ++i)
202  for(ordinal_type d=0; d <subcellDim; ++d)
203  refPtsSubcell(i,d) = subcellDofCoords(subcellBasis.getDofOrdinal(subcellDim, 0, i),d);
204 
205  // modified points with orientation
206  DynRankViewHostType ortPtsSubcell("ortPtsSubcell", ndofSubcell, subcellDim);
208  refPtsSubcell,
209  subcellTopo,
210  subcellOrt);
211 
212  // map to reference coordinates
213  DynRankViewHostType refPtsCell("refPtsCell", ndofSubcell, cellDim);
215  ortPtsSubcell,
216  subcellDim,
217  subcellId,
218  cellTopo);
219  //
220  // Basis evaluation on the collocation points
221  //
222  // evaluate values on the reference cell
223  DynRankViewHostType refValues("refValues", numCellBasis, ndofSubcell, cellDim);
224  cellBasis.getValues(refValues, refPtsCell, OPERATOR_VALUE);
225 
226 
227  DynRankViewHostType ortJacobian("ortJacobian", subcellDim, subcellDim);
228  Impl::OrientationTools::getJacobianOfOrientationMap(ortJacobian, subcellTopo, subcellOrt);
229  auto ortJacobianDet = (subcellDim == 1) ? ortJacobian(0,0) :
230  ortJacobian(0,0)*ortJacobian(1,1)-ortJacobian(1,0)*ortJacobian(0,1);
231 
232  //
233  // Restrict vector valued basis functions on the subcell dimensions
234  //
235  DynRankViewHostType sideNormal("sideNormal", cellDim);
236  CellTools<host_space_type>::getReferenceSideNormal(sideNormal, subcellId, cellTopo);
237 
238  if((cellBaseKey == shards::Triangle<>::key))
239  for(ordinal_type i=0; i<cellDim; ++i)
240  sideNormal(i) *= 2; //scale by reference edge length
241 
242  Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type>
243  refMat("refMat", ndofSubcell, ndofSubcell),
244  ortMat("ortMat", ndofSubcell, ndofSubcell);
245  DynRankViewHostType tmpValues("tmpValues", numCellBasis, ndofSubcell);
246  for (ordinal_type i=0;i<ndofSubcell;++i) {
247  const ordinal_type ii = cellBasis.getDofOrdinal(subcellDim, subcellId, i);
248  for (ordinal_type j=0;j<ndofSubcell;++j)
249  for (ordinal_type k=0;k<cellDim;++k) {
250  refMat(i,j) += refValues(ii,j,k)*sideNormal(k)* ortJacobianDet;
251  ortMat(i,j) = (i==j); //ortMat is the identity thanks to the kronecker property of the bases
252  }
253  }
254 
255  // solve the system
256  {
257  Teuchos::LAPACK<ordinal_type,value_type> lapack;
258  ordinal_type info = 0;
259  Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type> pivVec("pivVec", ndofSubcell, 1);
260 
261  lapack.GESV(ndofSubcell, ndofSubcell,
262  refMat.data(),
263  refMat.stride_1(),
264  (ordinal_type*)pivVec.data(),
265  ortMat.data(),
266  ortMat.stride_1(),
267  &info);
268 
269  if (info) {
270  std::stringstream ss;
271  ss << ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
272  << "LAPACK return with error code: "
273  << info;
274  INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
275  }
276 
277  // clean up numerical noise
278  {
279  const double eps = threshold();
280  for (ordinal_type i=0;i<ndofSubcell;++i)
281  for (ordinal_type j=0;j<ndofSubcell;++j) {
282  auto intOrtMat = std::round(ortMat(i,j));
283  ortMat(i,j) = (std::abs(ortMat(i,j) - std::round(ortMat(i,j))) < eps) ? intOrtMat : ortMat(i,j);
284  }
285  }
286  }
287 
288  {
289  // move the data to original device memory
290  const Kokkos::pair<ordinal_type,ordinal_type> range(0, ndofSubcell);
291  Kokkos::deep_copy(Kokkos::subview(output, range, range),
292  Kokkos::subview(ortMat, range, range));
293  }
294  }
295  }
296 
297 }
298 #endif
static void getCoeffMatrix_HDIV(OutputViewType &output, const subcellBasisType subcellBasis, const cellBasisType cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
Compute coefficient matrix for HDIV by collocating point values.
static void getReferenceSideNormal(Kokkos::DynRankView< refSideNormalValueType, refSideNormalProperties...> refSideNormal, const ordinal_type sideOrd, const shards::CellTopology parentCell)
Computes constant normal vectors to sides of 2D or 3D reference cells.
static void getJacobianOfOrientationMap(JacobianViewType jacobian, const shards::CellTopology cellTopo, const ordinal_type cellOrt)
Computes jacobian of the parameterization maps of 1- and 2-subcells with orientation.
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.
static void mapToModifiedReference(outPointViewType outPoints, const refPointViewType refPoints, const shards::CellTopology cellTopo, const ordinal_type cellOrt=0)
Computes modified parameterization maps of 1- and 2-subcells with orientation.