Intrepid2
Intrepid2_HGRAD_WEDGE_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_WEDGE_C1_FEM_HPP__
50 #define __INTREPID2_HGRAD_WEDGE_C1_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
55 
56 namespace Intrepid2 {
57 
88  namespace Impl {
89 
94  public:
95  typedef struct Wedge<6> cell_topology_type;
99  template<EOperator opType>
100  struct Serial {
101  template<typename OutputViewType,
102  typename inputViewType>
103  KOKKOS_INLINE_FUNCTION
104  static void
105  getValues( OutputViewType output,
106  const inputViewType input );
107 
108  };
109 
110  template<typename DeviceType,
111  typename outputValueValueType, class ...outputValueProperties,
112  typename inputPointValueType, class ...inputPointProperties>
113  static void
114  getValues( 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_D2 :
144  case OPERATOR_MAX : {
145  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
146  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
147  Serial<opType>::getValues( output, input );
148  break;
149  }
150  default: {
151  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
152  opType != OPERATOR_GRAD &&
153  opType != OPERATOR_D2 &&
154  opType != OPERATOR_MAX,
155  ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::Serial::getValues) operator is not supported");
156  }
157  }
158  }
159  };
160 
161 
162  };
163  }
164  template<typename DeviceType = void,
165  typename outputValueType = double,
166  typename pointValueType = double>
167  class Basis_HGRAD_WEDGE_C1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
168  public:
169  using OrdinalTypeArray1DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost;
170  using OrdinalTypeArray2DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost;
171  using OrdinalTypeArray3DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost;
172 
176 
180 
182 
183  virtual
184  void
185  getValues( 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_WEDGE_C1_FEM::
197  getValues<DeviceType>( outputValues,
198  inputPoints,
199  operatorType );
200  }
201 
202  virtual
203  void
204  getDofCoords( ScalarViewType dofCoords ) const override {
205 #ifdef HAVE_INTREPID2_DEBUG
206  // Verify rank of output array.
207  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
208  ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
209  // Verify 0th dimension of output array.
210  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
211  ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
212  // Verify 1st dimension of output array.
213  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
214  ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
215 #endif
216  Kokkos::deep_copy(dofCoords, this->dofCoords_);
217  }
218 
219  virtual
220  void
221  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
222 #ifdef HAVE_INTREPID2_DEBUG
223  // Verify rank of output array.
224  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
225  ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
226  // Verify 0th dimension of output array.
227  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
228  ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
229 #endif
230  Kokkos::deep_copy(dofCoeffs, 1.0);
231  }
232 
233  virtual
234  const char*
235  getName() const override {
236  return "Intrepid2_HGRAD_WEDGE_C1_FEM";
237  }
238 
247  BasisPtr<DeviceType,outputValueType,pointValueType>
248  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
249  if(subCellDim == 1) {
250  return Teuchos::rcp(new
252  } else if(subCellDim == 2) {
253  if(subCellOrd < 3) //lateral faces
255  else
257  }
258  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
259  }
260 
261  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
262  getHostBasis() const override{
264  }
265  };
266 }// namespace Intrepid2
267 
269 
270 #endif
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral cell...
Header file for the Intrepid2::Basis_HGRAD_TRI_C1_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge cell.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C1_FEM class.
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 (...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
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...
virtual const char * getName() const override
Returns basis name.
Definition file for FEM basis functions of degree 1 for H(grad) functions on WEDGE cells...
See Intrepid2::Basis_HGRAD_WEDGE_C1_FEM.
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.
Header file for the abstract base class Intrepid2::Basis.