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 
55 #ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HGRAD_HPP__
56 #define __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HGRAD_HPP__
57 
58 // disable clang warnings
59 #if defined (__clang__) && !defined (__INTEL_COMPILER)
60 #pragma clang system_header
61 #endif
62 
63 namespace Intrepid2 {
64 
65  namespace Impl {
66 
67  template<typename OutputViewType,
68  typename subcellBasisType,
69  typename cellBasisType>
70  inline
71  void
73  getCoeffMatrix_HGRAD(OutputViewType &output,
74  const subcellBasisType subcellBasis,
75  const cellBasisType cellBasis,
76  const ordinal_type subcellId,
77  const ordinal_type subcellOrt) {
78  typedef typename OutputViewType::execution_space space_type;
79  typedef typename OutputViewType::value_type value_type;
80 
81  // with shards, everything should be computed on host space
82  typedef typename
83  Kokkos::Impl::is_space<space_type>::host_mirror_space::execution_space host_space_type;
84 
85  typedef Kokkos::DynRankView<value_type,host_space_type> DynRankViewHostType;
86 
87 
88  //
89  // Topology
90  //
91 
92  // populate points on a subcell and map to subcell
93  const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
94  const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
95 
96  const ordinal_type cellDim = cellTopo.getDimension();
97  const ordinal_type subcellDim = subcellTopo.getDimension();
98 
99  INTREPID2_TEST_FOR_EXCEPTION( subcellDim >= cellDim,
100  std::logic_error,
101  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
102  "cellDim must be greater than subcellDim.");
103 
104  const auto subcellBaseKey = subcellTopo.getBaseKey();
105 
106  INTREPID2_TEST_FOR_EXCEPTION( subcellBaseKey != shards::Line<>::key &&
107  subcellBaseKey != shards::Quadrilateral<>::key &&
108  subcellBaseKey != shards::Triangle<>::key,
109  std::logic_error,
110  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
111  "subcellBasis must have line, quad, or triangle topology.");
112 
113 
114  //
115  // Function space
116  //
117 
118  {
119  const bool cellBasisIsHGRAD = cellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
120  const bool subcellBasisIsHGRAD = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
121  if (cellBasisIsHGRAD) {
122  INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHGRAD,
123  std::logic_error,
124  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
125  "subcellBasis function space is not consistent to cellBasis.");
126  }
127 
128  INTREPID2_TEST_FOR_EXCEPTION( subcellBasis.getDegree() != cellBasis.getDegree(),
129  std::logic_error,
130  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
131  "subcellBasis has a different polynomial degree from cellBasis' degree.");
132  }
133 
134 
135  //
136  // Collocation points
137  //
138 
139  const ordinal_type numCellBasis = cellBasis.getCardinality();
140  const ordinal_type numSubcellBasis = subcellBasis.getCardinality();
141 
142  const ordinal_type ordSubcell = cellBasis.getDofOrdinal(subcellDim, subcellId, 0);
143  INTREPID2_TEST_FOR_EXCEPTION( ordSubcell == -1,
144  std::logic_error,
145  ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
146  "Invalid subcellId returns -1 ordSubcell.");
147 
148  const ordinal_type ndofSubcell = cellBasis.getDofTag(ordSubcell)(3);
149 
150  // reference points on a subcell
151  DynRankViewHostType refPtsSubcell("refPtsSubcell", ndofSubcell, subcellDim);
152  DynRankViewHostType subcellDofCoords("subcellDofCoords", numSubcellBasis, subcellDim);
153  subcellBasis.getDofCoords(subcellDofCoords);
154  for(ordinal_type i=0; i<ndofSubcell; ++i)
155  for(ordinal_type d=0; d <subcellDim; ++d)
156  refPtsSubcell(i,d) = subcellDofCoords(subcellBasis.getDofOrdinal(subcellDim, 0, i),d);
157 
158  // modified points with orientation
159  DynRankViewHostType ortPtsSubcell("ortPtsSubcell", ndofSubcell, subcellDim);
161  refPtsSubcell,
162  subcellTopo,
163  subcellOrt);
164 
165  // map to reference coordinates
166  DynRankViewHostType refPtsCell("refPtsCell", ndofSubcell, cellDim);
168  ortPtsSubcell,
169  subcellDim,
170  subcellId,
171  cellTopo);
172 
173  //
174  // Basis evaluation on the collocation points
175  //
176 
177  DynRankViewHostType refValues("refValues", numCellBasis, ndofSubcell);
178  cellBasis.getValues(refValues, refPtsCell, OPERATOR_VALUE);
179 
180 
181  //
182  // Construct collocation matrix and solve problems
183  //
184 
185  // construct collocation matrix; using lapack, it should be left layout
186  Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type>
187  refMat("refMat", ndofSubcell, ndofSubcell),
188  ortMat("ortMat", ndofSubcell, ndofSubcell),
189  pivVec("pivVec", ndofSubcell, 1);
190 
191  for (ordinal_type i=0;i<ndofSubcell;++i) {
192  const ordinal_type iref = cellBasis.getDofOrdinal(subcellDim, subcellId, i);
193  for (ordinal_type j=0;j<ndofSubcell;++j) {
194  refMat(i,j) = refValues(iref,j);
195  ortMat(i,j) = (i==j); //thanks to the kronecher property of the basis functions
196  }
197  }
198 
199  // solve the system
200  {
201  Teuchos::LAPACK<ordinal_type,value_type> lapack;
202  ordinal_type info = 0;
203 
204  lapack.GESV(ndofSubcell, ndofSubcell,
205  refMat.data(),
206  refMat.stride_1(),
207  (ordinal_type*)pivVec.data(),
208  ortMat.data(),
209  ortMat.stride_1(),
210  &info);
211 
212  if (info) {
213  std::stringstream ss;
214  ss << ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): "
215  << "LAPACK return with error code: "
216  << info;
217  INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
218  }
219 
220  //Transpose matrix and clean up numerical noise
221  {
222  const double eps = threshold();
223  for (ordinal_type i=0;i<ndofSubcell;++i)
224  for (ordinal_type j=i;j<ndofSubcell;++j) {
225  auto intOrtMat = std::round(ortMat(i,j));
226  ortMat(i,j) = (std::abs(ortMat(i,j) - std::round(ortMat(i,j))) < eps) ? intOrtMat : ortMat(i,j);
227  }
228  }
229  }
230 
231  {
232  // move the data to original device memory
233  const Kokkos::pair<ordinal_type,ordinal_type> range(0, ndofSubcell);
234  Kokkos::deep_copy(Kokkos::subview(output, range, range), ortMat);
235  }
236  }
237  }
238 
239 }
240 #endif
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 getCoeffMatrix_HGRAD(OutputViewType &output, const subcellBasisType subcellBasis, const cellBasisType cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
Compute coefficient matrix for HGRAD by collocating point values.
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.