Intrepid2
Intrepid2_HGRAD_QUAD_C2_FEM.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 
49 #ifndef __INTREPID2_HGRAD_QUAD_C2_FEM_HPP__
50 #define __INTREPID2_HGRAD_QUAD_C2_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
54 
55 namespace Intrepid2 {
56 
102  namespace Impl {
103 
107  template<bool serendipity>
109  public:
110  typedef struct Quadrilateral<serendipity ? 8 : 9> cell_topology_type;
114  template<EOperator opType>
115  struct Serial {
116  template<typename OutputViewType,
117  typename inputViewType>
118  KOKKOS_INLINE_FUNCTION
119  static void
120  getValues( OutputViewType output,
121  const inputViewType input );
122 
123  };
124 
125  template<typename DeviceType,
126  typename outputValueValueType, class ...outputValueProperties,
127  typename inputPointValueType, class ...inputPointProperties>
128  static void
129  getValues( const typename DeviceType::execution_space& space,
130  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
131  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
132  const EOperator operatorType);
133 
137  template<typename outputValueViewType,
138  typename inputPointViewType,
139  EOperator opType>
140  struct Functor {
141  outputValueViewType _outputValues;
142  const inputPointViewType _inputPoints;
143 
144  KOKKOS_INLINE_FUNCTION
145  Functor( outputValueViewType outputValues_,
146  inputPointViewType inputPoints_ )
147  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
148 
149  KOKKOS_INLINE_FUNCTION
150  void operator()(const ordinal_type pt) const {
151  switch (opType) {
152  case OPERATOR_VALUE : {
153  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
154  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
155  Serial<opType>::getValues( output, input );
156  break;
157  }
158  case OPERATOR_GRAD :
159  case OPERATOR_CURL :
160  case OPERATOR_D1 :
161  case OPERATOR_D2 :
162  case OPERATOR_D3 :
163  case OPERATOR_D4 :
164  case OPERATOR_MAX : {
165  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
166  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
167  Serial<opType>::getValues( output, input );
168  break;
169  }
170  default: {
171  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
172  opType != OPERATOR_GRAD &&
173  opType != OPERATOR_CURL &&
174  opType != OPERATOR_D1 &&
175  opType != OPERATOR_D2 &&
176  opType != OPERATOR_D3 &&
177  opType != OPERATOR_D4 &&
178  opType != OPERATOR_MAX,
179  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::Serial::getValues) operator is not supported");
180  }
181  }
182  }
183  };
184  };
185  }
186 
187  template<bool serendipity,
188  typename DeviceType,
189  typename outputValueType,
190  typename pointValueType>
191  class Basis_HGRAD_QUAD_DEG2_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
192  public:
194  using typename BasisBase::ExecutionSpace;
195  using typename BasisBase::OrdinalTypeArray1DHost;
196  using typename BasisBase::OrdinalTypeArray2DHost;
197  using typename BasisBase::OrdinalTypeArray3DHost;
198 
202 
203  using typename BasisBase::OutputViewType;
204  using typename BasisBase::PointViewType;
205  using typename BasisBase::ScalarViewType;
206 
207  using BasisBase::getValues;
208 
209  virtual
210  void
211  getValues( const ExecutionSpace& space,
212  OutputViewType outputValues,
213  const PointViewType inputPoints,
214  const EOperator operatorType = OPERATOR_VALUE ) const override {
215 #ifdef HAVE_INTREPID2_DEBUG
216  // Verify arguments
217  Intrepid2::getValues_HGRAD_Args(outputValues,
218  inputPoints,
219  operatorType,
220  this->getBaseCellTopology(),
221  this->getCardinality() );
222 #endif
224  template getValues<DeviceType>(space,
225  outputValues,
226  inputPoints,
227  operatorType);
228  }
229 
230  virtual
231  void
232  getDofCoords( ScalarViewType dofCoords ) const override {
233 #ifdef HAVE_INTREPID2_DEBUG
234  // Verify rank of output array.
235  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
236  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getDofCoords) rank = 2 required for dofCoords array");
237  // Verify 0th dimension of output array.
238  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
239  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
240  // Verify 1st dimension of output array.
241  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
242  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
243 #endif
244  Kokkos::deep_copy(dofCoords, this->dofCoords_);
245  }
246 
247  virtual
248  void
249  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
250 #ifdef HAVE_INTREPID2_DEBUG
251  // Verify rank of output array.
252  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
253  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
254  // Verify 0th dimension of output array.
255  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
256  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
257 #endif
258  Kokkos::deep_copy(dofCoeffs, 1.0);
259  }
260 
261  virtual
262  const char*
263  getName() const override {
264  return serendipity ? "Intrepid2_HGRAD_QUAD_I2_FEM" : "Intrepid2_HGRAD_QUAD_C2_FEM";
265  }
266 
275  BasisPtr<DeviceType, outputValueType, pointValueType>
276  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
277  if(subCellDim == 1) {
278  return Teuchos::rcp(new
280  }
281  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
282  }
283 
284  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
285  getHostBasis() const override{
287  }
288  };
289 
290  template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
291  using Basis_HGRAD_QUAD_C2_FEM = Basis_HGRAD_QUAD_DEG2_FEM<false, DeviceType, outputValueType, pointValueType>;
292 
293  template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
294  using Basis_HGRAD_QUAD_I2_FEM = Basis_HGRAD_QUAD_DEG2_FEM<true, DeviceType, outputValueType, pointValueType>;
295 
296 }// namespace Intrepid2
297 
299 
300 #endif
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual void getValues(const ExecutionSpace &, OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Header file for the Intrepid2::Basis_HGRAD_LINE_C2_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell...
See Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM.
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
BasisPtr< typename Kokkos::HostSpace::device_type, outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Definition file for FEM basis functions of degree 2 for H(grad) functions on QUAD cells...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.