Intrepid2
Intrepid2_HGRAD_WEDGE_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
48
49#ifndef __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__
50#define __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
55
56namespace Intrepid2 {
57
129
130 namespace Impl {
131
135 template<bool serendipity>
137 public:
138 typedef struct Wedge<serendipity ? 15 : 18> cell_topology_type;
142 template<EOperator opType>
143 struct Serial {
144 template<typename OutputViewType,
145 typename inputViewType>
146 KOKKOS_INLINE_FUNCTION
147 static void
148 getValues( OutputViewType output,
149 const inputViewType input );
150
151 };
152
153 template<typename DeviceType,
154 typename outputValueValueType, class ...outputValueProperties,
155 typename inputPointValueType, class ...inputPointProperties>
156 static void
157 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
158 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
159 const EOperator operatorType);
160
164 template<typename outputValueViewType,
165 typename inputPointViewType,
166 EOperator opType>
167 struct Functor {
168 outputValueViewType _outputValues;
169 const inputPointViewType _inputPoints;
170
171 KOKKOS_INLINE_FUNCTION
172 Functor( outputValueViewType outputValues_,
173 inputPointViewType inputPoints_ )
174 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
175 KOKKOS_INLINE_FUNCTION
176 void operator()(const ordinal_type pt) const {
177 switch (opType) {
178 case OPERATOR_VALUE : {
179 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
180 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
181 Serial<opType>::getValues( output, input );
182 break;
183 }
184 case OPERATOR_GRAD :
185 case OPERATOR_D2 :
186 case OPERATOR_D3 :
187 case OPERATOR_D4 :
188 case OPERATOR_MAX : {
189 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
190 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
191 Serial<opType>::getValues( output, input );
192 break;
193 }
194 default: {
195 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
196 opType != OPERATOR_GRAD &&
197 opType != OPERATOR_D2 &&
198 opType != OPERATOR_MAX,
199 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::Serial::getValues) operator is not supported");
200 }
201 }
202 }
203 };
204 };
205 }
206
207 template<bool serendipity,
208 typename DeviceType = void,
209 typename outputValueType = double,
210 typename pointValueType = double>
211 class Basis_HGRAD_WEDGE_DEG2_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
212 public:
219
223
224 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
225
226 virtual
227 void
228 getValues( OutputViewType outputValues,
229 const PointViewType inputPoints,
230 const EOperator operatorType = OPERATOR_VALUE ) const override {
231#ifdef HAVE_INTREPID2_DEBUG
232 // Verify arguments
234 inputPoints,
235 operatorType,
236 this->getBaseCellTopology(),
237 this->getCardinality() );
238#endif
239 if constexpr (serendipity)
241 getValues<DeviceType>( outputValues,
242 inputPoints,
243 operatorType );
244 else
246 getValues<DeviceType>( outputValues,
247 inputPoints,
248 operatorType );;
249 }
250
251 virtual
252 void
253 getDofCoords( ScalarViewType dofCoords ) const override {
254#ifdef HAVE_INTREPID2_DEBUG
255 // Verify rank of output array.
256 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
257 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getDofCoords) rank = 2 required for dofCoords array");
258 // Verify 0th dimension of output array.
259 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
260 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
261 // Verify 1st dimension of output array.
262 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
263 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
264#endif
265 Kokkos::deep_copy(dofCoords, this->dofCoords_);
266 }
267
268 virtual
269 void
270 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
271#ifdef HAVE_INTREPID2_DEBUG
272 // Verify rank of output array.
273 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
274 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
275 // Verify 0th dimension of output array.
276 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
277 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
278#endif
279 Kokkos::deep_copy(dofCoeffs, 1.0);
280 }
281
282 virtual
283 const char*
284 getName() const override {
285 if constexpr (serendipity)
286 return "Intrepid2_HGRAD_WEDGE_I2_Serendipity_FEM";
287 else
288 return "Intrepid2_HGRAD_WEDGE_C2_FEM";
289 }
290
300 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
301 if(subCellDim == 1) {
302 return Teuchos::rcp(new
304 } else if(subCellDim == 2) {
305 if(subCellOrd < 3) //lateral faces
307 else
309 }
310 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
311 }
312
317 };
318
319 template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
320 using Basis_HGRAD_WEDGE_C2_FEM = Basis_HGRAD_WEDGE_DEG2_FEM<false, DeviceType, outputValueType, pointValueType>;
321
322 template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
323 using Basis_HGRAD_WEDGE_I2_Serendipity_FEM = Basis_HGRAD_WEDGE_DEG2_FEM<true, DeviceType, outputValueType, pointValueType>;
324
325}// namespace Intrepid2
326
328
329#endif
Header file for the abstract base class Intrepid2::Basis.
void getValues_HGRAD_Args(const outputValueViewType outputValues, const inputPointViewType inputPoints, const EOperator operatorType, const shards::CellTopology cellTopo, const ordinal_type basisCard)
Runtime check of the arguments for the getValues method in an HGRAD-conforming FEM basis....
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TRI_C2_FEM class.
Definition file for FEM basis functions of degree 2 for H(grad) functions on WEDGE cells.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
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 void getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
shards::CellTopology getBaseCellTopology() const
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM.