Nektar++
ExpList3DHomogeneous1D.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ExpList3DHomogeneous1D.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 ExpList 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
50 const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom,
51 const bool useFFT, const bool dealiasing)
52 : ExpListHomogeneous1D(e3DH1D, pSession, HomoBasis, lhom, useFFT,
53 dealiasing)
54{
55}
56
57// Constructor for ExpList3DHomogeneous1D to act as a Explist field
60 const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom,
61 const bool useFFT, const bool dealiasing,
62 const SpatialDomains::MeshGraphSharedPtr &graph2D, const std::string &var,
64 : ExpListHomogeneous1D(e3DH1D, pSession, HomoBasis, lhom, useFFT,
65 dealiasing)
66{
67 GenExpList3DHomogeneous1D(graph2D->GetExpansionInfo(var), ImpType);
69 m_graph = graph2D;
70}
71
72// Constructor for ExpList3DHomogeneous1D to act as a Explist field
75 const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom,
76 const bool useFFT, const bool dealiasing,
77 const SpatialDomains::ExpansionInfoMap &expansions,
79 : ExpListHomogeneous1D(e3DH1D, pSession, HomoBasis, lhom, useFFT,
80 dealiasing)
81{
82 GenExpList3DHomogeneous1D(expansions, ImpType);
83}
84
86 const SpatialDomains::ExpansionInfoMap &expansions,
88{
89 int n, j, nel;
90 bool False = false;
91 ExpListSharedPtr plane_zero;
92
93 // note that nzplanes can be larger than nzmodes
95 m_session, expansions, false, ImpType);
96
98 nel = m_planes[0]->GetExpSize();
99
100 for (j = 0; j < nel; ++j)
101 {
102 (*m_exp).push_back(m_planes[0]->GetExp(j));
103 }
104
105 for (n = 1; n < m_planes.size(); ++n)
106 {
107 m_planes[n] =
109 for (j = 0; j < nel; ++j)
110 {
111 (*m_exp).push_back((*m_exp)[j]);
112 }
113 }
114
115 SetCoeffPhys();
116}
117
118/**
119 * @param In ExpList3DHomogeneous1D object to copy.
120 */
122 bool DeclarePlanesSetCoeffPhys)
124{
125 if (DeclarePlanesSetCoeffPhys)
126 {
127 bool False = false;
128 ExpListSharedPtr zero_plane =
129 std::dynamic_pointer_cast<ExpList>(In.m_planes[0]);
130
131 for (int n = 0; n < m_planes.size(); ++n)
132 {
133 m_planes[n] =
135 }
136
137 SetCoeffPhys();
138 }
139}
140
141/**
142 * @param In ExpList3DHomogeneous1D object to copy.
143 * @param eIDs Id of elements that should be copied.
144 */
146 const ExpList3DHomogeneous1D &In, const std::vector<unsigned int> &eIDs,
147 bool DeclarePlanesSetCoeffPhys,
149 : ExpListHomogeneous1D(In, eIDs, ImpType)
150{
151 if (DeclarePlanesSetCoeffPhys)
152 {
153 bool False = false;
154 std::vector<unsigned int> eIDsPlane;
155 int nel = eIDs.size() / m_planes.size();
156
157 for (int i = 0; i < nel; ++i)
158 {
159 eIDsPlane.push_back(eIDs[i]);
160 }
161
162 ExpListSharedPtr zero_plane_old =
163 std::dynamic_pointer_cast<ExpList>(In.m_planes[0]);
164
166 *(zero_plane_old), eIDsPlane, DeclarePlanesSetCoeffPhys, ImpType);
167
168 for (int n = 0; n < m_planes.size(); ++n)
169 {
170 m_planes[n] =
172 }
173
174 SetCoeffPhys();
175 }
176}
177
178/**
179 * Destructor
180 */
182{
183}
184
186{
187 int i, n, cnt;
188 int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
189 int npoints_per_plane = m_planes[0]->GetTotPoints();
190
191 int nzplanes = m_planes.size();
192
193 // Set total coefficients and points
194 m_ncoeffs = ncoeffs_per_plane * nzplanes;
195 m_npoints = npoints_per_plane * nzplanes;
196
199
200 int nel = m_planes[0]->GetExpSize();
201 m_coeff_offset = Array<OneD, int>(nel * nzplanes);
202 m_phys_offset = Array<OneD, int>(nel * nzplanes);
203 Array<OneD, NekDouble> tmparray;
204
205 for (cnt = n = 0; n < nzplanes; ++n)
206 {
207 m_planes[n]->SetCoeffsArray(tmparray =
208 m_coeffs + ncoeffs_per_plane * n);
209 m_planes[n]->SetPhysArray(tmparray = m_phys + npoints_per_plane * n);
210
211 for (i = 0; i < nel; ++i)
212 {
213 m_coeff_offset[cnt] =
214 m_planes[n]->GetCoeff_Offset(i) + n * ncoeffs_per_plane;
215 m_phys_offset[cnt++] =
216 m_planes[n]->GetPhys_Offset(i) + n * npoints_per_plane;
217 }
218 }
219}
220
221// @TODO: Could potentially move the extra plane coordinate output here
226{
227 int n;
229 int nzplanes = m_planes.size();
230 int npoints = GetTotPoints(eid);
231
232 (*m_exp)[eid]->GetCoords(xc0, xc1);
233
234 // Fill z-direction
236 Array<OneD, NekDouble> local_pts(m_planes.size());
237
238 for (n = 0; n < m_planes.size(); n++)
239 {
240 local_pts[n] = pts[m_transposition->GetPlaneID(n)];
241 }
242
243 Array<OneD, NekDouble> z(nzplanes);
244
245 Vmath::Smul(nzplanes, m_lhom / 2.0, local_pts, 1, z, 1);
246 Vmath::Sadd(nzplanes, m_lhom / 2.0, z, 1, z, 1);
247
248 for (n = 0; n < nzplanes; ++n)
249 {
250 Vmath::Fill(npoints, z[n], tmp_xc = xc2 + npoints * n, 1);
251 if (n)
252 {
253 Vmath::Vcopy(npoints, xc0, 1, tmp_xc = xc0 + npoints * n, 1);
254 Vmath::Vcopy(npoints, xc1, 1, tmp_xc = xc1 + npoints * n, 1);
255 }
256 }
257}
258
259/**
260 * The operation calls the 2D plane coordinates through the
261 * function ExpList#GetCoords and then evaluated the third
262 * coordinate using the member \a m_lhom
263 *
264 * @param coord_0 After calculation, the \f$x_1\f$ coordinate
265 * will be stored in this array.
266 *
267 * @param coord_1 After calculation, the \f$x_2\f$ coordinate
268 * will be stored in this array.
269 *
270 * @param coord_2 After calculation, the \f$x_3\f$ coordinate
271 * will be stored in this array. This
272 * coordinate is evaluated using the
273 * predefined value \a m_lhom
274 */
278{
279 int n;
281 int nzplanes = m_planes.size();
282 int npoints = m_planes[0]->GetTotPoints();
283
284 m_planes[0]->GetCoords(xc0, xc1);
285
286 // Fill z-direction
288
289 Array<OneD, NekDouble> local_pts(m_planes.size());
290
291 for (n = 0; n < m_planes.size(); n++)
292 {
293 local_pts[n] = pts[m_transposition->GetPlaneID(n)];
294 }
295
296 Array<OneD, NekDouble> z(nzplanes);
297
298 Vmath::Smul(nzplanes, m_lhom / 2.0, local_pts, 1, z, 1);
299 Vmath::Sadd(nzplanes, m_lhom / 2.0, z, 1, z, 1);
300
301 for (n = 0; n < nzplanes; ++n)
302 {
303 Vmath::Fill(npoints, z[n], tmp_xc = xc2 + npoints * n, 1);
304 if (n)
305 {
306 Vmath::Vcopy(npoints, xc0, 1, tmp_xc = xc0 + npoints * n, 1);
307 Vmath::Vcopy(npoints, xc1, 1, tmp_xc = xc1 + npoints * n, 1);
308 }
309 }
310}
311
313 int expansion)
314{
315 ASSERTL0(expansion == -1,
316 "Multi-zone output not supported for homogeneous expansions.");
317
318 const int nPtsPlane = m_planes[0]->GetNpoints();
319 const int nElmt = m_planes[0]->GetExpSize();
320 const int nPlanes = m_planes.size();
321
322 for (int i = 0, cnt = 0; i < nElmt; ++i)
323 {
324 const int np0 = (*m_exp)[i]->GetNumPoints(0);
325 const int np1 = (*m_exp)[i]->GetNumPoints(1);
326
327 for (int n = 1; n < nPlanes; ++n)
328 {
329 const int o1 = (n - 1) * nPtsPlane;
330 const int o2 = n * nPtsPlane;
331 for (int j = 1; j < np1; ++j)
332 {
333 for (int k = 1; k < np0; ++k)
334 {
335 outfile << cnt + (j - 1) * np0 + (k - 1) + o1 + 1 << " ";
336 outfile << cnt + (j - 1) * np0 + (k - 1) + o2 + 1 << " ";
337 outfile << cnt + (j - 1) * np0 + k + o2 + 1 << " ";
338 outfile << cnt + (j - 1) * np0 + k + o1 + 1 << " ";
339 outfile << cnt + j * np0 + (k - 1) + o1 + 1 << " ";
340 outfile << cnt + j * np0 + (k - 1) + o2 + 1 << " ";
341 outfile << cnt + j * np0 + k + o2 + 1 << " ";
342 outfile << cnt + j * np0 + k + o1 + 1 << endl;
343 }
344 }
345 }
346
347 cnt += np0 * np1;
348 }
349}
350
352 int expansion, int istrip)
353{
354 // If there is only one plane (e.g. HalfMode), we write a 2D plane.
355 if (m_planes.size() == 1)
356 {
357 m_planes[0]->WriteVtkPieceHeader(outfile, expansion);
358 return;
359 }
360
361 // If we are using Fourier points, output extra plane to fill domain
362 int outputExtraPlane = 0;
363 if (m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier &&
364 m_homogeneousBasis->GetPointsType() ==
366 {
367 outputExtraPlane = 1;
368 }
369
370 int i, j, k;
371 int nq0 = (*m_exp)[expansion]->GetNumPoints(0);
372 int nq1 = (*m_exp)[expansion]->GetNumPoints(1);
373 int nq2 = m_planes.size() + outputExtraPlane;
374 int ntot = nq0 * nq1 * nq2;
375 int ntotminus = (nq0 - 1) * (nq1 - 1) * (nq2 - 1);
376
377 Array<OneD, NekDouble> coords[3];
378 coords[0] = Array<OneD, NekDouble>(ntot);
379 coords[1] = Array<OneD, NekDouble>(ntot);
380 coords[2] = Array<OneD, NekDouble>(ntot);
381 v_GetCoords(expansion, coords[0], coords[1], coords[2]);
382
383 if (outputExtraPlane)
384 {
385 // Copy coords[0] and coords[1] to extra plane
387 Vmath::Vcopy(nq0 * nq1, coords[0], 1,
388 tmp = coords[0] + (nq2 - 1) * nq0 * nq1, 1);
389 Vmath::Vcopy(nq0 * nq1, coords[1], 1,
390 tmp = coords[1] + (nq2 - 1) * nq0 * nq1, 1);
391 // Fill coords[2] for extra plane
392 NekDouble z = coords[2][nq0 * nq1 * m_planes.size() - 1] +
393 (coords[2][nq0 * nq1] - coords[2][0]);
394 Vmath::Fill(nq0 * nq1, z, tmp = coords[2] + (nq2 - 1) * nq0 * nq1, 1);
395 }
396
397 NekDouble DistStrip;
398 m_session->LoadParameter("DistStrip", DistStrip, 0);
399 // Reset the z-coords for homostrips
400 for (int i = 0; i < ntot; i++)
401 {
402 coords[2][i] += istrip * DistStrip;
403 }
404
405 outfile << " <Piece NumberOfPoints=\"" << ntot << "\" NumberOfCells=\""
406 << ntotminus << "\">" << endl;
407 outfile << " <Points>" << endl;
408 outfile << " <DataArray type=\"Float64\" "
409 << R"(NumberOfComponents="3" format="ascii">)" << endl;
410 outfile << " ";
411 for (i = 0; i < ntot; ++i)
412 {
413 for (j = 0; j < 3; ++j)
414 {
415 outfile << coords[j][i] << " ";
416 }
417 outfile << endl;
418 }
419 outfile << endl;
420 outfile << " </DataArray>" << endl;
421 outfile << " </Points>" << endl;
422 outfile << " <Cells>" << endl;
423 outfile << " <DataArray type=\"Int32\" "
424 << R"(Name="connectivity" format="ascii">)" << endl;
425 for (i = 0; i < nq0 - 1; ++i)
426 {
427 for (j = 0; j < nq1 - 1; ++j)
428 {
429 for (k = 0; k < nq2 - 1; ++k)
430 {
431 outfile << k * nq0 * nq1 + j * nq0 + i << " "
432 << k * nq0 * nq1 + j * nq0 + i + 1 << " "
433 << k * nq0 * nq1 + (j + 1) * nq0 + i + 1 << " "
434 << k * nq0 * nq1 + (j + 1) * nq0 + i << " "
435 << (k + 1) * nq0 * nq1 + j * nq0 + i << " "
436 << (k + 1) * nq0 * nq1 + j * nq0 + i + 1 << " "
437 << (k + 1) * nq0 * nq1 + (j + 1) * nq0 + i + 1 << " "
438 << (k + 1) * nq0 * nq1 + (j + 1) * nq0 + i << endl;
439 }
440 }
441 }
442 outfile << endl;
443 outfile << " </DataArray>" << endl;
444 outfile << " <DataArray type=\"Int32\" "
445 << R"(Name="offsets" format="ascii">)" << endl;
446 for (i = 0; i < ntotminus; ++i)
447 {
448 outfile << i * 8 + 8 << " ";
449 }
450 outfile << endl;
451 outfile << " </DataArray>" << endl;
452 outfile << " <DataArray type=\"UInt8\" "
453 << R"(Name="types" format="ascii">)" << endl;
454 for (i = 0; i < ntotminus; ++i)
455 {
456 outfile << "12 ";
457 }
458 outfile << endl;
459 outfile << " </DataArray>" << endl;
460 outfile << " </Cells>" << endl;
461 outfile << " <PointData>" << endl;
462}
463
465 const Array<OneD, const NekDouble> &inarray,
467{
468 int cnt = 0;
469 NekDouble errL2, err = 0.0;
471 Array<OneD, NekDouble> local_w(m_planes.size());
472
473 for (int n = 0; n < m_planes.size(); n++)
474 {
475 local_w[n] = w[m_transposition->GetPlaneID(n)];
476 }
477
478 if (soln == NullNekDouble1DArray)
479 {
480 for (int n = 0; n < m_planes.size(); ++n)
481 {
482 errL2 = m_planes[n]->L2(inarray + cnt);
483 cnt += m_planes[n]->GetTotPoints();
484 err += errL2 * errL2 * local_w[n] * m_lhom * 0.5;
485 }
486 }
487 else
488 {
489 for (int n = 0; n < m_planes.size(); ++n)
490 {
491 errL2 = m_planes[n]->L2(inarray + cnt, soln + cnt);
492 cnt += m_planes[n]->GetTotPoints();
493 err += errL2 * errL2 * local_w[n] * m_lhom * 0.5;
494 }
495 }
496 m_comm->GetColumnComm()->AllReduce(err, LibUtilities::ReduceSum);
497
498 return sqrt(err);
499}
500
502{
503 Array<OneD, NekDouble> energy(m_planes.size() / 2);
504 NekDouble area = 0.0;
505
506 // Calculate total area of elements.
507 for (int n = 0; n < m_planes[0]->GetExpSize(); ++n)
508 {
510 1.0);
511 area += m_planes[0]->GetExp(n)->Integral(inarray);
512 }
513
514 m_comm->GetRowComm()->AllReduce(area, LibUtilities::ReduceSum);
515
516 // Calculate L2 norm of real/imaginary planes.
517 for (int n = 0; n < m_planes.size(); n += 2)
518 {
519 NekDouble err;
520
521 energy[n / 2] = 0;
522
523 for (int i = 0; i < m_planes[n]->GetExpSize(); ++i)
524 {
525 LocalRegions::ExpansionSharedPtr exp = m_planes[n]->GetExp(i);
526 Array<OneD, NekDouble> phys(exp->GetTotPoints());
527 exp->BwdTrans(m_planes[n]->GetCoeffs() +
529 phys);
530 err = exp->L2(phys);
531 energy[n / 2] += err * err;
532
533 exp = m_planes[n + 1]->GetExp(i);
534 exp->BwdTrans(m_planes[n + 1]->GetCoeffs() +
535 m_planes[n + 1]->GetCoeff_Offset(i),
536 phys);
537 err = exp->L2(phys);
538 energy[n / 2] += err * err;
539 }
540
541 m_comm->GetRowComm()->AllReduce(energy[n / 2], LibUtilities::ReduceSum);
542 energy[n / 2] /= 2.0 * area;
543 }
544
545 return energy;
546}
547} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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 GenExpList3DHomogeneous1D(const SpatialDomains::ExpansionInfoMap &expansions, const Collections::ImplementationType ImpType)
NekDouble v_L2(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray) override
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...
Array< OneD, const NekDouble > v_HomogeneousEnergy(void) override
void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip) override
void v_WriteTecplotConnectivity(std::ostream &outfile, int expansion) override
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:1080
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:1952
int GetCoeff_Offset(int n) const
Get the start offset position for a local contiguous list of coeffs correspoinding to element n.
Definition: ExpList.h:2087
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1124
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1053
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1115
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1060
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2079
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:1057
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1055
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1126
ExpansionType m_expType
Exapnsion type.
Definition: ExpList.h:1048
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1096
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1557
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:74
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition: MeshGraph.h:141
std::vector< double > w(NPUPPER)
std::vector< double > z(NPUPPER)
static Array< OneD, NekDouble > NullNekDouble1DArray
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
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294