Intrepid2
Intrepid2_HDIV_TET_I1_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_HDIV_TET_I1_FEM_HPP__
50 #define __INTREPID2_HDIV_TET_I1_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
54 
55 namespace Intrepid2 {
56 
103  namespace Impl {
104 
109  public:
110  typedef struct Tetrahedron<4> 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( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
130  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
131  const EOperator operatorType);
132 
136  template<typename outputValueViewType,
137  typename inputPointViewType,
138  EOperator opType>
139  struct Functor {
140  outputValueViewType _outputValues;
141  const inputPointViewType _inputPoints;
142 
143  KOKKOS_INLINE_FUNCTION
144  Functor( outputValueViewType outputValues_,
145  inputPointViewType inputPoints_ )
146  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
147 
148  KOKKOS_INLINE_FUNCTION
149  void operator()(const ordinal_type pt) const {
150  switch (opType) {
151  case OPERATOR_VALUE : {
152  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
153  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
154  Serial<opType>::getValues( output, input );
155  break;
156  }
157  case OPERATOR_DIV : {
158  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
159  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
160  Serial<opType>::getValues( output, input );
161  break;
162  }
163  default: {
164  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
165  opType != OPERATOR_DIV,
166  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::Serial::getValues) operator is not supported");
167  }
168  }
169  }
170  };
171 
172  };
173  }
174 
175  template<typename DeviceType = void,
176  typename outputValueType = double,
177  typename pointValueType = double>
178  class Basis_HDIV_TET_I1_FEM: public Basis<DeviceType,outputValueType,pointValueType> {
179  public:
180  using OrdinalTypeArray1DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost;
181  using OrdinalTypeArray2DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost;
182  using OrdinalTypeArray3DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost;
183 
187 
191 
193 
194  virtual
195  void
196  getValues( OutputViewType outputValues,
197  const PointViewType inputPoints,
198  const EOperator operatorType = OPERATOR_VALUE ) const override {
199 #ifdef HAVE_INTREPID2_DEBUG
200  // Verify arguments
201  Intrepid2::getValues_HDIV_Args(outputValues,
202  inputPoints,
203  operatorType,
204  this->getBaseCellTopology(),
205  this->getCardinality() );
206 #endif
207  Impl::Basis_HDIV_TET_I1_FEM::
208  getValues<DeviceType>( outputValues,
209  inputPoints,
210  operatorType );
211  }
212 
213  virtual
214  void
215  getDofCoords( ScalarViewType dofCoords ) const override {
216 #ifdef HAVE_INTREPID2_DEBUG
217  // Verify rank of output array.
218  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
219  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
220  // Verify 0th dimension of output array.
221  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
222  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
223  // Verify 1st dimension of output array.
224  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
225  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
226 #endif
227  Kokkos::deep_copy(dofCoords, this->dofCoords_);
228  }
229 
230  virtual
231  void
232  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
233 #ifdef HAVE_INTREPID2_DEBUG
234  // Verify rank of output array.
235  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
236  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
237  // Verify 0th dimension of output array.
238  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
239  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
240  // Verify 1st dimension of output array.
241  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
242  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
243 #endif
244  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
245  }
246 
247  virtual
248  const char*
249  getName() const override {
250  return "Intrepid2_HDIV_TET_I1_FEM";
251  }
252 
253  virtual
254  bool
255  requireOrientation() const override {
256  return true;
257  }
258 
269  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
270 
271  if(subCellDim == 2) {
272  return Teuchos::rcp(new
273  Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Triangle<3> >()));
274  }
275  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
276  }
277 
279  getHostBasis() const override{
281  }
282  };
283 
284 }// namespace Intrepid2
285 
287 
288 #endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HDIV_Args(const outputValueViewType outputValues, const inputPointViewType inputPoints, const EOperator operatorType, const shards::CellTopology cellTopo, const ordinal_type basisCard)
Runtime check of the arguments for the getValues method in an HDIV-conforming FEM basis....
Definition file for FEM basis functions of degree 1 for H(div) functions on TET cells.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
Implementation of the default H(div)-compatible FEM basis of degree 1 on a Tetrahedron cell.
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 bool requireOrientation() const override
True if orientation is required.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
virtual const char * getName() const override
Returns basis name.
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...
Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral,...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
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 > OrdinalTypeArray2DHost
View type for 2d host array.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HDIV_TET_I1_FEM.