Intrepid2
Intrepid2_HGRAD_TET_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_TET_C2_FEM_HPP__
50 #define __INTREPID2_HGRAD_TET_C2_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
54 
55 namespace Intrepid2 {
56 
101  namespace Impl {
102 
107  public:
108  typedef struct Tetrahedron<10> cell_topology_type;
112  template<EOperator opType>
113  struct Serial {
114  template<typename OutputViewType,
115  typename inputViewType>
116  KOKKOS_INLINE_FUNCTION
117  static void
118  getValues( OutputViewType output,
119  const inputViewType input );
120 
121  };
122 
123  template<typename DeviceType,
124  typename outputValueValueType, class ...outputValueProperties,
125  typename inputPointValueType, class ...inputPointProperties>
126  static void
127  getValues( const typename DeviceType::execution_space& space,
128  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
129  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
130  const EOperator operatorType);
131 
135  template<typename outputValueViewType,
136  typename inputPointViewType,
137  EOperator opType>
138  struct Functor {
139  outputValueViewType _outputValues;
140  const inputPointViewType _inputPoints;
141 
142  KOKKOS_INLINE_FUNCTION
143  Functor( outputValueViewType outputValues_,
144  inputPointViewType inputPoints_ )
145  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
146 
147  KOKKOS_INLINE_FUNCTION
148  void operator()(const ordinal_type pt) const {
149  switch (opType) {
150  case OPERATOR_VALUE : {
151  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
152  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
153  Serial<opType>::getValues( output, input );
154  break;
155  }
156  case OPERATOR_GRAD :
157  case OPERATOR_D1 :
158  case OPERATOR_D2 :
159  case OPERATOR_MAX : {
160  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
161  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
162  Serial<opType>::getValues( output, input );
163  break;
164  }
165  default: {
166  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
167  opType != OPERATOR_GRAD &&
168  opType != OPERATOR_MAX,
169  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::Serial::getValues) operator is not supported");
170  }
171  }
172  }
173  };
174  };
175  }
176 
177  template<typename DeviceType = void,
178  typename outputValueType = double,
179  typename pointValueType = double>
180  class Basis_HGRAD_TET_C2_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
181  public:
183  using typename BasisBase::ExecutionSpace;
184  using typename BasisBase::OrdinalTypeArray1DHost;
185  using typename BasisBase::OrdinalTypeArray2DHost;
186  using typename BasisBase::OrdinalTypeArray3DHost;
187 
191 
192  using typename BasisBase::OutputViewType;
193  using typename BasisBase::PointViewType;
194  using typename BasisBase::ScalarViewType;
195 
196  using BasisBase::getValues;
197 
198  virtual
199  void
200  getValues( const ExecutionSpace& space,
201  OutputViewType outputValues,
202  const PointViewType inputPoints,
203  const EOperator operatorType = OPERATOR_VALUE ) const override {
204 #ifdef HAVE_INTREPID2_DEBUG
205  // Verify arguments
206  Intrepid2::getValues_HGRAD_Args(outputValues,
207  inputPoints,
208  operatorType,
209  this->getBaseCellTopology(),
210  this->getCardinality() );
211 #endif
212  Impl::Basis_HGRAD_TET_C2_FEM::
213  getValues<DeviceType>(space,
214  outputValues,
215  inputPoints,
216  operatorType);
217  }
218 
219  virtual
220  void
221  getDofCoords( ScalarViewType dofCoords ) const override {
222 #ifdef HAVE_INTREPID2_DEBUG
223  // Verify rank of output array.
224  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
225  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
226  // Verify 0th dimension of output array.
227  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
228  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
229  // Verify 1st dimension of output array.
230  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
231  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
232 #endif
233  Kokkos::deep_copy(dofCoords, this->dofCoords_);
234  }
235 
236 
237  virtual
238  void
239  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
240 #ifdef HAVE_INTREPID2_DEBUG
241  // Verify rank of output array.
242  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
243  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
244  // Verify 0th dimension of output array.
245  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
246  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
247 #endif
248  Kokkos::deep_copy(dofCoeffs, 1.0);
249  }
250 
251  virtual
252  const char*
253  getName() const override {
254  return "Intrepid2_HGRAD_TET_C2_FEM";
255  }
256 
265  BasisPtr<DeviceType,outputValueType,pointValueType>
266  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
267  if(subCellDim == 1) {
269  } else if(subCellDim == 2) {
271  }
272 
273  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
274  }
275 
276  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
277  getHostBasis() const override{
279  }
280  };
281 }// namespace Intrepid2
282 
284 
285 #endif
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...
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.
See Intrepid2::Basis_HGRAD_TET_C2_FEM.
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...
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
virtual const char * getName() const override
Returns basis name.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell...
Header file for the Intrepid2::Basis_HGRAD_TRI_C2_FEM class.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Definition file for FEM basis functions of degree 2 for H(grad) functions on TET cells.
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.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Tetrahedron cell...
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.
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.