Intrepid2
Intrepid2_HGRAD_QUAD_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
49#ifndef __INTREPID2_HGRAD_QUAD_C2_FEM_HPP__
50#define __INTREPID2_HGRAD_QUAD_C2_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
53
54namespace Intrepid2 {
55
93 namespace Impl {
94
99 public:
100 typedef struct Quadrilateral<9> cell_topology_type;
104 template<EOperator opType>
105 struct Serial {
106 template<typename OutputViewType,
107 typename inputViewType>
109 static void
110 getValues( OutputViewType output,
111 const inputViewType input );
112
113 };
114
115 template<typename DeviceType,
116 typename outputValueValueType, class ...outputValueProperties,
117 typename inputPointValueType, class ...inputPointProperties>
118 static void
119 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
120 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
121 const EOperator operatorType);
122
126 template<typename outputValueViewType,
127 typename inputPointViewType,
128 EOperator opType>
129 struct Functor {
130 outputValueViewType _outputValues;
131 const inputPointViewType _inputPoints;
132
135 inputPointViewType inputPoints_ )
136 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
137
139 void operator()(const ordinal_type pt) const {
140 switch (opType) {
141 case OPERATOR_VALUE : {
142 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
143 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
144 Serial<opType>::getValues( output, input );
145 break;
146 }
147 case OPERATOR_GRAD :
148 case OPERATOR_CURL :
149 case OPERATOR_D1 :
150 case OPERATOR_D2 :
151 case OPERATOR_D3 :
152 case OPERATOR_D4 :
153 case OPERATOR_MAX : {
154 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
155 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
156 Serial<opType>::getValues( output, input );
157 break;
158 }
159 default: {
160 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
161 opType != OPERATOR_GRAD &&
162 opType != OPERATOR_CURL &&
163 opType != OPERATOR_D1 &&
164 opType != OPERATOR_D2 &&
165 opType != OPERATOR_D3 &&
166 opType != OPERATOR_D4 &&
167 opType != OPERATOR_MAX,
168 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_I1_FEM::Serial::getValues) operator is not supported");
169 }
170 }
171 }
172 };
173 };
174 }
175
176 template<typename DeviceType = void,
177 typename outputValueType = double,
178 typename pointValueType = double>
179 class Basis_HGRAD_QUAD_C2_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
180 public:
184
188
192
193 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
194
195 virtual
196 void
197 getValues( OutputViewType outputValues,
198 const PointViewType inputPoints,
199 const EOperator operatorType = OPERATOR_VALUE ) const override {
200#ifdef HAVE_INTREPID2_DEBUG
201 // Verify arguments
203 inputPoints,
204 operatorType,
205 this->getBaseCellTopology(),
206 this->getCardinality() );
207#endif
208 Impl::Basis_HGRAD_QUAD_C2_FEM::
209 getValues<DeviceType>( outputValues,
210 inputPoints,
211 operatorType );
212 }
213
214 virtual
215 void
216 getDofCoords( ScalarViewType dofCoords ) const override {
217#ifdef HAVE_INTREPID2_DEBUG
218 // Verify rank of output array.
219 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
220 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
221 // Verify 0th dimension of output array.
222 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
224 // Verify 1st dimension of output array.
225 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
227#endif
228 Kokkos::deep_copy(dofCoords, this->dofCoords_);
229 }
230
231 virtual
232 void
233 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
234#ifdef HAVE_INTREPID2_DEBUG
235 // Verify rank of output array.
236 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
237 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
238 // Verify 0th dimension of output array.
239 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
240 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
241#endif
242 Kokkos::deep_copy(dofCoeffs, 1.0);
243 }
244
245 virtual
246 const char*
247 getName() const override {
248 return "Intrepid2_HGRAD_QUAD_C2_FEM";
249 }
250
255 };
256}// namespace Intrepid2
257
259
260#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 2 for H(grad) functions on QUAD cells.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell.
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 const char * getName() const override
Returns basis name.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference 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_QUAD_C2_FEM.