Nektar++
ExpList2DHomogeneous1D.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ExpList2DHomogeneous1D.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: An ExpList2D which is homogeneous in 1 direction and so
32// uses much of the functionality from a ExpList and its daughters
33//
34///////////////////////////////////////////////////////////////////////////////
35
38
39using namespace std;
40
42{
43// Forward declaration for typedefs
45{
46}
47
48// Constructor for ExpList2DHomogeneous1D to act as a Explist2D field
51 const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom,
52 const bool useFFT, const bool dealiasing,
55 : ExpListHomogeneous1D(e2DH1D, pSession, HomoBasis, lhom, useFFT,
56 dealiasing)
57{
58 int i, n, cnt;
59
60 if (comm)
61 {
62 m_comm = comm;
63 }
64
65 ASSERTL1(m_planes.size() == planes.size(),
66 "Size of basis number of points and number"
67 "of planes are not the same");
68
69 // Set up expansion list with elements from all planes.
71 planes.size() * planes[0]->GetExpSize());
72
73 for (cnt = n = 0; n < planes.size(); ++n)
74 {
75 m_planes[n] = planes[n];
76 for (i = 0; i < planes[n]->GetExpSize(); ++i)
77 {
78 (*m_exp)[cnt++] = planes[n]->GetExp(i);
79 }
80 }
81
83}
84
85// Constructor for ExpList2DHomogeneous1D to act as a Explist2D field
88 const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom,
89 const bool useFFT, const bool dealiasing,
92 : ExpListHomogeneous1D(e2DH1D, pSession, HomoBasis, lhom, useFFT,
93 dealiasing)
94{
95 int n, j, nel;
96 ExpListSharedPtr plane_zero;
97
98 // note that nzplanes can be larger than nzmodes
100 m_session, graph1D, false, "DefaultVar", ImpType);
101
103
104 nel = m_planes[0]->GetExpSize();
105
106 for (j = 0; j < nel; ++j)
107 {
108 (*m_exp).push_back(m_planes[0]->GetExp(j));
109 }
110
111 for (n = 1; n < m_planes.size(); ++n)
112 {
113 m_planes[n] =
115 for (j = 0; j < nel; ++j)
116 {
117 (*m_exp).push_back((*m_exp)[j]);
118 }
119 }
120
121 SetCoeffPhys();
122}
123
124/**
125 * @param In ExpList2DHomogeneous1D object to copy.
126 */
129{
130 ExpListSharedPtr zero_plane =
131 std::dynamic_pointer_cast<ExpList>(In.m_planes[0]);
132
133 for (int n = 0; n < m_planes.size(); ++n)
134 {
135 m_planes[n] =
137 }
138
139 SetCoeffPhys();
140}
141
142/**
143 * Destructor
144 */
146{
147}
148
150{
151 int i, n, cnt;
152 int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
153 int npoints_per_plane = m_planes[0]->GetTotPoints();
154
155 int nzplanes = m_planes.size();
156
157 // Set total coefficients and points
158 m_ncoeffs = ncoeffs_per_plane * nzplanes;
159 m_npoints = npoints_per_plane * nzplanes;
160
163
164 int nel = m_planes[0]->GetExpSize();
165 m_coeff_offset = Array<OneD, int>(nel * nzplanes);
166 m_phys_offset = Array<OneD, int>(nel * nzplanes);
167 Array<OneD, NekDouble> tmparray;
168
169 for (cnt = n = 0; n < nzplanes; ++n)
170 {
171 m_planes[n]->SetCoeffsArray(tmparray =
172 m_coeffs + ncoeffs_per_plane * n);
173 m_planes[n]->SetPhysArray(tmparray = m_phys + npoints_per_plane * n);
174
175 for (i = 0; i < nel; ++i)
176 {
177 m_coeff_offset[cnt] =
178 m_planes[n]->GetCoeff_Offset(i) + n * ncoeffs_per_plane;
179 m_phys_offset[cnt++] =
180 m_planes[n]->GetPhys_Offset(i) + n * npoints_per_plane;
181 }
182 }
183}
184
189{
190 int n, coordim;
191 Array<OneD, NekDouble> tmp_xc, xhom;
192 int nyplanes = m_planes.size();
193 int npoints = GetTotPoints(eid);
194
195 switch (coordim = GetCoordim(0))
196 {
197 case 1:
198 {
199 (*m_exp)[eid]->GetCoords(xc0);
200 xhom = xc1;
201 }
202 break;
203 case 2:
204 ASSERTL0(xc1.size() != 0, "output coord_1 is not defined");
205 {
206 (*m_exp)[eid]->GetCoords(xc0, xc1);
207 xhom = xc2;
208 }
209 break;
210 default:
211 ASSERTL0(false, "Cannot have coordim being three dimensions"
212 "in a homogeneous field");
213 break;
214 }
215
216 // Fill homogeneous-direction
218 Array<OneD, NekDouble> z(nyplanes);
219
220 Array<OneD, NekDouble> local_pts(m_planes.size());
221
222 for (n = 0; n < m_planes.size(); n++)
223 {
224 local_pts[n] = pts[m_transposition->GetPlaneID(n)];
225 }
226
227 Vmath::Smul(nyplanes, m_lhom / 2.0, local_pts, 1, z, 1);
228 Vmath::Sadd(nyplanes, m_lhom / 2.0, z, 1, z, 1);
229
230 for (n = 0; n < nyplanes; ++n)
231 {
232 Vmath::Fill(npoints, z[n], tmp_xc = xhom + npoints * n, 1);
233 if (n)
234 {
235 Vmath::Vcopy(npoints, xc0, 1, tmp_xc = xc0 + npoints * n, 1);
236 if (coordim == 2)
237 {
238 Vmath::Vcopy(npoints, xc1, 1, tmp_xc = xc1 + npoints * n, 1);
239 }
240 }
241 }
242}
243
244/**
245 * The operation calls the 2D plane coordinates through the
246 * function ExpList#GetCoords and then evaluated the third
247 * coordinate using the member \a m_lhom
248 *
249 * @param coord_0 After calculation, the \f$x_1\f$ coordinate
250 * will be stored in this array.
251 *
252 * @param coord_1 After calculation, the \f$x_2\f$ coordinate
253 * will be stored in this array. This
254 * coordinate might be evaluated using the
255 * predefined value \a m_lhom
256 *
257 * @param coord_2 After calculation, the \f$x_3\f$ coordinate
258 * will be stored in this array. This
259 * coordinate is evaluated using the
260 * predefined value \a m_lhom
261 */
265{
266 int n, coordim;
267 Array<OneD, NekDouble> tmp_xc, xhom;
268 int nyplanes = m_planes.size();
269 int npoints = m_planes[0]->GetTotPoints();
270
271 m_planes[0]->GetCoords(xc0, xc1);
272
273 if ((coordim = GetCoordim(0)) == 1)
274 {
275 xhom = xc1;
276 }
277 else
278 {
279 xhom = xc2;
280 }
281
282 // Fill z-direction
284 Array<OneD, NekDouble> z(nyplanes);
285 Array<OneD, NekDouble> local_pts(m_planes.size());
286
287 for (n = 0; n < m_planes.size(); n++)
288 {
289 local_pts[n] = pts[m_transposition->GetPlaneID(n)];
290 }
291
292 Vmath::Smul(nyplanes, m_lhom / 2.0, local_pts, 1, z, 1);
293 Vmath::Sadd(nyplanes, m_lhom / 2.0, z, 1, z, 1);
294
295 for (n = 0; n < nyplanes; ++n)
296 {
297 Vmath::Fill(npoints, z[n], tmp_xc = xhom + npoints * n, 1);
298 if (n)
299 {
300 Vmath::Vcopy(npoints, xc0, 1, tmp_xc = xc0 + npoints * n, 1);
301 if (coordim == 2)
302 {
303 Vmath::Vcopy(npoints, xc1, 1, tmp_xc = xc1 + npoints * n, 1);
304 }
305 }
306 }
307}
308
309/**
310 * Write Tecplot Files Zone
311 * @param outfile Output file name.
312 * @param expansion Expansion that is considered
313 */
315 int expansion)
316{
317 int i, j;
318
319 int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
320 int nquad1 = m_planes.size();
321
322 Array<OneD, NekDouble> coords[3];
323
324 coords[0] = Array<OneD, NekDouble>(3 * nquad0 * nquad1);
325 coords[1] = coords[0] + nquad0 * nquad1;
326 coords[2] = coords[1] + nquad0 * nquad1;
327
328 GetCoords(expansion, coords[0], coords[1], coords[2]);
329
330 outfile << "Zone, I=" << nquad0 << ", J=" << nquad1 << ", F=Block"
331 << std::endl;
332
333 for (j = 0; j < GetCoordim(0) + 1; ++j)
334 {
335 for (i = 0; i < nquad0 * nquad1; ++i)
336 {
337 outfile << coords[j][i] << " ";
338 }
339 outfile << std::endl;
340 }
341}
342
344 int expansion,
345 [[maybe_unused]] int istrip)
346{
347 // If there is only one plane (e.g. HalfMode), we write a 2D plane.
348 if (m_planes.size() == 1)
349 {
350 m_planes[0]->WriteVtkPieceHeader(outfile, expansion);
351 return;
352 }
353
354 // If we are using Fourier points, output extra plane to fill domain
355 int outputExtraPlane = 0;
356 if (m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier &&
357 m_homogeneousBasis->GetPointsType() ==
359 {
360 outputExtraPlane = 1;
361 }
362 int i, j;
363 int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
364 int nquad1 = m_planes.size() + outputExtraPlane;
365 int ntot = nquad0 * nquad1;
366 int ntotminus = (nquad0 - 1) * (nquad1 - 1);
367
368 Array<OneD, NekDouble> coords[3];
369 coords[0] = Array<OneD, NekDouble>(ntot);
370 coords[1] = Array<OneD, NekDouble>(ntot);
371 coords[2] = Array<OneD, NekDouble>(ntot);
372 GetCoords(expansion, coords[0], coords[1], coords[2]);
373
374 if (outputExtraPlane)
375 {
376 // Copy coords[0] and coords[1] to extra plane
378 Vmath::Vcopy(nquad0, coords[0], 1,
379 tmp = coords[0] + (nquad1 - 1) * nquad0, 1);
380 Vmath::Vcopy(nquad0, coords[1], 1,
381 tmp = coords[1] + (nquad1 - 1) * nquad0, 1);
382 // Fill coords[2] for extra plane
383 NekDouble z = coords[2][nquad0 * m_planes.size() - 1] +
384 (coords[2][nquad0] - coords[2][0]);
385 Vmath::Fill(nquad0, z, tmp = coords[2] + (nquad1 - 1) * nquad0, 1);
386 }
387
388 outfile << " <Piece NumberOfPoints=\"" << ntot << "\" NumberOfCells=\""
389 << ntotminus << "\">" << endl;
390 outfile << " <Points>" << endl;
391 outfile << " <DataArray type=\"Float64\" "
392 << R"(NumberOfComponents="3" format="ascii">)" << endl;
393 outfile << " ";
394 for (i = 0; i < ntot; ++i)
395 {
396 for (j = 0; j < 3; ++j)
397 {
398 outfile << coords[j][i] << " ";
399 }
400 outfile << endl;
401 }
402 outfile << endl;
403 outfile << " </DataArray>" << endl;
404 outfile << " </Points>" << endl;
405 outfile << " <Cells>" << endl;
406 outfile << " <DataArray type=\"Int32\" "
407 << R"(Name="connectivity" format="ascii">)" << endl;
408 for (i = 0; i < nquad0 - 1; ++i)
409 {
410 for (j = 0; j < nquad1 - 1; ++j)
411 {
412 outfile << j * nquad0 + i << " " << j * nquad0 + i + 1 << " "
413 << (j + 1) * nquad0 + i + 1 << " " << (j + 1) * nquad0 + i
414 << endl;
415 }
416 }
417 outfile << endl;
418 outfile << " </DataArray>" << endl;
419 outfile << " <DataArray type=\"Int32\" "
420 << R"(Name="offsets" format="ascii">)" << endl;
421 for (i = 0; i < ntotminus; ++i)
422 {
423 outfile << i * 4 + 4 << " ";
424 }
425 outfile << endl;
426 outfile << " </DataArray>" << endl;
427 outfile << " <DataArray type=\"UInt8\" "
428 << R"(Name="types" format="ascii">)" << endl;
429 for (i = 0; i < ntotminus; ++i)
430 {
431 outfile << "9 ";
432 }
433 outfile << endl;
434 outfile << " </DataArray>" << endl;
435 outfile << " </Cells>" << endl;
436 outfile << " <PointData>" << endl;
437}
438
441{
442 int nPlanes = m_planes.size();
443 int nPtsPlane = m_planes[0]->GetNpoints();
444 int nDim = GetCoordim(0) + 1;
445
446 ASSERTL1(normals.size() >= nDim,
447 "Output vector does not have sufficient dimensions to"
448 "match coordim");
449 ASSERTL1(normals[0].size() >= nPtsPlane,
450 "Output vector does not have sufficient dimensions to"
451 "match coordim");
452
453 // Calculate normals from plane 0.
454 m_planes[0]->GetNormals(normals);
455
456 // Copy remaining planes.
457 for (int i = 0; i < nDim; ++i)
458 {
459 for (int n = 1; n < nPlanes; ++n)
460 {
461 Vmath::Vcopy(nPtsPlane, &normals[i][0], 1,
462 &normals[i][n * nPtsPlane], 1);
463 }
464 }
465}
466} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
Describes the specification for a Basis.
Definition: Basis.h:45
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 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...
void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip) override
void v_WriteTecplotZone(std::ostream &outfile, int expansion) override
void v_GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals) override
Populate normals with the normals of all expansions.
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.
Definition: ExpList.h:1083
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1123
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 .
Definition: ExpList.h:1770
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1056
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1118
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1063
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2070
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1058
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1125
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1099
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1557
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1903
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:74
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::vector< double > z(NPUPPER)
double NekDouble
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.hpp:194
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
STL namespace.