Intrepid2
Intrepid2_HCURL_HEX_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_HEX_I1_FEM_HPP__
50 #define __INTREPID2_HCURL_HEX_I1_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
54 
55 namespace Intrepid2 {
56 
117  namespace Impl {
118 
123  public:
124  typedef struct Hexahedron<8> cell_topology_type;
125 
129  template<EOperator opType>
130  struct Serial {
131  template<typename OutputViewType,
132  typename inputViewType>
133  KOKKOS_INLINE_FUNCTION
134  static void
135  getValues( OutputViewType output,
136  const inputViewType input );
137  };
138 
139  template<typename DeviceType,
140  typename outputValueValueType, class ...outputValueProperties,
141  typename inputPointValueType, class ...inputPointProperties>
142  static void
143  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
144  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
145  const EOperator operatorType);
146 
150  template<typename outputValueViewType,
151  typename inputPointViewType,
152  EOperator opType>
153  struct Functor {
154  outputValueViewType _outputValues;
155  const inputPointViewType _inputPoints;
156 
157  KOKKOS_INLINE_FUNCTION
158  Functor( outputValueViewType outputValues_,
159  inputPointViewType inputPoints_ )
160  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
161 
162  KOKKOS_INLINE_FUNCTION
163  void operator()(const ordinal_type pt) const {
164  switch (opType) {
165  case OPERATOR_VALUE:
166  case OPERATOR_CURL: {
167  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
168  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
169  Serial<opType>::getValues( output, input );
170  break;
171  }
172  default: {
173  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
174  opType != OPERATOR_CURL,
175  ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_CI_FEM::Serial::getVAlues) operator is not supported");
176  break;
177  }
178  } //end switch
179  }
180  };
181  };
182  }
183 
184  template<typename DeviceType = void,
185  typename outputValueType = double,
186  typename pointValueType = double>
187  class Basis_HCURL_HEX_I1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
188  public:
189  using OrdinalTypeArray1DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost;
190  using OrdinalTypeArray2DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost;
191  using OrdinalTypeArray3DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost;
192 
196 
200 
202 
203  virtual
204  void
205  getValues( OutputViewType outputValues,
206  const PointViewType inputPoints,
207  const EOperator operatorType = OPERATOR_VALUE ) const override {
208 #ifdef HAVE_INTREPID2_DEBUG
209  Intrepid2::getValues_HCURL_Args(outputValues,
210  inputPoints,
211  operatorType,
212  this->getBaseCellTopology(),
213  this->getCardinality() );
214 #endif
215  Impl::Basis_HCURL_HEX_I1_FEM::
216  getValues<DeviceType>( outputValues,
217  inputPoints,
218  operatorType );
219  }
220 
221  virtual
222  void
223  getDofCoords( ScalarViewType dofCoords ) const override {
224 #ifdef HAVE_INTREPID2_DEBUG
225  // Verify rank of output array.
226  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
227  ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
228  // Verify 0th dimension of output array.
229  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
230  ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
231  // Verify 1st dimension of output array.
232  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
233  ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
234 #endif
235  Kokkos::deep_copy(dofCoords, this->dofCoords_);
236  }
237 
238 
239  virtual
240  void
241  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
242 #ifdef HAVE_INTREPID2_DEBUG
243  // Verify rank of output array.
244  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
245  ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
246  // Verify 0th dimension of output array.
247  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
248  ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
249  // Verify 1st dimension of output array.
250  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
251  ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
252 #endif
253  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
254  }
255 
256  virtual
257  const char*
258  getName() const override {
259  return "Intrepid2_HCURL_HEX_I1_FEM";
260  }
261 
262  virtual
263  bool
264  requireOrientation() const override {
265  return true;
266  }
267 
277  BasisPtr<DeviceType,outputValueType,pointValueType>
278  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
279 
280  if(subCellDim == 1)
281  return Teuchos::rcp( new
282  Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Line<2> >()));
283  else if(subCellDim == 2) {
284  return Teuchos::rcp(new
286  }
287 
288  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
289  }
290 
291  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
292  getHostBasis() const override{
294  }
295 
296  };
297 }// namespace Intrepid2
298 
300 
301 #endif
virtual const char * getName() const override
Returns basis name.
ordinal_type getCardinality() const
Returns cardinality of the basis.
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...
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Quadrilateral cell...
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 bool requireOrientation() const override
True if orientation is required.
Definition file for FEM basis functions of degree 1 for H(curl) functions on HEX cells.
Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral, hexahedron and tetrahedron cells.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
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 H(curl)-compatible FEM basis of degree 1 on Hexahedron cell...
See Intrepid2::Basis_HCURL_HEX_I1_FEM.
Header file for the Intrepid2::Basis_HCURL_QUAD_I1_FEM class.
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...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Header file for the abstract base class Intrepid2::Basis.