Nektar++
ExpListHomogeneous1D.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ExpListHomogeneous1D.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
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <LibUtilities/Foundations/ManagerAccess.h> // for PointsManager, etc
41
42using namespace std;
43
44namespace Nektar
45{
46namespace MultiRegions
47{
48// Forward declaration for typedefs
50 : ExpList(type), m_homogeneousBasis(LibUtilities::NullBasisSharedPtr),
51 m_lhom(1), m_homogeneous1DBlockMat(
52 MemoryManager<Homo1DBlockMatrixMap>::AllocateSharedPtr())
53{
54}
55
57 const ExpansionType type,
59 const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom,
60 const bool useFFT, const bool dealiasing)
61 : ExpList(type), m_useFFT(useFFT), m_lhom(lhom),
62 m_homogeneous1DBlockMat(
63 MemoryManager<Homo1DBlockMatrixMap>::AllocateSharedPtr()),
64 m_dealiasing(dealiasing)
65{
66 m_session = pSession;
67 m_comm = pSession->GetComm();
68
70 "Homogeneous Basis is a null basis");
71
73
74 m_StripZcomm = m_session->DefinesSolverInfo("HomoStrip")
75 ? m_comm->GetColumnComm()->GetColumnComm()
76 : m_comm->GetColumnComm();
77
78 ASSERTL0(m_homogeneousBasis->GetNumPoints() % m_StripZcomm->GetSize() == 0,
79 "HomModesZ should be a multiple of npz.");
80
81 if ((m_homogeneousBasis->GetBasisType() !=
83 (m_homogeneousBasis->GetBasisType() !=
85 {
87 (m_homogeneousBasis->GetNumPoints() / m_StripZcomm->GetSize()) %
88 2 ==
89 0,
90 "HomModesZ/npz should be an even integer.");
91 }
92
95 HomoBasis, m_comm, m_StripZcomm);
96
98 m_homogeneousBasis->GetNumPoints() / m_StripZcomm->GetSize());
99
100 if (m_useFFT)
101 {
103 "NekFFTW", m_homogeneousBasis->GetNumPoints());
104 }
105
106 if (m_dealiasing)
107 {
108 if (m_useFFT)
109 {
110 int size = m_homogeneousBasis->GetNumPoints() +
111 m_homogeneousBasis->GetNumPoints() / 2;
112 m_padsize = size + (size % 2);
114 "NekFFTW", m_padsize);
115 }
116 else
117 {
118 ASSERTL0(false, "Dealiasing available just in combination "
119 "with FFTW");
120 }
121 }
122}
123
124/**
125 * @param In ExpListHomogeneous1D object to copy.
126 */
128 : ExpList(In, false), m_transposition(In.m_transposition),
129 m_StripZcomm(In.m_StripZcomm), m_useFFT(In.m_useFFT), m_FFT(In.m_FFT),
130 m_tmpIN(In.m_tmpIN), m_tmpOUT(In.m_tmpOUT),
131 m_homogeneousBasis(In.m_homogeneousBasis), m_lhom(In.m_lhom),
132 m_homogeneous1DBlockMat(In.m_homogeneous1DBlockMat),
133 m_dealiasing(In.m_dealiasing), m_padsize(In.m_padsize)
134{
136}
137
139 const ExpListHomogeneous1D &In, const std::vector<unsigned int> &eIDs,
141 : ExpList(In, eIDs, false, ImpType), m_transposition(In.m_transposition),
142 m_useFFT(In.m_useFFT), m_FFT(In.m_FFT), m_tmpIN(In.m_tmpIN),
143 m_tmpOUT(In.m_tmpOUT), m_homogeneousBasis(In.m_homogeneousBasis),
144 m_lhom(In.m_lhom),
145 m_homogeneous1DBlockMat(
146 MemoryManager<Homo1DBlockMatrixMap>::AllocateSharedPtr()),
147 m_dealiasing(In.m_dealiasing), m_padsize(In.m_padsize)
148{
150}
151
152/**
153 * Destructor
154 */
156{
157}
158
160 const int npts, const Array<OneD, const NekDouble> &inarray,
161 Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
162{
163 // Forwards trans
164 Homogeneous1DTrans(npts, inarray, outarray, true, Shuff, UnShuff);
165}
166
168 const int npts, const Array<OneD, const NekDouble> &inarray,
169 Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
170{
171 // Backwards trans
172 Homogeneous1DTrans(npts, inarray, outarray, false, Shuff, UnShuff);
173}
174
175/**
176 * Dealiasing routine
177 *
178 * @param inarray1 First term of the product
179 * @param inarray2 Second term of the product
180 * @param outarray Dealiased product
181 */
183 const int num_dofs, const Array<OneD, NekDouble> &inarray1,
184 const Array<OneD, NekDouble> &inarray2, Array<OneD, NekDouble> &outarray)
185{
186 // int num_dofs = inarray1.size();
187 int N = m_homogeneousBasis->GetNumPoints();
188
189 Array<OneD, NekDouble> V1(num_dofs);
190 Array<OneD, NekDouble> V2(num_dofs);
191 Array<OneD, NekDouble> V1V2(num_dofs);
192
193 if (m_WaveSpace)
194 {
195 V1 = inarray1;
196 V2 = inarray2;
197 }
198 else
199 {
200 HomogeneousFwdTrans(num_dofs, inarray1, V1);
201 HomogeneousFwdTrans(num_dofs, inarray2, V2);
202 }
203
204 int num_points_per_plane = num_dofs / m_planes.size();
205 int num_proc;
206 if (!m_session->DefinesSolverInfo("HomoStrip"))
207 {
208 num_proc = m_comm->GetColumnComm()->GetSize();
209 }
210 else
211 {
212 num_proc = m_StripZcomm->GetSize();
213 }
214 int num_dfts_per_proc =
215 num_points_per_plane / num_proc + (num_points_per_plane % num_proc > 0);
216
217 Array<OneD, NekDouble> ShufV1(num_dfts_per_proc * N, 0.0);
218 Array<OneD, NekDouble> ShufV2(num_dfts_per_proc * N, 0.0);
219 Array<OneD, NekDouble> ShufV1V2(num_dfts_per_proc * N, 0.0);
220
221 Array<OneD, NekDouble> ShufV1_PAD_coef(m_padsize, 0.0);
222 Array<OneD, NekDouble> ShufV2_PAD_coef(m_padsize, 0.0);
223 Array<OneD, NekDouble> ShufV1_PAD_phys(m_padsize, 0.0);
224 Array<OneD, NekDouble> ShufV2_PAD_phys(m_padsize, 0.0);
225
226 Array<OneD, NekDouble> ShufV1V2_PAD_coef(m_padsize, 0.0);
227 Array<OneD, NekDouble> ShufV1V2_PAD_phys(m_padsize, 0.0);
228
229 m_transposition->Transpose(num_dofs, V1, ShufV1, false,
231 m_transposition->Transpose(num_dofs, V2, ShufV2, false,
233
234 // Looping on the pencils
235 for (int i = 0; i < num_dfts_per_proc; i++)
236 {
237 // Copying the i-th pencil pf lenght N into a bigger
238 // pencil of lenght 3/2N We are in Fourier space
239 Vmath::Vcopy(N, &(ShufV1[i * N]), 1, &(ShufV1_PAD_coef[0]), 1);
240 Vmath::Vcopy(N, &(ShufV2[i * N]), 1, &(ShufV2_PAD_coef[0]), 1);
241
242 // Moving to physical space using the padded system
243 m_FFT_deal->FFTBwdTrans(ShufV1_PAD_coef, ShufV1_PAD_phys);
244 m_FFT_deal->FFTBwdTrans(ShufV2_PAD_coef, ShufV2_PAD_phys);
245
246 // Perfroming the vectors multiplication in physical space on
247 // the padded system
248 Vmath::Vmul(m_padsize, ShufV1_PAD_phys, 1, ShufV2_PAD_phys, 1,
249 ShufV1V2_PAD_phys, 1);
250
251 // Moving back the result (V1*V2)_phys in Fourier space, padded
252 // system
253 m_FFT_deal->FFTFwdTrans(ShufV1V2_PAD_phys, ShufV1V2_PAD_coef);
254
255 // Copying the first half of the padded pencil in the full
256 // vector (Fourier space)
257 Vmath::Vcopy(N, &(ShufV1V2_PAD_coef[0]), 1, &(ShufV1V2[i * N]), 1);
258 }
259
260 // Moving the results to the output
261 if (m_WaveSpace)
262 {
263 m_transposition->Transpose(num_dofs, ShufV1V2, outarray, false,
265 }
266 else
267 {
268 m_transposition->Transpose(num_dofs, ShufV1V2, V1V2, false,
270 HomogeneousBwdTrans(num_dofs, V1V2, outarray);
271 }
272}
273
274/**
275 * Dealiasing routine for dot product
276 *
277 * @param inarray1 First term of the product with dimension ndim
278 * (e.g. U)
279 * @param inarray2 Second term of the product with dimension ndim*nvec
280 * (e.g. grad U)
281 * @param outarray Dealiased product with dimension nvec
282 */
284 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
285 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
287{
288 int ndim = inarray1.size();
289 ASSERTL1(inarray2.size() % ndim == 0,
290 "Wrong dimensions for DealiasedDotProd.");
291 int nvec = inarray2.size() / ndim;
292
293 // int num_dofs = inarray1[0].size();
294 int N = m_homogeneousBasis->GetNumPoints();
295
296 int num_points_per_plane = num_dofs / m_planes.size();
297 int num_proc;
298 if (!m_session->DefinesSolverInfo("HomoStrip"))
299 {
300 num_proc = m_comm->GetColumnComm()->GetSize();
301 }
302 else
303 {
304 num_proc = m_StripZcomm->GetSize();
305 }
306 int num_dfts_per_proc =
307 num_points_per_plane / num_proc + (num_points_per_plane % num_proc > 0);
308
309 // Get inputs in Fourier space
311 Array<OneD, Array<OneD, NekDouble>> V2(ndim * nvec);
312 if (m_WaveSpace)
313 {
314 for (int i = 0; i < ndim; i++)
315 {
316 V1[i] = inarray1[i];
317 }
318 for (int i = 0; i < ndim * nvec; i++)
319 {
320 V2[i] = inarray2[i];
321 }
322 }
323 else
324 {
325 for (int i = 0; i < ndim; i++)
326 {
327 V1[i] = Array<OneD, NekDouble>(num_dofs);
328 HomogeneousFwdTrans(num_dofs, inarray1[i], V1[i]);
329 }
330 for (int i = 0; i < ndim * nvec; i++)
331 {
332 V2[i] = Array<OneD, NekDouble>(num_dofs);
333 HomogeneousFwdTrans(num_dofs, inarray2[i], V2[i]);
334 }
335 }
336
337 // Allocate variables for ffts
339 Array<OneD, NekDouble> ShufV1_PAD_coef(m_padsize, 0.0);
340 Array<OneD, Array<OneD, NekDouble>> ShufV1_PAD_phys(ndim);
341 for (int i = 0; i < ndim; i++)
342 {
343 ShufV1[i] = Array<OneD, NekDouble>(num_dfts_per_proc * N, 0.0);
344 ShufV1_PAD_phys[i] = Array<OneD, NekDouble>(m_padsize, 0.0);
345 }
346
347 Array<OneD, Array<OneD, NekDouble>> ShufV2(ndim * nvec);
348 Array<OneD, NekDouble> ShufV2_PAD_coef(m_padsize, 0.0);
349 Array<OneD, Array<OneD, NekDouble>> ShufV2_PAD_phys(ndim * nvec);
350 for (int i = 0; i < ndim * nvec; i++)
351 {
352 ShufV2[i] = Array<OneD, NekDouble>(num_dfts_per_proc * N, 0.0);
353 ShufV2_PAD_phys[i] = Array<OneD, NekDouble>(m_padsize, 0.0);
354 }
355
357 Array<OneD, NekDouble> ShufV1V2_PAD_coef(m_padsize, 0.0);
358 Array<OneD, NekDouble> ShufV1V2_PAD_phys(m_padsize, 0.0);
359 for (int i = 0; i < nvec; i++)
360 {
361 ShufV1V2[i] = Array<OneD, NekDouble>(num_dfts_per_proc * N, 0.0);
362 }
363
364 for (int i = 0; i < ndim; i++)
365 {
366 m_transposition->Transpose(num_dofs, V1[i], ShufV1[i], false,
368 }
369 for (int i = 0; i < ndim * nvec; i++)
370 {
371 m_transposition->Transpose(num_dofs, V2[i], ShufV2[i], false,
373 }
374
375 // Looping on the pencils
376 for (int i = 0; i < num_dfts_per_proc; i++)
377 {
378 for (int j = 0; j < ndim; j++)
379 {
380 // Copying the i-th pencil pf lenght N into a bigger
381 // pencil of lenght 1.5N We are in Fourier space
382 Vmath::Vcopy(N, &(ShufV1[j][i * N]), 1, &(ShufV1_PAD_coef[0]), 1);
383 // Moving to physical space using the padded system
384 m_FFT_deal->FFTBwdTrans(ShufV1_PAD_coef, ShufV1_PAD_phys[j]);
385 }
386 for (int j = 0; j < ndim * nvec; j++)
387 {
388 Vmath::Vcopy(N, &(ShufV2[j][i * N]), 1, &(ShufV2_PAD_coef[0]), 1);
389 m_FFT_deal->FFTBwdTrans(ShufV2_PAD_coef, ShufV2_PAD_phys[j]);
390 }
391
392 // Performing the vectors multiplication in physical space on
393 // the padded system
394 for (int j = 0; j < nvec; j++)
395 {
396 Vmath::Zero(m_padsize, ShufV1V2_PAD_phys, 1);
397 for (int k = 0; k < ndim; k++)
398 {
399 Vmath::Vvtvp(m_padsize, ShufV1_PAD_phys[k], 1,
400 ShufV2_PAD_phys[j * ndim + k], 1,
401 ShufV1V2_PAD_phys, 1, ShufV1V2_PAD_phys, 1);
402 }
403 // Moving back the result (V1*V2)_phys in Fourier space,
404 // padded system
405 m_FFT_deal->FFTFwdTrans(ShufV1V2_PAD_phys, ShufV1V2_PAD_coef);
406 // Copying the first half of the padded pencil in the full
407 // vector (Fourier space)
408 Vmath::Vcopy(N, &(ShufV1V2_PAD_coef[0]), 1, &(ShufV1V2[j][i * N]),
409 1);
410 }
411 }
412
413 // Moving the results to the output
414 if (m_WaveSpace)
415 {
416 for (int j = 0; j < nvec; j++)
417 {
418 m_transposition->Transpose(num_dofs, ShufV1V2[j], outarray[j],
419 false, LibUtilities::eZtoXY);
420 }
421 }
422 else
423 {
424 Array<OneD, NekDouble> V1V2(num_dofs);
425 for (int j = 0; j < nvec; j++)
426 {
427 m_transposition->Transpose(num_dofs, ShufV1V2[j], V1V2, false,
429 HomogeneousBwdTrans(num_dofs, V1V2, outarray[j]);
430 }
431 }
432}
433
434/**
435 * Forward transform
436 */
438 const Array<OneD, const NekDouble> &inarray,
439 Array<OneD, NekDouble> &outarray)
440{
441 int cnt = 0, cnt1 = 0;
442 Array<OneD, NekDouble> tmparray;
443
444 for (int n = 0; n < m_planes.size(); ++n)
445 {
446 m_planes[n]->FwdTrans(inarray + cnt, tmparray = outarray + cnt1);
447 cnt += m_planes[n]->GetTotPoints();
448
449 cnt1 += m_planes[n]->GetNcoeffs(); // need to skip ncoeffs
450 }
451 if (!m_WaveSpace)
452 {
453 HomogeneousFwdTrans(cnt1, outarray, outarray);
454 }
455}
456
457/**
458 * Forward transform element by element
459 */
461 const Array<OneD, const NekDouble> &inarray,
462 Array<OneD, NekDouble> &outarray)
463{
464 int cnt = 0, cnt1 = 0;
465 Array<OneD, NekDouble> tmparray;
466
467 // spectral element FwdTrans plane by plane
468 for (int n = 0; n < m_planes.size(); ++n)
469 {
470 m_planes[n]->FwdTransLocalElmt(inarray + cnt,
471 tmparray = outarray + cnt1);
472 cnt += m_planes[n]->GetTotPoints();
473 cnt1 += m_planes[n]->GetNcoeffs();
474 }
475 if (!m_WaveSpace)
476 {
477 HomogeneousFwdTrans(cnt1, outarray, outarray);
478 }
479}
480
481/**
482 * Forward transform element by element with boundaries constrained
483 */
485 const Array<OneD, const NekDouble> &inarray,
486 Array<OneD, NekDouble> &outarray)
487{
488 int cnt = 0, cnt1 = 0;
489 Array<OneD, NekDouble> tmparray;
490
491 // spectral element FwdTrans plane by plane
492 for (int n = 0; n < m_planes.size(); ++n)
493 {
494 m_planes[n]->FwdTransBndConstrained(inarray + cnt,
495 tmparray = outarray + cnt1);
496 cnt += m_planes[n]->GetTotPoints();
497 cnt1 += m_planes[n]->GetNcoeffs();
498 }
499 if (!m_WaveSpace)
500 {
501 HomogeneousFwdTrans(cnt1, outarray, outarray);
502 }
503}
504
505/**
506 * Backward transform
507 */
509 const Array<OneD, const NekDouble> &inarray,
510 Array<OneD, NekDouble> &outarray)
511{
512 int cnt = 0, cnt1 = 0;
513 Array<OneD, NekDouble> tmparray;
514
515 for (int n = 0; n < m_planes.size(); ++n)
516 {
517 m_planes[n]->BwdTrans(inarray + cnt, tmparray = outarray + cnt1);
518 cnt += m_planes[n]->GetNcoeffs();
519 cnt1 += m_planes[n]->GetTotPoints();
520 }
521 if (!m_WaveSpace)
522 {
523 HomogeneousBwdTrans(cnt1, outarray, outarray);
524 }
525}
526
527/**
528 * Inner product
529 */
531 const Array<OneD, const NekDouble> &inarray,
532 Array<OneD, NekDouble> &outarray)
533{
534 ASSERTL1(inarray.size() >= m_npoints,
535 "Inarray is not of sufficient length");
536 int cnt = 0, cnt1 = 0;
537 Array<OneD, NekDouble> tmparray, tmpIn;
538
539 if (m_WaveSpace)
540 {
541 tmpIn = inarray;
542 }
543 else
544 {
545 tmpIn = Array<OneD, NekDouble>(m_npoints, 0.0);
546 HomogeneousFwdTrans(m_npoints, inarray, tmpIn);
547 }
548
549 for (int n = 0; n < m_planes.size(); ++n)
550 {
551 m_planes[n]->IProductWRTBase(tmpIn + cnt, tmparray = outarray + cnt1);
552
553 cnt1 += m_planes[n]->GetNcoeffs();
554 cnt += m_planes[n]->GetTotPoints();
555 }
556}
557
558/**
559 * Inner product with respect to the derivative of the Basis
560 * including one homogeneous direction.
561 *
562 * This function evaluates \f[\int_\Omega \partial_d \phi \mathbf{a}_d\f],
563 * where \f$\partial_d \phi\f$ is the derivative of the basis function
564 * (including a homogeneous basis) in direction \f$d\f$ and \f$\mathbf{a}_d\f$
565 * is the \f$d\f$-th component of a generic vector \f$\mathbf{a}\f$, which is
566 * fourier decomposed into modes of wavenumber \f$k\f$.
567 *
568 * Assuming a 3DH1D basis, the generic vector \f$\mathbf{a}\f$ is defined as
569 * \f[\mathbf{a}=\sum_k\begin{bmatrix}a_x^k(x,y)\\a_y^k(x,y)\\a_z^k(x,y)\end{bmatrix}e^{i
570 * \beta k z}\f].
571 *
572 */
574 const int dir, const Array<OneD, const NekDouble> &inarray,
575 Array<OneD, NekDouble> &outarray)
576{
577 int nT_pts = m_npoints; // number of total points = n. of Fourier
578 // points * n. of points per plane (nT_pts)
579 int nP_pts =
580 nT_pts / m_planes.size(); // number of points per plane = n of
581 // Fourier transform required (nP_pts)
582 int nT_coeffs = m_planes[0]->GetNcoeffs() * m_planes.size();
583 int nP_coeffs = nT_coeffs / m_planes.size();
584
585 Array<OneD, NekDouble> tmpIn(nT_pts);
586 Array<OneD, NekDouble> tmp1, tmp2;
587
588 // Check: input in Fourier space
589 if (m_WaveSpace)
590 {
591 tmpIn = inarray;
592 }
593 else
594 {
595 HomogeneousFwdTrans(m_npoints, inarray, tmpIn);
596 }
597
598 // Do 2D-IProductWRTDerivBase on each plane
599 if (dir == 0 || dir == 1)
600 {
601 for (int i = 0; i < m_planes.size(); i++)
602 {
603 // Call 2D routine on each plane
604 m_planes[i]->IProductWRTDerivBase(dir, tmpIn + i * nP_pts,
605 tmp1 = outarray + i * nP_coeffs);
606 }
607 }
608 // Do homogeneous derivative
609 else if (dir == 2)
610 {
611 if (m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier ||
612 m_homogeneousBasis->GetBasisType() ==
614 m_homogeneousBasis->GetBasisType() ==
616 m_homogeneousBasis->GetBasisType() ==
618 {
619 // TODO Check that nT_coeffs is correct for HalfMode
620 Array<OneD, NekDouble> tmpHom(nT_coeffs, 0.0);
621
622 NekDouble sign = -1.0;
624
625 // Half Modes
626 if (m_homogeneousBasis->GetBasisType() ==
628 {
629 // Do IProduct with 2D basis
630 m_planes[0]->IProductWRTBase(tmpIn, tmpHom);
631
632 // Add fourier coefficient
633 beta = sign * 2 * M_PI * (m_transposition->GetK(0)) / m_lhom;
634 Vmath::Smul(nT_coeffs, beta, tmpHom, 1, tmpHom, 1);
635 }
636 else if (m_homogeneousBasis->GetBasisType() ==
638 {
639 // Do IProduct with 2D basis
640 m_planes[0]->IProductWRTBase(tmpIn, tmpHom);
641
642 // Add fourier coefficient
643 beta = -sign * 2 * M_PI * (m_transposition->GetK(0)) / m_lhom;
644 Vmath::Smul(nT_coeffs, beta, tmpHom, 1, tmpHom, 1);
645 }
646 // Fully complex (Fourier or FourierSingleMode)
647 else
648 {
649 Array<OneD, NekDouble> tmpHomTwo(nT_coeffs, 0.0);
650 for (int i = 0; i < m_planes.size(); i++)
651 {
652 if (i != 1 || m_transposition->GetK(i) !=
653 0) // Mean + Null mode unchanged
654 {
655
656 // Do IProduct with 2D basis ie \int_\Omega \phi_{pq}
657 // inarray[ndim]
658 m_planes[i]->IProductWRTBase(tmpIn + i * nP_pts,
659 tmp1 = tmpHomTwo +
660 i * nP_coeffs);
661
662 // Add fourier coefficient (switch Real and Imaginary
663 // parts ie planes) Multiply by i changes real and
664 // imaginary parts
665 beta = sign * 2 * M_PI * (m_transposition->GetK(i)) /
666 m_lhom;
668 nP_coeffs, beta, tmp1 = tmpHomTwo + i * nP_coeffs,
669 1, tmp2 = tmpHom + (i - int(sign)) * nP_coeffs, 1);
670 }
671 sign = -1.0 * sign;
672 }
673 }
674
675 // Add last term to innerproduct
676 Vmath::Vcopy(nT_coeffs, tmpHom, 1, outarray, 1);
677 }
678 else
679 {
681 "The IProductWRTDerivBase routine with one homogeneous "
682 "direction is not implemented for the homogeneous basis "
683 "if it is is not of type Fourier "
684 "or FourierSingleMode/HalfModeRe/HalfModeIm");
685 }
686 }
687 else
688 {
689 NEKERROR(
691 "cannot handle direction give to IProductWRTDerivBase operator."
692 "dir != 0, 1 or 2.");
693 }
694}
695
696/**
697 * Inner product with respect to the derivative of the Basis
698 * including one homogeneous direction.
699 *
700 * This function evaluates \f[\int_\Omega \nabla \phi \cdot \mathbf{a}\f],
701 * where \f$\phi\f$ are the basis functions (including a homogeneous basis)
702 * and \f$\mathbf{a}\f$ is a generic vector \f$\mathbf{a}\f$, which is fourier
703 * decomposed into modes of wavenumber \f$k\f$.
704 *
705 * Assuming a 3DH1D basis the generic vector \f$\mathbf{a}\f$ is defined as
706 * \f[\mathbf{a}=\sum_k\begin{bmatrix}a_x^k(x,y)\\a_y^k(x,y)\\a_z^k(x,y)\end{bmatrix}e^{i
707 * \beta k z}\f].
708 *
709 */
711 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
712 Array<OneD, NekDouble> &outarray)
713{
714 int ndim = inarray.size(); // Dimension including homogeneous direction
715 // (e.g. 3DH1D, ndim=3D)
716 int nT_coeffs = m_planes[0]->GetNcoeffs() * m_planes.size();
717
718 // Write x-direction into outarray
719 v_IProductWRTDerivBase(0, inarray[0], outarray);
720
721 // Add other directions on top
722 // requires temporary array
723 Array<OneD, NekDouble> tmp(nT_coeffs, 0.0);
724 for (int dir = 1; dir < ndim; dir++)
725 {
726 v_IProductWRTDerivBase(dir, inarray[dir], tmp);
727 Vmath::Vadd(nT_coeffs, tmp, 1, outarray, 1, outarray, 1);
728 }
729}
730
731/**
732 * Homogeneous transform Bwd/Fwd (MVM and FFT)
733 */
735 const int num_dofs, const Array<OneD, const NekDouble> &inarray,
736 Array<OneD, NekDouble> &outarray, bool IsForwards, bool Shuff, bool UnShuff)
737{
738
739 if (m_useFFT)
740 {
741 int num_points_per_plane = num_dofs / m_planes.size();
742 int num_dfts_per_proc;
743 if (!m_session->DefinesSolverInfo("HomoStrip"))
744 {
745 int nP = m_comm->GetColumnComm()->GetSize();
746 num_dfts_per_proc =
747 num_points_per_plane / nP + (num_points_per_plane % nP > 0);
748 }
749 else
750 {
751 int nP = m_StripZcomm->GetSize();
752 num_dfts_per_proc =
753 num_points_per_plane / nP + (num_points_per_plane % nP > 0);
754 }
755
757 num_dfts_per_proc * m_homogeneousBasis->GetNumPoints(), 0.0);
759 num_dfts_per_proc * m_homogeneousBasis->GetNumPoints(), 0.0);
760
761 if (Shuff)
762 {
763 m_transposition->Transpose(num_dofs, inarray, fft_in, false,
765 }
766 else
767 {
768 Vmath::Vcopy(num_dfts_per_proc * m_homogeneousBasis->GetNumPoints(),
769 inarray, 1, fft_in, 1);
770 }
771
772 if (IsForwards)
773 {
774 for (int i = 0; i < num_dfts_per_proc; i++)
775 {
776 m_FFT->FFTFwdTrans(
777 m_tmpIN = fft_in + i * m_homogeneousBasis->GetNumPoints(),
778 m_tmpOUT =
779 fft_out + i * m_homogeneousBasis->GetNumPoints());
780 }
781 }
782 else
783 {
784 for (int i = 0; i < num_dfts_per_proc; i++)
785 {
786 m_FFT->FFTBwdTrans(
787 m_tmpIN = fft_in + i * m_homogeneousBasis->GetNumPoints(),
788 m_tmpOUT =
789 fft_out + i * m_homogeneousBasis->GetNumPoints());
790 }
791 }
792
793 if (UnShuff)
794 {
795 m_transposition->Transpose(num_dofs, fft_out, outarray, false,
797 }
798 else
799 {
800 Vmath::Vcopy(num_dfts_per_proc * m_homogeneousBasis->GetNumPoints(),
801 fft_out, 1, outarray, 1);
802 }
803 }
804 else
805 {
806 DNekBlkMatSharedPtr blkmat;
807
808 if (num_dofs == m_npoints) // transform phys space
809 {
810 if (IsForwards)
811 {
813 }
814 else
815 {
817 }
818 }
819 else
820 {
821 if (IsForwards)
822 {
824 }
825 else
826 {
828 }
829 }
830
831 int nrows = blkmat->GetRows();
832 int ncols = blkmat->GetColumns();
833
834 Array<OneD, NekDouble> sortedinarray(ncols, 0.0);
835 Array<OneD, NekDouble> sortedoutarray(nrows, 0.0);
836
837 if (Shuff)
838 {
839 m_transposition->Transpose(num_dofs, inarray, sortedinarray,
840 !IsForwards, LibUtilities::eXYtoZ);
841 }
842 else
843 {
844 Vmath::Vcopy(ncols, inarray, 1, sortedinarray, 1);
845 }
846
847 // Create NekVectors from the given data arrays
848 NekVector<NekDouble> in(ncols, sortedinarray, eWrapper);
849 NekVector<NekDouble> out(nrows, sortedoutarray, eWrapper);
850
851 // Perform matrix-vector multiply.
852 out = (*blkmat) * in;
853
854 if (UnShuff)
855 {
856 m_transposition->Transpose(num_dofs, sortedoutarray, outarray,
857 IsForwards, LibUtilities::eZtoXY);
858 }
859 else
860 {
861 Vmath::Vcopy(nrows, sortedoutarray, 1, outarray, 1);
862 }
863 }
864}
865
867 Homogeneous1DMatType mattype) const
868{
869 auto matrixIter = m_homogeneous1DBlockMat->find(mattype);
870
871 if (matrixIter == m_homogeneous1DBlockMat->end())
872 {
873 return ((*m_homogeneous1DBlockMat)[mattype] =
875 }
876 else
877 {
878 return matrixIter->second;
879 }
880}
881
883 Homogeneous1DMatType mattype) const
884{
885 DNekMatSharedPtr loc_mat;
886 DNekBlkMatSharedPtr BlkMatrix;
887 int n_exp = 0;
888 int num_trans_per_proc = 0;
889
890 if ((mattype == eForwardsCoeffSpace1D) ||
891 (mattype == eBackwardsCoeffSpace1D)) // will operate on coeffs
892 {
893 n_exp = m_planes[0]->GetNcoeffs();
894 }
895 else
896 {
897 n_exp = m_planes[0]->GetTotPoints(); // will operatore on phys
898 }
899
900 num_trans_per_proc = n_exp / m_comm->GetColumnComm()->GetSize() +
901 (n_exp % m_comm->GetColumnComm()->GetSize() > 0);
902
903 Array<OneD, unsigned int> nrows(num_trans_per_proc);
904 Array<OneD, unsigned int> ncols(num_trans_per_proc);
905
906 if ((mattype == eForwardsCoeffSpace1D) || (mattype == eForwardsPhysSpace1D))
907 {
908 nrows = Array<OneD, unsigned int>(num_trans_per_proc,
909 m_homogeneousBasis->GetNumModes());
910 ncols = Array<OneD, unsigned int>(num_trans_per_proc,
911 m_homogeneousBasis->GetNumPoints());
912 }
913 else
914 {
915 nrows = Array<OneD, unsigned int>(num_trans_per_proc,
916 m_homogeneousBasis->GetNumPoints());
917 ncols = Array<OneD, unsigned int>(num_trans_per_proc,
918 m_homogeneousBasis->GetNumModes());
919 }
920
921 MatrixStorage blkmatStorage = eDIAGONAL;
922 BlkMatrix = MemoryManager<DNekBlkMat>::AllocateSharedPtr(nrows, ncols,
923 blkmatStorage);
924
925 // Half Mode
926 if (m_homogeneousBasis->GetBasisType() ==
929 {
930 StdRegions::StdPointExp StdPoint(m_homogeneousBasis->GetBasisKey());
931
932 if ((mattype == eForwardsCoeffSpace1D) ||
933 (mattype == eForwardsPhysSpace1D))
934 {
936 StdPoint.DetShapeType(), StdPoint);
937
938 loc_mat = StdPoint.GetStdMatrix(matkey);
939 }
940 else
941 {
943 StdPoint.DetShapeType(), StdPoint);
944
945 loc_mat = StdPoint.GetStdMatrix(matkey);
946 }
947 }
948 // other cases
949 else
950 {
951 StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
952
953 if ((mattype == eForwardsCoeffSpace1D) ||
954 (mattype == eForwardsPhysSpace1D))
955 {
957 StdSeg.DetShapeType(), StdSeg);
958
959 loc_mat = StdSeg.GetStdMatrix(matkey);
960 }
961 else
962 {
964 StdSeg.DetShapeType(), StdSeg);
965
966 loc_mat = StdSeg.GetStdMatrix(matkey);
967 }
968 }
969
970 // set up array of block matrices.
971 for (int i = 0; i < num_trans_per_proc; ++i)
972 {
973 BlkMatrix->SetBlock(i, i, loc_mat);
974 }
975
976 return BlkMatrix;
977}
978
979std::vector<LibUtilities::FieldDefinitionsSharedPtr> ExpListHomogeneous1D::
981{
982 std::vector<LibUtilities::FieldDefinitionsSharedPtr> returnval;
983
984 // Set up Homogeneous length details.
986
987 std::vector<NekDouble> HomoLen;
988 HomoLen.push_back(m_lhom);
989
990 std::vector<unsigned int> StripsIDs;
991
992 bool strips;
993 m_session->MatchSolverInfo("HomoStrip", "True", strips, false);
994 if (strips)
995 {
996 StripsIDs.push_back(m_transposition->GetStripID());
997 }
998
999 std::vector<unsigned int> PlanesIDs;
1000 int IDoffset = 0;
1001
1002 // introduce a 2 plane offset for single mode case so can
1003 // be post-processed or used in MultiMode expansion.
1005 {
1006 IDoffset = 2;
1007 }
1008
1009 for (int i = 0; i < m_planes.size(); i++)
1010 {
1011 PlanesIDs.push_back(m_transposition->GetPlaneID(i) + IDoffset);
1012 }
1013
1014 m_planes[0]->GeneralGetFieldDefinitions(returnval, 1, HomoBasis, HomoLen,
1015 strips, StripsIDs, PlanesIDs);
1016 return returnval;
1017}
1018
1020 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
1021{
1022 // Set up Homogeneous length details.
1024
1025 std::vector<NekDouble> HomoLen;
1026 HomoLen.push_back(m_lhom);
1027
1028 std::vector<unsigned int> StripsIDs;
1029
1030 bool strips;
1031 m_session->MatchSolverInfo("HomoStrip", "True", strips, false);
1032 if (strips)
1033 {
1034 StripsIDs.push_back(m_transposition->GetStripID());
1035 }
1036
1037 std::vector<unsigned int> PlanesIDs;
1038 int IDoffset = 0;
1039
1041 {
1042 IDoffset = 2;
1043 }
1044
1045 for (int i = 0; i < m_planes.size(); i++)
1046 {
1047 PlanesIDs.push_back(m_transposition->GetPlaneID(i) + IDoffset);
1048 }
1049
1050 // enforce NumHomoDir == 1 by direct call
1051 m_planes[0]->GeneralGetFieldDefinitions(fielddef, 1, HomoBasis, HomoLen,
1052 strips, StripsIDs, PlanesIDs);
1053}
1054
1055/** This routine appends the data from the expansion list into
1056 the output format where each element is given by looping
1057 over its Fourier modes where as data in the expandion is
1058 stored with all consecutive elements and then the Fourier
1059 modes
1060 */
1063 std::vector<NekDouble> &fielddata, Array<OneD, NekDouble> &coeffs)
1064{
1065 int i, n;
1066 int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
1067
1068 // Determine mapping from element ids to location in
1069 // expansion list
1070 if (m_elmtToExpId.size() == 0)
1071 {
1072 for (i = 0; i < m_planes[0]->GetExpSize(); ++i)
1073 {
1074 m_elmtToExpId[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
1075 }
1076 }
1077
1078 for (i = 0; i < fielddef->m_elementIDs.size(); ++i)
1079 {
1080 int eid = m_elmtToExpId[fielddef->m_elementIDs[i]];
1081 int datalen = (*m_exp)[eid]->GetNcoeffs();
1082
1083 for (n = 0; n < m_planes.size(); ++n)
1084 {
1085 fielddata.insert(
1086 fielddata.end(),
1087 &coeffs[m_coeff_offset[eid] + n * ncoeffs_per_plane],
1088 &coeffs[m_coeff_offset[eid] + n * ncoeffs_per_plane] + datalen);
1089 }
1090 }
1091}
1092
1095 std::vector<NekDouble> &fielddata)
1096{
1097 v_AppendFieldData(fielddef, fielddata, m_coeffs);
1098}
1099
1100// Extract the data in fielddata into the m_coeff list
1103 std::vector<NekDouble> &fielddata, std::string &field,
1104 Array<OneD, NekDouble> &coeffs, std::unordered_map<int, int> zIdToPlane)
1105{
1106 int i, n;
1107 int offset = 0;
1108 int nzmodes = 1;
1109 int datalen = fielddata.size() / fielddef->m_fields.size();
1110 std::vector<unsigned int> fieldDefHomoZids;
1111
1112 // Find data location according to field definition
1113 for (i = 0; i < fielddef->m_fields.size(); ++i)
1114 {
1115 if (fielddef->m_fields[i] == field)
1116 {
1117 break;
1118 }
1119 offset += datalen;
1120 }
1121
1122 if (i == fielddef->m_fields.size())
1123 {
1124 cout << "Field " << field << "not found in data file. " << endl;
1125 }
1126 else
1127 {
1128
1129 int modes_offset = 0;
1130 int planes_offset = 0;
1131 Array<OneD, NekDouble> coeff_tmp;
1132
1133 // Build map of plane IDs lying on this processor and determine
1134 // mapping from element ids to location in expansion list.
1135 if (zIdToPlane.size() == 0)
1136 {
1137 int IDoffset = m_homogeneousBasis->GetBasisType() ==
1139 ? 2
1140 : 0;
1141 for (i = 0; i < m_planes.size(); ++i)
1142 {
1143 zIdToPlane[m_transposition->GetPlaneID(i) + IDoffset] = i;
1144 }
1145 }
1146
1147 if (m_elmtToExpId.size() == 0)
1148 {
1149 for (i = 0; i < m_planes[0]->GetExpSize(); ++i)
1150 {
1151 m_elmtToExpId[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
1152 }
1153 }
1154
1155 if (fielddef->m_numHomogeneousDir)
1156 {
1157 nzmodes = fielddef->m_homogeneousZIDs.size();
1158 fieldDefHomoZids = fielddef->m_homogeneousZIDs;
1159 }
1160 else // input file is 2D and so set nzmodes to 1
1161 {
1162 nzmodes = 1;
1163 fieldDefHomoZids.push_back(0);
1164 }
1165
1166 // calculate number of modes in the current partition
1167 int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
1168
1169 for (i = 0; i < fielddef->m_elementIDs.size(); ++i)
1170 {
1171 if (fielddef->m_uniOrder == true) // reset modes_offset to zero
1172 {
1173 modes_offset = 0;
1174 }
1175
1177 fielddef->m_shapeType, fielddef->m_numModes, modes_offset);
1178
1179 auto it = m_elmtToExpId.find(fielddef->m_elementIDs[i]);
1180
1181 // ensure element is on this partition for parallel case.
1182 if (it == m_elmtToExpId.end())
1183 {
1184 // increase offset for correct FieldData access
1185 offset += datalen * nzmodes;
1186 modes_offset +=
1187 (*m_exp)[0]->GetNumBases() + fielddef->m_numHomogeneousDir;
1188 continue;
1189 }
1190
1191 int eid = it->second;
1192 bool sameBasis = true;
1193 for (int j = 0; j < fielddef->m_basis.size() - 1; ++j)
1194 {
1195 if (fielddef->m_basis[j] != (*m_exp)[eid]->GetBasisType(j))
1196 {
1197 sameBasis = false;
1198 break;
1199 }
1200 }
1201
1202 for (n = 0; n < nzmodes; ++n, offset += datalen)
1203 {
1204
1205 it = zIdToPlane.find(fieldDefHomoZids[n]);
1206
1207 // Check to make sure this mode number lies in this field.
1208 if (it == zIdToPlane.end())
1209 {
1210 continue;
1211 }
1212
1213 planes_offset = it->second;
1214 if (datalen == (*m_exp)[eid]->GetNcoeffs() && sameBasis)
1215 {
1216 Vmath::Vcopy(datalen, &fielddata[offset], 1,
1217 &coeffs[m_coeff_offset[eid] +
1218 planes_offset * ncoeffs_per_plane],
1219 1);
1220 }
1221 else // unpack data to new order
1222 {
1223 (*m_exp)[eid]->ExtractDataToCoeffs(
1224 &fielddata[offset], fielddef->m_numModes, modes_offset,
1225 &coeffs[m_coeff_offset[eid] +
1226 planes_offset * ncoeffs_per_plane],
1227 fielddef->m_basis);
1228 }
1229 }
1230 modes_offset +=
1231 (*m_exp)[0]->GetNumBases() + fielddef->m_numHomogeneousDir;
1232 }
1233 }
1234}
1235
1236// Extract the data in fielddata into the m_coeff list
1238 const std::shared_ptr<ExpList> &fromExpList,
1239 const Array<OneD, const NekDouble> &fromCoeffs,
1240 Array<OneD, NekDouble> &toCoeffs)
1241{
1242 int i;
1243 int fromNcoeffs_per_plane = fromExpList->GetPlane(0)->GetNcoeffs();
1244 int toNcoeffs_per_plane = m_planes[0]->GetNcoeffs();
1245 Array<OneD, NekDouble> tocoeffs_tmp, fromcoeffs_tmp;
1246
1247 for (i = 0; i < m_planes.size(); ++i)
1248 {
1249 m_planes[i]->ExtractCoeffsToCoeffs(
1250 fromExpList->GetPlane(i),
1251 fromcoeffs_tmp = fromCoeffs + fromNcoeffs_per_plane * i,
1252 tocoeffs_tmp = toCoeffs + toNcoeffs_per_plane * i);
1253 }
1254}
1255
1257 int expansion, std::string var)
1258{
1259 // If there is only one plane (e.g. HalfMode), we write a 2D plane.
1260 if (m_planes.size() == 1)
1261 {
1262 m_planes[0]->WriteVtkPieceData(outfile, expansion, var);
1263 return;
1264 }
1265
1266 int i;
1267 int nq = (*m_exp)[expansion]->GetTotPoints();
1268 int npoints_per_plane = m_planes[0]->GetTotPoints();
1269
1270 // If we are using Fourier points, output extra plane to fill domain
1271 int outputExtraPlane = 0;
1272 Array<OneD, NekDouble> extraPlane;
1273 if (m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier &&
1274 m_homogeneousBasis->GetPointsType() ==
1276 {
1277 outputExtraPlane = 1;
1278 // Get extra plane data
1279 if (m_StripZcomm->GetSize() == 1)
1280 {
1281 extraPlane = m_phys + m_phys_offset[expansion];
1282 }
1283 else
1284 {
1285 // Determine to and from rank for communication
1286 int size = m_StripZcomm->GetSize();
1287 int rank = m_StripZcomm->GetRank();
1288 int fromRank = (rank + 1) % size;
1289 int toRank = (rank == 0) ? size - 1 : rank - 1;
1290 // Communicate using SendRecv
1291 extraPlane = Array<OneD, NekDouble>(nq);
1292 Array<OneD, NekDouble> send(nq, m_phys + m_phys_offset[expansion]);
1293
1294 m_StripZcomm->SendRecv(toRank, send, fromRank, extraPlane);
1295 }
1296 }
1297
1298 // printing the fields of that zone
1299 outfile << " <DataArray type=\"Float64\" Name=\"" << var << "\">"
1300 << endl;
1301 outfile << " ";
1302 for (int n = 0; n < m_planes.size(); ++n)
1303 {
1304 const Array<OneD, NekDouble> phys =
1305 m_phys + m_phys_offset[expansion] + n * npoints_per_plane;
1306
1307 for (i = 0; i < nq; ++i)
1308 {
1309 outfile << (fabs(phys[i]) < NekConstants::kNekZeroTol ? 0 : phys[i])
1310 << " ";
1311 }
1312 }
1313 if (outputExtraPlane)
1314 {
1315 for (i = 0; i < nq; ++i)
1316 {
1317 outfile << (fabs(extraPlane[i]) < NekConstants::kNekZeroTol
1318 ? 0
1319 : extraPlane[i])
1320 << " ";
1321 }
1322 }
1323 outfile << endl;
1324 outfile << " </DataArray>" << endl;
1325}
1326
1328 const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1329 Array<OneD, NekDouble> &outarray)
1330{
1331 int cnt, cnt1;
1332 Array<OneD, NekDouble> tmparray;
1333 cnt = m_planes[0]->GetTotPoints();
1334 cnt1 = m_planes[0]->Get1DScaledTotPoints(scale);
1335
1336 ASSERTL1(m_planes.size() * cnt1 <= outarray.size(),
1337 "size of outarray does not match internal estimage");
1338
1339 for (int i = 0; i < m_planes.size(); i++)
1340 {
1341
1342 m_planes[i]->PhysInterp1DScaled(scale, inarray + i * cnt,
1343 tmparray = outarray + i * cnt1);
1344 }
1345}
1346
1348 const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1349 Array<OneD, NekDouble> &outarray)
1350{
1351 int cnt, cnt1;
1352 Array<OneD, NekDouble> tmparray;
1353 cnt = m_planes[0]->Get1DScaledTotPoints(scale);
1354 cnt1 = m_planes[0]->GetTotPoints();
1355 ASSERTL1(m_planes.size() * cnt <= inarray.size(),
1356 "size of outarray does not match internal estimage");
1357
1358 for (int i = 0; i < m_planes.size(); i++)
1359 {
1360 m_planes[i]->PhysGalerkinProjection1DScaled(
1361 scale, inarray + i * cnt, tmparray = outarray + i * cnt1);
1362 }
1363}
1367{
1368 int nT_pts = m_npoints; // number of total points = n. of
1369 // Fourier points * n. of points per plane (nT_pts)
1370 int nP_pts =
1371 nT_pts / m_planes.size(); // number of points per plane = n of
1372 // Fourier transform required (nP_pts)
1373
1374 Array<OneD, NekDouble> temparray(nT_pts);
1375 Array<OneD, NekDouble> outarray(nT_pts);
1379
1380 for (int i = 0; i < m_planes.size(); i++)
1381 {
1382 m_planes[i]->PhysDeriv(inarray + i * nP_pts, tmp2 = out_d0 + i * nP_pts,
1383 tmp3 = out_d1 + i * nP_pts);
1384 }
1385
1386 if (out_d2 != NullNekDouble1DArray)
1387 {
1388 if (m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier ||
1389 m_homogeneousBasis->GetBasisType() ==
1391 m_homogeneousBasis->GetBasisType() ==
1393 m_homogeneousBasis->GetBasisType() ==
1395 {
1396 if (m_WaveSpace)
1397 {
1398 temparray = inarray;
1399 }
1400 else
1401 {
1402 HomogeneousFwdTrans(nT_pts, inarray, temparray);
1403 }
1404
1405 NekDouble sign = -1.0;
1407
1408 // Half Mode
1409 if (m_homogeneousBasis->GetBasisType() ==
1411 {
1412 beta = sign * 2 * M_PI * (m_transposition->GetK(0)) / m_lhom;
1413
1414 Vmath::Smul(nP_pts, beta, temparray, 1, outarray, 1);
1415 }
1416 else if (m_homogeneousBasis->GetBasisType() ==
1418 {
1419 beta = -sign * 2 * M_PI * (m_transposition->GetK(0)) / m_lhom;
1420
1421 Vmath::Smul(nP_pts, beta, temparray, 1, outarray, 1);
1422 }
1423
1424 // Fully complex
1425 else
1426 {
1427 for (int i = 0; i < m_planes.size(); i++)
1428 {
1429 beta =
1430 -sign * 2 * M_PI * (m_transposition->GetK(i)) / m_lhom;
1431
1432 Vmath::Smul(nP_pts, beta, tmp1 = temparray + i * nP_pts, 1,
1433 tmp2 = outarray + (i - int(sign)) * nP_pts, 1);
1434
1435 sign = -1.0 * sign;
1436 }
1437 }
1438
1439 if (m_WaveSpace)
1440 {
1441 out_d2 = outarray;
1442 }
1443 else
1444 {
1445 HomogeneousBwdTrans(nT_pts, outarray, out_d2);
1446 }
1447 }
1448 else
1449 {
1450 if (!m_session->DefinesSolverInfo("HomoStrip"))
1451 {
1452 ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,
1453 "Parallelisation in the homogeneous direction "
1454 "implemented just for Fourier basis");
1455 }
1456 else
1457 {
1458 ASSERTL0(m_StripZcomm->GetSize() == 1,
1459 "Parallelisation in the homogeneous direction "
1460 "implemented just for Fourier basis");
1461 }
1462
1463 if (m_WaveSpace)
1464 {
1465 ASSERTL0(false, "Semi-phyisical time-stepping not "
1466 "implemented yet for non-Fourier "
1467 "basis");
1468 }
1469 else
1470 {
1471 StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
1472
1473 m_transposition->Transpose(nT_pts, inarray, temparray, false,
1475
1476 for (int i = 0; i < nP_pts; i++)
1477 {
1478 StdSeg.PhysDeriv(temparray + i * m_planes.size(),
1479 tmp2 = outarray + i * m_planes.size());
1480 }
1481
1482 m_transposition->Transpose(nT_pts, outarray, out_d2, false,
1484
1485 Vmath::Smul(nT_pts, 2.0 / m_lhom, out_d2, 1, out_d2, 1);
1486 }
1487 }
1488 }
1489}
1490
1492 Direction edir, const Array<OneD, const NekDouble> &inarray,
1494
1495{
1496 int nT_pts = m_npoints; // number of total points = n. of Fourier
1497 // points * n. of points per plane (nT_pts)
1498 int nP_pts =
1499 nT_pts / m_planes.size(); // number of points per plane = n of
1500 // Fourier transform required (nP_pts)
1501
1502 int dir = (int)edir;
1503
1504 Array<OneD, NekDouble> temparray(nT_pts);
1505 Array<OneD, NekDouble> outarray(nT_pts);
1508
1509 if (dir < 2)
1510 {
1511 for (int i = 0; i < m_planes.size(); i++)
1512 {
1513 m_planes[i]->PhysDeriv(edir, inarray + i * nP_pts,
1514 tmp2 = out_d + i * nP_pts);
1515 }
1516 }
1517 else
1518 {
1519 if (m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier ||
1520 m_homogeneousBasis->GetBasisType() ==
1522 m_homogeneousBasis->GetBasisType() ==
1524 m_homogeneousBasis->GetBasisType() ==
1526 {
1527 if (m_WaveSpace)
1528 {
1529 temparray = inarray;
1530 }
1531 else
1532 {
1533 HomogeneousFwdTrans(nT_pts, inarray, temparray);
1534 }
1535
1536 NekDouble sign = -1.0;
1538
1539 // HalfMode
1540 if (m_homogeneousBasis->GetBasisType() ==
1542 {
1543 beta = 2 * sign * M_PI * (m_transposition->GetK(0)) / m_lhom;
1544
1545 Vmath::Smul(nP_pts, beta, temparray, 1, outarray, 1);
1546 }
1547 else if (m_homogeneousBasis->GetBasisType() ==
1549 {
1550 beta = -2 * sign * M_PI * (m_transposition->GetK(0)) / m_lhom;
1551
1552 Vmath::Smul(nP_pts, beta, temparray, 1, outarray, 1);
1553 }
1554 // Fully complex
1555 else
1556 {
1557 for (int i = 0; i < m_planes.size(); i++)
1558 {
1559 beta =
1560 -sign * 2 * M_PI * (m_transposition->GetK(i)) / m_lhom;
1561
1562 Vmath::Smul(nP_pts, beta, tmp1 = temparray + i * nP_pts, 1,
1563 tmp2 = outarray + (i - int(sign)) * nP_pts, 1);
1564
1565 sign = -1.0 * sign;
1566 }
1567 }
1568 if (m_WaveSpace)
1569 {
1570 out_d = outarray;
1571 }
1572 else
1573 {
1574 HomogeneousBwdTrans(nT_pts, outarray, out_d);
1575 }
1576 }
1577 else
1578 {
1579 if (!m_session->DefinesSolverInfo("HomoStrip"))
1580 {
1581 ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,
1582 "Parallelisation in the homogeneous direction "
1583 "implemented just for Fourier basis");
1584 }
1585 else
1586 {
1587 ASSERTL0(m_StripZcomm->GetSize() == 1,
1588 "Parallelisation in the homogeneous direction "
1589 "implemented just for Fourier basis");
1590 }
1591
1592 if (m_WaveSpace)
1593 {
1594 ASSERTL0(false, "Semi-phyisical time-stepping not implemented "
1595 "yet for non-Fourier basis");
1596 }
1597 else
1598 {
1599 StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
1600
1601 m_transposition->Transpose(nT_pts, inarray, temparray, false,
1603
1604 for (int i = 0; i < nP_pts; i++)
1605 {
1606 StdSeg.PhysDeriv(temparray + i * m_planes.size(),
1607 tmp2 = outarray + i * m_planes.size());
1608 }
1609
1610 m_transposition->Transpose(nT_pts, outarray, out_d, false,
1612
1613 Vmath::Smul(nT_pts, 2.0 / m_lhom, out_d, 1, out_d, 1);
1614 }
1615 }
1616 }
1617}
1618
1622
1623{
1624 v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
1625}
1626
1628 Direction edir, const Array<OneD, const NekDouble> &inarray,
1630{
1631 v_PhysDeriv(edir, inarray, out_d);
1632}
1633
1635 void)
1636{
1637 return m_transposition;
1638}
1639
1641{
1642 return m_lhom;
1643}
1644
1646{
1647 m_lhom = lhom;
1648}
1649
1651{
1652 return m_transposition->GetPlanesIDs();
1653}
1654
1656 const Array<OneD, const NekDouble> &inarray)
1657{
1658 NekDouble val = 0.0;
1659 int i = 0;
1660
1661 for (i = 0; i < (*m_exp).size(); ++i)
1662 {
1663 val += (*m_exp)[i]->Integral(inarray + m_phys_offset[i]);
1664 }
1665 val *= m_lhom / m_homogeneousBasis->GetNumModes();
1666
1667 m_comm->AllReduce(val, LibUtilities::ReduceSum);
1668
1669 return val;
1670}
1671} // namespace MultiRegions
1672} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:49
Describes the specification for a Basis.
Definition: Basis.h:47
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
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 PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void) override
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void) override
NekDouble m_lhom
Width of homogeneous direction.
LibUtilities::TranspositionSharedPtr m_transposition
virtual 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.
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var) override
virtual 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
virtual void v_FwdTransLocalElmt(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Array< OneD, ExpListSharedPtr > m_planes
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata) override
virtual void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
LibUtilities::NektarFFTSharedPtr m_FFT_deal
virtual void v_HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true) override
virtual void v_ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs) override
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
virtual 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
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
LibUtilities::NektarFFTSharedPtr m_FFT
virtual NekDouble v_GetHomoLen(void) override
void Homogeneous1DTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool IsForwards, bool Shuff=true, bool UnShuff=true)
virtual void v_HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true) override
virtual void v_SetHomoLen(const NekDouble lhom) override
DNekBlkMatSharedPtr GenHomogeneous1DBlockMatrix(Homogeneous1DMatType mattype) const
DNekBlkMatSharedPtr GetHomogeneous1DBlockMatrix(Homogeneous1DMatType mattype) const
virtual Array< OneD, const unsigned int > v_GetZIDs(void) override
virtual void v_DealiasedProd(const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray) override
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
ExpListHomogeneous1D(const ExpansionType type)
Default constructor.
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:102
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1092
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1527
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1136
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1065
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1127
void HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1872
std::unordered_map< int, int > m_elmtToExpId
Mapping from geometry ID of element to index inside m_exp.
Definition: ExpList.h:1148
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1067
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1138
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1108
void HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1863
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:609
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
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:855
Class representing a segment element in reference space All interface of this class sits in StdExpans...
Definition: StdSegExp.h:48
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:352
BasisManagerT & BasisManager(void)
NektarFFTFactory & GetNektarFFTFactory()
Definition: NektarFFT.cpp:69
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:186
int GetNumberOfCoefficients(ShapeType shape, std::vector< unsigned int > &modes, int offset=0)
Definition: ShapeType.hpp:310
std::shared_ptr< Transposition > TranspositionSharedPtr
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:61
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:76
@ eFourierSingleMode
Fourier ModifiedExpansion with just the first mode .
Definition: BasisType.h:66
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:70
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:68
@ eFourier
Fourier Expansion .
Definition: BasisType.h:57
std::map< Homogeneous1DMatType, DNekBlkMatSharedPtr > Homo1DBlockMatrixMap
A map between homo matrix keys and their associated block matrices.
static const NekDouble kNekZeroTol
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:77
static Array< OneD, NekDouble > NullNekDouble1DArray
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.cpp:207
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:569
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.cpp:354
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.cpp:245
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191