Intrepid2
Intrepid2_OrientationToolsDefCoeffMatrix_HVOL.hpp
Go to the documentation of this file.
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
23#ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HVOL_HPP__
24#define __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HVOL_HPP__
25
26// disable clang warnings
27#if defined (__clang__) && !defined (__INTEL_COMPILER)
28#pragma clang system_header
29#endif
30
31namespace Intrepid2 {
32
33namespace Impl {
34namespace Debug {
35
36#ifdef HAVE_INTREPID2_DEBUG
37template<typename cellBasisType>
38inline
39void
40check_getCoeffMatrix_HVOL(const cellBasisType& cellBasis,
41 const ordinal_type cellOrt) {
42
43 // populate points on a subcell and map to subcell
44 const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
45 const ordinal_type cellDim = cellTopo.getDimension();
46
47 INTREPID2_TEST_FOR_EXCEPTION( cellDim > 2,
48 std::logic_error,
49 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HVOL): " \
50 "HVOL orientation supported only for (side) cells with dimension less than 3.");
51
52 const auto cellBaseKey = cellTopo.getBaseKey();
53
54 INTREPID2_TEST_FOR_EXCEPTION( cellBaseKey != shards::Line<>::key &&
55 cellBaseKey != shards::Quadrilateral<>::key &&
56 cellBaseKey != shards::Triangle<>::key,
57 std::logic_error,
58 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HVOL): " \
59 "cellBasis must have line, quad, or triangle topology.");
60
61
62 //
63 // Function space
64 //
65
66 {
67 INTREPID2_TEST_FOR_EXCEPTION( cellBasis.getFunctionSpace() != FUNCTION_SPACE_HVOL,
68 std::logic_error,
69 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HVOL): " \
70 "cellBasis function space is not HVOL.");
71 }
72 }
73#endif
74} // Debug namespace
75
76template<typename OutputViewType,
77typename cellBasisHostType>
78inline
79void
81getCoeffMatrix_HVOL(OutputViewType &output,
82 const cellBasisHostType& cellBasis,
83 const ordinal_type cellOrt,
84 const bool inverse) {
85
86#ifdef HAVE_INTREPID2_DEBUG
87 Debug::check_getCoeffMatrix_HVOL(cellBasis,cellOrt);
88#endif
89
90 using host_device_type = typename Kokkos::HostSpace::device_type;
91 using value_type = typename OutputViewType::non_const_value_type;
92
93 //
94 // Topology
95 //
96
97 const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
98 const ordinal_type cellDim = cellTopo.getDimension();
99 const auto cellBaseKey = cellTopo.getBaseKey();
100 const ordinal_type cardinality = cellBasis.getCardinality();
101
102 //
103 // Reference points
104 //
105
106 // Reference points \xi_j on the subcell
107 Kokkos::DynRankView<value_type,host_device_type> refPtsCell("refPtsCell", cardinality, cellDim),refPtsCellNotOriented("refPtsCellNotOriented", cardinality, cellDim);
108
109 ordinal_type latticeOffset(1);
110
111 // this work for line and quadrilateral topologies
112 ordinal_type latticeOrder = (cellTopo.getBaseKey() == shards::Triangle<>::key) ?
113 cellBasis.getDegree() + 3 * latticeOffset : // triangle
114 cellBasis.getDegree() + 2 * latticeOffset; // line and quad
115
116 PointTools::getLattice(refPtsCellNotOriented, cellTopo, latticeOrder, 1, POINTTYPE_WARPBLEND);
117
118 // map the points into the parent, cell accounting for orientation
119 mapToModifiedReference(refPtsCell,refPtsCellNotOriented,cellBaseKey,cellOrt);
120
121 //
122 // Bases evaluation on the reference points
123 //
124
125 // cellBasisValues = \psi_k(\eta_o (\xi_j))
126 Kokkos::DynRankView<value_type,host_device_type> cellBasisValues("cellBasisValues", cardinality, cardinality);
127
128 // basisValues = \phi_i (\xi_j)
129 Kokkos::DynRankView<value_type,host_device_type> nonOrientedBasisValues("subcellBasisValues", cardinality, cardinality);
130
131 cellBasis.getValues(cellBasisValues, refPtsCell, OPERATOR_VALUE);
132 cellBasis.getValues(nonOrientedBasisValues, refPtsCellNotOriented, OPERATOR_VALUE);
133
134 //
135 // Compute Psi_jk = \psi_k(\eta_o (\xi_j)) det(J_\eta) and Phi_ji = \phi_i (\xi_j),
136 // and solve
137 // Psi A^T = Phi
138 //
139
140 // construct Psi and Phi matrices. LAPACK wants left layout
141 Kokkos::DynRankView<value_type,Kokkos::LayoutLeft,host_device_type>
142 PsiMat("PsiMat", cardinality, cardinality),
143 PhiMat("PhiMat", cardinality, cardinality),
144 RefMat,
145 OrtMat;
146
147 auto cellTagToOrdinal = cellBasis.getAllDofOrdinal();
148
149 Kokkos::DynRankView<value_type,host_device_type> jac("jacobian",cellDim,cellDim);
151 value_type jacDet(0);
152 if(cellDim == 2) {
153 jacDet = jac(0,0)*jac(1,1)-jac(0,1)*jac(1,0);
154 } else { //celldim == 1
155 jacDet = jac(0,0);
156 }
157
158 for (ordinal_type i=0;i<cardinality;++i) {
159 const ordinal_type ic = cellTagToOrdinal(cellDim, 0, i);
160 for (ordinal_type j=0;j<cardinality;++j) {
161 PsiMat(j, i) = cellBasisValues(ic,j)*jacDet;
162 PhiMat(j, i) = nonOrientedBasisValues(ic,j);
163 }
164 }
165
166 RefMat = inverse ? PhiMat : PsiMat;
167 OrtMat = inverse ? PsiMat : PhiMat;
168
169 // Solve the system
170 {
171 Teuchos::LAPACK<ordinal_type,value_type> lapack;
172 ordinal_type info = 0;
173
174 Kokkos::DynRankView<ordinal_type,Kokkos::LayoutLeft,host_device_type> pivVec("pivVec", cardinality);
175 lapack.GESV(cardinality, cardinality,
176 RefMat.data(),
177 RefMat.stride(1),
178 pivVec.data(),
179 OrtMat.data(),
180 OrtMat.stride(1),
181 &info);
182
183 if (info) {
184 std::stringstream ss;
185 ss << ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HVOL): "
186 << "LAPACK return with error code: "
187 << info;
188 INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
189 }
190
191 //After solving the system w/ LAPACK, Phi contains A^T
192
193 // transpose B and clean up numerical noise (for permutation matrices)
194 const double eps = tolerence();
195 for (ordinal_type i=0;i<cardinality;++i) {
196 auto intmatii = std::round(OrtMat(i,i));
197 OrtMat(i,i) = (std::abs(OrtMat(i,i) - intmatii) < eps) ? intmatii : OrtMat(i,i);
198 for (ordinal_type j=i+1;j<cardinality;++j) {
199 auto matij = OrtMat(i,j);
200
201 auto intmatji = std::round(OrtMat(j,i));
202 OrtMat(i,j) = (std::abs(OrtMat(j,i) - intmatji) < eps) ? intmatji : OrtMat(j,i);
203
204 auto intmatij = std::round(matij);
205 OrtMat(j,i) = (std::abs(matij - intmatij) < eps) ? intmatij : matij;
206 }
207 }
208
209 }
210
211 // Print A Matrix
212 /*
213 {
214 std::cout << "Ort: " << cellOrt << ": |";
215 for (ordinal_type i=0;i<cardinality;++i) {
216 for (ordinal_type j=0;j<cardinality;++j) {
217 std::cout << OrtMat(i,j) << " ";
218 }
219 std::cout << "| ";
220 }
221 std::cout <<std::endl;
222 }
223 */
224
225 {
226 // move the data to original device memory
227 const Kokkos::pair<ordinal_type,ordinal_type> range(0, cardinality);
228 auto suboutput = Kokkos::subview(output, range, range);
229 auto tmp = Kokkos::create_mirror_view_and_copy(typename OutputViewType::device_type::memory_space(), OrtMat);
230 Kokkos::deep_copy(suboutput, tmp);
231 }
232}
233}
234
235}
236#endif
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 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 void getCoeffMatrix_HVOL(OutputViewType &output, const cellBasisHostType &cellBasis, const ordinal_type cellOrt, const bool inverse=false)
Compute orientation matrix for HVOL basis for a given (2D or 1D) cell 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...