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"
53 
54 namespace Intrepid2 {
55 
102  namespace Impl {
103 
108  public:
109  typedef struct Tetrahedron<4> cell_topology_type;
113  template<EOperator opType>
114  struct Serial {
115  template<typename OutputViewType,
116  typename inputViewType>
117  KOKKOS_INLINE_FUNCTION
118  static void
119  getValues( OutputViewType output,
120  const inputViewType input );
121 
122  };
123 
124  template<typename ExecSpaceType,
125  typename outputValueValueType, class ...outputValueProperties,
126  typename inputPointValueType, class ...inputPointProperties>
127  static void
128  getValues( 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, Kokkos::ALL() );
152  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
153  Serial<opType>::getValues( output, input );
154  break;
155  }
156  case OPERATOR_DIV : {
157  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
158  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
159  Serial<opType>::getValues( output, input );
160  break;
161  }
162  default: {
163  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
164  opType != OPERATOR_DIV,
165  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::Serial::getValues) operator is not supported");
166  }
167  }
168  }
169  };
170 
171  };
172  }
173 
174  template<typename ExecSpaceType = void,
175  typename outputValueType = double,
176  typename pointValueType = double>
177  class Basis_HDIV_TET_I1_FEM: public Basis<ExecSpaceType,outputValueType,pointValueType> {
178  public:
182 
186 
190 
192 
193  virtual
194  void
195  getValues( OutputViewType outputValues,
196  const PointViewType inputPoints,
197  const EOperator operatorType = OPERATOR_VALUE ) const {
198 #ifdef HAVE_INTREPID2_DEBUG
199  // Verify arguments
200  Intrepid2::getValues_HDIV_Args(outputValues,
201  inputPoints,
202  operatorType,
203  this->getBaseCellTopology(),
204  this->getCardinality() );
205 #endif
206  Impl::Basis_HDIV_TET_I1_FEM::
207  getValues<ExecSpaceType>( outputValues,
208  inputPoints,
209  operatorType );
210  }
211 
212  virtual
213  void
214  getDofCoords( ScalarViewType dofCoords ) const {
215 #ifdef HAVE_INTREPID2_DEBUG
216  // Verify rank of output array.
217  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
218  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
219  // Verify 0th dimension of output array.
220  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
221  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
222  // Verify 1st dimension of output array.
223  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
224  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
225 #endif
226  Kokkos::deep_copy(dofCoords, this->dofCoords_);
227  }
228 
229  virtual
230  void
231  getDofCoeffs( ScalarViewType dofCoeffs ) const {
232 #ifdef HAVE_INTREPID2_DEBUG
233  // Verify rank of output array.
234  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
235  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
236  // Verify 0th dimension of output array.
237  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
238  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
239  // Verify 1st dimension of output array.
240  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
241  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
242 #endif
243  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
244  }
245 
246  virtual
247  const char*
248  getName() const {
249  return "Intrepid2_HDIV_TET_I1_FEM";
250  }
251 
252  virtual
253  bool
255  return true;
256  }
257 
258  };
259 
260 }// namespace Intrepid2
261 
263 
264 #endif
Definition file for FEM basis functions of degree 1 for H(div) functions on TET cells.
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
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, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
virtual void getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Implementation of the default H(div)-compatible FEM basis of degree 1 on a Tetrahedron cell...
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
virtual const char * getName() const
Returns basis name.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
virtual bool requireOrientation() const
True if orientation is required.
See Intrepid2::Basis_HDIV_TET_I1_FEM.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Header file for the abstract base class Intrepid2::Basis.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...