23 #ifndef Intrepid2_DerivedBasis_HGRAD_HEX_h
24 #define Intrepid2_DerivedBasis_HGRAD_HEX_h
33 template<
class HGRAD_LINE>
38 using ExecutionSpace =
typename HGRAD_LINE::ExecutionSpace;
39 using OutputValueType =
typename HGRAD_LINE::OutputValueType;
40 using PointValueType =
typename HGRAD_LINE::PointValueType;
42 using OutputViewType =
typename HGRAD_LINE::OutputViewType;
43 using PointViewType =
typename HGRAD_LINE::PointViewType ;
44 using ScalarViewType =
typename HGRAD_LINE::ScalarViewType;
46 using LineBasis = HGRAD_LINE;
48 using BasisBase =
typename HGRAD_LINE::BasisBase;
52 ordinal_type order_x_;
53 ordinal_type order_y_;
54 ordinal_type order_z_;
55 EPointType pointType_;
66 Teuchos::rcp( new LineBasis(polyOrder_z, pointType)))
68 this->functionSpace_ = FUNCTION_SPACE_HGRAD;
70 std::ostringstream basisName;
72 name_ = basisName.str();
74 order_x_ = polyOrder_x;
75 order_y_ = polyOrder_y;
76 order_z_ = polyOrder_z;
77 pointType_ = pointType;
79 this->setShardsTopologyAndTags();
94 return (this->getDofCount(1,0) > 1);
97 using BasisBase::getValues;
106 return name_.c_str();
118 Teuchos::RCP<BasisBase>
120 if(subCellDim == 1) {
126 return Teuchos::rcp(
new LineBasis(order_x_, pointType_) );
131 return Teuchos::rcp(
new LineBasis(order_y_, pointType_) );
136 return Teuchos::rcp(
new LineBasis(order_z_, pointType_) );
138 }
else if(subCellDim == 2) {
141 return Teuchos::rcp(
new QuadBasis(order_x_, order_z_, pointType_) );
143 return Teuchos::rcp(
new QuadBasis(order_y_,order_z_, pointType_) );
145 return Teuchos::rcp(
new QuadBasis(order_x_, order_z_, pointType_) );
147 return Teuchos::rcp(
new QuadBasis(order_z_, order_y_, pointType_) );
149 return Teuchos::rcp(
new QuadBasis(order_y_, order_x_, pointType_) );
151 return Teuchos::rcp(
new QuadBasis(order_x_, order_y_, pointType_) );
155 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
162 virtual HostBasisPtr<OutputValueType, PointValueType>
166 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.
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.
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.