51 namespace MultiRegions
65 const std::vector<unsigned int> &eIDs,
66 const bool DeclareCoeffPhysArrays,
68 ExpList(In, eIDs, DeclareCoeffPhysArrays)
79 if (DeclareCoeffPhysArrays)
116 SpatialDomains::ExpansionMap::const_iterator expIt;
117 for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
124 if((TetGeom = boost::dynamic_pointer_cast<SpatialDomains::TetGeom>(expIt->second->m_geomShPtr)))
134 (*m_exp).push_back(tet);
161 else if((HexGeom = boost::dynamic_pointer_cast<SpatialDomains::HexGeom>(expIt->second->m_geomShPtr)))
164 (*m_exp).push_back(hex);
171 ASSERTL0(
false,
"dynamic cast to a proper Geometry3D failed");
205 const std::string &variable,
218 = graph3D->GetExpansions(variable);
220 SpatialDomains::ExpansionMap::const_iterator expIt;
221 for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
228 if((TetGeom = boost::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 = boost::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 = boost::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 = boost::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 SpatialDomains::ExpansionMap::const_iterator expmap = expansions.find(i);
359 ASSERTL1(expmap != expansions.end(),
"Unable to find expansion.");
362 if((TetGeom = boost::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 = boost::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 = boost::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 = boost::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.");
483 int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
484 int nquad1 = (*m_exp)[expansion]->GetNumPoints(1);
485 int nquad2 = (*m_exp)[expansion]->GetNumPoints(2);
486 int ntot = nquad0*nquad1*nquad2;
487 int ntotminus = (nquad0-1)*(nquad1-1)*(nquad2-1);
493 (*m_exp)[expansion]->GetCoords(coords[0],coords[1],coords[2]);
495 outfile <<
" <Piece NumberOfPoints=\""
496 << ntot <<
"\" NumberOfCells=\""
497 << ntotminus <<
"\">" << endl;
498 outfile <<
" <Points>" << endl;
499 outfile <<
" <DataArray type=\"Float64\" "
500 <<
"NumberOfComponents=\"3\" format=\"ascii\">" << endl;
502 for (i = 0; i < ntot; ++i)
504 for (j = 0; j < 3; ++j)
506 outfile << setprecision(8) << scientific
507 << (float)coords[j][i] <<
" ";
513 outfile <<
" </DataArray>" << endl;
514 outfile <<
" </Points>" << endl;
515 outfile <<
" <Cells>" << endl;
516 outfile <<
" <DataArray type=\"Int32\" "
517 <<
"Name=\"connectivity\" format=\"ascii\">" << endl;
518 for (i = 0; i < nquad0-1; ++i)
520 for (j = 0; j < nquad1-1; ++j)
522 for (k = 0; k < nquad2-1; ++k)
524 outfile << k*nquad0*nquad1 + j*nquad0 + i <<
" "
525 << k*nquad0*nquad1 + j*nquad0 + i + 1 <<
" "
526 << k*nquad0*nquad1 + (j+1)*nquad0 + i + 1 <<
" "
527 << k*nquad0*nquad1 + (j+1)*nquad0 + i <<
" "
528 << (k+1)*nquad0*nquad1 + j*nquad0 + i <<
" "
529 << (k+1)*nquad0*nquad1 + j*nquad0 + i + 1 <<
" "
530 << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i + 1 <<
" "
531 << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i <<
" " << endl;
536 outfile <<
" </DataArray>" << endl;
537 outfile <<
" <DataArray type=\"Int32\" "
538 <<
"Name=\"offsets\" format=\"ascii\">" << endl;
539 for (i = 0; i < ntotminus; ++i)
541 outfile << i*8+8 <<
" ";
544 outfile <<
" </DataArray>" << endl;
545 outfile <<
" <DataArray type=\"UInt8\" "
546 <<
"Name=\"types\" format=\"ascii\">" << endl;
547 for (i = 0; i < ntotminus; ++i)
552 outfile <<
" </DataArray>" << endl;
553 outfile <<
" </Cells>" << endl;
554 outfile <<
" <PointData>" << endl;
560 for (i = 0; i <
m_exp->size(); ++i)
562 for (j = 0; j < (*m_exp)[i]->GetNfaces(); ++j)
564 (*m_exp)[i]->ComputeFaceNormal(j);
579 int pt0 = (*m_exp)[i]->GetNumPoints(0);
580 int pt1 = (*m_exp)[i]->GetNumPoints(1);
581 int pt2 = (*m_exp)[i]->GetNumPoints(2);
582 int npt0 = (int) pt0*scale;
583 int npt1 = (int) pt1*scale;
584 int npt2 = (int) pt2*scale;
592 (*
m_exp)[i]->GetBasis(1)->GetPointsKey(),
593 (*
m_exp)[i]->GetBasis(2)->GetPointsKey(),
594 &inarray[cnt], newPointsKey0,
595 newPointsKey1, newPointsKey2,
599 cnt1 += npt0*npt1*npt2;
613 int pt0 = (*m_exp)[i]->GetNumPoints(0);
614 int pt1 = (*m_exp)[i]->GetNumPoints(1);
615 int pt2 = (*m_exp)[i]->GetNumPoints(2);
616 int npt0 = (int) pt0*scale;
617 int npt1 = (int) pt1*scale;
618 int npt2 = (int) pt2*scale;
629 (*
m_exp)[i]->GetBasis(0)->GetPointsKey(),
630 (*
m_exp)[i]->GetBasis(1)->GetPointsKey(),
631 (*
m_exp)[i]->GetBasis(2)->GetPointsKey(),
634 cnt += npt0*npt1*npt2;
boost::shared_ptr< PyrGeom > PyrGeomSharedPtr
#define ASSERTL0(condition, msg)
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual void v_SetUpPhysNormals()
Set up the normals on each expansion.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
BasisType GetBasisType() const
Return type of expansion basis.
Lagrange Polynomials using the Gauss points .
boost::shared_ptr< HexExp > HexExpSharedPtr
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
boost::shared_ptr< HexGeom > HexGeomSharedPtr
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.
boost::shared_ptr< TetExp > TetExpSharedPtr
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
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.
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
boost::shared_ptr< PyrExp > PyrExpSharedPtr
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
int GetNumPoints() const
Return points order at which basis is defined.
void SetCoeffPhysOffsets()
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
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.
boost::shared_ptr< PrismExp > PrismExpSharedPtr
virtual void v_ReadGlobalOptimizationParameters()
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...
boost::shared_ptr< Expansion > ExpansionShPtr
boost::shared_ptr< PrismGeom > PrismGeomSharedPtr
void ReadGlobalOptimizationParameters()
boost::shared_ptr< TetGeom > TetGeomSharedPtr
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...
int GetNumModes() const
Returns the order of the basis.
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...
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Describes the specification for a Basis.
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)