43 namespace MultiRegions
54 const bool useFFT,
const bool dealiasing,
68 "Size of basis number of points and number"
69 "of planes are not the same");
73 planes.size() * planes[0]->GetExpSize());
75 for (cnt = n = 0; n < planes.size(); ++n)
78 for (i = 0; i < planes[n]->GetExpSize(); ++i)
80 (*m_exp)[cnt++] = planes[n]->GetExp(i);
91 const bool useFFT,
const bool dealiasing,
102 m_session, graph1D,
false,
"DefaultVar", ImpType);
108 for (j = 0; j < nel; ++j)
113 for (n = 1; n <
m_planes.size(); ++n)
117 for (j = 0; j < nel; ++j)
119 (*m_exp).push_back((*
m_exp)[j]);
133 std::dynamic_pointer_cast<ExpList>(In.
m_planes[0]);
135 for (
int n = 0; n <
m_planes.size(); ++n)
154 int ncoeffs_per_plane =
m_planes[0]->GetNcoeffs();
155 int npoints_per_plane =
m_planes[0]->GetTotPoints();
160 m_ncoeffs = ncoeffs_per_plane * nzplanes;
161 m_npoints = npoints_per_plane * nzplanes;
166 int nel =
m_planes[0]->GetExpSize();
171 for (cnt = n = 0; n < nzplanes; ++n)
173 m_planes[n]->SetCoeffsArray(tmparray =
175 m_planes[n]->SetPhysArray(tmparray =
m_phys + npoints_per_plane * n);
177 for (i = 0; i < nel; ++i)
180 m_planes[n]->GetCoeff_Offset(i) + n * ncoeffs_per_plane;
182 m_planes[n]->GetPhys_Offset(i) + n * npoints_per_plane;
201 (*m_exp)[eid]->GetCoords(xc0);
206 ASSERTL0(xc1.size() != 0,
"output coord_1 is not defined");
208 (*m_exp)[eid]->GetCoords(xc0, xc1);
213 ASSERTL0(
false,
"Cannot have coordim being three dimensions"
214 "in a homogeneous field");
224 for (n = 0; n <
m_planes.size(); n++)
232 for (n = 0; n < nyplanes; ++n)
234 Vmath::Fill(npoints, z[n], tmp_xc = xhom + npoints * n, 1);
237 Vmath::Vcopy(npoints, xc0, 1, tmp_xc = xc0 + npoints * n, 1);
240 Vmath::Vcopy(npoints, xc1, 1, tmp_xc = xc1 + npoints * n, 1);
271 int npoints =
m_planes[0]->GetTotPoints();
289 for (n = 0; n <
m_planes.size(); n++)
297 for (n = 0; n < nyplanes; ++n)
299 Vmath::Fill(npoints, z[n], tmp_xc = xhom + npoints * n, 1);
302 Vmath::Vcopy(npoints, xc0, 1, tmp_xc = xc0 + npoints * n, 1);
305 Vmath::Vcopy(npoints, xc1, 1, tmp_xc = xc1 + npoints * n, 1);
321 int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
327 coords[1] = coords[0] + nquad0 * nquad1;
328 coords[2] = coords[1] + nquad0 * nquad1;
330 GetCoords(expansion, coords[0], coords[1], coords[2]);
332 outfile <<
"Zone, I=" << nquad0 <<
", J=" << nquad1 <<
", F=Block"
337 for (i = 0; i < nquad0 * nquad1; ++i)
339 outfile << coords[j][i] <<
" ";
341 outfile << std::endl;
346 int expansion,
int istrip)
348 boost::ignore_unused(istrip);
353 m_planes[0]->WriteVtkPieceHeader(outfile, expansion);
358 int outputExtraPlane = 0;
363 outputExtraPlane = 1;
366 int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
367 int nquad1 =
m_planes.size() + outputExtraPlane;
368 int ntot = nquad0 * nquad1;
369 int ntotminus = (nquad0 - 1) * (nquad1 - 1);
375 GetCoords(expansion, coords[0], coords[1], coords[2]);
377 if (outputExtraPlane)
382 tmp = coords[0] + (nquad1 - 1) * nquad0, 1);
384 tmp = coords[1] + (nquad1 - 1) * nquad0, 1);
387 (coords[2][nquad0] - coords[2][0]);
388 Vmath::Fill(nquad0, z, tmp = coords[2] + (nquad1 - 1) * nquad0, 1);
391 outfile <<
" <Piece NumberOfPoints=\"" << ntot <<
"\" NumberOfCells=\""
392 << ntotminus <<
"\">" << endl;
393 outfile <<
" <Points>" << endl;
394 outfile <<
" <DataArray type=\"Float64\" "
395 <<
"NumberOfComponents=\"3\" format=\"ascii\">" << endl;
397 for (i = 0; i < ntot; ++i)
399 for (j = 0; j < 3; ++j)
401 outfile << coords[j][i] <<
" ";
406 outfile <<
" </DataArray>" << endl;
407 outfile <<
" </Points>" << endl;
408 outfile <<
" <Cells>" << endl;
409 outfile <<
" <DataArray type=\"Int32\" "
410 <<
"Name=\"connectivity\" format=\"ascii\">" << endl;
411 for (i = 0; i < nquad0 - 1; ++i)
413 for (j = 0; j < nquad1 - 1; ++j)
415 outfile << j * nquad0 + i <<
" " << j * nquad0 + i + 1 <<
" "
416 << (j + 1) * nquad0 + i + 1 <<
" " << (j + 1) * nquad0 + i
421 outfile <<
" </DataArray>" << endl;
422 outfile <<
" <DataArray type=\"Int32\" "
423 <<
"Name=\"offsets\" format=\"ascii\">" << endl;
424 for (i = 0; i < ntotminus; ++i)
426 outfile << i * 4 + 4 <<
" ";
429 outfile <<
" </DataArray>" << endl;
430 outfile <<
" <DataArray type=\"UInt8\" "
431 <<
"Name=\"types\" format=\"ascii\">" << endl;
432 for (i = 0; i < ntotminus; ++i)
437 outfile <<
" </DataArray>" << endl;
438 outfile <<
" </Cells>" << endl;
439 outfile <<
" <PointData>" << endl;
446 int nPtsPlane =
m_planes[0]->GetNpoints();
450 "Output vector does not have sufficient dimensions to"
452 ASSERTL1(normals[0].size() >= nPtsPlane,
453 "Output vector does not have sufficient dimensions to"
460 for (
int i = 0; i < nDim; ++i)
462 for (
int n = 1; n < nPlanes; ++n)
465 &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...
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2) override
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip) override
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble >> &normals) override
Populate normals with the normals of all expansions.
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion) override
virtual ~ExpList2DHomogeneous1D()
Destructor.
ExpList2DHomogeneous1D()
Default constructor.
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.
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 .
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 scalar y = alpha + x.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)