Intrepid2
Intrepid2_HCURL_TRI_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_TRI_I1_FEM_HPP
50#define INTREPID2_HCURL_TRI_I1_FEM_HPP
51
52#include "Intrepid2_Basis.hpp"
54
55namespace Intrepid2 {
56
97 namespace Impl {
98
103 public:
104 typedef struct Triangle<3> cell_topology_type;
108 template<EOperator opType>
109 struct Serial {
110 template<typename OutputViewType,
111 typename inputViewType>
113 static void
114 getValues( OutputViewType output,
115 const inputViewType input );
116
117 };
118
119 template<typename DeviceType,
120 typename outputValueValueType, class ...outputValueProperties,
121 typename inputPointValueType, class ...inputPointProperties>
122 static void
123 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
124 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
125 const EOperator operatorType);
126
130 template<typename outputValueViewType,
131 typename inputPointViewType,
132 EOperator opType>
133 struct Functor {
134 outputValueViewType _outputValues;
135 const inputPointViewType _inputPoints;
136
139 inputPointViewType inputPoints_ )
140 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
141
143 void operator()(const ordinal_type pt) const {
144 switch (opType) {
145 case OPERATOR_VALUE : {
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 case OPERATOR_CURL : {
152 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
153 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
154 Serial<opType>::getValues( output, input );
155 break;
156 }
157 case OPERATOR_DIV: {
158 INTREPID2_TEST_FOR_ABORT( (opType == OPERATOR_DIV),
159 ">>> ERROR (Basis_HCURL_TRI_I1_FEM): DIV is invalid operator for HCURL Basis Functions");
160 break;
161 }
162 case OPERATOR_GRAD: {
163 INTREPID2_TEST_FOR_ABORT( (opType == OPERATOR_GRAD),
164 ">>> ERROR (Basis_HCURL_TRI_I1_FEM): GRAD is invalid operator for HCURL Basis Functions");
165 break;
166 }
167 default: {
168 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
169 opType != OPERATOR_CURL,
170 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::Serial::getValues) operator is not supported");
171 }
172 }
173 }
174 };
175
176 };
177 }
178
179 template< typename DeviceType = void,
180 typename outputValueType = double,
181 typename pointValueType = double >
182 class Basis_HCURL_TRI_I1_FEM : public Basis<DeviceType, outputValueType, pointValueType> {
183 public:
187
191
195
196 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
197
198 virtual
199 void
200 getValues( OutputViewType outputValues,
201 const PointViewType inputPoints,
202 const EOperator operatorType = OPERATOR_VALUE ) const override {
203#ifdef HAVE_INTREPID2_DEBUG
204 // Verify arguments
206 inputPoints,
207 operatorType,
208 this->getBaseCellTopology(),
209 this->getCardinality() );
210#endif
211 Impl::Basis_HCURL_TRI_I1_FEM::
212 getValues<DeviceType>( outputValues,
213 inputPoints,
214 operatorType );
215 }
216
217 virtual
218 void
219 getDofCoords( ScalarViewType dofCoords ) const override {
220#ifdef HAVE_INTREPID2_DEBUG
221 // Verify rank of output array.
222 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
224 // Verify 0th dimension of output array.
225 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
227 // Verify 1st dimension of output array.
228 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
230#endif
231 Kokkos::deep_copy(dofCoords, this->dofCoords_);
232 }
233
234 virtual
235 void
236 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
237#ifdef HAVE_INTREPID2_DEBUG
238 // Verify rank of output array.
239 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
240 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
241 // Verify 0th dimension of output array.
242 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
243 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
244 // Verify 1st dimension of output array.
245 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
246 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
247#endif
248 Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
249 }
250
251
252 virtual
253 const char*
254 getName() const override {
255 return "Intrepid2_HCURL_TRI_I1_FEM";
256 }
257
258 virtual
259 bool
260 requireOrientation() const override {
261 return true;
262 }
263
274 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
275 if(subCellDim == 1)
276 return Teuchos::rcp( new
277 Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Line<2> >()));
278
279 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
280 }
281
286 };
287
288}// namespace Intrepid2
289
291
292#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 default FEM basis functions of degree 1 for H(curl) functions on Triangle cells.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Triangle cell.
virtual bool requireOrientation() const override
True if orientation is required.
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...
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 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...
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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_TRI_I1_FEM.