Intrepid2
Intrepid2_DerivedBasis_HVOL_QUAD.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 
18 #ifndef Intrepid2_DerivedBasis_HVOL_QUAD_h
19 #define Intrepid2_DerivedBasis_HVOL_QUAD_h
20 
23 
24 namespace Intrepid2
25 {
30  template<class HVOL_LINE>
32  :
33  public Basis_TensorBasis<typename HVOL_LINE::BasisBase>
34  {
35  protected:
36  std::string name_;
37  using LineBasis = HVOL_LINE;
39 
40  ordinal_type polyOrder_x_, polyOrder_y_;
41  EPointType pointType_;
42  public:
43  using ExecutionSpace = typename HVOL_LINE::ExecutionSpace;
44  using OutputValueType = typename HVOL_LINE::OutputValueType;
45  using PointValueType = typename HVOL_LINE::PointValueType;
46 
47  using OutputViewType = typename HVOL_LINE::OutputViewType;
48  using PointViewType = typename HVOL_LINE::PointViewType ;
49  using ScalarViewType = typename HVOL_LINE::ScalarViewType;
50 
51  using BasisBase = typename HVOL_LINE::BasisBase;
52 
58  Basis_Derived_HVOL_QUAD(int polyOrder_x, int polyOrder_y, const EPointType pointType=POINTTYPE_DEFAULT)
59  :
60  TensorBasis(Teuchos::rcp( new LineBasis(polyOrder_x, pointType)),
61  Teuchos::rcp( new LineBasis(polyOrder_y, pointType))),
62  polyOrder_x_(polyOrder_x),
63  polyOrder_y_(polyOrder_y),
64  pointType_(pointType)
65  {
66  this->functionSpace_ = FUNCTION_SPACE_HVOL;
67 
68  std::ostringstream basisName;
69  basisName << "HVOL_QUAD (" << this->TensorBasis::getName() << ")";
70  name_ = basisName.str();
71 
72  this->setShardsTopologyAndTags();
73  }
74 
79  Basis_Derived_HVOL_QUAD(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT) : Basis_Derived_HVOL_QUAD(polyOrder,polyOrder,pointType) {}
80 
85  virtual
86  const char*
87  getName() const override {
88  return name_.c_str();
89  }
90 
93  virtual bool requireOrientation() const override {
94  return false;
95  }
96 
97  virtual OperatorTensorDecomposition getSimpleOperatorDecomposition(const EOperator &operatorType) const override
98  {
99  const EOperator VALUE = Intrepid2::OPERATOR_VALUE;
100 
101  if (operatorType == VALUE)
102  {
103  return OperatorTensorDecomposition(VALUE, VALUE);
104  }
105  else
106  {
107  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"operator not yet supported");
108  }
109  }
110 
112 
117  virtual BasisPtr<typename Kokkos::HostSpace::device_type, typename BasisBase::OutputValueType, typename BasisBase::PointValueType>
118  getHostBasis() const override {
120  return Teuchos::rcp( new HostBasisType(polyOrder_x_, polyOrder_y_, pointType_) );
121  }
122  };
123 } // end namespace Intrepid2
124 
125 #endif /* Intrepid2_DerivedBasis_HVOL_QUAD_h */
Basis_Derived_HVOL_QUAD(int polyOrder_x, int polyOrder_y, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
virtual OperatorTensorDecomposition getSimpleOperatorDecomposition(const EOperator &operatorType) const override
Returns a simple decomposition of the specified operator: what operator(s) should be applied to basis...
Implementation of bases that are tensor products of two or three component bases. ...
virtual BasisPtr< typename Kokkos::HostSpace::device_type, typename BasisBase::OutputValueType, typename BasisBase::PointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Implementation of H(vol) basis on the quadrilateral that is templated on H(vol) on the line...
Basis_Derived_HVOL_QUAD(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
For a multi-component tensor basis, specifies the operators to be applied to the components to produc...
virtual void getValues(BasisValues< OutputValueType, DeviceType > outputValues, const TensorPoints< PointValueType, DeviceType > inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell, using point and output value containers that allow pre...
virtual const char * getName() const override
Returns basis name.
H(vol) basis on the line based on Legendre polynomials.
virtual const char * getName() const override
Returns basis name.
Basis defined as the tensor product of two component bases.
virtual bool requireOrientation() const override
True if orientation is required.