Intrepid2
Intrepid2_DerivedBasis_HGRAD_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
17#ifndef Intrepid2_DerivedBasis_HGRAD_QUAD_h
18#define Intrepid2_DerivedBasis_HGRAD_QUAD_h
19
21
22namespace Intrepid2
23{
24 template<class HGRAD_LINE>
26 : public Basis_TensorBasis<typename HGRAD_LINE::BasisBase>
27 {
28 protected:
29 std::string name_;
30 ordinal_type order_x_, order_y_;
31 EPointType pointType_;
32 public:
33 using ExecutionSpace = typename HGRAD_LINE::ExecutionSpace;
34 using OutputValueType = typename HGRAD_LINE::OutputValueType;
35 using PointValueType = typename HGRAD_LINE::PointValueType;
36
37 using OutputViewType = typename HGRAD_LINE::OutputViewType;
38 using PointViewType = typename HGRAD_LINE::PointViewType ;
39 using ScalarViewType = typename HGRAD_LINE::ScalarViewType;
40
41 using BasisBase = typename HGRAD_LINE::BasisBase;
42 using LineBasis = HGRAD_LINE;
44
50 Basis_Derived_HGRAD_QUAD(int polyOrder_x, int polyOrder_y, const EPointType pointType=POINTTYPE_DEFAULT)
51 :
52 TensorBasis(Teuchos::rcp( new LineBasis(polyOrder_x, pointType)),
53 Teuchos::rcp( new LineBasis(polyOrder_y, pointType)))
54 {
55 this->functionSpace_ = FUNCTION_SPACE_HGRAD;
56
57 std::ostringstream basisName;
58 basisName << "HGRAD_QUAD (" << this->TensorBasis::getName() << ")";
59 name_ = basisName.str();
60
61 order_x_= polyOrder_x;
62 order_y_ = polyOrder_y;
63 pointType_ = pointType;
64
65 this->setShardsTopologyAndTags();
66 }
67
72 Basis_Derived_HGRAD_QUAD(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT) : Basis_Derived_HGRAD_QUAD(polyOrder,polyOrder,pointType) {}
73
74
77 virtual bool requireOrientation() const override {
78 return (this->getDofCount(1,0) > 1); //if it has more than 1 DOF per edge, then it needs orientations
79 }
80
81 using BasisBase::getValues;
82
87 virtual
88 const char*
89 getName() const override {
90 return name_.c_str();
91 }
92
102 Teuchos::RCP<BasisBase>
103 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
104 if(subCellDim == 1) {
105 switch(subCellOrd) {
106 case 0:
107 case 2:
108 return Teuchos::rcp( new LineBasis(order_x_, pointType_) );
109 case 1:
110 case 3:
111 return Teuchos::rcp( new LineBasis(order_y_, pointType_) );
112 }
113 }
114
115 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
116 }
117
123 getHostBasis() const override {
125
126 auto hostBasis = Teuchos::rcp(new HostBasis(order_x_, order_y_, pointType_));
127
128 return hostBasis;
129 }
130 };
131} // end namespace Intrepid2
132
133#endif /* Intrepid2_DerivedBasis_HGRAD_QUAD_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 bases that are tensor products of two or three component bases.
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.
Basis_Derived_HGRAD_QUAD(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
Basis_Derived_HGRAD_QUAD(int polyOrder_x, int polyOrder_y, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
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.
virtual const char * getName() const override
Returns basis name.