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( const typename DeviceType::execution_space& space,
110  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
111  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
112  const EOperator operatorType);
113 
117  template<typename outputValueViewType,
118  typename inputPointViewType,
119  EOperator opType>
120  struct Functor {
121  outputValueViewType _outputValues;
122  const inputPointViewType _inputPoints;
123 
124  KOKKOS_INLINE_FUNCTION
125  Functor( outputValueViewType outputValues_,
126  inputPointViewType inputPoints_ )
127  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
128 
129  KOKKOS_INLINE_FUNCTION
130  void operator()(const ordinal_type pt) const {
131  switch (opType) {
132  case OPERATOR_VALUE : {
133  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
134  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
135  Serial<opType>::getValues( output, input );
136  break;
137  }
138  case OPERATOR_GRAD :
139  case OPERATOR_CURL :
140  case OPERATOR_D2 :
141  case OPERATOR_MAX : {
142  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
143  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
144  Serial<opType>::getValues( output, input );
145  break;
146  }
147  default: {
148  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
149  opType != OPERATOR_GRAD &&
150  opType != OPERATOR_CURL &&
151  opType != OPERATOR_D2 &&
152  opType != OPERATOR_MAX,
153  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::Serial::getValues) operator is not supported");
154  }
155  }
156  }
157  };
158  };
159  }
160 
161  template<typename DeviceType = void,
162  typename outputValueType = double,
163  typename pointValueType = double>
164  class Basis_HGRAD_QUAD_C1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
165  public:
167  using typename BasisBase::ExecutionSpace;
168  using typename BasisBase::OrdinalTypeArray1DHost;
169  using typename BasisBase::OrdinalTypeArray2DHost;
170  using typename BasisBase::OrdinalTypeArray3DHost;
171 
175 
176  using typename BasisBase::OutputViewType;
177  using typename BasisBase::PointViewType;
178  using typename BasisBase::ScalarViewType;
179 
180  using BasisBase::getValues;
181 
182  virtual
183  void
184  getValues( const ExecutionSpace& space,
185  OutputViewType outputValues,
186  const PointViewType inputPoints,
187  const EOperator operatorType = OPERATOR_VALUE ) const override {
188 #ifdef HAVE_INTREPID2_DEBUG
189  // Verify arguments
190  Intrepid2::getValues_HGRAD_Args(outputValues,
191  inputPoints,
192  operatorType,
193  this->getBaseCellTopology(),
194  this->getCardinality() );
195 #endif
196  Impl::Basis_HGRAD_QUAD_C1_FEM::
197  getValues<DeviceType>(space,
198  outputValues,
199  inputPoints,
200  operatorType);
201  }
202 
203  virtual
204  void
205  getDofCoords( ScalarViewType dofCoords ) const override {
206 #ifdef HAVE_INTREPID2_DEBUG
207  // Verify rank of output array.
208  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
209  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
210  // Verify 0th dimension of output array.
211  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
212  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
213  // Verify 1st dimension of output array.
214  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
215  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
216 #endif
217  Kokkos::deep_copy(dofCoords, this->dofCoords_);
218  }
219 
220  virtual
221  void
222  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
223 #ifdef HAVE_INTREPID2_DEBUG
224  // Verify rank of output array.
225  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
226  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
227  // Verify 0th dimension of output array.
228  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
229  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
230 #endif
231  Kokkos::deep_copy(dofCoeffs, 1.0);
232  }
233 
234  virtual
235  const char*
236  getName() const override {
237  return "Intrepid2_HGRAD_QUAD_C1_FEM";
238  }
239 
248  BasisPtr<DeviceType, outputValueType, pointValueType>
249  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
250  if(subCellDim == 1) {
251  return Teuchos::rcp(new
253  }
254  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
255  }
256 
257  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
258  getHostBasis() const override {
260  }
261  };
262 }// namespace Intrepid2
263 
265 
266 #endif
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral cell...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Definition file for FEM basis functions of degree 1 for H(grad) functions on QUAD cells...
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 const char * getName() const override
Returns basis name.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
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.
See Intrepid2::Basis_HGRAD_QUAD_C1_FEM.
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< 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::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Header file for the Intrepid2::Basis_HGRAD_LINE_C1_FEM class.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
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.
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.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.