Intrepid2
Intrepid2_HDIV_WEDGE_I1_FEM.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Intrepid2 Package
4//
5// Copyright 2007 NTESS and the Intrepid2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
15#ifndef __INTREPID2_HDIV_WEDGE_I1_FEM_HPP__
16#define __INTREPID2_HDIV_WEDGE_I1_FEM_HPP__
17
18#include "Intrepid2_Basis.hpp"
20
21namespace Intrepid2 {
22
53 namespace Impl {
54
59 public:
60 typedef struct Wedge<6> cell_topology_type;
64 template<EOperator opType>
65 struct Serial {
66 template<typename OutputViewType,
67 typename inputViewType>
68 KOKKOS_INLINE_FUNCTION
69 static void
70 getValues( OutputViewType output,
71 const inputViewType input );
72
73 };
74
75 template<typename DeviceType,
76 typename outputValueValueType, class ...outputValueProperties,
77 typename inputPointValueType, class ...inputPointProperties>
78 static void
79 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
80 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
81 const EOperator operatorType);
82
86 template<typename outputValueViewType,
87 typename inputPointViewType,
88 EOperator opType>
89 struct Functor {
90 outputValueViewType _outputValues;
91 const inputPointViewType _inputPoints;
92
93 KOKKOS_INLINE_FUNCTION
94 Functor( outputValueViewType outputValues_,
95 inputPointViewType inputPoints_ )
96 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
97
98 KOKKOS_INLINE_FUNCTION
99 void operator()(const ordinal_type pt) const {
100 switch (opType) {
101 case OPERATOR_VALUE : {
102 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
103 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
104 Serial<opType>::getValues( output, input );
105 break;
106 }
107 case OPERATOR_DIV : {
108 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
109 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
110 Serial<opType>::getValues( output, input );
111 break;
112 }
113 default: {
114 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
115 opType != OPERATOR_DIV ,
116 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::Serial::getValues) operator is not supported");
117 }
118 }
119 }
120 };
121 };
122 }
123
124 template<typename DeviceType = void,
125 typename outputValueType = double,
126 typename pointValueType = double>
127 class Basis_HDIV_WEDGE_I1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
128 public:
132
136
140
141 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
142
143 virtual
144 void
145 getValues( OutputViewType outputValues,
146 const PointViewType inputPoints,
147 const EOperator operatorType = OPERATOR_VALUE ) const override {
148#ifdef HAVE_INTREPID2_DEBUG
149 // Verify arguments
151 inputPoints,
152 operatorType,
153 this->getBaseCellTopology(),
154 this->getCardinality() );
155#endif
156 Impl::Basis_HDIV_WEDGE_I1_FEM::
157 getValues<DeviceType>( outputValues,
158 inputPoints,
159 operatorType );
160 }
161
162 virtual void
163 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
164 ordinal_type& perThreadSpaceSize,
165 const PointViewType inputPointsconst,
166 const EOperator operatorType = OPERATOR_VALUE) const override;
167
168 KOKKOS_INLINE_FUNCTION
169 virtual void
170 getValues(
171 OutputViewType outputValues,
172 const PointViewType inputPoints,
173 const EOperator operatorType,
174 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
175 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
176 const ordinal_type subcellDim = -1,
177 const ordinal_type subcellOrdinal = -1) const override;
178
179 virtual
180 void
181 getDofCoords( ScalarViewType dofCoords ) const override {
182#ifdef HAVE_INTREPID2_DEBUG
183 // Verify rank of output array.
184 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
185 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
186 // Verify 0th dimension of output array.
187 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
188 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
189 // Verify 1st dimension of output array.
190 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
191 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
192#endif
193 Kokkos::deep_copy(dofCoords, this->dofCoords_);
194 }
195
196
197 virtual
198 void
199 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
200 #ifdef HAVE_INTREPID2_DEBUG
201 // Verify rank of output array.
202 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
203 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
204 // Verify 0th dimension of output array.
205 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
206 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
207 // Verify 1st dimension of output array.
208 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
209 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
210 #endif
211 Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
212 }
213
214 virtual
215 const char*
216 getName() const override {
217 return "Intrepid2_HDIV_WEDGE_I1_FEM";
218 }
219
220 virtual
221 bool
222 requireOrientation() const override {
223 return true;
224 }
225
236 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
237
238 if(subCellDim == 2) {
239 if((subCellOrd >=0) && (subCellOrd<3))
240 return Teuchos::rcp(new
241 Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Quadrilateral<4> >()));
242 if((subCellOrd==3)||(subCellOrd==4))
243 return Teuchos::rcp(new
244 Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Triangle<3> >()));
245 }
246 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
247 }
248
249
254 };
255}// namespace Intrepid2
256
258
259#endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HDIV_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 HDIV-conforming FEM basis....
Definition file for FEM basis functions of degree 1 for H(div) functions on WEDGE cells.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge cell.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual bool requireOrientation() const override
True if orientation is required.
virtual const char * getName() const override
Returns basis name.
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 HVOL-compatible FEM contstant basis on triangle, quadrilateral,...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
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.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HDIV_WEDGE_I1_FEM.