Intrepid2
Intrepid2_HGRAD_QUAD_C1_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_HGRAD_QUAD_C1_FEM_HPP__
50 #define __INTREPID2_HGRAD_QUAD_C1_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
54 
55 namespace Intrepid2 {
56 
83  namespace Impl {
84 
89  public:
90  typedef struct Quadrilateral<4> cell_topology_type;
94  template<EOperator opType>
95  struct Serial {
96  template<typename OutputViewType,
97  typename inputViewType>
98  KOKKOS_INLINE_FUNCTION
99  static void
100  getValues( OutputViewType output,
101  const inputViewType input );
102 
103  };
104 
105  template<typename DeviceType,
106  typename outputValueValueType, class ...outputValueProperties,
107  typename inputPointValueType, class ...inputPointProperties>
108  static void
109  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
110  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
111  const EOperator operatorType);
112 
116  template<typename outputValueViewType,
117  typename inputPointViewType,
118  EOperator opType>
119  struct Functor {
120  outputValueViewType _outputValues;
121  const inputPointViewType _inputPoints;
122 
123  KOKKOS_INLINE_FUNCTION
124  Functor( outputValueViewType outputValues_,
125  inputPointViewType inputPoints_ )
126  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
127 
128  KOKKOS_INLINE_FUNCTION
129  void operator()(const ordinal_type pt) const {
130  switch (opType) {
131  case OPERATOR_VALUE : {
132  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
133  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
134  Serial<opType>::getValues( output, input );
135  break;
136  }
137  case OPERATOR_GRAD :
138  case OPERATOR_CURL :
139  case OPERATOR_D2 :
140  case OPERATOR_MAX : {
141  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
142  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
143  Serial<opType>::getValues( output, input );
144  break;
145  }
146  default: {
147  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
148  opType != OPERATOR_GRAD &&
149  opType != OPERATOR_CURL &&
150  opType != OPERATOR_D2 &&
151  opType != OPERATOR_MAX,
152  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::Serial::getValues) operator is not supported");
153  }
154  }
155  }
156  };
157  };
158  }
159 
160  template<typename DeviceType = void,
161  typename outputValueType = double,
162  typename pointValueType = double>
163  class Basis_HGRAD_QUAD_C1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
164  public:
165  using OrdinalTypeArray1DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost;
166  using OrdinalTypeArray2DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost;
167  using OrdinalTypeArray3DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost;
168 
172 
176 
178 
179  virtual
180  void
181  getValues( OutputViewType outputValues,
182  const PointViewType inputPoints,
183  const EOperator operatorType = OPERATOR_VALUE ) const override {
184 #ifdef HAVE_INTREPID2_DEBUG
185  // Verify arguments
186  Intrepid2::getValues_HGRAD_Args(outputValues,
187  inputPoints,
188  operatorType,
189  this->getBaseCellTopology(),
190  this->getCardinality() );
191 #endif
192  Impl::Basis_HGRAD_QUAD_C1_FEM::
193  getValues<DeviceType>( outputValues,
194  inputPoints,
195  operatorType );
196  }
197 
198  virtual
199  void
200  getDofCoords( ScalarViewType dofCoords ) const override {
201 #ifdef HAVE_INTREPID2_DEBUG
202  // Verify rank of output array.
203  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
204  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
205  // Verify 0th dimension of output array.
206  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
207  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
208  // Verify 1st dimension of output array.
209  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
210  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
211 #endif
212  Kokkos::deep_copy(dofCoords, this->dofCoords_);
213  }
214 
215  virtual
216  void
217  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
218 #ifdef HAVE_INTREPID2_DEBUG
219  // Verify rank of output array.
220  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
221  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
222  // Verify 0th dimension of output array.
223  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
224  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
225 #endif
226  Kokkos::deep_copy(dofCoeffs, 1.0);
227  }
228 
229  virtual
230  const char*
231  getName() const override {
232  return "Intrepid2_HGRAD_QUAD_C1_FEM";
233  }
234 
244  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
245  if(subCellDim == 1) {
246  return Teuchos::rcp(new
248  }
249  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
250  }
251 
253  getHostBasis() const override {
255  }
256  };
257 }// namespace Intrepid2
258 
260 
261 #endif
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral cell...
Definition file for FEM basis functions of degree 1 for H(grad) functions on QUAD cells...
virtual const char * getName() const override
Returns basis name.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
ordinal_type getCardinality() const
Returns cardinality of the basis.
See Intrepid2::Basis_HGRAD_QUAD_C1_FEM.
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...
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
void getValues_HGRAD_Args(const outputValueViewType outputValues, const inputPointViewType inputPoints, const EOperator operatorType, const shards::CellTopology cellTopo, const ordinal_type basisCard)
Runtime check of the arguments for the getValues method in an HGRAD-conforming FEM basis...
Header file for the Intrepid2::Basis_HGRAD_LINE_C1_FEM class.
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(grad)-compatible FEM basis of degree 1 on Line cell.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
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.