Intrepid2
Intrepid2_DerivedBasis_HGRAD_WEDGE.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 
17 #ifndef Intrepid2_DerivedBasis_HGRAD_WEDGE_h
18 #define Intrepid2_DerivedBasis_HGRAD_WEDGE_h
19 
22 
23 namespace Intrepid2
24 {
25  template<class HGRAD_TRI, class HGRAD_LINE>
27  : public Basis_TensorBasis<typename HGRAD_LINE::BasisBase>
28  {
29  protected:
30  std::string name_;
31  ordinal_type order_xy_, order_z_;
32  EPointType pointType_;
33  public:
34  using ExecutionSpace = typename HGRAD_LINE::ExecutionSpace;
35  using OutputValueType = typename HGRAD_LINE::OutputValueType;
36  using PointValueType = typename HGRAD_LINE::PointValueType;
37 
38  using OutputViewType = typename HGRAD_LINE::OutputViewType;
39  using PointViewType = typename HGRAD_LINE::PointViewType ;
40  using ScalarViewType = typename HGRAD_LINE::ScalarViewType;
41 
42  using BasisBase = typename HGRAD_LINE::BasisBase;
43  using TriBasis = HGRAD_TRI;
44  using LineBasis = HGRAD_LINE;
46 
52  Basis_Derived_HGRAD_WEDGE(int polyOrder_xy, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT)
53  :
54  TensorBasis(Teuchos::rcp( new TriBasis(polyOrder_xy, pointType)),
55  Teuchos::rcp( new LineBasis(polyOrder_z, pointType)))
56  {
57  this->functionSpace_ = FUNCTION_SPACE_HGRAD;
58 
59  std::ostringstream basisName;
60  basisName << "HGRAD_WEDGE (" << this->TensorBasis::getName() << ")";
61  name_ = basisName.str();
62 
63  order_xy_= polyOrder_xy;
64  order_z_ = polyOrder_z;
65  pointType_ = pointType;
66 
67  this->setShardsTopologyAndTags();
68  }
69 
74  Basis_Derived_HGRAD_WEDGE(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT) : Basis_Derived_HGRAD_WEDGE(polyOrder,polyOrder,pointType) {}
75 
76 
79  virtual bool requireOrientation() const override {
80  return (this->getDofCount(1,0) > 1); //if it has more than 1 DOF per edge, than it needs orientations
81  }
82 
83  using BasisBase::getValues;
84 
89  virtual
90  const char*
91  getName() const override {
92  return name_.c_str();
93  }
94 
104  Teuchos::RCP<BasisBase>
105  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
106  {
107  if(subCellDim == 1) {
108  if(subCellOrd < 6) //edges associated to basal and upper triagular faces
109  return Teuchos::rcp( new LineBasis(order_xy_, pointType_) );
110  else
111  return Teuchos::rcp( new LineBasis(order_z_, pointType_) );
112  }
113  else if(subCellDim == 2) {
114  switch(subCellOrd) {
115  case 0:
116  return Teuchos::rcp( new Intrepid2::Basis_Derived_HGRAD_QUAD<LineBasis>(order_xy_, order_z_, pointType_) );
117  case 1:
118  return Teuchos::rcp( new Intrepid2::Basis_Derived_HGRAD_QUAD<LineBasis>(order_xy_, order_z_, pointType_) );
119  case 2:
120  return Teuchos::rcp( new Intrepid2::Basis_Derived_HGRAD_QUAD<LineBasis>(order_z_, order_xy_, pointType_) );
121  case 3:
122  return Teuchos::rcp( new TriBasis(order_xy_, pointType_) );
123  case 4:
124  return Teuchos::rcp( new TriBasis(order_xy_, pointType_) );
125  default:
126  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"subCellOrd is out of bounds");
127  }
128  }
129  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"subCellDim is out of bounds");
130  }
131 
136  virtual HostBasisPtr<OutputValueType, PointValueType>
137  getHostBasis() const override {
139 
140  auto hostBasis = Teuchos::rcp(new HostBasis(order_xy_, order_z_, pointType_));
141 
142  return hostBasis;
143  }
144  };
145 } // end namespace Intrepid2
146 
147 #endif /* Intrepid2_DerivedBasis_HGRAD_WEDGE_h */
Basis_Derived_HGRAD_WEDGE(int polyOrder_xy, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
Implementation of bases that are tensor products of two or three component bases. ...
Basis_Derived_HGRAD_WEDGE(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
virtual const char * getName() const override
Returns basis name.
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.
Teuchos::RCP< BasisBase > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual HostBasisPtr< OutputValueType, PointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Basis defined as the tensor product of two component bases.