Intrepid2
Intrepid2_HDIV_TRI_I1_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_HDIV_TRI_I1_FEM_HPP__
50 #define __INTREPID2_HDIV_TRI_I1_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
53 
54 namespace Intrepid2 {
55 
101  namespace Impl {
102 
107  public:
108  typedef struct Triangle<3> cell_topology_type;
112  template<EOperator opType>
113  struct Serial {
114  template<typename OutputViewType,
115  typename inputViewType>
116  KOKKOS_INLINE_FUNCTION
117  static void
118  getValues( OutputViewType output,
119  const inputViewType input );
120 
121  };
122 
123  template<typename ExecSpaceType,
124  typename outputValueValueType, class ...outputValueProperties,
125  typename inputPointValueType, class ...inputPointProperties>
126  static void
127  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
128  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
129  const EOperator operatorType);
130 
134  template<typename outputValueViewType,
135  typename inputPointViewType,
136  EOperator opType>
137  struct Functor {
138  outputValueViewType _outputValues;
139  const inputPointViewType _inputPoints;
140 
141  KOKKOS_INLINE_FUNCTION
142  Functor( outputValueViewType outputValues_,
143  inputPointViewType inputPoints_ )
144  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
145 
146  KOKKOS_INLINE_FUNCTION
147  void operator()(const ordinal_type pt) const {
148  switch (opType) {
149  case OPERATOR_VALUE : {
150  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
151  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
152  Serial<opType>::getValues( output, input );
153  break;
154  }
155  case OPERATOR_DIV : {
156  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
157  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
158  Serial<opType>::getValues( output, input );
159  break;
160  }
161  default: {
162  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
163  opType != OPERATOR_DIV,
164  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::Serial::getValues) operator is not supported");
165  }
166  }
167  }
168  };
169  };
170  }
171 
172  template<typename ExecSpaceType = void,
173  typename outputValueType = double,
174  typename pointValueType = double>
175  class Basis_HDIV_TRI_I1_FEM: public Basis<ExecSpaceType,outputValueType,pointValueType> {
176  public:
180 
181  using ordinal_type_array_1d_host INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use OrdinalTypeArray1DHost instead","OrdinalTypeArray1DHost") = OrdinalTypeArray1DHost INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use OrdinalTypeArray1DHost instead");
182  using ordinal_type_array_2d_host INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use OrdinalTypeArray2DHost instead","OrdinalTypeArray2DHost") = OrdinalTypeArray2DHost INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use OrdinalTypeArray2DHost instead");
183  using ordinal_type_array_3d_host INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use OrdinalTypeArray3DHost instead","OrdinalTypeArray3DHost") = OrdinalTypeArray3DHost INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use OrdinalTypeArray3DHost instead");
184 
188 
192 
193  using outputViewType INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use OutputViewType instead","OutputViewType") = OutputViewType INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use OutputViewType instead");
194  using pointViewType INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use PointViewType instead","PointViewType") = PointViewType INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use PointViewType instead");
195  using scalarViewType INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use ScalarViewType instead","ScalarViewType") = ScalarViewType INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use ScalarViewType instead");
196 
198 
199  virtual
200  void
201  getValues( OutputViewType outputValues,
202  const PointViewType inputPoints,
203  const EOperator operatorType = OPERATOR_VALUE ) const {
204 #ifdef HAVE_INTREPID2_DEBUG
205  // Verify arguments
206  Intrepid2::getValues_HDIV_Args(outputValues,
207  inputPoints,
208  operatorType,
209  this->getBaseCellTopology(),
210  this->getCardinality() );
211 #endif
212  Impl::Basis_HDIV_TRI_I1_FEM::
213  getValues<ExecSpaceType>( outputValues,
214  inputPoints,
215  operatorType );
216  }
217 
218  virtual
219  void
220  getDofCoords( ScalarViewType dofCoords ) const {
221 #ifdef HAVE_INTREPID2_DEBUG
222  // Verify rank of output array.
223  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
224  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
225  // Verify 0th dimension of output array.
226  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
227  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
228  // Verify 1st dimension of output array.
229  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
230  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
231 #endif
232  Kokkos::deep_copy(dofCoords, this->dofCoords_);
233  }
234 
235  virtual
236  void
237  getDofCoeffs( ScalarViewType dofCoeffs ) const {
238 #ifdef HAVE_INTREPID2_DEBUG
239  // Verify rank of output array.
240  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
241  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
242  // Verify 0th dimension of output array.
243  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
244  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
245  // Verify 1st dimension of output array.
246  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
247  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
248 #endif
249  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
250  }
251 
252  virtual
253  const char*
254  getName() const {
255  return "Intrepid2_HDIV_TRI_I1_FEM";
256  }
257 
258  virtual
259  bool
261  return true;
262  }
263 
264  };
265 
266 }// namespace Intrepid2
267 
269 
270 #endif
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.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
See Intrepid2::Basis_HDIV_TRI_I1_FEM.
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 getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Implementation of the default H(div)-compatible FEM basis of degree 1 on a Triangle cell...
virtual bool requireOrientation() const
True if orientation is required.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
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.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Definition file for FEM basis functions of degree 1 for H(div) functions on TRI cells.
Header file for the abstract base class Intrepid2::Basis.
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...