Intrepid2
Intrepid2_HGRAD_TRI_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_TRI_C1_FEM_HPP__
50 #define __INTREPID2_HGRAD_TRI_C1_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
54 
55 namespace Intrepid2 {
56 
80  namespace Impl {
81 
86  public:
87  typedef struct Triangle<3> cell_topology_type;
91  template<EOperator opType>
92  struct Serial {
93  template<typename OutputViewType,
94  typename inputViewType>
95  KOKKOS_INLINE_FUNCTION
96  static void
97  getValues( OutputViewType output,
98  const inputViewType input );
99 
100  };
101 
102  template<typename DeviceType,
103  typename outputValueValueType, class ...outputValueProperties,
104  typename inputPointValueType, class ...inputPointProperties>
105  static void
106  getValues( const typename DeviceType::execution_space& space,
107  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
108  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
109  const EOperator operatorType);
110 
114  template<typename outputValueViewType,
115  typename inputPointViewType,
116  EOperator opType>
117  struct Functor {
118  outputValueViewType _outputValues;
119  const inputPointViewType _inputPoints;
120 
121  KOKKOS_INLINE_FUNCTION
122  Functor( outputValueViewType outputValues_,
123  inputPointViewType inputPoints_ )
124  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
125 
126  KOKKOS_INLINE_FUNCTION
127  void operator()(const ordinal_type pt) const {
128  switch (opType) {
129  case OPERATOR_VALUE : {
130  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
131  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
132  Serial<opType>::getValues( output, input );
133  break;
134  }
135  case OPERATOR_GRAD :
136  case OPERATOR_CURL :
137  case OPERATOR_MAX : {
138  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
139  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
140  Serial<opType>::getValues( output, input );
141  break;
142  }
143  default: {
144  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
145  opType != OPERATOR_GRAD &&
146  opType != OPERATOR_CURL &&
147  opType != OPERATOR_MAX,
148  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C1_FEM::Serial::getValues) operator is not supported");
149  }
150  }
151  }
152  };
153  };
154  }
155 
156  template<typename DeviceType = void,
157  typename outputValueType = double,
158  typename pointValueType = double>
159  class Basis_HGRAD_TRI_C1_FEM: public Basis<DeviceType,outputValueType,pointValueType> {
160  public:
162  using typename BasisBase::ExecutionSpace;
163  using typename BasisBase::OrdinalTypeArray1DHost;
164  using typename BasisBase::OrdinalTypeArray2DHost;
165  using typename BasisBase::OrdinalTypeArray3DHost;
166 
170 
171  using typename BasisBase::OutputViewType;
172  using typename BasisBase::PointViewType;
173  using typename BasisBase::ScalarViewType;
174 
175  using BasisBase::getValues;
176 
177  virtual
178  void
179  getValues( const ExecutionSpace& space,
180  OutputViewType outputValues,
181  const PointViewType inputPoints,
182  const EOperator operatorType = OPERATOR_VALUE ) const override {
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_TRI_C1_FEM::
192  getValues<DeviceType>(space,
193  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( rank(dofCoords) != 2, std::invalid_argument,
204  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_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_TRI_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_TRI_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( rank(dofCoeffs) != 1, std::invalid_argument,
221  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_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_TRI_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_TRI_C1_FEM";
233  }
234 
243  BasisPtr<DeviceType, outputValueType, pointValueType>
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 
252  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
253  getHostBasis() const override{
255  }
256  };
257 
258 }// namespace Intrepid2
259 
261 
262 #endif
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...
See Intrepid2::Basis_HGRAD_TRI_C1_FEM.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Triangle cell...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
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.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
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.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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 > OrdinalTypeArray2DHost
View type for 2d host array.
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::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Line cell.
virtual const char * getName() const override
Returns basis name.
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.
Definition file for FEM basis functions of degree 1 for H(grad) functions on TRI cells.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.