Intrepid2
Intrepid2_HCURL_HEX_I1_FEM.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Intrepid2 Package
4//
5// Copyright 2007 NTESS and the Intrepid2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
16#ifndef __INTREPID2_HCURL_HEX_I1_FEM_HPP__
17#define __INTREPID2_HCURL_HEX_I1_FEM_HPP__
18
19#include "Intrepid2_Basis.hpp"
21
22namespace Intrepid2 {
23
84 namespace Impl {
85
90 public:
91 typedef struct Hexahedron<8> cell_topology_type;
92
96 template<EOperator opType>
97 struct Serial {
98 template<typename OutputViewType,
99 typename inputViewType>
100 KOKKOS_INLINE_FUNCTION
101 static void
102 getValues( OutputViewType output,
103 const inputViewType input );
104 };
105
106 template<typename DeviceType,
107 typename outputValueValueType, class ...outputValueProperties,
108 typename inputPointValueType, class ...inputPointProperties>
109 static void
110 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
111 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
112 const EOperator operatorType);
113
117 template<typename outputValueViewType,
118 typename inputPointViewType,
119 EOperator opType>
120 struct Functor {
121 outputValueViewType _outputValues;
122 const inputPointViewType _inputPoints;
123
124 KOKKOS_INLINE_FUNCTION
125 Functor( outputValueViewType outputValues_,
126 inputPointViewType inputPoints_ )
127 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
128
129 KOKKOS_INLINE_FUNCTION
130 void operator()(const ordinal_type pt) const {
131 switch (opType) {
132 case OPERATOR_VALUE:
133 case OPERATOR_CURL: {
134 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
135 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
136 Serial<opType>::getValues( output, input );
137 break;
138 }
139 default: {
140 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
141 opType != OPERATOR_CURL,
142 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_CI_FEM::Serial::getVAlues) operator is not supported");
143 break;
144 }
145 } //end switch
146 }
147 };
148 };
149 }
150
151 template<typename DeviceType = void,
152 typename outputValueType = double,
153 typename pointValueType = double>
154 class Basis_HCURL_HEX_I1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
155 public:
159
163
167
168 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
169
170 virtual
171 void
172 getValues( OutputViewType outputValues,
173 const PointViewType inputPoints,
174 const EOperator operatorType = OPERATOR_VALUE ) const override {
175#ifdef HAVE_INTREPID2_DEBUG
177 inputPoints,
178 operatorType,
179 this->getBaseCellTopology(),
180 this->getCardinality() );
181#endif
182 Impl::Basis_HCURL_HEX_I1_FEM::
183 getValues<DeviceType>( outputValues,
184 inputPoints,
185 operatorType );
186 }
187
188 virtual void
189 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
190 ordinal_type& perThreadSpaceSize,
191 const PointViewType inputPointsconst,
192 const EOperator operatorType = OPERATOR_VALUE) const override;
193
194 KOKKOS_INLINE_FUNCTION
195 virtual void
196 getValues(
197 OutputViewType outputValues,
198 const PointViewType inputPoints,
199 const EOperator operatorType,
200 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
201 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
202 const ordinal_type subcellDim = -1,
203 const ordinal_type subcellOrdinal = -1) const override;
204
205 virtual
206 void
207 getDofCoords( ScalarViewType dofCoords ) const override {
208#ifdef HAVE_INTREPID2_DEBUG
209 // Verify rank of output array.
210 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
211 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
212 // Verify 0th dimension of output array.
213 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
214 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
215 // Verify 1st dimension of output array.
216 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
217 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
218#endif
219 Kokkos::deep_copy(dofCoords, this->dofCoords_);
220 }
221
222
223 virtual
224 void
225 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
226#ifdef HAVE_INTREPID2_DEBUG
227 // Verify rank of output array.
228 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
230 // Verify 0th dimension of output array.
231 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
232 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
233 // Verify 1st dimension of output array.
234 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
236#endif
237 Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
238 }
239
240 virtual
241 const char*
242 getName() const override {
243 return "Intrepid2_HCURL_HEX_I1_FEM";
244 }
245
246 virtual
247 bool
248 requireOrientation() const override {
249 return true;
250 }
251
262 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
263
264 if(subCellDim == 1)
265 return Teuchos::rcp( new
266 Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Line<2> >()));
267 else if(subCellDim == 2) {
268 return Teuchos::rcp(new
270 }
271
272 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
273 }
274
279
280 };
281}// namespace Intrepid2
282
284
285#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 HEX cells.
Header file for the Intrepid2::Basis_HCURL_QUAD_I1_FEM class.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Hexahedron cell.
virtual bool requireOrientation() const override
True if orientation is required.
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 const char * getName() const override
Returns basis name.
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 H(curl)-compatible FEM basis of degree 1 on Quadrilateral 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_HEX_I1_FEM.