Intrepid2
Intrepid2_HCURL_TRI_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_TRI_I1_FEM_HPP
50 #define INTREPID2_HCURL_TRI_I1_FEM_HPP
51 
52 #include "Intrepid2_Basis.hpp"
54 
55 namespace Intrepid2 {
56 
97  namespace Impl {
98 
103  public:
104  typedef struct Triangle<3> cell_topology_type;
108  template<EOperator opType>
109  struct Serial {
110  template<typename OutputViewType,
111  typename inputViewType>
112  KOKKOS_INLINE_FUNCTION
113  static void
114  getValues( OutputViewType output,
115  const inputViewType input );
116 
117  };
118 
119  template<typename DeviceType,
120  typename outputValueValueType, class ...outputValueProperties,
121  typename inputPointValueType, class ...inputPointProperties>
122  static void
123  getValues( const typename DeviceType::execution_space& space,
124  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
125  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
126  const EOperator operatorType);
127 
131  template<typename outputValueViewType,
132  typename inputPointViewType,
133  EOperator opType>
134  struct Functor {
135  outputValueViewType _outputValues;
136  const inputPointViewType _inputPoints;
137 
138  KOKKOS_INLINE_FUNCTION
139  Functor( outputValueViewType outputValues_,
140  inputPointViewType inputPoints_ )
141  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
142 
143  KOKKOS_INLINE_FUNCTION
144  void operator()(const ordinal_type pt) const {
145  switch (opType) {
146  case OPERATOR_VALUE : {
147  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
148  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
149  Serial<opType>::getValues( output, input );
150  break;
151  }
152  case OPERATOR_CURL : {
153  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
154  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
155  Serial<opType>::getValues( output, input );
156  break;
157  }
158  case OPERATOR_DIV: {
159  INTREPID2_TEST_FOR_ABORT( (opType == OPERATOR_DIV),
160  ">>> ERROR (Basis_HCURL_TRI_I1_FEM): DIV is invalid operator for HCURL Basis Functions");
161  break;
162  }
163  case OPERATOR_GRAD: {
164  INTREPID2_TEST_FOR_ABORT( (opType == OPERATOR_GRAD),
165  ">>> ERROR (Basis_HCURL_TRI_I1_FEM): GRAD is invalid operator for HCURL Basis Functions");
166  break;
167  }
168  default: {
169  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
170  opType != OPERATOR_CURL,
171  ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::Serial::getValues) operator is not supported");
172  }
173  }
174  }
175  };
176 
177  };
178  }
179 
180  template< typename DeviceType = void,
181  typename outputValueType = double,
182  typename pointValueType = double >
183  class Basis_HCURL_TRI_I1_FEM : public Basis<DeviceType, outputValueType, pointValueType> {
184  public:
186  using typename BasisBase::ExecutionSpace;
187 
188  using typename BasisBase::OrdinalTypeArray1DHost;
189  using typename BasisBase::OrdinalTypeArray2DHost;
190  using typename BasisBase::OrdinalTypeArray3DHost;
191 
192  using typename BasisBase::OutputViewType;
193  using typename BasisBase::PointViewType ;
194  using typename BasisBase::ScalarViewType;
195 
199 
200  using BasisBase::getValues;
201 
202  virtual
203  void
204  getValues( const ExecutionSpace& space,
205  OutputViewType outputValues,
206  const PointViewType inputPoints,
207  const EOperator operatorType = OPERATOR_VALUE ) const override {
208 #ifdef HAVE_INTREPID2_DEBUG
209  // Verify arguments
210  Intrepid2::getValues_HCURL_Args(outputValues,
211  inputPoints,
212  operatorType,
213  this->getBaseCellTopology(),
214  this->getCardinality() );
215 #endif
216  Impl::Basis_HCURL_TRI_I1_FEM::
217  getValues<DeviceType>(space,
218  outputValues,
219  inputPoints,
220  operatorType);
221  }
222 
223  virtual
224  void
225  getDofCoords( ScalarViewType dofCoords ) const override {
226 #ifdef HAVE_INTREPID2_DEBUG
227  // Verify rank of output array.
228  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
229  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
230  // Verify 0th dimension of output array.
231  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
232  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
233  // Verify 1st dimension of output array.
234  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
235  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
236 #endif
237  Kokkos::deep_copy(dofCoords, this->dofCoords_);
238  }
239 
240  virtual
241  void
242  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
243 #ifdef HAVE_INTREPID2_DEBUG
244  // Verify rank of output array.
245  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 2, std::invalid_argument,
246  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
247  // Verify 0th dimension of output array.
248  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
249  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
250  // Verify 1st dimension of output array.
251  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
252  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
253 #endif
254  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
255  }
256 
257 
258  virtual
259  const char*
260  getName() const override {
261  return "Intrepid2_HCURL_TRI_I1_FEM";
262  }
263 
264  virtual
265  bool
266  requireOrientation() const override {
267  return true;
268  }
269 
279  BasisPtr<DeviceType,outputValueType,pointValueType>
280  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
281  if(subCellDim == 1)
282  return Teuchos::rcp( new
283  Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Line<2> >()));
284 
285  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
286  }
287 
288  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
289  getHostBasis() const override{
291  }
292  };
293 
294 }// namespace Intrepid2
295 
297 
298 #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.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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 bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Definition file for default FEM basis functions of degree 1 for H(curl) functions on Triangle cells...
See Intrepid2::Basis_HCURL_TRI_I1_FEM.
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 > OrdinalTypeArray2DHost
View type for 2d host array.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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. ...
virtual const char * getName() const override
Returns basis name.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
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...
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...
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.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.