Intrepid2
Intrepid2_DerivedBasis_HGRAD_HEX.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
23#ifndef Intrepid2_DerivedBasis_HGRAD_HEX_h
24#define Intrepid2_DerivedBasis_HGRAD_HEX_h
25
27
29
30namespace Intrepid2
31{
32 // TODO: make this a subclass of TensorBasis3 instead, following what we've done for H(curl) and H(div)
33 template<class HGRAD_LINE>
35 : public Basis_TensorBasis<typename HGRAD_LINE::BasisBase>
36 {
37 public:
38 using ExecutionSpace = typename HGRAD_LINE::ExecutionSpace;
39 using OutputValueType = typename HGRAD_LINE::OutputValueType;
40 using PointValueType = typename HGRAD_LINE::PointValueType;
41
42 using OutputViewType = typename HGRAD_LINE::OutputViewType;
43 using PointViewType = typename HGRAD_LINE::PointViewType ;
44 using ScalarViewType = typename HGRAD_LINE::ScalarViewType;
45
46 using LineBasis = HGRAD_LINE;
48 using BasisBase = typename HGRAD_LINE::BasisBase;
50
51 std::string name_;
52 ordinal_type order_x_;
53 ordinal_type order_y_;
54 ordinal_type order_z_;
55 EPointType pointType_;
56
63 Basis_Derived_HGRAD_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT)
64 :
65 TensorBasis(Teuchos::rcp( new QuadBasis(polyOrder_x,polyOrder_y, pointType)),
66 Teuchos::rcp( new LineBasis(polyOrder_z, pointType)))
67 {
68 this->functionSpace_ = FUNCTION_SPACE_HGRAD;
69
70 std::ostringstream basisName;
71 basisName << "HGRAD_HEX (" << this->TensorBasis::getName() << ")";
72 name_ = basisName.str();
73
74 order_x_ = polyOrder_x;
75 order_y_ = polyOrder_y;
76 order_z_ = polyOrder_z;
77 pointType_ = pointType;
78
79 this->setShardsTopologyAndTags();
80 }
81
82
87 Basis_Derived_HGRAD_HEX(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT) :
88 Basis_Derived_HGRAD_HEX(polyOrder, polyOrder, polyOrder, pointType) {}
89
90
93 virtual bool requireOrientation() const override {
94 return (this->getDofCount(1,0) > 1); //if it has more than 1 DOF per edge, then it needs orientations
95 }
96
97 using BasisBase::getValues;
98
103 virtual
104 const char*
105 getName() const override {
106 return name_.c_str();
107 }
108
118 Teuchos::RCP<BasisBase>
119 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
120 if(subCellDim == 1) {
121 switch(subCellOrd) {
122 case 0:
123 case 2:
124 case 4:
125 case 6:
126 return Teuchos::rcp( new LineBasis(order_x_, pointType_) );
127 case 1:
128 case 3:
129 case 5:
130 case 7:
131 return Teuchos::rcp( new LineBasis(order_y_, pointType_) );
132 case 8:
133 case 9:
134 case 10:
135 case 11:
136 return Teuchos::rcp( new LineBasis(order_z_, pointType_) );
137 }
138 } else if(subCellDim == 2) {
139 switch(subCellOrd) {
140 case 0:
141 return Teuchos::rcp( new QuadBasis(order_x_, order_z_, pointType_) );
142 case 1:
143 return Teuchos::rcp( new QuadBasis(order_y_,order_z_, pointType_) );
144 case 2:
145 return Teuchos::rcp( new QuadBasis(order_x_, order_z_, pointType_) );
146 case 3:
147 return Teuchos::rcp( new QuadBasis(order_z_, order_y_, pointType_) );
148 case 4:
149 return Teuchos::rcp( new QuadBasis(order_y_, order_x_, pointType_) );
150 case 5:
151 return Teuchos::rcp( new QuadBasis(order_x_, order_y_, pointType_) );
152 }
153 }
154
155 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
156 }
157
163 getHostBasis() const override {
165
166 auto hostBasis = Teuchos::rcp(new HostBasis(order_x_, order_y_, order_z_, pointType_));
167
168 return hostBasis;
169 }
170 };
171} // end namespace Intrepid2
172
173#endif /* Intrepid2_DerivedBasis_HGRAD_HEX_h */
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type),...
Implementation of H(grad) basis on the quadrilateral that is templated on H(grad) on the line.
Implementation of bases that are tensor products of two or three component bases.
virtual HostBasisPtr< 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.
virtual bool requireOrientation() const override
True if orientation is required.
Teuchos::RCP< BasisBase > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Basis_Derived_HGRAD_HEX(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
Basis_Derived_HGRAD_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
Basis defined as the tensor product of two component bases.
virtual const char * getName() const override
Returns basis name.