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"
53 
54 namespace Intrepid2 {
55 
82  namespace Impl {
83 
88  public:
89  typedef struct Quadrilateral<4> cell_topology_type;
93  template<EOperator opType>
94  struct Serial {
95  template<typename OutputViewType,
96  typename inputViewType>
97  KOKKOS_INLINE_FUNCTION
98  static void
99  getValues( OutputViewType output,
100  const inputViewType input );
101 
102  };
103 
104  template<typename ExecSpaceType,
105  typename outputValueValueType, class ...outputValueProperties,
106  typename inputPointValueType, class ...inputPointProperties>
107  static void
108  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
109  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
110  const EOperator operatorType);
111 
115  template<typename outputValueViewType,
116  typename inputPointViewType,
117  EOperator opType>
118  struct Functor {
119  outputValueViewType _outputValues;
120  const inputPointViewType _inputPoints;
121 
122  KOKKOS_INLINE_FUNCTION
123  Functor( outputValueViewType outputValues_,
124  inputPointViewType inputPoints_ )
125  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
126 
127  KOKKOS_INLINE_FUNCTION
128  void operator()(const ordinal_type pt) const {
129  switch (opType) {
130  case OPERATOR_VALUE : {
131  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
132  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
133  Serial<opType>::getValues( output, input );
134  break;
135  }
136  case OPERATOR_GRAD :
137  case OPERATOR_CURL :
138  case OPERATOR_D2 :
139  case OPERATOR_MAX : {
140  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
141  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
142  Serial<opType>::getValues( output, input );
143  break;
144  }
145  default: {
146  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
147  opType != OPERATOR_GRAD &&
148  opType != OPERATOR_CURL &&
149  opType != OPERATOR_D2 &&
150  opType != OPERATOR_MAX,
151  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::Serial::getValues) operator is not supported");
152  }
153  }
154  }
155  };
156  };
157  }
158 
159  template<typename ExecSpaceType = void,
160  typename outputValueType = double,
161  typename pointValueType = double>
162  class Basis_HGRAD_QUAD_C1_FEM : public Basis<ExecSpaceType,outputValueType,pointValueType> {
163  public:
167 
171 
175 
177 
178  virtual
179  void
180  getValues( OutputViewType outputValues,
181  const PointViewType inputPoints,
182  const EOperator operatorType = OPERATOR_VALUE ) const {
183 #ifdef HAVE_INTREPID2_DEBUG
184  // Verify arguments
185  Intrepid2::getValues_HGRAD_Args(outputValues,
186  inputPoints,
187  operatorType,
188  this->getBaseCellTopology(),
189  this->getCardinality() );
190 #endif
191  Impl::Basis_HGRAD_QUAD_C1_FEM::
192  getValues<ExecSpaceType>( outputValues,
193  inputPoints,
194  operatorType );
195  }
196 
197  virtual
198  void
199  getDofCoords( ScalarViewType dofCoords ) const {
200 #ifdef HAVE_INTREPID2_DEBUG
201  // Verify rank of output array.
202  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
203  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
204  // Verify 0th dimension of output array.
205  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
206  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
207  // Verify 1st dimension of output array.
208  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
209  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
210 #endif
211  Kokkos::deep_copy(dofCoords, this->dofCoords_);
212  }
213 
214  virtual
215  void
216  getDofCoeffs( ScalarViewType dofCoeffs ) const {
217 #ifdef HAVE_INTREPID2_DEBUG
218  // Verify rank of output array.
219  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
220  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
221  // Verify 0th dimension of output array.
222  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
223  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
224 #endif
225  Kokkos::deep_copy(dofCoeffs, 1.0);
226  }
227 
228  virtual
229  const char*
230  getName() const {
231  return "Intrepid2_HGRAD_QUAD_C1_FEM";
232  }
233 
234  };
235 }// namespace Intrepid2
236 
238 
239 #endif
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral cell...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Definition file for FEM basis functions of degree 1 for H(grad) functions on QUAD cells...
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
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, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
See Intrepid2::Basis_HGRAD_QUAD_C1_FEM.
virtual const char * getName() const
Returns basis name.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
virtual void getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
Header file for the abstract base class Intrepid2::Basis.