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 
59 #ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HDIV_HPP__
60 #define __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HDIV_HPP__
61 
62 // disable clang warnings
63 #if defined (__clang__) && !defined (__INTEL_COMPILER)
64 #pragma clang system_header
65 #endif
66 
67 namespace Intrepid2 {
68 
69 namespace Impl {
70 
71 namespace Debug {
72 
73 #ifdef HAVE_INTREPID2_DEBUG
74 template<typename subcellBasisType,
75 typename cellBasisType>
76 inline
77 void
78 check_getCoeffMatrix_HDIV(const subcellBasisType& subcellBasis,
79  const cellBasisType& cellBasis,
80  const ordinal_type subcellId,
81  const ordinal_type subcellOrt) {
82  const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
83  const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
84 
85  const ordinal_type cellDim = cellTopo.getDimension();
86  const ordinal_type subcellDim = subcellTopo.getDimension();
87 
88  INTREPID2_TEST_FOR_EXCEPTION( subcellDim >= cellDim,
89  std::logic_error,
90  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
91  "cellDim must be greater than subcellDim.");
92 
93  const auto subcellBaseKey = subcellTopo.getBaseKey();
94  const auto cellBaseKey = cellTopo.getBaseKey();
95 
96  INTREPID2_TEST_FOR_EXCEPTION( subcellBaseKey != shards::Line<>::key &&
97  subcellBaseKey != shards::Quadrilateral<>::key &&
98  subcellBaseKey != shards::Triangle<>::key,
99  std::logic_error,
100  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
101  "subcellBasis must have line, quad, or triangle topology.");
102 
103  INTREPID2_TEST_FOR_EXCEPTION( cellBaseKey != shards::Quadrilateral<>::key &&
104  cellBaseKey != shards::Triangle<>::key &&
105  cellBaseKey != shards::Hexahedron<>::key &&
106  cellBaseKey != shards::Wedge<>::key &&
107  cellBaseKey != shards::Tetrahedron<>::key &&
108  cellBaseKey != shards::Pyramid<>::key,
109  std::logic_error,
110  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
111  "cellBasis must have quad, triangle, hexhedron or tetrahedron topology.");
112 
113  //
114  // Function space
115  //
116  {
117  const bool isHDIV = cellBasis.getFunctionSpace() == FUNCTION_SPACE_HDIV;
118  INTREPID2_TEST_FOR_EXCEPTION( !isHDIV,
119  std::logic_error,
120  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
121  "cellBasis is not HDIV.");
122  {
123  const bool subcellBasisIsHGRAD = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
124  const bool subcellBasisIsHVOL = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HVOL;
125  const bool cellIsTri = cellBaseKey == shards::Triangle<>::key;
126  const bool cellIsTet = cellBaseKey == shards::Tetrahedron<>::key;
127  const bool cellIsHex = cellBaseKey == shards::Hexahedron<>::key;
128  const bool cellIsQuad = cellBaseKey == shards::Quadrilateral<>::key;
129 
130  switch (subcellDim) {
131  case 1: {
132  //TODO: Hex, QUAD, TET and TRI element should have the same 1d basis
133  if (cellIsHex || cellIsQuad) {
134  INTREPID2_TEST_FOR_EXCEPTION( !(subcellBasisIsHGRAD||subcellBasisIsHVOL),
135  std::logic_error,
136  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_DIV): "
137  "subcellBasis function space (1d) is not consistent to cellBasis, which should be open line hgrad, order -1.");
138  } else if (cellIsTet || cellIsTri) {
139  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHVOL,
140  std::logic_error,
141  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_DIV): "
142  "subcellBasis function space (1d) is not consistent to cellBasis, which should be HVOL line, order -1.");
143  }
144  break;
145  }
146  case 2: {
147  if (subcellBaseKey == shards::Quadrilateral<>::key) {
148  // quad face basis is tensor product of open line basis functions
149  INTREPID2_TEST_FOR_EXCEPTION( !(subcellBasisIsHGRAD||subcellBasisIsHVOL),
150  std::logic_error,
151  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
152  "subcellBasis function space is not compatible, which should be open line hgrad, order -1.");
153  } else if (subcellBaseKey == shards::Triangle<>::key) {
154  // triangle face basis comes from HVOL basis
155  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHVOL,
156  std::logic_error,
157  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
158  "subcellBasis function space is not compatible, which should HVOL, order-1.");
159  }
160  break;
161  }
162  }
163  }
164  }
165 }
166 #endif
167 } //Debug Namespace
168 
169 template<typename OutputViewType,
170 typename subcellBasisHostType,
171 typename cellBasisHostType>
172 inline
173 void
175 getCoeffMatrix_HDIV(OutputViewType &output,
176  const subcellBasisHostType& subcellBasis,
177  const cellBasisHostType& cellBasis,
178  const ordinal_type subcellId,
179  const ordinal_type subcellOrt,
180  const bool inverse) {
181 
182 #ifdef HAVE_INTREPID2_DEBUG
183  Debug::check_getCoeffMatrix_HDIV(subcellBasis,cellBasis,subcellId,subcellOrt);
184 #endif
185 
186  using value_type = typename OutputViewType::non_const_value_type;
187  using host_device_type = Kokkos::HostSpace::device_type;
188 
189  //
190  // Collocation points
191  //
192  const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
193  const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
194  const ordinal_type cellDim = cellTopo.getDimension();
195  const ordinal_type subcellDim = subcellTopo.getDimension();
196  const auto subcellBaseKey = subcellTopo.getBaseKey();
197  const ordinal_type numCellBasis = cellBasis.getCardinality();
198  const ordinal_type numSubcellBasis = subcellBasis.getCardinality();
199  const ordinal_type ndofSubcell = cellBasis.getDofCount(subcellDim,subcellId);
200 
201  //
202  // Reference points
203  //
204 
205  //use lattice to compute reference subcell points \xi_j
206  auto latticeDegree = (subcellBaseKey == shards::Triangle<>::key) ?
207  cellBasis.getDegree()+2 : cellBasis.getDegree()+1;
208 
209  // Reference points \xi_j on the subcell
210  Kokkos::DynRankView<value_type,host_device_type> refPtsSubcell("refPtsSubcell", ndofSubcell, subcellDim);
211  auto latticeSize=PointTools::getLatticeSize(subcellTopo, latticeDegree, 1);
212  INTREPID2_TEST_FOR_EXCEPTION( latticeSize != ndofSubcell,
213  std::logic_error,
214  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
215  "Lattice size should be equal to the onber of subcell internal DoFs");
216  PointTools::getLattice(refPtsSubcell, subcellTopo, latticeDegree, 1);//, POINTTYPE_WARPBLEND);
217 
218  // evaluate values on the modified cell
219  auto subcellParam = Intrepid2::RefSubcellParametrization<host_device_type>::get(subcellDim, cellTopo.getKey());
220 
221  // refPtsCell = F_s (\eta_o (refPtsSubcell))
222  Kokkos::DynRankView<value_type,host_device_type> refPtsCell("refPtsCell", ndofSubcell, cellDim);
223  // map points from the subcell manifold into the cell one
224  mapSubcellCoordsToRefCell(refPtsCell,refPtsSubcell, subcellParam, subcellBaseKey, subcellId, subcellOrt);
225 
226  //computing normal to the subcell accounting for orientation
227  Kokkos::DynRankView<value_type,host_device_type> tangentsAndNormal("trJacobianF", cellDim, cellDim );
228  OrientationTools::getRefSideTangentsAndNormal(tangentsAndNormal, subcellParam, subcellBaseKey, subcellId, subcellOrt);
229  auto sideNormal = Kokkos::subview(tangentsAndNormal, cellDim-1, Kokkos::ALL());
230 
231 
232  //
233  // Basis evaluation on the collocation points
234  //
235 
236  // cellBasisValues = \psi_k(F_s (\eta_o (\xi_j)))
237  Kokkos::DynRankView<value_type,host_device_type> cellBasisValues("cellBasisValues", numCellBasis, ndofSubcell, cellDim);
238  cellBasis.getValues(cellBasisValues, refPtsCell, OPERATOR_VALUE);
239 
240  // subcellBasisValues = \phi_i (\xi_j)
241  Kokkos::DynRankView<value_type,host_device_type> subCellValues("subCellValues", numSubcellBasis, ndofSubcell);
242  subcellBasis.getValues(subCellValues, refPtsSubcell, OPERATOR_VALUE);
243 
244  //
245  // Compute Psi_jk = \psi_k(F_s (\eta_o (\xi_j))) \cdot (n_s det(J_\eta))
246  // and Phi_ji = \phi_i (\xi_j), and solve
247  // Psi A^T = Phi
248  //
249 
250  // construct Psi and Phi matrices. LAPACK wants left layout
251  Kokkos::DynRankView<value_type,Kokkos::LayoutLeft,host_device_type>
252  PsiMat("PsiMat", ndofSubcell, ndofSubcell),
253  PhiMat("PhiMat", ndofSubcell, ndofSubcell),
254  RefMat,
255  OrtMat;
256 
257  auto cellTagToOrdinal = cellBasis.getAllDofOrdinal();
258  auto subcellTagToOrdinal = subcellBasis.getAllDofOrdinal();
259 
260  for (ordinal_type i=0;i<ndofSubcell;++i) {
261  const ordinal_type ic = cellTagToOrdinal(subcellDim, subcellId, i);
262  const ordinal_type isc = subcellTagToOrdinal(subcellDim, 0, i);
263  for (ordinal_type j=0;j<ndofSubcell;++j) {
264  PhiMat(j,i) = get_scalar_value(subCellValues(isc,j));
265  for (ordinal_type k=0;k<cellDim;++k)
266  PsiMat(j,i) += get_scalar_value(cellBasisValues(ic,j,k))*sideNormal(k);
267  }
268  }
269 
270  RefMat = inverse ? PhiMat : PsiMat;
271  OrtMat = inverse ? PsiMat : PhiMat;
272 
273  // solve the system
274  {
275  Teuchos::LAPACK<ordinal_type,value_type> lapack;
276  ordinal_type info = 0;
277  Kokkos::DynRankView<ordinal_type,host_device_type> pivVec("pivVec", ndofSubcell);
278 
279  lapack.GESV(ndofSubcell, ndofSubcell,
280  RefMat.data(),
281  RefMat.stride_1(),
282  pivVec.data(),
283  OrtMat.data(),
284  OrtMat.stride_1(),
285  &info);
286 
287  if (info) {
288  std::stringstream ss;
289  ss << ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
290  << "LAPACK return with error code: "
291  << info;
292  INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
293  }
294 
295  //After solving the system w/ LAPACK, Phi contains A^T
296 
297  // transpose and clean up numerical noise (for permutation matrices)
298  const double eps = tolerence();
299  for (ordinal_type i=0;i<ndofSubcell;++i) {
300  auto intmatii = std::round(OrtMat(i,i));
301  OrtMat(i,i) = (std::abs(OrtMat(i,i) - intmatii) < eps) ? intmatii : OrtMat(i,i);
302  for (ordinal_type j=i+1;j<ndofSubcell;++j) {
303  auto matij = OrtMat(i,j);
304 
305  auto intmatji = std::round(OrtMat(j,i));
306  OrtMat(i,j) = (std::abs(OrtMat(j,i) - intmatji) < eps) ? intmatji : OrtMat(j,i);
307 
308  auto intmatij = std::round(matij);
309  OrtMat(j,i) = (std::abs(matij - intmatij) < eps) ? intmatij : matij;
310  }
311  }
312 
313  // Print Matrix A
314  /*
315  {
316  std::cout << "|";
317  for (ordinal_type i=0;i<ndofSubcell;++i) {
318  for (ordinal_type j=0;j<ndofSubcell;++j) {
319  std::cout << OrtMat(i,j) << " ";
320  }
321  std::cout << "| ";
322  }
323  std::cout <<std::endl;
324  }
325  */
326 
327  }
328 
329  {
330  // move the data to original device memory
331  const Kokkos::pair<ordinal_type,ordinal_type> range(0, ndofSubcell);
332  auto suboutput = Kokkos::subview(output, range, range);
333  auto tmp = Kokkos::create_mirror_view_and_copy(typename OutputViewType::device_type::memory_space(), OrtMat);
334  Kokkos::deep_copy(suboutput, tmp);
335  }
336 }
337 }
338 
339 }
340 #endif
static KOKKOS_INLINE_FUNCTION void mapSubcellCoordsToRefCell(coordsViewType cellCoords, const subcellCoordsViewType subCellCoords, const ParamViewType subcellParametrization, const unsigned subcellTopoKey, const ordinal_type subCellOrd, const ordinal_type ort)
Maps points defined on the subCell manifold into the parent Cell accounting for orientation.
static ConstViewType get(const ordinal_type subcellDim, const unsigned parentCellKey)
Returns a Kokkos view with the coefficients of the parametrization maps for the edges or faces of a r...
static void getLattice(Kokkos::DynRankView< pointValueType, pointProperties...> points, const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0, const EPointType pointType=POINTTYPE_EQUISPACED)
Computes a lattice of points of a given order on a reference simplex, quadrilateral or hexahedron (cu...
static void getCoeffMatrix_HDIV(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt, const bool inverse=false)
Compute orientation matrix for HDIV basis for a given subcell and its reference basis.
static KOKKOS_INLINE_FUNCTION void getRefSideTangentsAndNormal(TanNormViewType tangentsAndNormal, const ParamViewType subCellParametrization, const unsigned subcellTopoKey, const ordinal_type subCellOrd, const ordinal_type ort)
Computes the (oriented) tangents and normal of the side subCell The normals are only defined for side...
static ordinal_type getLatticeSize(const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0)
Computes the number of points in a lattice of a given order on a simplex (currently disabled for othe...