Intrepid2
Intrepid2_HGRAD_TRI_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
48
49#ifndef __INTREPID2_HGRAD_TRI_C2_FEM_HPP__
50#define __INTREPID2_HGRAD_TRI_C2_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
54
55namespace Intrepid2 {
56
86
87 namespace Impl {
88
93 public:
94 typedef struct Triangle<6> cell_topology_type;
98 template<EOperator opType>
99 struct Serial {
100 template<typename OutputViewType,
101 typename inputViewType>
102 KOKKOS_INLINE_FUNCTION
103 static void
104 getValues( OutputViewType output,
105 const inputViewType input );
106
107 };
108
109 template<typename DeviceType,
110 typename outputValueValueType, class ...outputValueProperties,
111 typename inputPointValueType, class ...inputPointProperties>
112 static void
113 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
114 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
115 const EOperator operatorType);
116
120 template<typename outputValueViewType,
121 typename inputPointViewType,
122 EOperator opType>
123 struct Functor {
124 outputValueViewType _outputValues;
125 const inputPointViewType _inputPoints;
126
127 KOKKOS_INLINE_FUNCTION
128 Functor( outputValueViewType outputValues_,
129 inputPointViewType inputPoints_ )
130 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
131
132 KOKKOS_INLINE_FUNCTION
133 void operator()(const ordinal_type pt) const {
134 switch (opType) {
135 case OPERATOR_VALUE : {
136 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
137 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
138 Serial<opType>::getValues( output, input );
139 break;
140 }
141 case OPERATOR_GRAD :
142 case OPERATOR_CURL :
143 case OPERATOR_D1 :
144 case OPERATOR_D2 :
145 case OPERATOR_MAX : {
146 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
147 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
148 Serial<opType>::getValues( output, input );
149 break;
150 }
151 default: {
152 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
153 opType != OPERATOR_GRAD &&
154 opType != OPERATOR_CURL &&
155 opType != OPERATOR_D1 &&
156 opType != OPERATOR_D2 &&
157 opType != OPERATOR_MAX,
158 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::Serial::getValues) operator is not supported");
159 }
160 }
161 }
162 };
163 };
164 }
165
166 template<typename DeviceType = void,
167 typename outputValueType = double,
168 typename pointValueType = double>
169 class Basis_HGRAD_TRI_C2_FEM: public Basis<DeviceType,outputValueType,pointValueType> {
170 public:
174
178
182
183 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
184
185 virtual
186 void
187 getValues( OutputViewType outputValues,
188 const PointViewType inputPoints,
189 const EOperator operatorType = OPERATOR_VALUE ) const override {
190#ifdef HAVE_INTREPID2_DEBUG
191 // Verify arguments
193 inputPoints,
194 operatorType,
195 this->getBaseCellTopology(),
196 this->getCardinality() );
197#endif
198 Impl::Basis_HGRAD_TRI_C2_FEM::
199 getValues<DeviceType>( outputValues,
200 inputPoints,
201 operatorType );
202 }
203
204 virtual
205 void
206 getDofCoords( ScalarViewType dofCoords ) const override {
207#ifdef HAVE_INTREPID2_DEBUG
208 // Verify rank of output array.
209 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
210 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
211 // Verify 0th dimension of output array.
212 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
213 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
214 // Verify 1st dimension of output array.
215 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
216 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
217#endif
218 Kokkos::deep_copy(dofCoords, this->dofCoords_);
219 }
220
221 virtual
222 void
223 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
224#ifdef HAVE_INTREPID2_DEBUG
225 // Verify rank of output array.
226 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
227 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
228 // Verify 0th dimension of output array.
229 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
230 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
231#endif
232 Kokkos::deep_copy(dofCoeffs, 1.0);
233 }
234
235 virtual
236 const char*
237 getName() const override {
238 return "Intrepid2_HGRAD_TRI_C2_FEM";
239 }
240
241
251 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
252 if(subCellDim==1)
254
255 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
256 }
257
262 };
263
264}// namespace Intrepid2
265
267
268#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.
Header file for the Intrepid2::Basis_HGRAD_LINE_C2_FEM class.
Definition file for FEM basis functions of degree 2 for H(grad) functions on TRI cells.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
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 getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the 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< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
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.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
shards::CellTopology getBaseCellTopology() const
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HGRAD_TRI_C2_FEM.