49 #ifndef __INTREPID2_HCURL_TRI_IN_FEM_HPP__ 50 #define __INTREPID2_HCURL_TRI_IN_FEM_HPP__ 57 #include "Teuchos_LAPACK.hpp" 81 #define CardinalityHCurlTri(order) (order*(order+2)) 95 template<EOperator opType>
97 template<
typename outputValueViewType,
98 typename inputPointViewType,
99 typename workViewType,
100 typename vinvViewType>
101 KOKKOS_INLINE_FUNCTION
103 getValues( outputValueViewType outputValues,
104 const inputPointViewType inputPoints,
106 const vinvViewType vinv );
108 KOKKOS_INLINE_FUNCTION
110 getWorkSizePerPoint(ordinal_type order) {
111 auto cardinality = CardinalityHCurlTri(order);
116 return 5*cardinality;
118 return getDkCardinality<opType,2>()*cardinality;
123 template<
typename DeviceType, ordinal_type numPtsPerEval,
124 typename outputValueValueType,
class ...outputValueProperties,
125 typename inputPointValueType,
class ...inputPointProperties,
126 typename vinvValueType,
class ...vinvProperties>
128 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
129 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
130 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
136 template<
typename outputValueViewType,
137 typename inputPointViewType,
138 typename vinvViewType,
139 typename workViewType,
141 ordinal_type numPtsEval>
143 outputValueViewType _outputValues;
144 const inputPointViewType _inputPoints;
145 const vinvViewType _coeffs;
148 KOKKOS_INLINE_FUNCTION
149 Functor( outputValueViewType outputValues_,
150 inputPointViewType inputPoints_,
151 vinvViewType coeffs_,
153 : _outputValues(outputValues_), _inputPoints(inputPoints_),
154 _coeffs(coeffs_), _work(work_) {}
156 KOKKOS_INLINE_FUNCTION
157 void operator()(
const size_type iter)
const {
161 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
162 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
165 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
167 auto vcprop = Kokkos::common_view_alloc_prop(_work);
168 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
171 case OPERATOR_VALUE : {
172 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
176 case OPERATOR_CURL: {
177 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
182 INTREPID2_TEST_FOR_ABORT(
true,
183 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::Functor) operator is not supported");
192 template<
typename DeviceType = void,
193 typename outputValueType = double,
194 typename pointValueType =
double>
196 :
public Basis<DeviceType,outputValueType,pointValueType> {
205 const EPointType pointType = POINTTYPE_EQUISPACED);
218 getValues( OutputViewType outputValues,
219 const PointViewType inputPoints,
220 const EOperator operatorType = OPERATOR_VALUE)
const override {
221 #ifdef HAVE_INTREPID2_DEBUG 229 Impl::Basis_HCURL_TRI_In_FEM::
230 getValues<DeviceType,numPtsPerEval>( outputValues,
238 getDofCoords( ScalarViewType dofCoords )
const override {
239 #ifdef HAVE_INTREPID2_DEBUG 241 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
242 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
244 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
245 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
247 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
248 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
250 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
255 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
256 #ifdef HAVE_INTREPID2_DEBUG 258 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
259 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
261 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
262 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
264 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
265 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
267 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
271 getExpansionCoeffs( ScalarViewType coeffs )
const {
273 Kokkos::deep_copy(coeffs, this->
coeffs_);
279 return "Intrepid2_HCURL_TRI_In_FEM";
299 if(subCellDim == 1) {
300 return Teuchos::rcp(
new 304 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
315 Kokkos::DynRankView<scalarType,DeviceType>
coeffs_;
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
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...
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...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Basis_HCURL_TRI_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual bool requireOrientation() const override
True if orientation is required.
ordinal_type getCardinality() const
Returns cardinality of the basis.
EPointType pointType_
type of lattice used for creating the DoF coordinates
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
See Intrepid2::Basis_HCURL_TRI_In_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
See Intrepid2::Basis_HCURL_TRI_In_FEM.
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
EPointType
Enumeration of types of point distributions in Intrepid.
Triangle topology, 3 nodes.
See Intrepid2::Basis_HCURL_TRI_In_FEM.
Definition file for FEM basis functions of degree n for H(curl) functions on TRI. ...
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Tr...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Header file for the abstract base class Intrepid2::Basis.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...