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