43 namespace MultiRegions
46 ExpList2DHomogeneous1D::ExpList2DHomogeneous1D():
57 const bool dealiasing,
61 lhom,useFFT,dealiasing)
71 "Size of basis number of points and number"
72 "of planes are not the same");
77 planes.size() * planes[0]->GetExpSize());
79 for(cnt = n = 0; n < planes.size(); ++n)
82 for (i = 0; i < planes[n]->GetExpSize(); ++i)
84 (*m_exp)[cnt++] = planes[n]->GetExp(i);
97 const bool dealiasing,
101 lhom,useFFT,dealiasing)
114 for (j = 0; j < nel; ++j)
119 for (n = 1; n <
m_planes.size(); ++n)
123 for(j = 0; j < nel; ++j)
125 (*m_exp).push_back((*
m_exp)[j]);
141 std::dynamic_pointer_cast<ExpList> (In.
m_planes[0]);
143 for (
int n = 0; n <
m_planes.size(); ++n)
162 int ncoeffs_per_plane =
m_planes[0]->GetNcoeffs();
163 int npoints_per_plane =
m_planes[0]->GetTotPoints();
174 int nel =
m_planes[0]->GetExpSize();
179 for (cnt = n = 0; n < nzplanes; ++n)
182 tmparray=
m_coeffs + ncoeffs_per_plane*n);
184 tmparray =
m_phys + npoints_per_plane*n);
186 for(i = 0; i < nel; ++i)
189 + n*ncoeffs_per_plane;
191 + n*npoints_per_plane;
211 (*m_exp)[eid]->GetCoords(xc0);
217 "output coord_1 is not defined");
219 (*m_exp)[eid]->GetCoords(xc0,xc1);
224 ASSERTL0(
false,
"Cannot have coordim being three dimensions"
225 "in a homogeneous field");
235 for(n = 0; n <
m_planes.size(); n++)
243 for (n = 0; n < nyplanes; ++n)
245 Vmath::Fill(npoints,z[n],tmp_xc = xhom + npoints*n,1);
283 int npoints =
m_planes[0]->GetTotPoints();
301 for(n = 0; n <
m_planes.size(); n++)
309 for(n = 0; n < nyplanes; ++n)
311 Vmath::Fill(npoints,z[n],tmp_xc = xhom + npoints*n,1);
330 std::ostream &outfile,
335 int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
341 coords[1] = coords[0] + nquad0*nquad1;
342 coords[2] = coords[1] + nquad0*nquad1;
344 GetCoords(expansion,coords[0],coords[1],coords[2]);
346 outfile <<
"Zone, I=" << nquad0 <<
", J=" << nquad1
347 <<
", F=Block" << std::endl;
351 for(i = 0; i < nquad0*nquad1; ++i)
353 outfile << coords[j][i] <<
" ";
355 outfile << std::endl;
361 std::ostream &outfile,
365 boost::ignore_unused(istrip);
370 m_planes[0]->WriteVtkPieceHeader(outfile, expansion);
375 int outputExtraPlane = 0;
380 outputExtraPlane = 1;
383 int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
384 int nquad1 =
m_planes.size() + outputExtraPlane;
385 int ntot = nquad0*nquad1;
386 int ntotminus = (nquad0-1)*(nquad1-1);
392 GetCoords(expansion,coords[0],coords[1],coords[2]);
394 if (outputExtraPlane)
399 tmp = coords[0] + (nquad1-1)*nquad0, 1);
401 tmp = coords[1] + (nquad1-1)*nquad0, 1);
404 (coords[2][nquad0] - coords[2][0]);
405 Vmath::Fill(nquad0, z, tmp = coords[2] + (nquad1-1)*nquad0, 1);
408 outfile <<
" <Piece NumberOfPoints=\""
409 << ntot <<
"\" NumberOfCells=\""
410 << ntotminus <<
"\">" << endl;
411 outfile <<
" <Points>" << endl;
412 outfile <<
" <DataArray type=\"Float64\" "
413 <<
"NumberOfComponents=\"3\" format=\"ascii\">" << endl;
415 for (i = 0; i < ntot; ++i)
417 for (j = 0; j < 3; ++j)
419 outfile << coords[j][i] <<
" ";
424 outfile <<
" </DataArray>" << endl;
425 outfile <<
" </Points>" << endl;
426 outfile <<
" <Cells>" << endl;
427 outfile <<
" <DataArray type=\"Int32\" "
428 <<
"Name=\"connectivity\" format=\"ascii\">" << endl;
429 for (i = 0; i < nquad0-1; ++i)
431 for (j = 0; j < nquad1-1; ++j)
433 outfile << j*nquad0 + i <<
" "
434 << j*nquad0 + i + 1 <<
" "
435 << (j+1)*nquad0 + i + 1 <<
" "
436 << (j+1)*nquad0 + i << endl;
440 outfile <<
" </DataArray>" << endl;
441 outfile <<
" <DataArray type=\"Int32\" "
442 <<
"Name=\"offsets\" format=\"ascii\">" << endl;
443 for (i = 0; i < ntotminus; ++i)
445 outfile << i*4+4 <<
" ";
448 outfile <<
" </DataArray>" << endl;
449 outfile <<
" <DataArray type=\"UInt8\" "
450 <<
"Name=\"types\" format=\"ascii\">" << endl;
451 for (i = 0; i < ntotminus; ++i)
456 outfile <<
" </DataArray>" << endl;
457 outfile <<
" </Cells>" << endl;
458 outfile <<
" <PointData>" << endl;
465 int nPtsPlane =
m_planes[0]->GetNpoints();
469 "Output vector does not have sufficient dimensions to"
471 ASSERTL1(normals[0].size() >= nPtsPlane,
472 "Output vector does not have sufficient dimensions to"
479 for (
int i = 0; i < nDim; ++i)
481 for (
int n = 1; n < nPlanes; ++n)
484 &normals[i][n*nPtsPlane], 1);
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Describes the specification for a Basis.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Populate normals with the normals of all expansions.
virtual ~ExpList2DHomogeneous1D()
Destructor.
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion)
ExpList2DHomogeneous1D()
Default constructor.
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
void GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
This function calculates the coordinates of all the elemental quadrature points .
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2)
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
NekDouble m_lhom
Width of homogeneous direction.
LibUtilities::TranspositionSharedPtr m_transposition
Array< OneD, ExpListSharedPtr > m_planes
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
LibUtilities::CommSharedPtr m_comm
Communicator.
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
LibUtilities::SessionReaderSharedPtr m_session
Session.
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
@ eFourier
Fourier Expansion .
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
The above copyright notice and this permission notice shall be included.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha - x.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)