Intrepid2
Intrepid2_HGRAD_HEX_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_HEX_C1_FEM_HPP__
50#define __INTREPID2_HGRAD_HEX_C1_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
54
55namespace Intrepid2 {
56
92 namespace Impl {
93
98 public:
99 typedef struct Hexahedron<8> cell_topology_type;
103 template<EOperator opType>
104 struct Serial {
105 template<typename OutputViewType,
106 typename inputViewType>
107 KOKKOS_INLINE_FUNCTION
108 static void
109 getValues( OutputViewType output,
110 const inputViewType input );
111
112 };
113
114 template<typename DeviceType,
115 typename outputValueValueType, class ...outputValueProperties,
116 typename inputPointValueType, class ...inputPointProperties>
117 static void
118 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
119 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
120 const EOperator operatorType);
121
125 template<typename outputValueViewType,
126 typename inputPointViewType,
127 EOperator opType>
128 struct Functor {
129 outputValueViewType _outputValues;
130 const inputPointViewType _inputPoints;
131
132 KOKKOS_INLINE_FUNCTION
133 Functor( outputValueViewType outputValues_,
134 inputPointViewType inputPoints_ )
135 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
136
137 KOKKOS_INLINE_FUNCTION
138 void operator()(const ordinal_type pt) const {
139 switch (opType) {
140 case OPERATOR_VALUE : {
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 case OPERATOR_GRAD :
147 case OPERATOR_D2 :
148 case OPERATOR_D3 :
149 case OPERATOR_MAX : {
150 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
151 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
152 Serial<opType>::getValues( output, input );
153 break;
154 }
155 default: {
156 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
157 opType != OPERATOR_GRAD &&
158 opType != OPERATOR_D2 &&
159 opType != OPERATOR_D3 &&
160 opType != OPERATOR_MAX,
161 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_C1_FEM::Serial::getValues) operator is not supported");
162 }
163 }
164 }
165 };
166 };
167 }
168
169 template<typename DeviceType = void,
170 typename outputValueType = double,
171 typename pointValueType = double>
172 class Basis_HGRAD_HEX_C1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
173 public:
177
181
185
186 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
187
188 virtual
189 void
190 getValues( OutputViewType outputValues,
191 const PointViewType inputPoints,
192 const EOperator operatorType = OPERATOR_VALUE ) const override {
193#ifdef HAVE_INTREPID2_DEBUG
194 // Verify arguments
196 inputPoints,
197 operatorType,
198 this->getBaseCellTopology(),
199 this->getCardinality() );
200#endif
201 Impl::Basis_HGRAD_HEX_C1_FEM::
202 getValues<DeviceType>( outputValues,
203 inputPoints,
204 operatorType );
205 }
206
207 virtual
208 void
209 getDofCoords( ScalarViewType dofCoords ) const override {
210#ifdef HAVE_INTREPID2_DEBUG
211 // Verify rank of output array.
212 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
213 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
214 // Verify 0th dimension of output array.
215 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
216 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
217 // Verify 1st dimension of output array.
218 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
219 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
220#endif
221 Kokkos::deep_copy(dofCoords, this->dofCoords_);
222 }
223
224 virtual
225 void
226 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
227#ifdef HAVE_INTREPID2_DEBUG
228 // Verify rank of output array.
229 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
230 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
231 // Verify 0th dimension of output array.
232 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
233 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
234#endif
235 Kokkos::deep_copy(dofCoeffs, 1.0);
236 }
237
238 virtual
239 const char*
240 getName() const override {
241 return "Intrepid2_HGRAD_HEX_C1_FEM";
242 }
243
253 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
254 if(subCellDim == 1) {
255 return Teuchos::rcp(new
257 } else if(subCellDim == 2) {
258 return Teuchos::rcp(new
260 }
261 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
262 }
263
265 getHostBasis() const override{
267 }
268 };
269}// namespace Intrepid2
270
272
273#endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
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....
Definition file for FEM basis functions of degree 1 for H(grad) functions on HEX cells.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C1_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Hexahedron 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 const char * getName() const override
Returns basis name.
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 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...
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral cell.
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_HEX_C1_FEM.