Intrepid2
Intrepid2_OrientationToolsDefCoeffMatrix_HGRAD.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 
56 #ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HGRAD_HPP__
57 #define __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HGRAD_HPP__
58 
59 // disable clang warnings
60 #if defined (__clang__) && !defined (__INTEL_COMPILER)
61 #pragma clang system_header
62 #endif
63 
64 namespace Intrepid2 {
65 
66 namespace Impl {
67 namespace Debug {
68 
69 #ifdef HAVE_INTREPID2_DEBUG
70 template<typename subcellBasisType,
71 typename cellBasisType>
72 inline
73 void
74 check_getCoeffMatrix_HGRAD(const subcellBasisType& subcellBasis,
75  const cellBasisType& cellBasis,
76  const ordinal_type subcellId,
77  const ordinal_type subcellOrt) {
78 
79  // populate points on a subcell and map to subcell
80  const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
81  const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
82 
83  const ordinal_type cellDim = cellTopo.getDimension();
84  const ordinal_type subcellDim = subcellTopo.getDimension();
85 
86  INTREPID2_TEST_FOR_EXCEPTION( subcellDim >= cellDim,
87  std::logic_error,
88  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
89  "cellDim must be greater than subcellDim.");
90 
91  const auto subcellBaseKey = subcellTopo.getBaseKey();
92 
93  INTREPID2_TEST_FOR_EXCEPTION( subcellBaseKey != shards::Line<>::key &&
94  subcellBaseKey != shards::Quadrilateral<>::key &&
95  subcellBaseKey != shards::Triangle<>::key,
96  std::logic_error,
97  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
98  "subcellBasis must have line, quad, or triangle topology.");
99 
100 
101  //
102  // Function space
103  //
104 
105  {
106  const bool cellBasisIsHGRAD = cellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
107  const bool subcellBasisIsHGRAD = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
108  if (cellBasisIsHGRAD) {
109  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHGRAD,
110  std::logic_error,
111  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
112  "subcellBasis function space is not consistent to cellBasis.");
113  }
114 
115  INTREPID2_TEST_FOR_EXCEPTION( subcellBasis.getDegree() != cellBasis.getDegree(),
116  std::logic_error,
117  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
118  "subcellBasis has a different polynomial degree from cellBasis' degree.");
119  }
120 }
121 #endif
122 } // Debug namespace
123 
124 template<typename OutputViewType,
125 typename subcellBasisHostType,
126 typename cellBasisHostType>
127 inline
128 void
130 getCoeffMatrix_HGRAD(OutputViewType &output,
131  const subcellBasisHostType& subcellBasis,
132  const cellBasisHostType& cellBasis,
133  const ordinal_type subcellId,
134  const ordinal_type subcellOrt) {
135 
136 #ifdef HAVE_INTREPID2_DEBUG
137  Debug::check_getCoeffMatrix_HGRAD(subcellBasis,cellBasis,subcellId,subcellOrt);
138 #endif
139 
140  using host_device_type = typename Kokkos::HostSpace::device_type;
141  using value_type = typename OutputViewType::non_const_value_type;
142 
143  //
144  // Topology
145  //
146 
147  const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
148  const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
149  const ordinal_type cellDim = cellTopo.getDimension();
150  const ordinal_type subcellDim = subcellTopo.getDimension();
151  const auto subcellBaseKey = subcellTopo.getBaseKey();
152  const ordinal_type numCellBasis = cellBasis.getCardinality();
153  const ordinal_type numSubcellBasis = subcellBasis.getCardinality();
154  const ordinal_type ndofSubcell = cellBasis.getDofCount(subcellDim,subcellId);
155 
156  //
157  // Reference points
158  //
159 
160  // Reference points \xi_j on the subcell
161  Kokkos::DynRankView<value_type,host_device_type> refPtsSubcell("refPtsSubcell", ndofSubcell, subcellDim);
162  auto latticeSize=PointTools::getLatticeSize(subcellTopo, subcellBasis.getDegree(), 1);
163 
164  INTREPID2_TEST_FOR_EXCEPTION( latticeSize != ndofSubcell,
165  std::logic_error,
166  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
167  "Lattice size should be equal to the onber of subcell internal DoFs");
168  PointTools::getLattice(refPtsSubcell, subcellTopo, subcellBasis.getDegree(), 1, POINTTYPE_WARPBLEND);
169 
170  // map the points into the parent, cell accounting for orientation
171  auto subcellParam = Intrepid2::RefSubcellParametrization<host_device_type>::get(subcellDim, cellTopo.getKey());
172  Kokkos::DynRankView<value_type,host_device_type> refPtsCell("refPtsCell", ndofSubcell, cellDim);
173  // refPtsCell = F_s (\eta_o (refPtsSubcell))
174  mapSubcellCoordsToRefCell(refPtsCell,refPtsSubcell, subcellParam, subcellBaseKey, subcellId, subcellOrt);
175 
176  //
177  // Bases evaluation on the reference points
178  //
179 
180  // cellBasisValues = \psi_k(F_s (\eta_o (\xi_j)))
181  Kokkos::DynRankView<value_type,host_device_type> cellBasisValues("cellBasisValues", numCellBasis, ndofSubcell);
182 
183  // subcellBasisValues = \phi_i (\xi_j)
184  Kokkos::DynRankView<value_type,host_device_type> subcellBasisValues("subcellBasisValues", numSubcellBasis, ndofSubcell);
185 
186  cellBasis.getValues(cellBasisValues, refPtsCell, OPERATOR_VALUE);
187  subcellBasis.getValues(subcellBasisValues, refPtsSubcell, OPERATOR_VALUE);
188 
189  //
190  // Compute Psi_jk = \psi_k(F_s (\eta_o (\xi_j))) and Phi_ji = \phi_i (\xi_j),
191  // and solve
192  // Psi A^T = Phi
193  //
194 
195  // construct Psi and Phi matrices. LAPACK wants left layout
196  Kokkos::DynRankView<value_type,Kokkos::LayoutLeft,host_device_type>
197  PsiMat("PsiMat", ndofSubcell, ndofSubcell),
198  PhiMat("PhiMat", ndofSubcell, ndofSubcell);
199 
200  auto cellTagToOrdinal = cellBasis.getAllDofOrdinal();
201  auto subcellTagToOrdinal = subcellBasis.getAllDofOrdinal();
202 
203  for (ordinal_type i=0;i<ndofSubcell;++i) {
204  const ordinal_type ic = cellTagToOrdinal(subcellDim, subcellId, i);
205  const ordinal_type isc = subcellTagToOrdinal(subcellDim, 0, i);
206  for (ordinal_type j=0;j<ndofSubcell;++j) {
207  PsiMat(j, i) = cellBasisValues(ic,j);
208  PhiMat(j, i) = subcellBasisValues(isc,j);
209  }
210  }
211 
212  // Solve the system
213  {
214  Teuchos::LAPACK<ordinal_type,value_type> lapack;
215  ordinal_type info = 0;
216 
217  Kokkos::DynRankView<ordinal_type,Kokkos::LayoutLeft,host_device_type> pivVec("pivVec", ndofSubcell);
218  lapack.GESV(ndofSubcell, ndofSubcell,
219  PsiMat.data(),
220  PsiMat.stride_1(),
221  pivVec.data(),
222  PhiMat.data(),
223  PhiMat.stride_1(),
224  &info);
225 
226  if (info) {
227  std::stringstream ss;
228  ss << ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): "
229  << "LAPACK return with error code: "
230  << info;
231  INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
232  }
233 
234  //After solving the system w/ LAPACK, Phi contains A^T
235 
236  // transpose B and clean up numerical noise (for permutation matrices)
237  const double eps = tolerence();
238  for (ordinal_type i=0;i<ndofSubcell;++i) {
239  auto intmatii = std::round(PhiMat(i,i));
240  PhiMat(i,i) = (std::abs(PhiMat(i,i) - intmatii) < eps) ? intmatii : PhiMat(i,i);
241  for (ordinal_type j=i+1;j<ndofSubcell;++j) {
242  auto matij = PhiMat(i,j);
243 
244  auto intmatji = std::round(PhiMat(j,i));
245  PhiMat(i,j) = (std::abs(PhiMat(j,i) - intmatji) < eps) ? intmatji : PhiMat(j,i);
246 
247  auto intmatij = std::round(matij);
248  PhiMat(j,i) = (std::abs(matij - intmatij) < eps) ? intmatij : matij;
249  }
250  }
251 
252  }
253 
254  // Print A Matrix
255  /*
256  {
257  std::cout << "|";
258  for (ordinal_type i=0;i<ndofSubcell;++i) {
259  for (ordinal_type j=0;j<ndofSubcell;++j) {
260  std::cout << PhiMat(i,j) << " ";
261  }
262  std::cout << "| ";
263  }
264  std::cout <<std::endl;
265  }
266  */
267 
268  {
269  // move the data to original device memory
270  const Kokkos::pair<ordinal_type,ordinal_type> range(0, ndofSubcell);
271  auto suboutput = Kokkos::subview(output, range, range);
272  auto tmp = Kokkos::create_mirror_view_and_copy(typename OutputViewType::device_type::memory_space(), PhiMat);
273  Kokkos::deep_copy(suboutput, tmp);
274  }
275 }
276 }
277 
278 }
279 #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.
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 getCoeffMatrix_HGRAD(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
Compute orientation matrix for HGRAD basis for a given subcell and its reference basis.
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 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...