16 #ifndef Intrepid2_SerendipityBasis_h
17 #define Intrepid2_SerendipityBasis_h
27 template<
typename BasisBaseClass =
void>
33 using BasisBase = BasisBaseClass;
34 using BasisPtr = Teuchos::RCP<BasisBase>;
41 int numTensorialExtrusions_;
43 using DeviceType =
typename BasisBase::DeviceType;
44 using ExecutionSpace =
typename BasisBase::ExecutionSpace;
45 using OutputValueType =
typename BasisBase::OutputValueType;
46 using PointValueType =
typename BasisBase::PointValueType;
48 using OrdinalTypeArray1D =
typename BasisBase::OrdinalTypeArray1D;
49 using OrdinalTypeArray1DHost =
typename BasisBase::OrdinalTypeArray1DHost;
50 using OrdinalTypeArray2DHost =
typename BasisBase::OrdinalTypeArray2DHost;
51 using OutputViewType =
typename BasisBase::OutputViewType;
52 using PointViewType =
typename BasisBase::PointViewType;
54 OrdinalTypeArray1D ordinalMap_;
63 INTREPID2_TEST_FOR_EXCEPTION(fullBasis_->getBasisType() != BASIS_FEM_HIERARCHICAL, std::invalid_argument,
"SerendipityBasis only supports full bases whose type is BASIS_FEM_HIERARCHICAL");
64 this->basisType_ = fullBasis_->getBasisType();
66 this->functionSpace_ = fullBasis_->getFunctionSpace();
67 this->basisDegree_ = fullBasis_->getDegree();
70 std::ostringstream basisName;
71 basisName <<
"Serendipity(" << fullBasis_->getName() <<
")";
72 name_ = basisName.str();
76 vector<ordinal_type> fullBasisOrdinals;
77 vector<vector<ordinal_type>> polynomialDegreeOfField;
78 vector<vector<ordinal_type>> polynomialH1DegreeOfField;
79 ordinal_type fullCardinality = fullBasis_->getCardinality();
80 ordinal_type basisH1Degree = (this->functionSpace_ == FUNCTION_SPACE_HVOL) ? this->basisDegree_ + 1 : this->basisDegree_;
81 for (ordinal_type i=0; i<fullCardinality; i++)
83 vector<ordinal_type> fieldDegree = fullBasis_->getPolynomialDegreeOfFieldAsVector(i);
84 vector<ordinal_type> fieldH1Degree = fullBasis_->getH1PolynomialDegreeOfFieldAsVector(i);
85 ordinal_type superlinearDegree = 0;
86 for (
const ordinal_type & p : fieldH1Degree)
91 superlinearDegree += p;
94 if (superlinearDegree <= basisH1Degree)
97 fullBasisOrdinals.push_back(i);
98 polynomialDegreeOfField.push_back(fieldDegree);
99 polynomialH1DegreeOfField.push_back(fieldH1Degree);
102 this->basisCardinality_ = fullBasisOrdinals.size();
103 ordinalMap_ = OrdinalTypeArray1D(
"serendipity ordinal map",fullBasisOrdinals.size());
105 auto ordinalMapHost = Kokkos::create_mirror_view(ordinalMap_);
106 const ordinal_type fullBasisCardinality = fullBasisOrdinals.size();
107 for (ordinal_type i=0; i<fullBasisCardinality; i++)
109 ordinalMapHost(i) = fullBasisOrdinals[i];
111 Kokkos::deep_copy(ordinalMap_, ordinalMapHost);
114 this->basisCellTopologyKey_ = fullBasis_->getBaseCellTopology().getKey();
115 this->numTensorialExtrusions_ = fullBasis_->getNumTensorialExtrusions();
116 const ordinal_type spaceDim = fullBasis_->getBaseCellTopology().getDimension() + numTensorialExtrusions_;
118 this->basisCoordinates_ = COORDINATES_CARTESIAN;
121 int degreeSize = fullBasis_->getPolynomialDegreeLength();
122 this->fieldOrdinalPolynomialDegree_ = OrdinalTypeArray2DHost(
"TensorBasis - field ordinal polynomial degree", this->basisCardinality_, degreeSize);
123 this->fieldOrdinalH1PolynomialDegree_ = OrdinalTypeArray2DHost(
"TensorBasis - field ordinal H^1 degree", this->basisCardinality_, degreeSize);
125 for (ordinal_type fieldOrdinal1 = 0; fieldOrdinal1 < this->basisCardinality_; fieldOrdinal1++)
127 for (
int d=0; d<degreeSize; d++)
129 this->fieldOrdinalPolynomialDegree_ (fieldOrdinal1,d) = polynomialDegreeOfField [fieldOrdinal1][d];
130 this->fieldOrdinalH1PolynomialDegree_(fieldOrdinal1,d) = polynomialH1DegreeOfField[fieldOrdinal1][d];
135 const auto & cardinality = this->basisCardinality_;
137 const ordinal_type tagSize = 4;
138 const ordinal_type posScDim = 0;
139 const ordinal_type posScOrd = 1;
140 const ordinal_type posDfOrd = 2;
141 const ordinal_type posDfCnt = 3;
143 OrdinalTypeArray1DHost tagView(
"tag view", cardinality*tagSize);
145 const OrdinalTypeArray2DHost fullBasisOrdinalToTag = fullBasis->getAllDofTags();
148 vector<vector<ordinal_type>> subcellDofCount(spaceDim+1);
149 vector<vector<ordinal_type>> subcellDofOrdinal(spaceDim+1);
150 for (ordinal_type d=0; d<=spaceDim; d++)
152 const ordinal_type numSubcells = cellTopo->getSubcellCount(d);
153 subcellDofCount[d] = vector<ordinal_type>(numSubcells,0);
154 subcellDofOrdinal[d] = vector<ordinal_type>(numSubcells,0);
156 for (ordinal_type fieldOrdinal=0; fieldOrdinal<fullBasisCardinality; fieldOrdinal++)
158 const ordinal_type fullFieldOrdinal = fullBasisOrdinals[fieldOrdinal];
159 const ordinal_type subcellDim = fullBasisOrdinalToTag(fullFieldOrdinal,posScDim);
160 const ordinal_type subcellOrd = fullBasisOrdinalToTag(fullFieldOrdinal,posScOrd);
161 subcellDofCount[subcellDim][subcellOrd]++;
163 for (ordinal_type fieldOrdinal=0; fieldOrdinal<fullBasisCardinality; fieldOrdinal++)
165 const ordinal_type fullFieldOrdinal = fullBasisOrdinals[fieldOrdinal];
166 const ordinal_type subcellDim = fullBasisOrdinalToTag(fullFieldOrdinal,posScDim);
167 const ordinal_type subcellOrd = fullBasisOrdinalToTag(fullFieldOrdinal,posScOrd);
168 const ordinal_type subcellDfCnt = subcellDofCount[subcellDim][subcellOrd];
169 const ordinal_type subcellDfOrd = subcellDofOrdinal[subcellDim][subcellOrd]++;
171 tagView(tagSize*fieldOrdinal + posScDim) = subcellDim;
172 tagView(tagSize*fieldOrdinal + posScOrd) = subcellOrd;
173 tagView(tagSize*fieldOrdinal + posDfOrd) = subcellDfOrd;
174 tagView(tagSize*fieldOrdinal + posDfCnt) = subcellDfCnt;
178 this->setOrdinalTagData(this->tagToOrdinal_,
181 this->basisCardinality_,
188 virtual int getNumTensorialExtrusions()
const override
190 return numTensorialExtrusions_;
199 auto basisValues = fullBasis_->allocateBasisValues(points,operatorType);
200 basisValues.setOrdinalFilter(this->ordinalMap_);
207 using BasisBase::getValues;
216 return name_.c_str();
234 const EOperator operatorType = OPERATOR_VALUE )
const override
236 fullBasis_->getValues(outputValues, inputPoints, operatorType);
237 outputValues.setOrdinalFilter(this->ordinalMap_);
259 void getValues( OutputViewType outputValues,
const PointViewType inputPoints,
260 const EOperator operatorType = OPERATOR_VALUE )
const override
262 INTREPID2_TEST_FOR_EXCEPTION(
true, std::logic_error,
263 ">>> ERROR (Basis::getValues): this method is not supported by SerendipityBasis (use the getValues() method that accepts a BasisValues object instead).");
270 virtual HostBasisPtr<OutputValueType, PointValueType>
274 auto hostBasis = Teuchos::rcp(
new HostBasis(fullBasis_->getHostBasis()));
Implements arbitrary-dimensional extrusion of a base shards::CellTopology.
Serendipity Basis, defined as the sub-basis of a provided basis, consisting of basis elements for whi...
View-like interface to tensor points; point components are stored separately; the appropriate coordin...
virtual const char * getName() const override
Returns basis name.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
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...
The data containers in Intrepid2 that support sum factorization and other reduced-data optimizations ...
virtual BasisValues< OutputValueType, DeviceType > allocateBasisValues(TensorPoints< PointValueType, DeviceType > points, const EOperator operatorType=OPERATOR_VALUE) const override
Allocate BasisValues container suitable for passing to the getValues() variant that takes a TensorPoi...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
BasisPtr getUnderlyingBasis() const
Returns a pointer to the underlying full basis.
OrdinalTypeArray1D ordinalMap() const
Returns the ordinal map from the Serendipity basis ordinal to the ordinal in the underlying full basi...
static CellTopoPtr cellTopology(const shards::CellTopology &shardsCellTopo, ordinal_type tensorialDegree=0)
static accessor that returns a CellTopoPtr; these are lazily constructed and cached.
virtual HostBasisPtr< OutputValueType, PointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
SerendipityBasis(BasisPtr fullBasis)
Constructor.