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"
54 
55 namespace Intrepid2 {
56 
104  namespace Impl {
105 
110  public:
111  typedef struct Tetrahedron<4> cell_topology_type;
112 
116  template<EOperator opType>
117  struct Serial {
118  template<typename OutputViewType,
119  typename inputViewType>
120  KOKKOS_INLINE_FUNCTION
121  static void
122  getValues( OutputViewType output,
123  const inputViewType input );
124 
125  };
126 
127  template<typename DeviceType,
128  typename outputValueValueType, class ...outputValueProperties,
129  typename inputPointValueType, class ...inputPointProperties>
130  static void
131  getValues( const typename DeviceType::execution_space& space,
132  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
133  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
134  const EOperator operatorType);
135 
139  template<typename outputValueViewType,
140  typename inputPointViewType,
141  EOperator opType>
142  struct Functor {
143  outputValueViewType _outputValues;
144  const inputPointViewType _inputPoints;
145 
146  KOKKOS_INLINE_FUNCTION
147  Functor( outputValueViewType outputValues_,
148  inputPointViewType inputPoints_ )
149  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
150 
151  KOKKOS_INLINE_FUNCTION
152  void operator()(const ordinal_type pt) const {
153  switch (opType) {
154  case OPERATOR_VALUE : {
155  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
156  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
157  Serial<opType>::getValues( output, input );
158  break;
159  }
160  case OPERATOR_CURL : {
161  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
162  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
163  Serial<opType>::getValues( output, input );
164  break;
165  }
166  default: {
167  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
168  opType != OPERATOR_CURL,
169  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_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_HCURL_TET_I1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
181  public:
183  using typename BasisBase::ExecutionSpace;
184 
185  using typename BasisBase::OrdinalTypeArray1DHost;
186  using typename BasisBase::OrdinalTypeArray2DHost;
187  using typename BasisBase::OrdinalTypeArray3DHost;
188 
189  using typename BasisBase::OutputViewType;
190  using typename BasisBase::PointViewType ;
191  using typename BasisBase::ScalarViewType;
192 
196 
197  using BasisBase::getValues;
198 
199  virtual
200  void
201  getValues( const ExecutionSpace& space,
202  OutputViewType outputValues,
203  const PointViewType inputPoints,
204  const EOperator operatorType = OPERATOR_VALUE ) const override {
205 #ifdef HAVE_INTREPID2_DEBUG
206  // Verify arguments
207  Intrepid2::getValues_HCURL_Args(outputValues,
208  inputPoints,
209  operatorType,
210  this->getBaseCellTopology(),
211  this->getCardinality() );
212 #endif
213  Impl::Basis_HCURL_TET_I1_FEM::
214  getValues<DeviceType>(space,
215  outputValues,
216  inputPoints,
217  operatorType);
218  }
219 
220  virtual
221  void
222  getDofCoords( ScalarViewType dofCoords ) const override {
223 #ifdef HAVE_INTREPID2_DEBUG
224  // Verify rank of output array.
225  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
226  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
227  // Verify 0th dimension of output array.
228  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
229  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
230  // Verify 1st dimension of output array.
231  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
232  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
233 #endif
234  Kokkos::deep_copy(dofCoords, this->dofCoords_);
235  }
236 
237 
238  virtual
239  void
240  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
241 #ifdef HAVE_INTREPID2_DEBUG
242  // Verify rank of output array.
243  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 2, std::invalid_argument,
244  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
245  // Verify 0th dimension of output array.
246  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
247  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
248  // Verify 1st dimension of output array.
249  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
250  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
251 #endif
252  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
253  }
254 
255  virtual
256  const char*
257  getName() const override {
258  return "Intrepid2_HCURL_TET_I1_FEM";
259  }
260 
261  virtual
262  bool
263  requireOrientation() const override {
264  return true;
265  }
266 
276  BasisPtr<DeviceType,outputValueType,pointValueType>
277  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
278 
279  if(subCellDim == 1)
280  return Teuchos::rcp( new
281  Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Line<2> >()));
282  else if(subCellDim == 2) {
283  return Teuchos::rcp(new
285  }
286 
287  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
288  }
289 
290  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
291  getHostBasis() const override{
293  }
294 
295 
296  };
297 }// namespace Intrepid2
298 
300 
301 #endif
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_HCURL_TET_I1_FEM.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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.
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, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
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...
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< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Definition file for FEM basis functions of degree 1 for H(curl) functions on TET cells.
Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral, hexahedron and tetrahedron cells.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Triangle cell...
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
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 > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
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...
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Tetrahedron cell...
Header file for the Intrepid2::Basis_HCURL_TRI_I1_FEM class.
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.
virtual const char * getName() const override
Returns basis name.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.