51 namespace MultiRegions
65 const std::vector<unsigned int> &eIDs,
66 const bool DeclareCoeffPhysArrays)
67 :
ExpList(In, eIDs, DeclareCoeffPhysArrays)
107 SpatialDomains::ExpansionMap::const_iterator expIt;
108 for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
115 if((TetGeom = boost::dynamic_pointer_cast<SpatialDomains::TetGeom>(expIt->second->m_geomShPtr)))
125 (*m_exp).push_back(tet);
152 else if((HexGeom = boost::dynamic_pointer_cast<SpatialDomains::HexGeom>(expIt->second->m_geomShPtr)))
155 (*m_exp).push_back(hex);
162 ASSERTL0(
false,
"dynamic cast to a proper Geometry3D failed");
192 const std::string &variable) :
204 = graph3D->GetExpansions(variable);
206 SpatialDomains::ExpansionMap::const_iterator expIt;
207 for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
214 if((TetGeom = boost::dynamic_pointer_cast<
218 = expIt->second->m_basisKeyVector[0];
220 = expIt->second->m_basisKeyVector[1];
222 = expIt->second->m_basisKeyVector[2];
227 ASSERTL0(
false,
"LocalRegions::NodalTetExp is not "
235 tet->SetElmtId(elmtid++);
236 (*m_exp).push_back(tet);
239 else if((PrismGeom = boost::dynamic_pointer_cast<SpatialDomains
240 ::PrismGeom>(expIt->second->m_geomShPtr)))
243 = expIt->second->m_basisKeyVector[0];
245 = expIt->second->m_basisKeyVector[1];
247 = expIt->second->m_basisKeyVector[2];
252 prism->SetElmtId(elmtid++);
253 (*m_exp).push_back(prism);
255 else if((PyrGeom = boost::dynamic_pointer_cast<
259 = expIt->second->m_basisKeyVector[0];
261 = expIt->second->m_basisKeyVector[1];
263 = expIt->second->m_basisKeyVector[2];
268 pyramid->SetElmtId(elmtid++);
269 (*m_exp).push_back(pyramid);
271 else if((HexGeom = boost::dynamic_pointer_cast<
275 = expIt->second->m_basisKeyVector[0];
277 = expIt->second->m_basisKeyVector[1];
279 = expIt->second->m_basisKeyVector[2];
284 hex->SetElmtId(elmtid++);
285 (*m_exp).push_back(hex);
289 ASSERTL0(
false,
"dynamic cast to a proper Geometry3D "
331 for(
int i = 0; i < expansions.size(); ++i)
338 SpatialDomains::ExpansionMap::const_iterator expmap = expansions.find(i);
339 ASSERTL1(expmap != expansions.end(),
"Unable to find expansion.");
342 if((TetGeom = boost::dynamic_pointer_cast<
346 = exp->m_basisKeyVector[0];
348 = exp->m_basisKeyVector[1];
350 = exp->m_basisKeyVector[2];
355 ASSERTL0(
false,
"LocalRegions::NodalTetExp is not "
363 tet->SetElmtId(elmtid++);
364 (*m_exp).push_back(tet);
367 else if((PrismGeom = boost::dynamic_pointer_cast<
371 = exp->m_basisKeyVector[0];
373 = exp->m_basisKeyVector[1];
375 = exp->m_basisKeyVector[2];
380 prism->SetElmtId(elmtid++);
381 (*m_exp).push_back(prism);
383 else if((PyrGeom = boost::dynamic_pointer_cast<
387 = exp->m_basisKeyVector[0];
389 = exp->m_basisKeyVector[1];
391 = exp->m_basisKeyVector[2];
396 pyramid->SetElmtId(elmtid++);
397 (*m_exp).push_back(pyramid);
399 else if((HexGeom = boost::dynamic_pointer_cast<
403 = exp->m_basisKeyVector[0];
405 = exp->m_basisKeyVector[1];
407 = exp->m_basisKeyVector[2];
412 hex->SetElmtId(elmtid++);
413 (*m_exp).push_back(hex);
417 ASSERTL0(
false,
"dynamic cast to a proper Geometry3D "
455 for(i = 0; i <
m_exp->size(); ++i)
459 m_offset_elmt_id[i] = i;
460 m_ncoeffs += (*m_exp)[i]->GetNcoeffs();
461 m_npoints += (*m_exp)[i]->GetTotPoints();
474 switch ((*
m_exp)[i]->DetShapeType())
481 ASSERTL0(
false,
"Unknown expansion type.");
494 int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
495 int nquad1 = (*m_exp)[expansion]->GetNumPoints(1);
496 int nquad2 = (*m_exp)[expansion]->GetNumPoints(2);
497 int ntot = nquad0*nquad1*nquad2;
498 int ntotminus = (nquad0-1)*(nquad1-1)*(nquad2-1);
504 (*m_exp)[expansion]->GetCoords(coords[0],coords[1],coords[2]);
506 outfile <<
" <Piece NumberOfPoints=\""
507 << ntot <<
"\" NumberOfCells=\""
508 << ntotminus <<
"\">" << endl;
509 outfile <<
" <Points>" << endl;
510 outfile <<
" <DataArray type=\"Float64\" "
511 <<
"NumberOfComponents=\"3\" format=\"ascii\">" << endl;
513 for (i = 0; i < ntot; ++i)
515 for (j = 0; j < 3; ++j)
517 outfile << setprecision(8) << scientific
518 << (float)coords[j][i] <<
" ";
524 outfile <<
" </DataArray>" << endl;
525 outfile <<
" </Points>" << endl;
526 outfile <<
" <Cells>" << endl;
527 outfile <<
" <DataArray type=\"Int32\" "
528 <<
"Name=\"connectivity\" format=\"ascii\">" << endl;
529 for (i = 0; i < nquad0-1; ++i)
531 for (j = 0; j < nquad1-1; ++j)
533 for (k = 0; k < nquad2-1; ++k)
535 outfile << k*nquad0*nquad1 + j*nquad0 + i <<
" "
536 << k*nquad0*nquad1 + j*nquad0 + i + 1 <<
" "
537 << k*nquad0*nquad1 + (j+1)*nquad0 + i + 1 <<
" "
538 << k*nquad0*nquad1 + (j+1)*nquad0 + i <<
" "
539 << (k+1)*nquad0*nquad1 + j*nquad0 + i <<
" "
540 << (k+1)*nquad0*nquad1 + j*nquad0 + i + 1 <<
" "
541 << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i + 1 <<
" "
542 << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i <<
" " << endl;
547 outfile <<
" </DataArray>" << endl;
548 outfile <<
" <DataArray type=\"Int32\" "
549 <<
"Name=\"offsets\" format=\"ascii\">" << endl;
550 for (i = 0; i < ntotminus; ++i)
552 outfile << i*8+8 <<
" ";
555 outfile <<
" </DataArray>" << endl;
556 outfile <<
" <DataArray type=\"UInt8\" "
557 <<
"Name=\"types\" format=\"ascii\">" << endl;
558 for (i = 0; i < ntotminus; ++i)
563 outfile <<
" </DataArray>" << endl;
564 outfile <<
" </Cells>" << endl;
565 outfile <<
" <PointData>" << endl;
571 for (i = 0; i <
m_exp->size(); ++i)
573 for (j = 0; j < (*m_exp)[i]->GetNfaces(); ++j)
575 (*m_exp)[i]->ComputeFaceNormal(j);
590 int pt0 = (*m_exp)[i]->GetNumPoints(0);
591 int pt1 = (*m_exp)[i]->GetNumPoints(1);
592 int pt2 = (*m_exp)[i]->GetNumPoints(2);
593 int npt0 = (int) pt0*scale;
594 int npt1 = (int) pt1*scale;
595 int npt2 = (int) pt2*scale;
603 (*
m_exp)[i]->GetBasis(1)->GetPointsKey(),
604 (*
m_exp)[i]->GetBasis(2)->GetPointsKey(),
605 &inarray[cnt], newPointsKey0,
606 newPointsKey1, newPointsKey2,
610 cnt1 += npt0*npt1*npt2;
624 int pt0 = (*m_exp)[i]->GetNumPoints(0);
625 int pt1 = (*m_exp)[i]->GetNumPoints(1);
626 int pt2 = (*m_exp)[i]->GetNumPoints(2);
627 int npt0 = (int) pt0*scale;
628 int npt1 = (int) pt1*scale;
629 int npt2 = (int) pt2*scale;
640 (*
m_exp)[i]->GetBasis(0)->GetPointsKey(),
641 (*
m_exp)[i]->GetBasis(1)->GetPointsKey(),
642 (*
m_exp)[i]->GetBasis(2)->GetPointsKey(),
645 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)
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Base class for all multi-elemental spectral/hp expansions.
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
boost::shared_ptr< PyrExp > PyrExpSharedPtr
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Array< OneD, int > m_offset_elmt_id
Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coef...
int GetNumPoints() const
Return points order at which basis is defined.
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.
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
#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)