Intrepid2
Intrepid2_OrientationToolsDefCoeffMatrix_HCURL.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_HCURL_HPP__
60 #define __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HCURL_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  template<typename OutputViewType,
72  typename subcellBasisType,
73  typename cellBasisType>
74  inline
75  void
77  getCoeffMatrix_HCURL(OutputViewType &output,
78  const subcellBasisType subcellBasis,
79  const cellBasisType cellBasis,
80  const ordinal_type subcellId,
81  const ordinal_type subcellOrt) {
82  typedef typename OutputViewType::execution_space space_type;
83  typedef typename OutputViewType::value_type value_type;
84 
85  // with shards, everything should be computed on host space
86  typedef typename
87  Kokkos::Impl::is_space<space_type>::host_mirror_space::execution_space host_space_type;
88 
89  typedef Kokkos::DynRankView<value_type,host_space_type> DynRankViewHostType;
90 
91  //
92  // Topology
93  //
94  // populate points on a subcell and map to subcell
95  const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
96  const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
97 
98  const ordinal_type cellDim = cellTopo.getDimension();
99  const ordinal_type subcellDim = subcellTopo.getDimension();
100 
101  INTREPID2_TEST_FOR_EXCEPTION( subcellDim >= cellDim,
102  std::logic_error,
103  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
104  "cellDim must be greater than subcellDim.");
105 
106  const auto subcellBaseKey = subcellTopo.getBaseKey();
107  const auto cellBaseKey = cellTopo.getBaseKey();
108 
109  INTREPID2_TEST_FOR_EXCEPTION( subcellBaseKey != shards::Line<>::key &&
110  subcellBaseKey != shards::Quadrilateral<>::key &&
111  subcellBaseKey != shards::Triangle<>::key,
112  std::logic_error,
113  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
114  "subcellBasis must have line, quad, or triangle topology.");
115 
116  INTREPID2_TEST_FOR_EXCEPTION( cellBaseKey != shards::Quadrilateral<>::key &&
117  cellBaseKey != shards::Triangle<>::key &&
118  cellBaseKey != shards::Hexahedron<>::key &&
119  cellBaseKey != shards::Tetrahedron<>::key,
120  std::logic_error,
121  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
122  "cellBasis must have quad, triangle, hexhedron or tetrahedron topology.");
123 
124  //
125  // Function space
126  //
127  {
128  const bool cellBasisIsHCURL = cellBasis.getFunctionSpace() == FUNCTION_SPACE_HCURL;
129  if (cellBasisIsHCURL) {
130  const bool subcellBasisIsHGRAD = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
131  const bool subcellBasisIsHVOL = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HVOL;
132  const bool subcellBasisIsHCURL = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HCURL;
133  const bool cellIsTri = cellBaseKey == shards::Triangle<>::key;
134  const bool cellIsTet = cellBaseKey == shards::Tetrahedron<>::key;
135  const bool cellIsHex = cellBaseKey == shards::Hexahedron<>::key;
136  const bool cellIsQuad = cellBaseKey == shards::Quadrilateral<>::key;
137  // edge hcurl is hgrad with gauss legendre points
138  switch (subcellDim) {
139  case 1: {
140  //TODO: Hex, QUAD, TET and TRI element should have the same 1d basis
141  if (cellIsHex || cellIsQuad) {
142  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHGRAD,
143  std::logic_error,
144  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "
145  "subcellBasis function space (1d) is not consistent to cellBasis.");
146  } else if (cellIsTet || cellIsTri) {
147  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHVOL,
148  std::logic_error,
149  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "
150  "subcellBasis function space (1d) is not consistent to cellBasis.");
151  }
152  break;
153  }
154  case 2: {
155  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHCURL,
156  std::logic_error,
157  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "
158  "subcellBasis function space (2d) is not consistent to cellBasis.");
159  break;
160  }
161  }
162  }
163  }
164 
165  const ordinal_type numCellBasis = cellBasis.getCardinality();
166  const ordinal_type numSubcellBasis = subcellBasis.getCardinality();
167 
168  const ordinal_type ordSubcell = cellBasis.getDofOrdinal(subcellDim, subcellId, 0);
169  INTREPID2_TEST_FOR_EXCEPTION( ordSubcell == -1,
170  std::logic_error,
171  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
172  "Invalid subcellId returns -1 ordSubcell.");
173 
174  const ordinal_type ndofSubcell = cellBasis.getDofTag(ordSubcell)(3);
175 
176  // reference points on a subcell
177  DynRankViewHostType refPtsSubcell("refPtsSubcell", ndofSubcell, subcellDim);
178  DynRankViewHostType subcellDofCoords("subcellDofCoords", numSubcellBasis, subcellDim);
179  subcellBasis.getDofCoords(subcellDofCoords);
180  for(ordinal_type i=0; i<ndofSubcell; ++i)
181  for(ordinal_type d=0; d <subcellDim; ++d)
182  refPtsSubcell(i,d) = subcellDofCoords(subcellBasis.getDofOrdinal(subcellDim, 0, i),d);
183 
184  // modified points with orientation
185  DynRankViewHostType ortPtsSubcell("ortPtsSubcell", ndofSubcell, subcellDim);
187  refPtsSubcell,
188  subcellTopo,
189  subcellOrt);
190 
191  // map to reference coordinates
192  DynRankViewHostType refPtsCell("refPtsCell", ndofSubcell, cellDim);
194  ortPtsSubcell,
195  subcellDim,
196  subcellId,
197  cellTopo);
198 
199  //
200  // Basis evaluation on the reference points
201  //
202 
203  // evaluate values on the reference cell
204  DynRankViewHostType refValues("refValues", numCellBasis, ndofSubcell, cellDim);
205  cellBasis.getValues(refValues, refPtsCell, OPERATOR_VALUE);
206 
207  // evaluate values on the modified cell
208  DynRankViewHostType subcellTangents("subcellTangents", numSubcellBasis, subcellDim);
209  DynRankViewHostType ortJacobian("ortJacobian", subcellDim, subcellDim);
210  Impl::OrientationTools::getJacobianOfOrientationMap(ortJacobian, subcellTopo, subcellOrt);
211 
212  //
213  // Compute jacobianF
214  //
215 
216  DynRankViewHostType jacobianF("jacobianF", cellDim, subcellDim );
217  switch (subcellBaseKey) {
218  case shards::Line<>::key: {
219  auto lineDofCoeffs = Kokkos::subview(subcellTangents, Kokkos::ALL(),0);
220  subcellBasis.getDofCoeffs(lineDofCoeffs);
221  auto edgeTan = Kokkos::subview(jacobianF, Kokkos::ALL(),0);
222  CellTools<host_space_type>::getReferenceEdgeTangent(edgeTan, subcellId, cellTopo);
223  if((cellBaseKey == shards::Triangle<>::key) || (cellBaseKey == shards::Tetrahedron<>::key))
224  for(ordinal_type i=0; i<cellDim; ++i)
225  edgeTan(i) *= 2.0; //scale by reference tangent
226  break;
227  }
228  case shards::Quadrilateral<>::key:
229  case shards::Triangle<>::key: {
230  subcellBasis.getDofCoeffs(subcellTangents);
231 
232  auto faceTanU = Kokkos::subview(jacobianF, Kokkos::ALL(), 0);
233  auto faceTanV = Kokkos::subview(jacobianF, Kokkos::ALL(), 1);
234  CellTools<host_space_type>::getReferenceFaceTangents(faceTanU, faceTanV,subcellId, cellTopo);
235  break;
236  }
237  default: {
238  INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, "Should not come here" );
239  }
240  }
241 
242  Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type> // left layout for lapack
243  refMat("refMat", ndofSubcell, ndofSubcell),
244  ortMat("ortMat", ndofSubcell, ndofSubcell);
245  for (ordinal_type i=0;i<ndofSubcell;++i) {
246  const ordinal_type iout = cellBasis.getDofOrdinal(subcellDim, subcellId, i);
247  for (ordinal_type j=0;j<ndofSubcell;++j) {
248  value_type tmp = 0;
249  const ordinal_type jsc = subcellBasis.getDofOrdinal(subcellDim, 0, j);
250  for (ordinal_type k=0;k<subcellDim;++k)
251  for (ordinal_type l=0;l<subcellDim;++l)
252  for (ordinal_type d=0; d<cellDim; ++d)
253  tmp += refValues(iout,j,d)*jacobianF(d,l)*ortJacobian(l,k)*subcellTangents(jsc,k);
254  refMat(i,j) = tmp;
255  ortMat(j,i) = (i==j); //identity because of the basis Kronecher property
256  }
257  }
258 
259  //
260  // Construct collocation matrix and solve problems
261  //
262 
263  // solve the system using Lapack
264  {
265  Teuchos::LAPACK<ordinal_type,value_type> lapack;
266  ordinal_type info = 0;
267  Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type> pivVec("pivVec", ndofSubcell, 1);
268 
269  lapack.GESV(ndofSubcell, ndofSubcell,
270  refMat.data(),
271  refMat.stride_1(),
272  (ordinal_type*)pivVec.data(),
273  ortMat.data(),
274  ortMat.stride_1(),
275  &info);
276 
277  if (info) {
278  std::stringstream ss;
279  ss << ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "
280  << "LAPACK return with error code: "
281  << info;
282  INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
283  }
284 
285 
286  // Clean up numerical noise
287  {
288  const double eps = threshold();
289  for (ordinal_type i=0;i<ndofSubcell;++i)
290  for (ordinal_type j=i;j<ndofSubcell;++j) {
291  auto intOrtMat = std::round(ortMat(i,j));
292  ortMat(i,j) = (std::abs(ortMat(i,j) - std::round(ortMat(i,j))) < eps) ? intOrtMat : ortMat(i,j);
293  }
294  }
295  }
296 
297  {
298  // move the data to original device memory
299  const Kokkos::pair<ordinal_type,ordinal_type> range(0, ndofSubcell);
300  Kokkos::deep_copy(Kokkos::subview(output, range, range),
301  Kokkos::subview(ortMat, range, range));
302  }
303  }
304  }
305 
306 }
307 #endif
static void getReferenceEdgeTangent(Kokkos::DynRankView< refEdgeTangentValueType, refEdgeTangentProperties...> refEdgeTangent, const ordinal_type edgeOrd, const shards::CellTopology parentCell)
Computes constant tangent vectors to edges of 2D or 3D reference cells.
static void getCoeffMatrix_HCURL(OutputViewType &output, const subcellBasisType subcellBasis, const cellBasisType cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
Compute coefficient matrix for HCURL by collocating point values.
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 mapToReferenceSubcell(Kokkos::DynRankView< refSubcellPointValueType, refSubcellPointProperties...> refSubcellPoints, const Kokkos::DynRankView< paramPointValueType, paramPointProperties...> paramPoints, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell)
Computes parameterization maps of 1- and 2-subcells of reference cells.
static void getReferenceFaceTangents(Kokkos::DynRankView< refFaceTanUValueType, refFaceTanUProperties...> refFaceTanU, Kokkos::DynRankView< refFaceTanVValueType, refFaceTanVProperties...> refFaceTanV, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes pairs of constant tangent vectors to faces of a 3D reference cells.
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.