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"
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>
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
128 inputPointViewType inputPoints_ )
129 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
130
132 void operator()(const ordinal_type pt) const {
133 switch (opType) {
134 case OPERATOR_VALUE : {
135 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
136 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
137 Serial<opType>::getValues( output, input );
138 break;
139 }
140 case OPERATOR_GRAD :
141 case OPERATOR_D2 :
142 case OPERATOR_MAX : {
143 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
144 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
145 Serial<opType>::getValues( output, input );
146 break;
147 }
148 default: {
149 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
150 opType != OPERATOR_GRAD &&
151 opType != OPERATOR_D2 &&
152 opType != OPERATOR_MAX,
153 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::Serial::getValues) operator is not supported");
154 }
155 }
156 }
157 };
158
159
160 };
161 }
162 template<typename DeviceType = void,
163 typename outputValueType = double,
164 typename pointValueType = double>
165 class Basis_HGRAD_WEDGE_C1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
166 public:
170
174
178
179 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
180
181 virtual
182 void
183 getValues( OutputViewType outputValues,
184 const PointViewType inputPoints,
185 const EOperator operatorType = OPERATOR_VALUE ) const override {
186#ifdef HAVE_INTREPID2_DEBUG
187 // Verify arguments
189 inputPoints,
190 operatorType,
191 this->getBaseCellTopology(),
192 this->getCardinality() );
193#endif
194 Impl::Basis_HGRAD_WEDGE_C1_FEM::
195 getValues<DeviceType>( outputValues,
196 inputPoints,
197 operatorType );
198 }
199
200 virtual
201 void
202 getDofCoords( ScalarViewType dofCoords ) const override {
203#ifdef HAVE_INTREPID2_DEBUG
204 // Verify rank of output array.
205 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
206 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
207 // Verify 0th dimension of output array.
208 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
209 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
210 // Verify 1st dimension of output array.
211 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
212 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
213#endif
214 Kokkos::deep_copy(dofCoords, this->dofCoords_);
215 }
216
217 virtual
218 void
219 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
220#ifdef HAVE_INTREPID2_DEBUG
221 // Verify rank of output array.
222 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
224 // Verify 0th dimension of output array.
225 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
227#endif
228 Kokkos::deep_copy(dofCoeffs, 1.0);
229 }
230
231 virtual
232 const char*
233 getName() const override {
234 return "Intrepid2_HGRAD_WEDGE_C1_FEM";
235 }
236
241 };
242}// namespace Intrepid2
243
245
246#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.
Definition file for FEM basis functions of degree 1 for H(grad) functions on WEDGE cells.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge cell.
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...
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 getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual const char * getName() const override
Returns basis name.
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::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_HGRAD_WEDGE_C1_FEM.