Intrepid2
Intrepid2_HCURL_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_HCURL_TET_I1_FEM_HPP__
50 #define __INTREPID2_HCURL_TET_I1_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
53 
54 namespace Intrepid2 {
55 
103  namespace Impl {
104 
109  public:
110  typedef struct Tetrahedron<4> cell_topology_type;
111 
115  template<EOperator opType>
116  struct Serial {
117  template<typename OutputViewType,
118  typename inputViewType>
119  KOKKOS_INLINE_FUNCTION
120  static void
121  getValues( OutputViewType output,
122  const inputViewType input );
123 
124  };
125 
126  template<typename ExecSpaceType,
127  typename outputValueValueType, class ...outputValueProperties,
128  typename inputPointValueType, class ...inputPointProperties>
129  static void
130  getValues( 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, Kokkos::ALL() );
154  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
155  Serial<opType>::getValues( output, input );
156  break;
157  }
158  case OPERATOR_CURL : {
159  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
160  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
161  Serial<opType>::getValues( output, input );
162  break;
163  }
164  default: {
165  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
166  opType != OPERATOR_CURL,
167  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::Serial::getValues) operator is not supported");
168  }
169  }
170  }
171  };
172  };
173  }
174 
175  template<typename ExecSpaceType = void,
176  typename outputValueType = double,
177  typename pointValueType = double>
178  class Basis_HCURL_TET_I1_FEM : public Basis<ExecSpaceType,outputValueType,pointValueType> {
179  public:
183 
187 
191 
193 
194  virtual
195  void
196  getValues( OutputViewType outputValues,
197  const PointViewType inputPoints,
198  const EOperator operatorType = OPERATOR_VALUE ) const {
199 #ifdef HAVE_INTREPID2_DEBUG
200  // Verify arguments
201  Intrepid2::getValues_HCURL_Args(outputValues,
202  inputPoints,
203  operatorType,
204  this->getBaseCellTopology(),
205  this->getCardinality() );
206 #endif
207  Impl::Basis_HCURL_TET_I1_FEM::
208  getValues<ExecSpaceType>( outputValues,
209  inputPoints,
210  operatorType );
211  }
212 
213  virtual
214  void
215  getDofCoords( ScalarViewType dofCoords ) const {
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_HCURL_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_HCURL_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_HCURL_TET_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
226 #endif
227  Kokkos::deep_copy(dofCoords, this->dofCoords_);
228  }
229 
230 
231  virtual
232  void
233  getDofCoeffs( ScalarViewType dofCoeffs ) const {
234 #ifdef HAVE_INTREPID2_DEBUG
235  // Verify rank of output array.
236  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
237  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
238  // Verify 0th dimension of output array.
239  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
240  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
241  // Verify 1st dimension of output array.
242  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
243  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
244 #endif
245  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
246  }
247 
248  virtual
249  const char*
250  getName() const {
251  return "Intrepid2_HCURL_TET_I1_FEM";
252  }
253 
254  virtual
255  bool
257  return true;
258  }
259 
260  };
261 }// namespace Intrepid2
262 
264 
265 #endif
virtual void getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HCURL_TET_I1_FEM.
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 bool requireOrientation() const
True if orientation is required.
ordinal_type getCardinality() const
Returns cardinality of the 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...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
virtual const char * getName() const
Returns basis name.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
Definition file for FEM basis functions of degree 1 for H(curl) functions on TET cells.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
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...
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Tetrahedron cell...
Header file for the abstract base class Intrepid2::Basis.
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...