Intrepid2
Intrepid2_HGRAD_TET_C1_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_C1_FEM_HPP__
50 #define __INTREPID2_HGRAD_TET_C1_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
54 
55 namespace Intrepid2 {
56 
82  namespace Impl {
83 
88  public:
89  typedef struct Tetrahedron<4> cell_topology_type;
93  template<EOperator opType>
94  struct Serial {
95  template<typename OutputViewType,
96  typename inputViewType>
97  KOKKOS_INLINE_FUNCTION
98  static void
99  getValues( OutputViewType output,
100  const inputViewType input );
101 
102  };
103 
104  template<typename DeviceType,
105  typename outputValueValueType, class ...outputValueProperties,
106  typename inputPointValueType, class ...inputPointProperties>
107  static void
108  getValues( const typename DeviceType::execution_space& space,
109  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
110  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
111  const EOperator operatorType);
112 
116  template<typename outputValueViewType,
117  typename inputPointViewType,
118  EOperator opType>
119  struct Functor {
120  outputValueViewType _outputValues;
121  const inputPointViewType _inputPoints;
122 
123  KOKKOS_INLINE_FUNCTION
124  Functor( outputValueViewType outputValues_,
125  inputPointViewType inputPoints_ )
126  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
127 
128  KOKKOS_INLINE_FUNCTION
129  void operator()(const ordinal_type pt) const {
130  switch (opType) {
131  case OPERATOR_VALUE : {
132  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
133  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
134  Serial<opType>::getValues( output, input );
135  break;
136  }
137  case OPERATOR_GRAD :
138  case OPERATOR_MAX : {
139  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
140  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
141  Serial<opType>::getValues( output, input );
142  break;
143  }
144  default: {
145  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
146  opType != OPERATOR_GRAD &&
147  opType != OPERATOR_MAX,
148  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C1_FEM::Serial::getValues) operator is not supported");
149  }
150  }
151  }
152  };
153 
154 
155  };
156  }
157 
158  template<typename DeviceType = void,
159  typename outputValueType = double,
160  typename pointValueType = double>
161  class Basis_HGRAD_TET_C1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
162  public:
164  using typename BasisBase::ExecutionSpace;
165  using typename BasisBase::OrdinalTypeArray1DHost;
166  using typename BasisBase::OrdinalTypeArray2DHost;
167  using typename BasisBase::OrdinalTypeArray3DHost;
168 
172 
173  using typename BasisBase::OutputViewType;
174  using typename BasisBase::PointViewType;
175  using typename BasisBase::ScalarViewType;
176 
177  using BasisBase::getValues;
178 
179  virtual
180  void
181  getValues( const ExecutionSpace& space,
182  OutputViewType outputValues,
183  const PointViewType inputPoints,
184  const EOperator operatorType = OPERATOR_VALUE ) const override {
185 #ifdef HAVE_INTREPID2_DEBUG
186  // Verify arguments
187  Intrepid2::getValues_HGRAD_Args(outputValues,
188  inputPoints,
189  operatorType,
190  this->getBaseCellTopology(),
191  this->getCardinality() );
192 #endif
193  Impl::Basis_HGRAD_TET_C1_FEM::
194  getValues<DeviceType>(space,
195  outputValues,
196  inputPoints,
197  operatorType);
198  }
199 
200  virtual
201  void
202  getDofCoords( ScalarViewType dofCoords ) const override {
203 #ifdef HAVE_INTREPID2_DEBUG
204  // Verify rank of output array.
205  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
206  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
207  // Verify 0th dimension of output array.
208  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
209  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
210  // Verify 1st dimension of output array.
211  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
212  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
213 #endif
214  Kokkos::deep_copy(dofCoords, this->dofCoords_);
215  }
216 
217  virtual
218  void
219  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
220 #ifdef HAVE_INTREPID2_DEBUG
221  // Verify rank of output array.
222  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
223  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
224  // Verify 0th dimension of output array.
225  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
226  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
227 #endif
228  Kokkos::deep_copy(dofCoeffs, 1.0);
229  }
230 
231  virtual
232  const char*
233  getName() const override {
234  return "Intrepid2_HGRAD_TET_C1_FEM";
235  }
236 
245  BasisPtr<DeviceType, outputValueType, pointValueType>
246  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
247  if(subCellDim == 1) {
248  return Teuchos::rcp(new
250  } else if(subCellDim == 1) {
251  return Teuchos::rcp(new
253  }
254  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
255  }
256 
257  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
258  getHostBasis() const override{
260  }
261  };
262 }// namespace Intrepid2
263 
265 
266 #endif
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::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Header file for the Intrepid2::Basis_HGRAD_TRI_C1_FEM class.
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< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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 1 for H(grad) functions on TET cells.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Triangle cell...
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...
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::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
See Intrepid2::Basis_HGRAD_TET_C1_FEM.
virtual const char * getName() const override
Returns basis name.
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.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Tetrahedron cell...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference 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.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.