Intrepid2
Intrepid2_OrientationToolsDefCoeffMatrix_HCURL.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Intrepid2 Package
4 //
5 // Copyright 2007 NTESS and the Intrepid2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 
29 #ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HCURL_HPP__
30 #define __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HCURL_HPP__
31 
35 
36 // disable clang warnings
37 #if defined (__clang__) && !defined (__INTEL_COMPILER)
38 #pragma clang system_header
39 #endif
40 
41 namespace Intrepid2 {
42 
43 namespace Impl {
44 
45 namespace Debug {
46 
47 #ifdef HAVE_INTREPID2_DEBUG
48 template<typename subcellBasisType,
49 typename cellBasisType>
50 inline
51 void
52 check_getCoeffMatrix_HCURL(const subcellBasisType& subcellBasis,
53  const cellBasisType& cellBasis,
54  const ordinal_type subcellId,
55  const ordinal_type subcellOrt) {
56  const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
57  const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
58 
59  const ordinal_type cellDim = cellTopo.getDimension();
60  const ordinal_type subcellDim = subcellTopo.getDimension();
61 
62 
63 
64  INTREPID2_TEST_FOR_EXCEPTION( subcellDim > cellDim,
65  std::logic_error,
66  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
67  "cellDim cannot be smaller than subcellDim.");
68 
69  const auto subcellBaseKey = subcellTopo.getBaseKey();
70  const auto cellBaseKey = cellTopo.getBaseKey();
71 
72  INTREPID2_TEST_FOR_EXCEPTION( subcellBaseKey != shards::Line<>::key &&
73  subcellBaseKey != shards::Quadrilateral<>::key &&
74  subcellBaseKey != shards::Triangle<>::key,
75  std::logic_error,
76  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
77  "subcellBasis must have line, quad, or triangle topology.");
78 
79  INTREPID2_TEST_FOR_EXCEPTION( cellBaseKey != shards::Quadrilateral<>::key &&
80  cellBaseKey != shards::Triangle<>::key &&
81  cellBaseKey != shards::Hexahedron<>::key &&
82  cellBaseKey != shards::Wedge<>::key &&
83  cellBaseKey != shards::Tetrahedron<>::key,
84  std::logic_error,
85  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
86  "cellBasis must have quad, triangle, hexhedron or tetrahedron topology.");
87 
88  //
89  // Function space
90  //
91  {
92  const bool cellBasisIsHCURL = cellBasis.getFunctionSpace() == FUNCTION_SPACE_HCURL;
93 
94 
95  if (cellBasisIsHCURL) {
96  const bool subcellBasisIsHGRAD = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
97  const bool subcellBasisIsHVOL = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HVOL;
98  const bool subcellBasisIsHCURL = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HCURL;
99  const bool cellIsTri = cellBaseKey == shards::Triangle<>::key;
100  const bool cellIsTet = cellBaseKey == shards::Tetrahedron<>::key;
101  const bool cellIsHex = cellBaseKey == shards::Hexahedron<>::key;
102  const bool cellIsQuad = cellBaseKey == shards::Quadrilateral<>::key;
103 
104 
105  // edge hcurl is hgrad with gauss legendre points
106  switch (subcellDim) {
107  case 1: {
108  //TODO: Hex, QUAD, TET and TRI element should have the same 1d basis
109  if (cellIsHex || cellIsQuad) {
110  INTREPID2_TEST_FOR_EXCEPTION( !(subcellBasisIsHGRAD||subcellBasisIsHVOL),
111  std::logic_error,
112  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "
113  "subcellBasis function space (1d) is not consistent to cellBasis.");
114  } else if (cellIsTet || cellIsTri) {
115  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHVOL,
116  std::logic_error,
117  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "
118  "subcellBasis function space (1d) is not consistent to cellBasis.");
119  }
120  break;
121  }
122  case 2: {
123  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHCURL,
124  std::logic_error,
125  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "
126  "subcellBasis function space (2d) is not consistent to cellBasis.");
127  break;
128  }
129  }
130  }
131  }
132 }
133 #endif
134 }
135 
136 template<typename OutputViewType,
137 typename subcellBasisHostType,
138 typename cellBasisHostType>
139 inline
140 void
142 getCoeffMatrix_HCURL(OutputViewType &output,
143  const subcellBasisHostType& subcellBasis,
144  const cellBasisHostType& cellBasis,
145  const ordinal_type subcellId,
146  const ordinal_type subcellOrt,
147  const bool inverse) {
148 
149 #ifdef HAVE_INTREPID2_DEBUG
150  Debug::check_getCoeffMatrix_HCURL(subcellBasis,cellBasis,subcellId,subcellOrt);
151 #endif
152 
153  using value_type = typename OutputViewType::non_const_value_type;
154  using host_device_type = typename Kokkos::HostSpace::device_type;
155 
156  //
157  // Topology
158  //
159  // populate points on a subcell and map to subcell
160  const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
161  const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
162  const ordinal_type cellDim = cellTopo.getDimension();
163  const ordinal_type subcellDim = subcellTopo.getDimension();
164  const auto subcellBaseKey = subcellTopo.getBaseKey();
165  const ordinal_type numCellBasis = cellBasis.getCardinality();
166  const ordinal_type numSubcellBasis = subcellBasis.getCardinality();
167  const ordinal_type ndofSubcell = cellBasis.getDofCount(subcellDim,subcellId);
168 
169  // Compute reference points \xi_j and tangents t_j on the subcell
170  // To do so we use the DoF coordinates and DoF coefficients of a Lagrangian HCURL basis
171  // on the subcell spanning the same space as the bases \phi_j
172 
173  // Tangents t_j
174  Kokkos::DynRankView<value_type,host_device_type> subcellTangents("subcellTangents", numSubcellBasis, subcellDim);
175  auto degree = subcellBasis.getDegree();
176  BasisPtr<host_device_type, value_type, value_type> basisPtr;
177  if(subcellBaseKey == shards::Line<>::key) {
179  basisPtr->getDofCoeffs(Kokkos::subview(subcellTangents, Kokkos::ALL(),0));
180  } else if (subcellBaseKey == shards::Triangle<>::key) {
182  basisPtr->getDofCoeffs(subcellTangents);
183  } else if (subcellBaseKey == shards::Quadrilateral<>::key) {
185  basisPtr->getDofCoeffs(subcellTangents);
186  }
187 
188  // coordinates \xi_j
189  Kokkos::DynRankView<value_type,host_device_type> subcellDofCoords("subcellDofCoords", basisPtr->getCardinality(), subcellDim);
190  basisPtr->getDofCoords(subcellDofCoords);
191  INTREPID2_TEST_FOR_EXCEPTION( basisPtr->getDofCount(subcellDim,0) != ndofSubcell,
192  std::logic_error,
193  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
194  "the number of basisPtr internal DoFs should equate those of the subcell");
195 
196  // restrict \xi_j (and corresponding t_j) to points internal to the HCURL basis
197  Kokkos::DynRankView<value_type,host_device_type> refPtsSubcell("refPtsSubcell", ndofSubcell, subcellDim);
198  Kokkos::DynRankView<value_type,host_device_type> refSubcellTangents("subcellTangents", ndofSubcell, subcellDim);
199  auto tagToOrdinal = basisPtr->getAllDofOrdinal();
200  for (ordinal_type i=0;i<ndofSubcell;++i) {
201  ordinal_type isc = tagToOrdinal(subcellDim, 0, i);
202  for(ordinal_type d=0; d <subcellDim; ++d){
203  refPtsSubcell(i,d) = subcellDofCoords(isc,d);
204  for(ordinal_type k=0; k <subcellDim; ++k)
205  refSubcellTangents(i,d) = subcellTangents(isc,d);
206  }
207  }
208 
209  //
210  // Bases evaluation on the reference points
211  //
212 
213  // subcellBasisValues = \phi_i (\xi_j)
214  Kokkos::DynRankView<value_type,host_device_type> subCellValues("subCellValues", numSubcellBasis, ndofSubcell, subcellDim);
215  if(subcellDim==1) {
216  auto lineValues = Kokkos::subview(subCellValues, Kokkos::ALL(), Kokkos::ALL(), 0);
217  subcellBasis.getValues(lineValues, refPtsSubcell, OPERATOR_VALUE);
218  } else {
219  subcellBasis.getValues(subCellValues, refPtsSubcell, OPERATOR_VALUE);
220  }
221 
222  //
223  // Basis evaluation on the reference points
224  //
225 
226  Kokkos::DynRankView<value_type,host_device_type> refPtsCell("refPtsCell", ndofSubcell, cellDim);
227  Kokkos::DynRankView<value_type,host_device_type> trJacobianF("trJacobianF", subcellDim, cellDim );
228 
229  if(cellDim == subcellDim) {
230  // refPtsCell = \eta_o (refPtsSubcell)
231  mapToModifiedReference(refPtsCell,refPtsSubcell,subcellBaseKey,subcellOrt);
232 
233  //mapping tangents t_j into parent cell, i.e. computing J_F J_\eta t_j
234  Kokkos::DynRankView<value_type,host_device_type> jac("data", subcellDim, subcellDim);
235  getJacobianOfOrientationMap(jac,subcellBaseKey,subcellOrt);
236  for(ordinal_type d=0; d<subcellDim; ++d)
237  for(ordinal_type j=0; j<cellDim; ++j) {
238  trJacobianF(j,d) = jac(d,j);
239  }
240  }
241  else {
242  auto subcellParam = Intrepid2::RefSubcellParametrization<host_device_type>::get(subcellDim, cellTopo.getKey());
243  // refPtsCell = F_s (\eta_o (refPtsSubcell))
244  mapSubcellCoordsToRefCell(refPtsCell,refPtsSubcell, subcellParam, subcellBaseKey, subcellId, subcellOrt);
245 
246  //mapping tangents t_j into parent cell, i.e. computing J_F J_\eta t_j
247  OrientationTools::getRefSubcellTangents(trJacobianF, subcellParam, subcellBaseKey, subcellId, subcellOrt);
248  }
249 
250  // cellBasisValues = \psi_k(F_s (\eta_o (\xi_j)))
251  Kokkos::DynRankView<value_type,host_device_type> cellBasisValues("cellBasisValues", numCellBasis, ndofSubcell, cellDim);
252  cellBasis.getValues(cellBasisValues, refPtsCell, OPERATOR_VALUE);
253 
254  //
255  // Compute Psi_jk = \psi_k(F_s (\eta_o (\xi_j))) \cdot (J_F J_\eta t_j)
256  // and Phi_ji = \phi_i (\xi_j) \ctot t_j, and solve
257  // Psi A^T = Phi
258  //
259 
260  // construct Psi and Phi matrices. LAPACK wants left layout
261  Kokkos::DynRankView<value_type,Kokkos::LayoutLeft,host_device_type> // left layout for lapack
262  PsiMat("PsiMat", ndofSubcell, ndofSubcell),
263  PhiMat("PhiMat", ndofSubcell, ndofSubcell),
264  RefMat,
265  OrtMat;
266 
267  auto cellTagToOrdinal = cellBasis.getAllDofOrdinal();
268  auto subcellTagToOrdinal = subcellBasis.getAllDofOrdinal();
269 
270  for (ordinal_type i=0;i<ndofSubcell;++i) {
271  const ordinal_type ic = cellTagToOrdinal(subcellDim, subcellId, i);
272  for (ordinal_type j=0;j<ndofSubcell;++j) {
273  const ordinal_type isc = subcellTagToOrdinal(subcellDim, 0, i);
274  value_type refEntry = 0, ortEntry =0;
275  for (ordinal_type k=0;k<subcellDim;++k) {
276  ortEntry += subCellValues(isc,j,k)*refSubcellTangents(j,k);
277  for (ordinal_type d=0; d<cellDim; ++d)
278  refEntry += cellBasisValues(ic,j,d)*trJacobianF(k,d)*refSubcellTangents(j,k);
279  }
280  PsiMat(j,i) = refEntry;
281  PhiMat(j,i) = ortEntry;
282  }
283  }
284 
285  RefMat = inverse ? PhiMat : PsiMat;
286  OrtMat = inverse ? PsiMat : PhiMat;
287 
288  // Solve the system using Lapack
289  {
290  Teuchos::LAPACK<ordinal_type,value_type> lapack;
291  ordinal_type info = 0;
292 
293 
294  /*
295  Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type> work("pivVec", 2*ndofSubcell, 1);
296  lapack.GELS('N', ndofSubcell*card, ndofSubcell, ndofSubcell,
297  RefMat.data(),
298  RefMat.stride_1(),
299  OrtMat.data(),
300  OrtMat.stride_1(),
301  work.data(), work.extent(0),
302  &info);
303 
304  */
305  Kokkos::DynRankView<ordinal_type,host_device_type> pivVec("pivVec", ndofSubcell);
306  lapack.GESV(ndofSubcell, ndofSubcell,
307  RefMat.data(),
308  RefMat.stride_1(),
309  pivVec.data(),
310  OrtMat.data(),
311  OrtMat.stride_1(),
312  &info);
313  //*/
314  if (info) {
315  std::stringstream ss;
316  ss << ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "
317  << "LAPACK return with error code: "
318  << info;
319  INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
320  }
321 
322  {
323  //After solving the system w/ LAPACK, Phi contains A^T (or A^-T)
324  // transpose and clean up numerical noise (for permutation matrices)
325  const double eps = tolerence();
326  for (ordinal_type i=0;i<ndofSubcell;++i) {
327  auto intmatii = std::round(OrtMat(i,i));
328  OrtMat(i,i) = (std::abs(OrtMat(i,i) - intmatii) < eps) ? intmatii : OrtMat(i,i);
329  for (ordinal_type j=i+1;j<ndofSubcell;++j) {
330  auto matij = OrtMat(i,j);
331 
332  auto intmatji = std::round(OrtMat(j,i));
333  OrtMat(i,j) = (std::abs(OrtMat(j,i) - intmatji) < eps) ? intmatji : OrtMat(j,i);
334 
335  auto intmatij = std::round(matij);
336  OrtMat(j,i) = (std::abs(matij - intmatij) < eps) ? intmatij : matij;
337  }
338  }
339  }
340 
341 
342 
343  // Print A matrix
344  /*
345  {
346  std::cout << "|";
347  for (ordinal_type i=0;i<ndofSubcell;++i) {
348  for (ordinal_type j=0;j<ndofSubcell;++j) {
349  std::cout << OrtMat(i,j) << " ";
350  }
351  std::cout << "| ";
352  }
353  std::cout <<std::endl;
354  }
355  */
356 
357 
358  }
359 
360  {
361  // move the data to original device memory
362  const Kokkos::pair<ordinal_type,ordinal_type> range(0, ndofSubcell);
363  auto suboutput = Kokkos::subview(output, range, range);
364  auto tmp = Kokkos::create_mirror_view_and_copy(typename OutputViewType::device_type::memory_space(), OrtMat);
365  Kokkos::deep_copy(suboutput, tmp);
366  }
367 }
368 }
369 
370 }
371 #endif
Implementation of the default H(curl)-compatible FEM basis on Quadrilateral cell. ...
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 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_HCURL(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt, const bool inverse=false)
Compute orientation matrix for HCURL basis for a given subcell and its reference basis.
static KOKKOS_INLINE_FUNCTION void getRefSubcellTangents(TanViewType tangents, const ParamViewType subCellParametrization, const unsigned subcellTopoKey, const ordinal_type subCellOrd, const ordinal_type ort)
Computes the (oriented) subCell tangents.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Header file for the Intrepid2::Basis_HCURL_QUAD_In_FEM class.
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Tr...
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.
Header file for the Intrepid2::Basis_HCURL_TRI_In_FEM class.