Intrepid2
Intrepid2_HCURL_QUAD_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_QUAD_I1_FEM_HPP__
50#define __INTREPID2_HCURL_QUAD_I1_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
54
55namespace Intrepid2 {
56
100 namespace Impl {
101
106 public:
107 typedef struct Quadrilateral<4> cell_topology_type;
108
112 template<EOperator opType>
113 struct Serial {
114 template<typename OutputViewType,
115 typename inputViewType>
117 static void
118 getValues( OutputViewType output,
119 const inputViewType input );
120
121 };
122
123 template<typename DeviceType,
124 typename outputValueValueType, class ...outputValueProperties,
125 typename inputPointValueType, class ...inputPointProperties>
126 static void
127 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
128 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
129 const EOperator operatorType);
130
134 template<typename outputValueViewType,
135 typename inputPointViewType,
136 EOperator opType>
137 struct Functor {
138 outputValueViewType _outputValues;
139 const inputPointViewType _inputPoints;
140
143 inputPointViewType inputPoints_)
144 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
145
147 void operator()(const ordinal_type pt) const {
148 switch (opType) {
149 case OPERATOR_VALUE : {
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
156 case OPERATOR_CURL: {
157 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
158 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
159 Serial<opType>::getValues( output, input );
160 break;
161 }
162 default: {
163 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
164 opType != OPERATOR_CURL,
165 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_CI_FEM::Serial::getVAlues) operator is not supported");
166 break;
167 }
168 } //end switch
169 }
170
171 };
172 };
173 }
174
175 template<typename DeviceType = void,
176 typename outputValueType = double,
177 typename pointValueType = double>
178 class Basis_HCURL_QUAD_I1_FEM : public Basis<DeviceType,outputValueType, pointValueType> {
179 public:
183
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_QUAD_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_QUAD_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_QUAD_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_QUAD_I1_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() != 2, std::invalid_argument,
237 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::getDofCoeffs) rank = 2 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_HCURL_QUAD_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
241 // Verify 1st dimension of output array.
242 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
243 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
244#endif
245 Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
246 }
247
248 virtual
249 const char*
250 getName() const override {
251 return "Intrepid2_HCURL_QUAD_I1_FEM";
252 }
253
254 virtual
255 bool
256 requireOrientation() const override {
257 return true;
258 }
259
270 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
271 if(subCellDim == 1)
272 return Teuchos::rcp( new
273 Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Line<2> >()));
274
275 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
276 }
277
282
283 };
284}// namespace Intrepid2
285
287
288#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 Qadrilateral cells.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Quadrilateral cell.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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 bool requireOrientation() const override
True if orientation is required.
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.
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...
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_QUAD_I1_FEM.