Nektar++
ExpListHomogeneous2D.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ExpListHomogeneous2D.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 2-directions
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <LibUtilities/Foundations/ManagerAccess.h> // for PointsManager, etc
40
41using namespace std;
42
44{
45// Forward declaration for typedefs
47 : ExpList(type), m_homogeneousBasis_y(LibUtilities::NullBasisSharedPtr),
48 m_homogeneousBasis_z(LibUtilities::NullBasisSharedPtr), m_lhom_y(1),
49 m_lhom_z(1), m_homogeneous2DBlockMat(
50 MemoryManager<Homo2DBlockMatrixMap>::AllocateSharedPtr())
51{
52}
53
55 const ExpansionType type,
57 const LibUtilities::BasisKey &HomoBasis_y,
58 const LibUtilities::BasisKey &HomoBasis_z, const NekDouble lhom_y,
59 const NekDouble lhom_z, const bool useFFT, const bool dealiasing)
60 : ExpList(type), m_useFFT(useFFT), m_lhom_y(lhom_y), m_lhom_z(lhom_z),
61 m_homogeneous2DBlockMat(
62 MemoryManager<Homo2DBlockMatrixMap>::AllocateSharedPtr()),
63 m_dealiasing(dealiasing)
64{
65 m_session = pSession;
66 m_comm = pSession->GetComm();
67
69 "Homogeneous Basis in y direction is a null basis");
71 "Homogeneous Basis in z direction is a null basis");
72
75
78 HomoBasis_y, HomoBasis_z, m_comm->GetColumnComm());
79
80 m_Ycomm = m_comm->GetColumnComm()->GetRowComm();
81 m_Zcomm = m_comm->GetColumnComm()->GetRowComm();
82
83 m_ny = m_homogeneousBasis_y->GetNumPoints() / m_Ycomm->GetSize();
84 m_nz = m_homogeneousBasis_z->GetNumPoints() / m_Zcomm->GetSize();
85
87
88 if (m_useFFT)
89 {
90 m_FFT_y =
92 m_FFT_z =
94 }
95
96 if (m_dealiasing)
97 {
98 ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,
99 "Remove dealiasing if you want to run in parallel");
101 }
102}
103
104/**
105 * @param In ExpListHomogeneous2D object to copy.
106 */
108 : ExpList(In, false), m_useFFT(In.m_useFFT), m_FFT_y(In.m_FFT_y),
109 m_FFT_z(In.m_FFT_z), m_transposition(In.m_transposition),
110 m_Ycomm(In.m_Ycomm), m_Zcomm(In.m_Ycomm),
111 m_homogeneousBasis_y(In.m_homogeneousBasis_y),
112 m_homogeneousBasis_z(In.m_homogeneousBasis_z), m_lhom_y(In.m_lhom_y),
113 m_lhom_z(In.m_lhom_z),
114 m_homogeneous2DBlockMat(In.m_homogeneous2DBlockMat), m_ny(In.m_ny),
115 m_nz(In.m_nz), m_dealiasing(In.m_dealiasing), m_padsize_y(In.m_padsize_y),
116 m_padsize_z(In.m_padsize_z), MatFwdPAD(In.MatFwdPAD),
117 MatBwdPAD(In.MatBwdPAD)
118{
120}
121
123 const ExpListHomogeneous2D &In, const std::vector<unsigned int> &eIDs)
124 : ExpList(In, eIDs, false), m_useFFT(In.m_useFFT), m_FFT_y(In.m_FFT_y),
125 m_FFT_z(In.m_FFT_z), m_transposition(In.m_transposition),
126 m_Ycomm(In.m_Ycomm), m_Zcomm(In.m_Ycomm),
127 m_homogeneousBasis_y(In.m_homogeneousBasis_y),
128 m_homogeneousBasis_z(In.m_homogeneousBasis_z), m_lhom_y(In.m_lhom_y),
129 m_lhom_z(In.m_lhom_z),
130 m_homogeneous2DBlockMat(
131 MemoryManager<Homo2DBlockMatrixMap>::AllocateSharedPtr()),
132 m_ny(In.m_ny), m_nz(In.m_nz), m_dealiasing(In.m_dealiasing),
133 m_padsize_y(In.m_padsize_y), m_padsize_z(In.m_padsize_z),
134 MatFwdPAD(In.MatFwdPAD), MatBwdPAD(In.MatBwdPAD)
135{
137}
138
139/**
140 * Destructor
141 */
143{
144}
145
147 const int npts, const Array<OneD, const NekDouble> &inarray,
148 Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
149{
150 // Forwards trans
151 Homogeneous2DTrans(npts, inarray, outarray, true, Shuff, UnShuff);
152}
153
155 const int npts, const Array<OneD, const NekDouble> &inarray,
156 Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
157{
158 // Backwards trans
159 Homogeneous2DTrans(npts, inarray, outarray, false, Shuff, UnShuff);
160}
161
163 const int npoints, const Array<OneD, NekDouble> &inarray1,
164 const Array<OneD, NekDouble> &inarray2, Array<OneD, NekDouble> &outarray)
165{
166 // int npoints = outarray.size(); // number of total physical points
167 int nlines =
168 m_lines.size(); // number of lines == number of Fourier modes = number
169 // of Fourier coeff = number of points per slab
170 int nslabs = npoints /
171 nlines; // number of slabs = numebr of physical points per line
172
173 Array<OneD, NekDouble> V1(npoints);
174 Array<OneD, NekDouble> V2(npoints);
175 Array<OneD, NekDouble> V1V2(npoints);
176 Array<OneD, NekDouble> ShufV1(npoints);
177 Array<OneD, NekDouble> ShufV2(npoints);
178 Array<OneD, NekDouble> ShufV1V2(npoints);
179
180 if (m_WaveSpace)
181 {
182 V1 = inarray1;
183 V2 = inarray2;
184 }
185 else
186 {
187 HomogeneousFwdTrans(npoints, inarray1, V1);
188 HomogeneousFwdTrans(npoints, inarray2, V2);
189 }
190
191 m_transposition->Transpose(npoints, V1, ShufV1, false,
193 m_transposition->Transpose(npoints, V2, ShufV2, false,
195
196 Array<OneD, NekDouble> PadV1_slab_coeff(m_padsize_y * m_padsize_z, 0.0);
197 Array<OneD, NekDouble> PadV2_slab_coeff(m_padsize_y * m_padsize_z, 0.0);
198 Array<OneD, NekDouble> PadRe_slab_coeff(m_padsize_y * m_padsize_z, 0.0);
199
200 Array<OneD, NekDouble> PadV1_slab_phys(m_padsize_y * m_padsize_z, 0.0);
201 Array<OneD, NekDouble> PadV2_slab_phys(m_padsize_y * m_padsize_z, 0.0);
202 Array<OneD, NekDouble> PadRe_slab_phys(m_padsize_y * m_padsize_z, 0.0);
203
204 NekVector<NekDouble> PadIN_V1(m_padsize_y * m_padsize_z, PadV1_slab_coeff,
205 eWrapper);
206 NekVector<NekDouble> PadOUT_V1(m_padsize_y * m_padsize_z, PadV1_slab_phys,
207 eWrapper);
208
209 NekVector<NekDouble> PadIN_V2(m_padsize_y * m_padsize_z, PadV2_slab_coeff,
210 eWrapper);
211 NekVector<NekDouble> PadOUT_V2(m_padsize_y * m_padsize_z, PadV2_slab_phys,
212 eWrapper);
213
214 NekVector<NekDouble> PadIN_Re(m_padsize_y * m_padsize_z, PadRe_slab_phys,
215 eWrapper);
216 NekVector<NekDouble> PadOUT_Re(m_padsize_y * m_padsize_z, PadRe_slab_coeff,
217 eWrapper);
218
219 // Looping on the slabs
220 for (int j = 0; j < nslabs; j++)
221 {
222 // Copying the j-th slab of size N*M into a bigger slab of lenght 2*N*M
223 // We are in Fourier space
224 for (int i = 0; i < m_nz; i++)
225 {
226 Vmath::Vcopy(m_ny, &(ShufV1[i * m_ny + j * nlines]), 1,
227 &(PadV1_slab_coeff[i * 2 * m_ny]), 1);
228 Vmath::Vcopy(m_ny, &(ShufV2[i * m_ny + j * nlines]), 1,
229 &(PadV2_slab_coeff[i * 2 * m_ny]), 1);
230 }
231
232 // Moving to physical space using the padded system
233 PadOUT_V1 = (*MatBwdPAD) * PadIN_V1;
234 PadOUT_V2 = (*MatBwdPAD) * PadIN_V2;
235
236 // Perfroming the vectors multiplication in physical
237 // space on the padded system
238 Vmath::Vmul(m_padsize_y * m_padsize_z, PadV1_slab_phys, 1,
239 PadV2_slab_phys, 1, PadRe_slab_phys, 1);
240
241 // Moving back the result (V1*V2)_phys in Fourier
242 // space, padded system
243 PadOUT_Re = (*MatFwdPAD) * PadIN_Re;
244
245 // Copying the first half of the padded pencil in the
246 // full vector (Fourier space)
247 for (int i = 0; i < m_nz; i++)
248 {
249 Vmath::Vcopy(m_ny, &(PadRe_slab_coeff[i * 2 * m_ny]), 1,
250 &(ShufV1V2[i * m_ny + j * nlines]), 1);
251 }
252 }
253
254 if (m_WaveSpace)
255 {
256 m_transposition->Transpose(npoints, ShufV1V2, outarray, false,
258 }
259 else
260 {
261 m_transposition->Transpose(npoints, ShufV1V2, V1V2, false,
263
264 // Moving the results in physical space for the output
265 HomogeneousBwdTrans(npoints, V1V2, outarray);
266 }
267}
268
270 const int npts, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
271 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
273{
274 // TODO Proper implementation of this
275 size_t ndim = inarray1.size();
276 ASSERTL1(inarray2.size() % ndim == 0,
277 "Wrong dimensions for DealiasedDotProd.");
278
279 size_t nvec = inarray2.size() % ndim;
280
281 Array<OneD, NekDouble> out(npts);
282 for (size_t i = 0; i < nvec; i++)
283 {
284 Vmath::Zero(npts, outarray[i], 1);
285 for (size_t j = 0; j < ndim; j++)
286 {
287 DealiasedProd(npts, inarray1[j], inarray2[i * ndim + j], out);
288 Vmath::Vadd(npts, outarray[i], 1, out, 1, outarray[i], 1);
289 }
290 }
291}
292
294 const Array<OneD, const NekDouble> &inarray,
295 Array<OneD, NekDouble> &outarray)
296{
297 size_t cnt = 0, cnt1 = 0;
298 Array<OneD, NekDouble> tmparray;
299 size_t nlines = m_lines.size();
300
301 for (size_t n = 0; n < nlines; ++n)
302 {
303 m_lines[n]->FwdTrans(inarray + cnt, tmparray = outarray + cnt1);
304 cnt += m_lines[n]->GetTotPoints();
305 cnt1 += m_lines[n]->GetNcoeffs();
306 }
307 if (!m_WaveSpace)
308 {
309 HomogeneousFwdTrans(cnt1, outarray, outarray);
310 }
311}
312
314 const Array<OneD, const NekDouble> &inarray,
315 Array<OneD, NekDouble> &outarray)
316{
317 size_t cnt = 0, cnt1 = 0;
318 Array<OneD, NekDouble> tmparray;
319 size_t nlines = m_lines.size();
320
321 for (size_t n = 0; n < nlines; ++n)
322 {
323 m_lines[n]->FwdTransLocalElmt(inarray + cnt,
324 tmparray = outarray + cnt1);
325
326 cnt += m_lines[n]->GetTotPoints();
327 cnt1 += m_lines[n]->GetNcoeffs();
328 }
329 if (!m_WaveSpace)
330 {
331 HomogeneousFwdTrans(cnt1, outarray, outarray);
332 }
333}
334
336 const Array<OneD, const NekDouble> &inarray,
337 Array<OneD, NekDouble> &outarray)
338{
339 int cnt = 0, cnt1 = 0;
340 Array<OneD, NekDouble> tmparray;
341 int nlines = m_lines.size();
342
343 for (int n = 0; n < nlines; ++n)
344 {
345 m_lines[n]->FwdTransBndConstrained(inarray + cnt,
346 tmparray = outarray + cnt1);
347
348 cnt += m_lines[n]->GetTotPoints();
349 cnt1 += m_lines[n]->GetNcoeffs();
350 }
351 if (!m_WaveSpace)
352 {
353 HomogeneousFwdTrans(cnt1, outarray, outarray);
354 }
355}
356
358 const Array<OneD, const NekDouble> &inarray,
359 Array<OneD, NekDouble> &outarray)
360{
361 size_t cnt = 0, cnt1 = 0;
362 Array<OneD, NekDouble> tmparray;
363 size_t nlines = m_lines.size();
364
365 for (size_t n = 0; n < nlines; ++n)
366 {
367 m_lines[n]->BwdTrans(inarray + cnt, tmparray = outarray + cnt1);
368 cnt += m_lines[n]->GetNcoeffs();
369 cnt1 += m_lines[n]->GetTotPoints();
370 }
371 if (!m_WaveSpace)
372 {
373 HomogeneousBwdTrans(cnt1, outarray, outarray);
374 }
375}
376
378 const Array<OneD, const NekDouble> &inarray,
379 Array<OneD, NekDouble> &outarray)
380{
381 size_t cnt = 0, cnt1 = 0;
382 Array<OneD, NekDouble> tmparray;
383 size_t nlines = m_lines.size();
384
385 for (size_t n = 0; n < nlines; ++n)
386 {
387 m_lines[n]->IProductWRTBase(inarray + cnt, tmparray = outarray + cnt1);
388
389 cnt += m_lines[n]->GetNcoeffs();
390 cnt1 += m_lines[n]->GetTotPoints();
391 }
392}
393
395 const int num_dofs, const Array<OneD, const NekDouble> &inarray,
396 Array<OneD, NekDouble> &outarray, bool IsForwards,
397 [[maybe_unused]] bool Shuff, [[maybe_unused]] bool UnShuff)
398{
399 if (m_useFFT)
400 {
401
402 int n = m_lines.size(); // number of Fourier points in the Fourier
403 // directions (x-z grid)
404 int p =
405 num_dofs /
406 n; // number of points per line = n of Fourier transform required
407
408 Array<OneD, NekDouble> fft_in(num_dofs);
409 Array<OneD, NekDouble> fft_out(num_dofs);
410
411 m_transposition->Transpose(num_dofs, inarray, fft_in, false,
413
414 if (IsForwards)
415 {
416 for (int i = 0; i < (p * m_nz); i++)
417 {
418 m_FFT_y->FFTFwdTrans(m_tmpIN = fft_in + i * m_ny,
419 m_tmpOUT = fft_out + i * m_ny);
420 }
421 }
422 else
423 {
424 for (int i = 0; i < (p * m_nz); i++)
425 {
426 m_FFT_y->FFTBwdTrans(m_tmpIN = fft_in + i * m_ny,
427 m_tmpOUT = fft_out + i * m_ny);
428 }
429 }
430
431 m_transposition->Transpose(num_dofs, fft_out, fft_in, false,
433
434 if (IsForwards)
435 {
436 for (int i = 0; i < (p * m_ny); i++)
437 {
438 m_FFT_z->FFTFwdTrans(m_tmpIN = fft_in + i * m_nz,
439 m_tmpOUT = fft_out + i * m_nz);
440 }
441 }
442 else
443 {
444 for (int i = 0; i < (p * m_ny); i++)
445 {
446 m_FFT_z->FFTBwdTrans(m_tmpIN = fft_in + i * m_nz,
447 m_tmpOUT = fft_out + i * m_nz);
448 }
449 }
450
451 // TODO: required ZYtoX routine
452 m_transposition->Transpose(num_dofs, fft_out, fft_in, false,
454
455 m_transposition->Transpose(num_dofs, fft_in, outarray, false,
457 }
458 else
459 {
460 DNekBlkMatSharedPtr blkmatY;
461 DNekBlkMatSharedPtr blkmatZ;
462
463 if (num_dofs == m_npoints) // transform phys space
464 {
465 if (IsForwards)
466 {
469 }
470 else
471 {
474 }
475 }
476 else
477 {
478 if (IsForwards)
479 {
482 }
483 else
484 {
487 }
488 }
489
490 int nrowsY = blkmatY->GetRows();
491 int ncolsY = blkmatY->GetColumns();
492
493 Array<OneD, NekDouble> sortedinarrayY(ncolsY);
494 Array<OneD, NekDouble> sortedoutarrayY(nrowsY);
495
496 int nrowsZ = blkmatZ->GetRows();
497 int ncolsZ = blkmatZ->GetColumns();
498
499 Array<OneD, NekDouble> sortedinarrayZ(ncolsZ);
500 Array<OneD, NekDouble> sortedoutarrayZ(nrowsZ);
501
502 NekVector<NekDouble> inY(ncolsY, sortedinarrayY, eWrapper);
503 NekVector<NekDouble> outY(nrowsY, sortedoutarrayY, eWrapper);
504
505 NekVector<NekDouble> inZ(ncolsZ, sortedinarrayZ, eWrapper);
506 NekVector<NekDouble> outZ(nrowsZ, sortedoutarrayZ, eWrapper);
507
508 m_transposition->Transpose(num_dofs, inarray, sortedinarrayY,
509 !IsForwards, LibUtilities::eXtoYZ);
510
511 outY = (*blkmatY) * inY;
512
513 m_transposition->Transpose(sortedoutarrayY.size(), sortedoutarrayY,
514 sortedinarrayZ, false,
516
517 outZ = (*blkmatZ) * inZ;
518
519 m_transposition->Transpose(sortedoutarrayZ.size(), sortedoutarrayZ,
520 sortedoutarrayY, false,
522
523 m_transposition->Transpose(sortedoutarrayY.size(), sortedoutarrayY,
524 outarray, false, LibUtilities::eYZtoX);
525 }
526}
527
529 Homogeneous2DMatType mattype) const
530{
531 auto matrixIter = m_homogeneous2DBlockMat->find(mattype);
532
533 if (matrixIter == m_homogeneous2DBlockMat->end())
534 {
535 return ((*m_homogeneous2DBlockMat)[mattype] =
537 }
538 else
539 {
540 return matrixIter->second;
541 }
542}
543
545 Homogeneous2DMatType mattype) const
546{
547 int i;
548 int n_exp = 0;
549
550 DNekMatSharedPtr loc_mat;
551 DNekBlkMatSharedPtr BlkMatrix;
552
554
555 int NumPoints = 0;
556 int NumModes = 0;
557 int NumPencils = 0;
558
559 if ((mattype == eForwardsCoeffSpaceY1D) ||
560 (mattype == eBackwardsCoeffSpaceY1D) ||
561 (mattype == eForwardsPhysSpaceY1D) ||
562 (mattype == eBackwardsPhysSpaceY1D))
563 {
565 NumPoints = m_homogeneousBasis_y->GetNumModes();
566 NumModes = m_homogeneousBasis_y->GetNumPoints();
567 NumPencils = m_homogeneousBasis_z->GetNumPoints();
568 }
569 else
570 {
572 NumPoints = m_homogeneousBasis_z->GetNumModes();
573 NumModes = m_homogeneousBasis_z->GetNumPoints();
574 NumPencils = m_homogeneousBasis_y->GetNumPoints();
575 }
576
577 if ((mattype == eForwardsCoeffSpaceY1D) ||
578 (mattype == eForwardsCoeffSpaceZ1D) ||
579 (mattype == eBackwardsCoeffSpaceY1D) ||
580 (mattype == eBackwardsCoeffSpaceZ1D))
581 {
582 n_exp = m_lines[0]->GetNcoeffs();
583 }
584 else
585 {
586 n_exp = m_lines[0]->GetTotPoints(); // will operatore on m_phys
587 }
588
589 Array<OneD, unsigned int> nrows(n_exp);
590 Array<OneD, unsigned int> ncols(n_exp);
591
592 if ((mattype == eForwardsCoeffSpaceY1D) ||
593 (mattype == eForwardsPhysSpaceY1D) ||
594 (mattype == eForwardsCoeffSpaceZ1D) ||
595 (mattype == eForwardsPhysSpaceZ1D))
596 {
597 nrows = Array<OneD, unsigned int>(n_exp * NumPencils, NumModes);
598 ncols = Array<OneD, unsigned int>(n_exp * NumPencils, NumPoints);
599 }
600 else
601 {
602 nrows = Array<OneD, unsigned int>(n_exp * NumPencils, NumPoints);
603 ncols = Array<OneD, unsigned int>(n_exp * NumPencils, NumModes);
604 }
605
606 MatrixStorage blkmatStorage = eDIAGONAL;
607 BlkMatrix = MemoryManager<DNekBlkMat>::AllocateSharedPtr(nrows, ncols,
608 blkmatStorage);
609
611
612 if ((mattype == eForwardsCoeffSpaceY1D) ||
613 (mattype == eForwardsPhysSpaceY1D) ||
614 (mattype == eForwardsCoeffSpaceZ1D) ||
615 (mattype == eForwardsPhysSpaceZ1D))
616 {
618 StdSeg.DetShapeType(), StdSeg);
619
620 loc_mat = StdSeg.GetStdMatrix(matkey);
621 }
622 else
623 {
625 StdSeg.DetShapeType(), StdSeg);
626
627 loc_mat = StdSeg.GetStdMatrix(matkey);
628 }
629
630 // set up array of block matrices.
631 for (i = 0; i < (n_exp * NumPencils); ++i)
632 {
633 BlkMatrix->SetBlock(i, i, loc_mat);
634 }
635
636 return BlkMatrix;
637}
638
639std::vector<LibUtilities::FieldDefinitionsSharedPtr> ExpListHomogeneous2D::
641{
642 std::vector<LibUtilities::FieldDefinitionsSharedPtr> returnval;
643 // Set up Homogeneous length details.
645 HomoBasis[0] = m_homogeneousBasis_y;
646 HomoBasis[1] = m_homogeneousBasis_z;
647
648 std::vector<NekDouble> HomoLen(2);
649 HomoLen[0] = m_lhom_y;
650 HomoLen[1] = m_lhom_z;
651
652 int nhom_modes_y = m_homogeneousBasis_y->GetNumModes();
653 int nhom_modes_z = m_homogeneousBasis_z->GetNumModes();
654
655 std::vector<unsigned int> sIDs = LibUtilities::NullUnsignedIntVector;
656
657 std::vector<unsigned int> yIDs;
658 std::vector<unsigned int> zIDs;
659
660 for (int n = 0; n < nhom_modes_z; ++n)
661 {
662 for (int m = 0; m < nhom_modes_y; ++m)
663 {
664 zIDs.push_back(n);
665 yIDs.push_back(m);
666 }
667 }
668
669 m_lines[0]->GeneralGetFieldDefinitions(returnval, 2, HomoBasis, HomoLen,
670 false, sIDs, zIDs, yIDs);
671 return returnval;
672}
673
675 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
676{
677 // Set up Homogeneous length details.
679 HomoBasis[0] = m_homogeneousBasis_y;
680 HomoBasis[1] = m_homogeneousBasis_z;
681 std::vector<NekDouble> HomoLen(2);
682 HomoLen[0] = m_lhom_y;
683 HomoLen[1] = m_lhom_z;
684
685 int nhom_modes_y = m_homogeneousBasis_y->GetNumModes();
686 int nhom_modes_z = m_homogeneousBasis_z->GetNumModes();
687
688 std::vector<unsigned int> sIDs = LibUtilities::NullUnsignedIntVector;
689
690 std::vector<unsigned int> yIDs;
691 std::vector<unsigned int> zIDs;
692
693 for (int n = 0; n < nhom_modes_z; ++n)
694 {
695 for (int m = 0; m < nhom_modes_y; ++m)
696 {
697 zIDs.push_back(n);
698 yIDs.push_back(m);
699 }
700 }
701
702 // enforce NumHomoDir == 1 by direct call
703 m_lines[0]->GeneralGetFieldDefinitions(fielddef, 2, HomoBasis, HomoLen,
704 false, sIDs, zIDs, yIDs);
705}
706
709 std::vector<NekDouble> &fielddata, Array<OneD, NekDouble> &coeffs)
710{
711 int i, k;
712
713 int NumMod_y = m_homogeneousBasis_y->GetNumModes();
714 int NumMod_z = m_homogeneousBasis_z->GetNumModes();
715
716 int ncoeffs_per_line = m_lines[0]->GetNcoeffs();
717
718 // Determine mapping from element ids to location in
719 // expansion list
720 map<int, int> ElmtID_to_ExpID;
721 for (i = 0; i < m_lines[0]->GetExpSize(); ++i)
722 {
723 ElmtID_to_ExpID[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
724 }
725
726 for (i = 0; i < fielddef->m_elementIDs.size(); ++i)
727 {
728 int eid = ElmtID_to_ExpID[fielddef->m_elementIDs[i]];
729 int datalen = (*m_exp)[eid]->GetNcoeffs();
730
731 for (k = 0; k < (NumMod_y * NumMod_z); ++k)
732 {
733 fielddata.insert(
734 fielddata.end(),
735 &coeffs[m_coeff_offset[eid] + k * ncoeffs_per_line],
736 &coeffs[m_coeff_offset[eid] + k * ncoeffs_per_line] + datalen);
737 }
738 }
739}
740
743 std::vector<NekDouble> &fielddata)
744{
745 v_AppendFieldData(fielddef, fielddata, m_coeffs);
746}
747// Extract the data in fielddata into the m_coeff list
750 std::vector<NekDouble> &fielddata, std::string &field,
752 [[maybe_unused]] std::unordered_map<int, int> zIdToPlane)
753{
754 int i, k;
755 int offset = 0;
756 int datalen = fielddata.size() / fielddef->m_fields.size();
757 int ncoeffs_per_line = m_lines[0]->GetNcoeffs();
758 int NumMod_y = m_homogeneousBasis_y->GetNumModes();
759 int NumMod_z = m_homogeneousBasis_z->GetNumModes();
760 // Find data location according to field definition
761 for (i = 0; i < fielddef->m_fields.size(); ++i)
762 {
763 if (fielddef->m_fields[i] == field)
764 {
765 break;
766 }
767 offset += datalen;
768 }
769 ASSERTL0(i != fielddef->m_fields.size(), "Field not found in data file");
770 // Determine mapping from element ids to location in
771 // expansion list
772 map<int, int> ElmtID_to_ExpID;
773 for (i = 0; i < m_lines[0]->GetExpSize(); ++i)
774 {
775 ElmtID_to_ExpID[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
776 }
777 for (i = 0; i < fielddef->m_elementIDs.size(); ++i)
778 {
779 int eid = ElmtID_to_ExpID[fielddef->m_elementIDs[i]];
780 int datalen = (*m_exp)[eid]->GetNcoeffs();
781 for (k = 0; k < (NumMod_y * NumMod_z); ++k)
782 {
783 Vmath::Vcopy(datalen, &fielddata[offset], 1,
784 &coeffs[m_coeff_offset[eid] + k * ncoeffs_per_line],
785 1);
786 offset += datalen;
787 }
788 }
789}
791 int expansion, std::string var)
792{
793 int i;
794 int nq = (*m_exp)[expansion]->GetTotPoints();
795 int npoints_per_line = m_lines[0]->GetTotPoints();
796
797 // printing the fields of that zone
798 outfile << R"( <DataArray type="Float64" Name=")" << var << "\">"
799 << endl;
800 outfile << " ";
801 for (int n = 0; n < m_lines.size(); ++n)
802 {
803 const Array<OneD, NekDouble> phys =
804 m_phys + m_phys_offset[expansion] + n * npoints_per_line;
805
806 for (i = 0; i < nq; ++i)
807 {
808 outfile << (fabs(phys[i]) < NekConstants::kNekZeroTol ? 0 : phys[i])
809 << " ";
810 }
811 }
812 outfile << endl;
813 outfile << " </DataArray>" << endl;
814}
815
819
820{
821 int nyzlines = m_lines.size(); // number of Fourier points in the Fourier
822 // directions (nF_pts)
823 int n_points_line = m_npoints / nyzlines; // number of points per line
824
831
832 for (int i = 0; i < nyzlines; i++)
833 {
834 m_lines[i]->PhysDeriv(tmp1 = inarray + i * n_points_line,
835 tmp2 = out_d0 + i * n_points_line);
836 }
837
838 if (m_homogeneousBasis_y->GetBasisType() == LibUtilities::eFourier &&
840 {
841 if (m_WaveSpace)
842 {
843 temparray = inarray;
844 }
845 else
846 {
847 HomogeneousFwdTrans(m_npoints, inarray, temparray);
848 }
849 NekDouble sign = -1.0;
851
852 // along y
853 for (int i = 0; i < m_ny; i++)
854 {
855 beta = -sign * 2 * M_PI * (i / 2) / m_lhom_y;
856
857 for (int j = 0; j < m_nz; j++)
858 {
859 Vmath::Smul(n_points_line, beta,
860 tmp1 = temparray + n_points_line * (i + j * m_ny),
861 1,
862 tmp2 = temparray1 +
863 n_points_line * ((i - int(sign)) + j * m_ny),
864 1);
865 }
866
867 sign = -1.0 * sign;
868 }
869
870 // along z
871 sign = -1.0;
872 for (int i = 0; i < m_nz; i++)
873 {
874 beta = -sign * 2 * M_PI * (i / 2) / m_lhom_z;
876 m_ny * n_points_line, beta,
877 tmp1 = temparray + i * m_ny * n_points_line, 1,
878 tmp2 = temparray2 + (i - int(sign)) * m_ny * n_points_line, 1);
879 sign = -1.0 * sign;
880 }
881 if (m_WaveSpace)
882 {
883 out_d1 = temparray1;
884 out_d2 = temparray2;
885 }
886 else
887 {
888 HomogeneousBwdTrans(m_npoints, temparray1, out_d1);
889 HomogeneousBwdTrans(m_npoints, temparray2, out_d2);
890 }
891 }
892 else
893 {
894 if (m_WaveSpace)
895 {
896 ASSERTL0(false, "Semi-phyisical time-stepping not implemented yet "
897 "for non-Fourier basis")
898 }
899 else
900 {
901 StdRegions::StdQuadExp StdQuad(m_homogeneousBasis_y->GetBasisKey(),
902 m_homogeneousBasis_z->GetBasisKey());
903
904 m_transposition->Transpose(m_npoints, inarray, temparray, false,
906
907 for (int i = 0; i < n_points_line; i++)
908 {
909 StdQuad.PhysDeriv(tmp1 = temparray + i * nyzlines,
910 tmp2 = temparray1 + i * nyzlines,
911 tmp3 = temparray2 + i * nyzlines);
912 }
913
914 m_transposition->Transpose(m_npoints, temparray1, out_d1, false,
916 m_transposition->Transpose(m_npoints, temparray2, out_d2, false,
918 Vmath::Smul(m_npoints, 2.0 / m_lhom_y, out_d1, 1, out_d1, 1);
919 Vmath::Smul(m_npoints, 2.0 / m_lhom_z, out_d2, 1, out_d2, 1);
920 }
921 }
922}
923
925 Direction edir, const Array<OneD, const NekDouble> &inarray,
927
928{
929 int nyzlines = m_lines.size(); // number of Fourier points in the Fourier
930 // directions (nF_pts)
931 int n_points_line = m_npoints / nyzlines; // number of points per line
932 // convert enum into int
933 int dir = (int)edir;
934
941
942 if (dir < 1)
943 {
944 for (int i = 0; i < nyzlines; i++)
945 {
946 m_lines[i]->PhysDeriv(tmp1 = inarray + i * n_points_line,
947 tmp2 = out_d + i * n_points_line);
948 }
949 }
950 else
951 {
952 if (m_homogeneousBasis_y->GetBasisType() == LibUtilities::eFourier &&
954 {
955 if (m_WaveSpace)
956 {
957 temparray = inarray;
958 }
959 else
960 {
961 HomogeneousFwdTrans(m_npoints, inarray, temparray);
962 }
963 NekDouble sign = -1.0;
965
966 if (dir == 1)
967 {
968 // along y
969 for (int i = 0; i < m_ny; i++)
970 {
971 beta = -sign * 2 * M_PI * (i / 2) / m_lhom_y;
972
973 for (int j = 0; j < m_nz; j++)
974 {
976 n_points_line, beta,
977 tmp1 = temparray + n_points_line * (i + j * m_ny),
978 1,
979 tmp2 = temparray1 +
980 n_points_line * ((i - int(sign)) + j * m_ny),
981 1);
982 }
983 sign = -1.0 * sign;
984 }
985 if (m_WaveSpace)
986 {
987 out_d = temparray1;
988 }
989 else
990 {
991 HomogeneousBwdTrans(m_npoints, temparray1, out_d);
992 }
993 }
994 else
995 {
996 // along z
997 for (int i = 0; i < m_nz; i++)
998 {
999 beta = -sign * 2 * M_PI * (i / 2) / m_lhom_z;
1000 Vmath::Smul(m_ny * n_points_line, beta,
1001 tmp1 = temparray + i * m_ny * n_points_line, 1,
1002 tmp2 = temparray2 +
1003 (i - int(sign)) * m_ny * n_points_line,
1004 1);
1005 sign = -1.0 * sign;
1006 }
1007 if (m_WaveSpace)
1008 {
1009 out_d = temparray2;
1010 }
1011 else
1012 {
1013 HomogeneousBwdTrans(m_npoints, temparray2, out_d);
1014 }
1015 }
1016 }
1017 else
1018 {
1019 if (m_WaveSpace)
1020 {
1021 ASSERTL0(false, "Semi-phyisical time-stepping not implemented "
1022 "yet for non-Fourier basis")
1023 }
1024 else
1025 {
1026 StdRegions::StdQuadExp StdQuad(
1027 m_homogeneousBasis_y->GetBasisKey(),
1028 m_homogeneousBasis_z->GetBasisKey());
1029
1030 m_transposition->Transpose(m_npoints, inarray, temparray, false,
1032
1033 for (int i = 0; i < n_points_line; i++)
1034 {
1035 StdQuad.PhysDeriv(tmp1 = temparray + i * nyzlines,
1036 tmp2 = temparray1 + i * nyzlines,
1037 tmp3 = temparray2 + i * nyzlines);
1038 }
1039
1040 if (dir == 1)
1041 {
1042 m_transposition->Transpose(m_npoints, temparray1, out_d,
1043 false, LibUtilities::eYZtoX);
1044 Vmath::Smul(m_npoints, 2.0 / m_lhom_y, out_d, 1, out_d, 1);
1045 }
1046 else
1047 {
1048 m_transposition->Transpose(m_npoints, temparray2, out_d,
1049 false, LibUtilities::eYZtoX);
1050 Vmath::Smul(m_npoints, 2.0 / m_lhom_z, out_d, 1, out_d, 1);
1051 }
1052 }
1053 }
1054 }
1055}
1056
1060
1061{
1062 v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
1063}
1064
1066 Direction edir, const Array<OneD, const NekDouble> &inarray,
1068{
1069 // convert int into enum
1070 v_PhysDeriv(edir, inarray, out_d);
1071}
1072
1074{
1075 NekDouble size_y = 1.5 * m_ny;
1076 NekDouble size_z = 1.5 * m_nz;
1077 m_padsize_y = int(size_y);
1078 m_padsize_z = int(size_z);
1079
1083 Ppad_y);
1084
1088 Ppad_z);
1089
1092
1093 StdRegions::StdQuadExp StdQuad(m_paddingBasis_y->GetBasisKey(),
1094 m_paddingBasis_z->GetBasisKey());
1095
1097 StdQuad.DetShapeType(), StdQuad);
1099 StdQuad.DetShapeType(), StdQuad);
1100
1101 MatFwdPAD = StdQuad.GetStdMatrix(matkey1);
1102 MatBwdPAD = StdQuad.GetStdMatrix(matkey2);
1103}
1104} // 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
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:47
Represents a basis of a given type.
Definition: Basis.h:198
const BasisKey GetBasisKey() const
Returns the specification for the Basis.
Definition: Basis.h:313
Describes the specification for a Basis.
Definition: Basis.h:45
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
Defines a specification for a set of points.
Definition: Points.h:50
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
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_DealiasedProd(const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray) override
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_DealiasedDotProd(const int num_dofs, const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray) override
void v_FwdTransLocalElmt(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void) override
DNekBlkMatSharedPtr GenHomogeneous2DBlockMatrix(Homogeneous2DMatType mattype) const
void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
LibUtilities::NektarFFTSharedPtr m_FFT_y
ExpListHomogeneous2D(const ExpansionType type)
Default constructor.
int m_nz
Number of modes = number of poitns in z direction.
LibUtilities::NektarFFTSharedPtr m_FFT_z
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
NekDouble m_lhom_z
Width of homogeneous direction z.
int m_ny
Number of modes = number of poitns in y direction.
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
void v_HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true) override
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var) override
LibUtilities::TranspositionSharedPtr m_transposition
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
LibUtilities::BasisSharedPtr m_paddingBasis_z
Base expansion in z direction.
void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata) override
void v_ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs, std::unordered_map< int, int > zIdToPlane) override
Extract data from raw field data into expansion list.
void Homogeneous2DTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool IsForwards, bool Shuff=true, bool UnShuff=true)
DNekBlkMatSharedPtr GetHomogeneous2DBlockMatrix(Homogeneous2DMatType mattype) const
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
LibUtilities::BasisSharedPtr m_paddingBasis_y
Base expansion in y direction.
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
NekDouble m_lhom_y
Width of homogeneous direction y.
void v_HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true) override
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:98
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1080
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1124
void DealiasedProd(const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1866
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1053
void HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1857
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
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1096
void HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1848
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:367
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:849
Class representing a segment element in reference space All interface of this class sits in StdExpans...
Definition: StdSegExp.h:45
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:350
static std::vector< unsigned int > NullUnsignedIntVector
BasisManagerT & BasisManager(void)
std::shared_ptr< Basis > BasisSharedPtr
NektarFFTFactory & GetNektarFFTFactory()
Definition: NektarFFT.cpp:65
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:184
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:74
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
std::map< Homogeneous2DMatType, DNekBlkMatSharedPtr > Homo2DBlockMatrixMap
A map between homo matrix keys and their associated block matrices.
static const NekDouble kNekZeroTol
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:77
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
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 Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825