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 std::string cellBasisName(cellBasis.getName());
134  INTREPID2_TEST_FOR_EXCEPTION( cellBasisName.find("HDIV") == std::string::npos,
135  std::logic_error,
136  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
137  "cellBasis is not HDIV.");
138  {
139  const std::string subcellBasisName(subcellBasis.getName());
140  switch (subcellDim) {
141  case 1: {
142  //TODO: Hex, QUAD, TET and TRI element should have the same 1d basis
143  if ((cellBasisName.find("HEX") != std::string::npos) || (cellBasisName.find("QUAD") != std::string::npos)) {
144  INTREPID2_TEST_FOR_EXCEPTION( subcellBasisName.find("HGRAD") == std::string::npos,
145  std::logic_error,
146  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_DIV): "
147  "subcellBasis function space (1d) is not consistent to cellBasis, which should be open line hgrad, order -1.");
148  } else if ((cellBasisName.find("TET") != std::string::npos) || (cellBasisName.find("TRI") != std::string::npos)) {
149  INTREPID2_TEST_FOR_EXCEPTION( subcellBasisName.find("HVOL") == std::string::npos,
150  std::logic_error,
151  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_DIV): "
152  "subcellBasis function space (1d) is not consistent to cellBasis, which should be HVOL line, order -1.");
153  }
154  break;
155  }
156  case 2: {
157  if (subcellBaseKey == shards::Quadrilateral<>::key) {
158  // quad face basis is tensor product of open line basis functions
159  INTREPID2_TEST_FOR_EXCEPTION( subcellBasisName.find("HGRAD") == std::string::npos,
160  std::logic_error,
161  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
162  "subcellBasis function space is not compatible, which should be open line hgrad, order -1.");
163  } else if (subcellBaseKey == shards::Triangle<>::key) {
164  // triangle face basis comes from HVOL basis
165  INTREPID2_TEST_FOR_EXCEPTION( subcellBasisName.find("HVOL") == std::string::npos,
166  std::logic_error,
167  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
168  "subcellBasis function space is not compatible, which should HVOL, order-1.");
169  }
170  break;
171  }
172  }
173  }
174  }
175 
176  //
177  // Collocation points
178  //
179  const ordinal_type numCellBasis = cellBasis.getCardinality();
180  const ordinal_type numSubcellBasis = subcellBasis.getCardinality();
181 
182  const ordinal_type ordSubcell = cellBasis.getDofOrdinal(subcellDim, subcellId, 0);
183  INTREPID2_TEST_FOR_EXCEPTION( ordSubcell == -1,
184  std::logic_error,
185  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
186  "Invalid subcellId returns -1 ordSubcell.");
187 
188  const ordinal_type ndofSubcell = cellBasis.getDofTag(ordSubcell)(3);
189 
190  // reference points on a subcell
191  DynRankViewHostType refPtsSubcell("refPtsSubcell", ndofSubcell, subcellDim);
192  DynRankViewHostType subcellDofCoords("subcellDofCoords", numSubcellBasis, subcellDim);
193  subcellBasis.getDofCoords(subcellDofCoords);
194  for(ordinal_type i=0; i<ndofSubcell; ++i)
195  for(ordinal_type d=0; d <subcellDim; ++d)
196  refPtsSubcell(i,d) = subcellDofCoords(subcellBasis.getDofOrdinal(subcellDim, 0, i),d);
197 
198  // modified points with orientation
199  DynRankViewHostType ortPtsSubcell("ortPtsSubcell", ndofSubcell, subcellDim);
201  refPtsSubcell,
202  subcellTopo,
203  subcellOrt);
204 
205  // map to reference coordinates
206  DynRankViewHostType refPtsCell("refPtsCell", ndofSubcell, cellDim);
208  ortPtsSubcell,
209  subcellDim,
210  subcellId,
211  cellTopo);
212  //
213  // Basis evaluation on the collocation points
214  //
215  // evaluate values on the reference cell
216  DynRankViewHostType refValues("refValues", numCellBasis, ndofSubcell, cellDim);
217  cellBasis.getValues(refValues, refPtsCell, OPERATOR_VALUE);
218 
219 
220  DynRankViewHostType ortJacobian("ortJacobian", subcellDim, subcellDim);
221  Impl::OrientationTools::getJacobianOfOrientationMap(ortJacobian, subcellTopo, subcellOrt);
222  auto ortJacobianDet = (subcellDim == 1) ? ortJacobian(0,0) :
223  ortJacobian(0,0)*ortJacobian(1,1)-ortJacobian(1,0)*ortJacobian(0,1);
224 
225  //
226  // Restrict vector valued basis functions on the subcell dimensions
227  //
228  DynRankViewHostType sideNormal("sideNormal", cellDim);
229  CellTools<host_space_type>::getReferenceSideNormal(sideNormal, subcellId, cellTopo);
230 
231  if((cellBaseKey == shards::Triangle<>::key))
232  for(ordinal_type i=0; i<cellDim; ++i)
233  sideNormal(i) *= 2; //scale by reference edge length
234 
235  Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type>
236  refMat("refMat", ndofSubcell, ndofSubcell),
237  ortMat("ortMat", ndofSubcell, ndofSubcell);
238  DynRankViewHostType tmpValues("tmpValues", numCellBasis, ndofSubcell);
239  for (ordinal_type i=0;i<ndofSubcell;++i) {
240  const ordinal_type ii = cellBasis.getDofOrdinal(subcellDim, subcellId, i);
241  for (ordinal_type j=0;j<ndofSubcell;++j)
242  for (ordinal_type k=0;k<cellDim;++k) {
243  refMat(i,j) += refValues(ii,j,k)*sideNormal(k)* ortJacobianDet;
244  ortMat(i,j) = (i==j); //ortMat is the identity thanks to the kronecker property of the bases
245  }
246  }
247 
248  // solve the system
249  {
250  Teuchos::LAPACK<ordinal_type,value_type> lapack;
251  ordinal_type info = 0;
252  Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type> pivVec("pivVec", ndofSubcell, 1);
253 
254  lapack.GESV(ndofSubcell, ndofSubcell,
255  refMat.data(),
256  refMat.stride_1(),
257  (ordinal_type*)pivVec.data(),
258  ortMat.data(),
259  ortMat.stride_1(),
260  &info);
261 
262  if (info) {
263  std::stringstream ss;
264  ss << ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
265  << "LAPACK return with error code: "
266  << info;
267  INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
268  }
269 
270  // clean up numerical noise
271  {
272  const double eps = threshold();
273  for (ordinal_type i=0;i<ndofSubcell;++i)
274  for (ordinal_type j=0;j<ndofSubcell;++j) {
275  auto intOrtMat = std::round(ortMat(i,j));
276  ortMat(i,j) = (std::abs(ortMat(i,j) - std::round(ortMat(i,j))) < eps) ? intOrtMat : ortMat(i,j);
277  }
278  }
279  }
280 
281  {
282  // move the data to original device memory
283  const Kokkos::pair<ordinal_type,ordinal_type> range(0, ndofSubcell);
284  Kokkos::deep_copy(Kokkos::subview(output, range, range),
285  Kokkos::subview(ortMat, range, range));
286  }
287  }
288  }
289 
290 }
291 #endif
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 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 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.