Intrepid2
Intrepid2_HDIV_WEDGE_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
48#ifndef __INTREPID2_HDIV_WEDGE_I1_FEM_HPP__
49#define __INTREPID2_HDIV_WEDGE_I1_FEM_HPP__
50
51#include "Intrepid2_Basis.hpp"
53
54namespace Intrepid2 {
55
86 namespace Impl {
87
92 public:
93 typedef struct Wedge<6> cell_topology_type;
97 template<EOperator opType>
98 struct Serial {
99 template<typename OutputViewType,
100 typename inputViewType>
101 KOKKOS_INLINE_FUNCTION
102 static void
103 getValues( OutputViewType output,
104 const inputViewType input );
105
106 };
107
108 template<typename DeviceType,
109 typename outputValueValueType, class ...outputValueProperties,
110 typename inputPointValueType, class ...inputPointProperties>
111 static void
112 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
113 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
114 const EOperator operatorType);
115
119 template<typename outputValueViewType,
120 typename inputPointViewType,
121 EOperator opType>
122 struct Functor {
123 outputValueViewType _outputValues;
124 const inputPointViewType _inputPoints;
125
126 KOKKOS_INLINE_FUNCTION
127 Functor( outputValueViewType outputValues_,
128 inputPointViewType inputPoints_ )
129 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
130
131 KOKKOS_INLINE_FUNCTION
132 void operator()(const ordinal_type pt) const {
133 switch (opType) {
134 case OPERATOR_VALUE : {
135 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
136 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
137 Serial<opType>::getValues( output, input );
138 break;
139 }
140 case OPERATOR_DIV : {
141 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
142 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
143 Serial<opType>::getValues( output, input );
144 break;
145 }
146 default: {
147 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
148 opType != OPERATOR_DIV ,
149 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::Serial::getValues) operator is not supported");
150 }
151 }
152 }
153 };
154 };
155 }
156
157 template<typename DeviceType = void,
158 typename outputValueType = double,
159 typename pointValueType = double>
160 class Basis_HDIV_WEDGE_I1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
161 public:
165
169
173
174 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
175
176 virtual
177 void
178 getValues( OutputViewType outputValues,
179 const PointViewType inputPoints,
180 const EOperator operatorType = OPERATOR_VALUE ) const override {
181#ifdef HAVE_INTREPID2_DEBUG
182 // Verify arguments
184 inputPoints,
185 operatorType,
186 this->getBaseCellTopology(),
187 this->getCardinality() );
188#endif
189 Impl::Basis_HDIV_WEDGE_I1_FEM::
190 getValues<DeviceType>( outputValues,
191 inputPoints,
192 operatorType );
193 }
194
195 virtual
196 void
197 getDofCoords( ScalarViewType dofCoords ) const override {
198#ifdef HAVE_INTREPID2_DEBUG
199 // Verify rank of output array.
200 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
201 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
202 // Verify 0th dimension of output array.
203 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
204 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
205 // Verify 1st dimension of output array.
206 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
207 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
208#endif
209 Kokkos::deep_copy(dofCoords, this->dofCoords_);
210 }
211
212
213 virtual
214 void
215 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
216 #ifdef HAVE_INTREPID2_DEBUG
217 // Verify rank of output array.
218 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
219 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
220 // Verify 0th dimension of output array.
221 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
222 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
223 // Verify 1st dimension of output array.
224 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
225 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
226 #endif
227 Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
228 }
229
230 virtual
231 const char*
232 getName() const override {
233 return "Intrepid2_HDIV_WEDGE_I1_FEM";
234 }
235
236 virtual
237 bool
238 requireOrientation() const override {
239 return true;
240 }
241
252 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
253
254 if(subCellDim == 2) {
255 if((subCellOrd >=0) && (subCellOrd<3))
256 return Teuchos::rcp(new
257 Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Quadrilateral<4> >()));
258 if((subCellOrd==3)||(subCellOrd==4))
259 return Teuchos::rcp(new
260 Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Triangle<3> >()));
261 }
262 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
263 }
264
265
267 getHostBasis() const override{
269 }
270 };
271}// namespace Intrepid2
272
274
275#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.
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 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...
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 getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
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.