Intrepid2
Intrepid2_HGRAD_TRI_C2_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_TRI_C2_FEM_HPP__
50 #define __INTREPID2_HGRAD_TRI_C2_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
54 
55 namespace Intrepid2 {
56 
87  namespace Impl {
88 
93  public:
94  typedef struct Triangle<6> cell_topology_type;
98  template<EOperator opType>
99  struct Serial {
100  template<typename OutputViewType,
101  typename inputViewType>
102  KOKKOS_INLINE_FUNCTION
103  static void
104  getValues( OutputViewType output,
105  const inputViewType input );
106 
107  };
108 
109  template<typename DeviceType,
110  typename outputValueValueType, class ...outputValueProperties,
111  typename inputPointValueType, class ...inputPointProperties>
112  static void
113  getValues( const typename DeviceType::execution_space& space,
114  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
115  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
116  const EOperator operatorType);
117 
121  template<typename outputValueViewType,
122  typename inputPointViewType,
123  EOperator opType>
124  struct Functor {
125  outputValueViewType _outputValues;
126  const inputPointViewType _inputPoints;
127 
128  KOKKOS_INLINE_FUNCTION
129  Functor( outputValueViewType outputValues_,
130  inputPointViewType inputPoints_ )
131  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
132 
133  KOKKOS_INLINE_FUNCTION
134  void operator()(const ordinal_type pt) const {
135  switch (opType) {
136  case OPERATOR_VALUE : {
137  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
138  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
139  Serial<opType>::getValues( output, input );
140  break;
141  }
142  case OPERATOR_GRAD :
143  case OPERATOR_CURL :
144  case OPERATOR_D1 :
145  case OPERATOR_D2 :
146  case OPERATOR_MAX : {
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  default: {
153  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
154  opType != OPERATOR_GRAD &&
155  opType != OPERATOR_CURL &&
156  opType != OPERATOR_D1 &&
157  opType != OPERATOR_D2 &&
158  opType != OPERATOR_MAX,
159  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::Serial::getValues) operator is not supported");
160  }
161  }
162  }
163  };
164  };
165  }
166 
167  template<typename DeviceType = void,
168  typename outputValueType = double,
169  typename pointValueType = double>
170  class Basis_HGRAD_TRI_C2_FEM: public Basis<DeviceType,outputValueType,pointValueType> {
171  public:
173  using typename BasisBase::ExecutionSpace;
174  using typename BasisBase::OrdinalTypeArray1DHost;
175  using typename BasisBase::OrdinalTypeArray2DHost;
176  using typename BasisBase::OrdinalTypeArray3DHost;
177 
181 
182  using typename BasisBase::OutputViewType;
183  using typename BasisBase::PointViewType;
184  using typename BasisBase::ScalarViewType;
185 
186  using BasisBase::getValues;
187 
188  virtual
189  void
190  getValues( const ExecutionSpace& space,
191  OutputViewType outputValues,
192  const PointViewType inputPoints,
193  const EOperator operatorType = OPERATOR_VALUE ) const override {
194 #ifdef HAVE_INTREPID2_DEBUG
195  // Verify arguments
196  Intrepid2::getValues_HGRAD_Args(outputValues,
197  inputPoints,
198  operatorType,
199  this->getBaseCellTopology(),
200  this->getCardinality() );
201 #endif
202  Impl::Basis_HGRAD_TRI_C2_FEM::
203  getValues<DeviceType>(space,
204  outputValues,
205  inputPoints,
206  operatorType);
207  }
208 
209  virtual
210  void
211  getDofCoords( ScalarViewType dofCoords ) const override {
212 #ifdef HAVE_INTREPID2_DEBUG
213  // Verify rank of output array.
214  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
215  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
216  // Verify 0th dimension of output array.
217  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
218  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
219  // Verify 1st dimension of output array.
220  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
221  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
222 #endif
223  Kokkos::deep_copy(dofCoords, this->dofCoords_);
224  }
225 
226  virtual
227  void
228  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
229 #ifdef HAVE_INTREPID2_DEBUG
230  // Verify rank of output array.
231  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
232  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
233  // Verify 0th dimension of output array.
234  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
235  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
236 #endif
237  Kokkos::deep_copy(dofCoeffs, 1.0);
238  }
239 
240  virtual
241  const char*
242  getName() const override {
243  return "Intrepid2_HGRAD_TRI_C2_FEM";
244  }
245 
246 
255  BasisPtr<DeviceType,outputValueType,pointValueType>
256  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
257  if(subCellDim==1)
259 
260  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
261  }
262 
263  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
264  getHostBasis() const override{
266  }
267  };
268 
269 }// namespace Intrepid2
270 
272 
273 #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< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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_TRI_C2_FEM.
Header file for the Intrepid2::Basis_HGRAD_LINE_C2_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle 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...
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Definition file for FEM basis functions of degree 2 for H(grad) functions on TRI cells.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
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.
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::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.
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.