57 #ifndef Intrepid2_DerivedBasis_HGRAD_HEX_h 58 #define Intrepid2_DerivedBasis_HGRAD_HEX_h 67 template<
class HGRAD_LINE>
72 using ExecutionSpace =
typename HGRAD_LINE::ExecutionSpace;
73 using OutputValueType =
typename HGRAD_LINE::OutputValueType;
74 using PointValueType =
typename HGRAD_LINE::PointValueType;
76 using OutputViewType =
typename HGRAD_LINE::OutputViewType;
77 using PointViewType =
typename HGRAD_LINE::PointViewType ;
78 using ScalarViewType =
typename HGRAD_LINE::ScalarViewType;
80 using LineBasis = HGRAD_LINE;
82 using BasisBase =
typename HGRAD_LINE::BasisBase;
86 ordinal_type order_x_;
87 ordinal_type order_y_;
88 ordinal_type order_z_;
100 Teuchos::rcp( new LineBasis(polyOrder_z, pointType)))
102 this->functionSpace_ = FUNCTION_SPACE_HGRAD;
104 std::ostringstream basisName;
106 name_ = basisName.str();
108 order_x_ = polyOrder_x;
109 order_y_ = polyOrder_y;
110 order_z_ = polyOrder_z;
111 pointType_ = pointType;
113 this->setShardsTopologyAndTags();
128 return (this->getDofCount(1,0) > 1);
131 using BasisBase::getValues;
140 return name_.c_str();
152 Teuchos::RCP<BasisBase>
154 if(subCellDim == 1) {
160 return Teuchos::rcp(
new LineBasis(order_x_, pointType_) );
165 return Teuchos::rcp(
new LineBasis(order_y_, pointType_) );
170 return Teuchos::rcp(
new LineBasis(order_z_, pointType_) );
172 }
else if(subCellDim == 2) {
175 return Teuchos::rcp(
new QuadBasis(order_x_, order_z_, pointType_) );
177 return Teuchos::rcp(
new QuadBasis(order_y_,order_z_, pointType_) );
179 return Teuchos::rcp(
new QuadBasis(order_x_, order_z_, pointType_) );
181 return Teuchos::rcp(
new QuadBasis(order_z_, order_y_, pointType_) );
183 return Teuchos::rcp(
new QuadBasis(order_y_, order_x_, pointType_) );
185 return Teuchos::rcp(
new QuadBasis(order_x_, order_y_, pointType_) );
189 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
200 auto hostBasis = Teuchos::rcp(
new HostBasis(order_x_, order_y_, order_z_, pointType_));
Implementation of bases that are tensor products of two or three component bases. ...
Teuchos::RCP< BasisBase > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type), allowing host access to input and output views, and ensuring host execution of basis evaluation.
Implementation of H(grad) basis on the quadrilateral that is templated on H(grad) on the line...
virtual bool requireOrientation() const override
True if orientation is required.
virtual const char * getName() const override
Returns basis name.
EPointType
Enumeration of types of point distributions in Intrepid.
virtual const char * getName() const override
Returns basis name.
Basis_Derived_HGRAD_HEX(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
Basis defined as the tensor product of two component bases.
virtual HostBasisPtr< OutputValueType, PointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Basis_Derived_HGRAD_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.