37 #include <boost/core/ignore_unused.hpp> 53 namespace MultiRegions
67 const std::vector<unsigned int> &eIDs,
68 const bool DeclareCoeffPhysArrays,
70 ExpList(In, eIDs, DeclareCoeffPhysArrays)
81 if (DeclareCoeffPhysArrays)
118 for (
auto &expIt : expansions)
125 if((TetGeom = std::dynamic_pointer_cast<SpatialDomains::TetGeom>(expIt.second->m_geomShPtr)))
135 (*m_exp).push_back(tet);
162 else if((HexGeom = std::dynamic_pointer_cast<SpatialDomains::HexGeom>(expIt.second->m_geomShPtr)))
165 (*m_exp).push_back(hex);
172 ASSERTL0(
false,
"dynamic cast to a proper Geometry3D failed");
206 const std::string &variable,
219 = graph3D->GetExpansions(variable);
221 for (
auto &expIt : expansions)
228 if((TetGeom = std::dynamic_pointer_cast<
232 = expIt.second->m_basisKeyVector[0];
234 = expIt.second->m_basisKeyVector[1];
236 = expIt.second->m_basisKeyVector[2];
241 ASSERTL0(
false,
"LocalRegions::NodalTetExp is not " 249 tet->SetElmtId(elmtid++);
250 (*m_exp).push_back(tet);
253 else if((PrismGeom = std::dynamic_pointer_cast<SpatialDomains
254 ::PrismGeom>(expIt.second->m_geomShPtr)))
257 = expIt.second->m_basisKeyVector[0];
259 = expIt.second->m_basisKeyVector[1];
261 = expIt.second->m_basisKeyVector[2];
266 prism->SetElmtId(elmtid++);
267 (*m_exp).push_back(prism);
269 else if((PyrGeom = std::dynamic_pointer_cast<
273 = expIt.second->m_basisKeyVector[0];
275 = expIt.second->m_basisKeyVector[1];
277 = expIt.second->m_basisKeyVector[2];
282 pyramid->SetElmtId(elmtid++);
283 (*m_exp).push_back(pyramid);
285 else if((HexGeom = std::dynamic_pointer_cast<
289 = expIt.second->m_basisKeyVector[0];
291 = expIt.second->m_basisKeyVector[1];
293 = expIt.second->m_basisKeyVector[2];
298 hex->SetElmtId(elmtid++);
299 (*m_exp).push_back(hex);
303 ASSERTL0(
false,
"dynamic cast to a proper Geometry3D " 351 for(
int i = 0; i < expansions.size(); ++i)
358 auto expmap = expansions.find(i);
359 ASSERTL1(expmap != expansions.end(),
"Unable to find expansion.");
362 if((TetGeom = std::dynamic_pointer_cast<
366 = exp->m_basisKeyVector[0];
368 = exp->m_basisKeyVector[1];
370 = exp->m_basisKeyVector[2];
375 ASSERTL0(
false,
"LocalRegions::NodalTetExp is not " 383 tet->SetElmtId(elmtid++);
384 (*m_exp).push_back(tet);
387 else if((PrismGeom = std::dynamic_pointer_cast<
391 = exp->m_basisKeyVector[0];
393 = exp->m_basisKeyVector[1];
395 = exp->m_basisKeyVector[2];
400 prism->SetElmtId(elmtid++);
401 (*m_exp).push_back(prism);
403 else if((PyrGeom = std::dynamic_pointer_cast<
407 = exp->m_basisKeyVector[0];
409 = exp->m_basisKeyVector[1];
411 = exp->m_basisKeyVector[2];
416 pyramid->SetElmtId(elmtid++);
417 (*m_exp).push_back(pyramid);
419 else if((HexGeom = std::dynamic_pointer_cast<
423 = exp->m_basisKeyVector[0];
425 = exp->m_basisKeyVector[1];
427 = exp->m_basisKeyVector[2];
432 hex->SetElmtId(elmtid++);
433 (*m_exp).push_back(hex);
437 ASSERTL0(
false,
"dynamic cast to a proper Geometry3D " 463 switch ((*
m_exp)[i]->DetShapeType())
470 ASSERTL0(
false,
"Unknown expansion type.");
482 boost::ignore_unused(istrip);
485 int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
486 int nquad1 = (*m_exp)[expansion]->GetNumPoints(1);
487 int nquad2 = (*m_exp)[expansion]->GetNumPoints(2);
488 int ntot = nquad0*nquad1*nquad2;
489 int ntotminus = (nquad0-1)*(nquad1-1)*(nquad2-1);
495 (*m_exp)[expansion]->GetCoords(coords[0],coords[1],coords[2]);
497 outfile <<
" <Piece NumberOfPoints=\"" 498 << ntot <<
"\" NumberOfCells=\"" 499 << ntotminus <<
"\">" << endl;
500 outfile <<
" <Points>" << endl;
501 outfile <<
" <DataArray type=\"Float64\" " 502 <<
"NumberOfComponents=\"3\" format=\"ascii\">" << endl;
504 for (i = 0; i < ntot; ++i)
506 for (j = 0; j < 3; ++j)
508 outfile << setprecision(8) << scientific
509 << (float)coords[j][i] <<
" ";
515 outfile <<
" </DataArray>" << endl;
516 outfile <<
" </Points>" << endl;
517 outfile <<
" <Cells>" << endl;
518 outfile <<
" <DataArray type=\"Int32\" " 519 <<
"Name=\"connectivity\" format=\"ascii\">" << endl;
520 for (i = 0; i < nquad0-1; ++i)
522 for (j = 0; j < nquad1-1; ++j)
524 for (k = 0; k < nquad2-1; ++k)
526 outfile << k*nquad0*nquad1 + j*nquad0 + i <<
" " 527 << k*nquad0*nquad1 + j*nquad0 + i + 1 <<
" " 528 << k*nquad0*nquad1 + (j+1)*nquad0 + i + 1 <<
" " 529 << k*nquad0*nquad1 + (j+1)*nquad0 + i <<
" " 530 << (k+1)*nquad0*nquad1 + j*nquad0 + i <<
" " 531 << (k+1)*nquad0*nquad1 + j*nquad0 + i + 1 <<
" " 532 << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i + 1 <<
" " 533 << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i <<
" " << endl;
538 outfile <<
" </DataArray>" << endl;
539 outfile <<
" <DataArray type=\"Int32\" " 540 <<
"Name=\"offsets\" format=\"ascii\">" << endl;
541 for (i = 0; i < ntotminus; ++i)
543 outfile << i*8+8 <<
" ";
546 outfile <<
" </DataArray>" << endl;
547 outfile <<
" <DataArray type=\"UInt8\" " 548 <<
"Name=\"types\" format=\"ascii\">" << endl;
549 for (i = 0; i < ntotminus; ++i)
554 outfile <<
" </DataArray>" << endl;
555 outfile <<
" </Cells>" << endl;
556 outfile <<
" <PointData>" << endl;
562 for (i = 0; i <
m_exp->size(); ++i)
564 for (j = 0; j < (*m_exp)[i]->GetNfaces(); ++j)
566 (*m_exp)[i]->ComputeFaceNormal(j);
581 int pt0 = (*m_exp)[i]->GetNumPoints(0);
582 int pt1 = (*m_exp)[i]->GetNumPoints(1);
583 int pt2 = (*m_exp)[i]->GetNumPoints(2);
584 int npt0 = (int) pt0*scale;
585 int npt1 = (int) pt1*scale;
586 int npt2 = (int) pt2*scale;
594 (*
m_exp)[i]->GetBasis(1)->GetPointsKey(),
595 (*
m_exp)[i]->GetBasis(2)->GetPointsKey(),
596 &inarray[cnt], newPointsKey0,
597 newPointsKey1, newPointsKey2,
601 cnt1 += npt0*npt1*npt2;
615 int pt0 = (*m_exp)[i]->GetNumPoints(0);
616 int pt1 = (*m_exp)[i]->GetNumPoints(1);
617 int pt2 = (*m_exp)[i]->GetNumPoints(2);
618 int npt0 = (int) pt0*scale;
619 int npt1 = (int) pt1*scale;
620 int npt2 = (int) pt2*scale;
631 (*
m_exp)[i]->GetBasis(0)->GetPointsKey(),
632 (*
m_exp)[i]->GetBasis(1)->GetPointsKey(),
633 (*
m_exp)[i]->GetBasis(2)->GetPointsKey(),
636 cnt += npt0*npt1*npt2;
int GetNumPoints() const
Return points order at which basis is defined.
#define ASSERTL0(condition, msg)
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
virtual void v_SetUpPhysNormals()
Set up the normals on each expansion.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
std::shared_ptr< TetGeom > TetGeomSharedPtr
std::shared_ptr< HexExp > HexExpSharedPtr
BasisType GetBasisType() const
Return type of expansion basis.
Lagrange Polynomials using the Gauss points .
std::shared_ptr< TetExp > TetExpSharedPtr
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
int GetExpSize(void)
This function returns the number of elements in the expansion.
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Base class for all multi-elemental spectral/hp expansions.
int GetNumModes() const
Returns the order of the basis.
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
void SetCoeffPhysOffsets()
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
LibUtilities::SessionReaderSharedPtr m_session
Session.
void PhysGalerkinProject3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
Defines a specification for a set of points.
std::shared_ptr< Expansion > ExpansionShPtr
virtual void v_ReadGlobalOptimizationParameters()
std::shared_ptr< PyrExp > PyrExpSharedPtr
std::shared_ptr< PyrGeom > PyrGeomSharedPtr
virtual ~ExpList3D()
Destructor.
void Interp3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
this function interpolates a 3D function evaluated at the quadrature points of the 3D basis...
std::shared_ptr< PrismExp > PrismExpSharedPtr
std::shared_ptr< HexGeom > HexGeomSharedPtr
void ReadGlobalOptimizationParameters()
Abstraction of a three-dimensional multi-elemental expansion which is merely a collection of local ex...
ExpList3D()
Default constructor.
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
std::shared_ptr< PrismGeom > PrismGeomSharedPtr
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Describes the specification for a Basis.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
std::map< int, ExpansionShPtr > ExpansionMap
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)