Intrepid2
Intrepid2_HCURL_TET_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
49#ifndef __INTREPID2_HCURL_TET_I1_FEM_HPP__
50#define __INTREPID2_HCURL_TET_I1_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
54
55namespace Intrepid2 {
56
104 namespace Impl {
105
110 public:
111 typedef struct Tetrahedron<4> cell_topology_type;
112
116 template<EOperator opType>
117 struct Serial {
118 template<typename OutputViewType,
119 typename inputViewType>
120 KOKKOS_INLINE_FUNCTION
121 static void
122 getValues( OutputViewType output,
123 const inputViewType input );
124
125 };
126
127 template<typename DeviceType,
128 typename outputValueValueType, class ...outputValueProperties,
129 typename inputPointValueType, class ...inputPointProperties>
130 static void
131 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
132 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
133 const EOperator operatorType);
134
138 template<typename outputValueViewType,
139 typename inputPointViewType,
140 EOperator opType>
141 struct Functor {
142 outputValueViewType _outputValues;
143 const inputPointViewType _inputPoints;
144
145 KOKKOS_INLINE_FUNCTION
146 Functor( outputValueViewType outputValues_,
147 inputPointViewType inputPoints_ )
148 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
149
150 KOKKOS_INLINE_FUNCTION
151 void operator()(const ordinal_type pt) const {
152 switch (opType) {
153 case OPERATOR_VALUE : {
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 case OPERATOR_CURL : {
160 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
161 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
162 Serial<opType>::getValues( output, input );
163 break;
164 }
165 default: {
166 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
167 opType != OPERATOR_CURL,
168 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_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_HCURL_TET_I1_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_HCURL_TET_I1_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_HCURL_TET_I1_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->basisCardinality_, std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_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->basisCellTopology_.getDimension(), std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
227#endif
228 Kokkos::deep_copy(dofCoords, this->dofCoords_);
229 }
230
231
232 virtual
233 void
234 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
235#ifdef HAVE_INTREPID2_DEBUG
236 // Verify rank of output array.
237 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
239 // Verify 0th dimension of output array.
240 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
242 // Verify 1st dimension of output array.
243 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
244 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
245#endif
246 Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
247 }
248
249 virtual
250 const char*
251 getName() const override {
252 return "Intrepid2_HCURL_TET_I1_FEM";
253 }
254
255 virtual
256 bool
257 requireOrientation() const override {
258 return true;
259 }
260
271 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
272
273 if(subCellDim == 1)
274 return Teuchos::rcp( new
275 Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Line<2> >()));
276 else if(subCellDim == 2) {
277 return Teuchos::rcp(new
279 }
280
281 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
282 }
283
285 getHostBasis() const override{
287 }
288
289
290 };
291}// namespace Intrepid2
292
294
295#endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HCURL_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 HCURL-conforming FEM basis....
Definition file for FEM basis functions of degree 1 for H(curl) functions on TET cells.
Header file for the Intrepid2::Basis_HCURL_TRI_I1_FEM class.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Tetrahedron cell.
virtual bool requireOrientation() const override
True if orientation is required.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference 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...
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 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 getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Triangle 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_HCURL_TET_I1_FEM.