Intrepid2
Intrepid2_DerivedBasis_HVOL_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_HVOL_HEX_h
24#define Intrepid2_DerivedBasis_HVOL_HEX_h
25
28
29namespace Intrepid2
30{
35 template<class HVOL_LINE>
37 :
38 public Basis_TensorBasis<typename HVOL_LINE::BasisBase>
39 // TODO: make this a subclass of TensorBasis3 instead, following what we've done for H(curl) and H(div)
40 {
41 std::string name_;
42 ordinal_type polyOrder_x_, polyOrder_y_, polyOrder_z_;
43 EPointType pointType_;
44 public:
45 using ExecutionSpace = typename HVOL_LINE::ExecutionSpace;
46 using OutputValueType = typename HVOL_LINE::OutputValueType;
47 using PointValueType = typename HVOL_LINE::PointValueType;
48
49 using OutputViewType = typename HVOL_LINE::OutputViewType;
50 using PointViewType = typename HVOL_LINE::PointViewType ;
51 using ScalarViewType = typename HVOL_LINE::ScalarViewType;
52
53 using BasisBase = typename HVOL_LINE::BasisBase;
54
55 using LineBasis = HVOL_LINE;
58
65 Basis_Derived_HVOL_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT)
66 :
67 TensorBasis(Teuchos::rcp(new QuadBasis(polyOrder_x,polyOrder_y,pointType)),
68 Teuchos::rcp(new LineBasis(polyOrder_z,pointType))),
69 polyOrder_x_(polyOrder_x),
70 polyOrder_y_(polyOrder_y),
71 polyOrder_z_(polyOrder_z),
72 pointType_(pointType)
73 {
74 this->functionSpace_ = FUNCTION_SPACE_HVOL;
75
76 std::ostringstream basisName;
77 basisName << "HVOL_HEX (" << this->TensorBasis::getName() << ")";
78 name_ = basisName.str();
79
80 this->setShardsTopologyAndTags();
81 }
82
87 Basis_Derived_HVOL_HEX(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT) : Basis_Derived_HVOL_HEX(polyOrder, polyOrder, polyOrder,pointType) {}
88
93 virtual
94 const char*
95 getName() const override {
96 return name_.c_str();
97 }
98
101 virtual bool requireOrientation() const override {
102 return false;
103 }
104
105 virtual OperatorTensorDecomposition getSimpleOperatorDecomposition(const EOperator &operatorType) const override
106 {
107 const EOperator VALUE = Intrepid2::OPERATOR_VALUE;
108
109 if (operatorType == VALUE)
110 {
111 return OperatorTensorDecomposition(VALUE, VALUE);
112 }
113 else
114 {
115 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"operator not yet supported");
116 }
117 }
118
120
126 getHostBasis() const override {
128 return Teuchos::rcp( new HostBasisType(polyOrder_x_, polyOrder_y_, polyOrder_z_, pointType_) );
129 }
130 };
131} // end namespace Intrepid2
132
133#endif /* Intrepid2_DerivedBasis_HVOL_HEX_h */
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Implementation of H(vol) basis on the quadrilateral that is templated on H(vol) on the line.
Implementation of bases that are tensor products of two or three component bases.
Implementation of H(vol) basis on the quadrilateral that is templated on H(vol) on the line.
virtual bool requireOrientation() const override
True if orientation is required.
Basis_Derived_HVOL_HEX(int polyOrder_x, int polyOrder_y, int polyOrder_z, 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...
virtual const char * getName() const override
Returns basis name.
Basis_Derived_HVOL_HEX(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
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 defined as the tensor product of two component bases.
virtual const char * getName() const override
Returns basis name.
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...
For a multi-component tensor basis, specifies the operators to be applied to the components to produc...