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
36 #include <LocalRegions/Expansion.h>
39 #include <StdRegions/StdPointExp.h>
40 #include <StdRegions/StdSegExp.h>
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46 namespace MultiRegions
47 {
48 // Forward declaration for typedefs
49 ExpListHomogeneous1D::ExpListHomogeneous1D(const ExpansionType type)
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  {
86  ASSERTL0(
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,
140  const Collections::ImplementationType ImpType)
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,
286  Array<OneD, Array<OneD, NekDouble>> &outarray)
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 
356  Array<OneD, Array<OneD, NekDouble>> ShufV1V2(nvec);
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  * Homogeneous transform Bwd/Fwd (MVM and FFT)
560  */
562  const int num_dofs, const Array<OneD, const NekDouble> &inarray,
563  Array<OneD, NekDouble> &outarray, bool IsForwards, bool Shuff, bool UnShuff)
564 {
565 
566  if (m_useFFT)
567  {
568  int num_points_per_plane = num_dofs / m_planes.size();
569  int num_dfts_per_proc;
570  if (!m_session->DefinesSolverInfo("HomoStrip"))
571  {
572  int nP = m_comm->GetColumnComm()->GetSize();
573  num_dfts_per_proc =
574  num_points_per_plane / nP + (num_points_per_plane % nP > 0);
575  }
576  else
577  {
578  int nP = m_StripZcomm->GetSize();
579  num_dfts_per_proc =
580  num_points_per_plane / nP + (num_points_per_plane % nP > 0);
581  }
582 
583  Array<OneD, NekDouble> fft_in(
584  num_dfts_per_proc * m_homogeneousBasis->GetNumPoints(), 0.0);
585  Array<OneD, NekDouble> fft_out(
586  num_dfts_per_proc * m_homogeneousBasis->GetNumPoints(), 0.0);
587 
588  if (Shuff)
589  {
590  m_transposition->Transpose(num_dofs, inarray, fft_in, false,
592  }
593  else
594  {
595  Vmath::Vcopy(num_dfts_per_proc * m_homogeneousBasis->GetNumPoints(),
596  inarray, 1, fft_in, 1);
597  }
598 
599  if (IsForwards)
600  {
601  for (int i = 0; i < num_dfts_per_proc; i++)
602  {
603  m_FFT->FFTFwdTrans(
604  m_tmpIN = fft_in + i * m_homogeneousBasis->GetNumPoints(),
605  m_tmpOUT =
606  fft_out + i * m_homogeneousBasis->GetNumPoints());
607  }
608  }
609  else
610  {
611  for (int i = 0; i < num_dfts_per_proc; i++)
612  {
613  m_FFT->FFTBwdTrans(
614  m_tmpIN = fft_in + i * m_homogeneousBasis->GetNumPoints(),
615  m_tmpOUT =
616  fft_out + i * m_homogeneousBasis->GetNumPoints());
617  }
618  }
619 
620  if (UnShuff)
621  {
622  m_transposition->Transpose(num_dofs, fft_out, outarray, false,
624  }
625  else
626  {
627  Vmath::Vcopy(num_dfts_per_proc * m_homogeneousBasis->GetNumPoints(),
628  fft_out, 1, outarray, 1);
629  }
630  }
631  else
632  {
633  DNekBlkMatSharedPtr blkmat;
634 
635  if (num_dofs == m_npoints) // transform phys space
636  {
637  if (IsForwards)
638  {
640  }
641  else
642  {
644  }
645  }
646  else
647  {
648  if (IsForwards)
649  {
651  }
652  else
653  {
655  }
656  }
657 
658  int nrows = blkmat->GetRows();
659  int ncols = blkmat->GetColumns();
660 
661  Array<OneD, NekDouble> sortedinarray(ncols, 0.0);
662  Array<OneD, NekDouble> sortedoutarray(nrows, 0.0);
663 
664  if (Shuff)
665  {
666  m_transposition->Transpose(num_dofs, inarray, sortedinarray,
667  !IsForwards, LibUtilities::eXYtoZ);
668  }
669  else
670  {
671  Vmath::Vcopy(ncols, inarray, 1, sortedinarray, 1);
672  }
673 
674  // Create NekVectors from the given data arrays
675  NekVector<NekDouble> in(ncols, sortedinarray, eWrapper);
676  NekVector<NekDouble> out(nrows, sortedoutarray, eWrapper);
677 
678  // Perform matrix-vector multiply.
679  out = (*blkmat) * in;
680 
681  if (UnShuff)
682  {
683  m_transposition->Transpose(num_dofs, sortedoutarray, outarray,
684  IsForwards, LibUtilities::eZtoXY);
685  }
686  else
687  {
688  Vmath::Vcopy(nrows, sortedoutarray, 1, outarray, 1);
689  }
690  }
691 }
692 
694  Homogeneous1DMatType mattype) const
695 {
696  auto matrixIter = m_homogeneous1DBlockMat->find(mattype);
697 
698  if (matrixIter == m_homogeneous1DBlockMat->end())
699  {
700  return ((*m_homogeneous1DBlockMat)[mattype] =
701  GenHomogeneous1DBlockMatrix(mattype));
702  }
703  else
704  {
705  return matrixIter->second;
706  }
707 }
708 
710  Homogeneous1DMatType mattype) const
711 {
712  DNekMatSharedPtr loc_mat;
713  DNekBlkMatSharedPtr BlkMatrix;
714  int n_exp = 0;
715  int num_trans_per_proc = 0;
716 
717  if ((mattype == eForwardsCoeffSpace1D) ||
718  (mattype == eBackwardsCoeffSpace1D)) // will operate on coeffs
719  {
720  n_exp = m_planes[0]->GetNcoeffs();
721  }
722  else
723  {
724  n_exp = m_planes[0]->GetTotPoints(); // will operatore on phys
725  }
726 
727  num_trans_per_proc = n_exp / m_comm->GetColumnComm()->GetSize() +
728  (n_exp % m_comm->GetColumnComm()->GetSize() > 0);
729 
730  Array<OneD, unsigned int> nrows(num_trans_per_proc);
731  Array<OneD, unsigned int> ncols(num_trans_per_proc);
732 
733  if ((mattype == eForwardsCoeffSpace1D) || (mattype == eForwardsPhysSpace1D))
734  {
735  nrows = Array<OneD, unsigned int>(num_trans_per_proc,
736  m_homogeneousBasis->GetNumModes());
737  ncols = Array<OneD, unsigned int>(num_trans_per_proc,
738  m_homogeneousBasis->GetNumPoints());
739  }
740  else
741  {
742  nrows = Array<OneD, unsigned int>(num_trans_per_proc,
743  m_homogeneousBasis->GetNumPoints());
744  ncols = Array<OneD, unsigned int>(num_trans_per_proc,
745  m_homogeneousBasis->GetNumModes());
746  }
747 
748  MatrixStorage blkmatStorage = eDIAGONAL;
749  BlkMatrix = MemoryManager<DNekBlkMat>::AllocateSharedPtr(nrows, ncols,
750  blkmatStorage);
751 
752  // Half Mode
753  if (m_homogeneousBasis->GetBasisType() ==
756  {
757  StdRegions::StdPointExp StdPoint(m_homogeneousBasis->GetBasisKey());
758 
759  if ((mattype == eForwardsCoeffSpace1D) ||
760  (mattype == eForwardsPhysSpace1D))
761  {
763  StdPoint.DetShapeType(), StdPoint);
764 
765  loc_mat = StdPoint.GetStdMatrix(matkey);
766  }
767  else
768  {
770  StdPoint.DetShapeType(), StdPoint);
771 
772  loc_mat = StdPoint.GetStdMatrix(matkey);
773  }
774  }
775  // other cases
776  else
777  {
778  StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
779 
780  if ((mattype == eForwardsCoeffSpace1D) ||
781  (mattype == eForwardsPhysSpace1D))
782  {
784  StdSeg.DetShapeType(), StdSeg);
785 
786  loc_mat = StdSeg.GetStdMatrix(matkey);
787  }
788  else
789  {
791  StdSeg.DetShapeType(), StdSeg);
792 
793  loc_mat = StdSeg.GetStdMatrix(matkey);
794  }
795  }
796 
797  // set up array of block matrices.
798  for (int i = 0; i < num_trans_per_proc; ++i)
799  {
800  BlkMatrix->SetBlock(i, i, loc_mat);
801  }
802 
803  return BlkMatrix;
804 }
805 
806 std::vector<LibUtilities::FieldDefinitionsSharedPtr> ExpListHomogeneous1D::
808 {
809  std::vector<LibUtilities::FieldDefinitionsSharedPtr> returnval;
810 
811  // Set up Homogeneous length details.
813 
814  std::vector<NekDouble> HomoLen;
815  HomoLen.push_back(m_lhom);
816 
817  std::vector<unsigned int> StripsIDs;
818 
819  bool strips;
820  m_session->MatchSolverInfo("HomoStrip", "True", strips, false);
821  if (strips)
822  {
823  StripsIDs.push_back(m_transposition->GetStripID());
824  }
825 
826  std::vector<unsigned int> PlanesIDs;
827  int IDoffset = 0;
828 
829  // introduce a 2 plane offset for single mode case so can
830  // be post-processed or used in MultiMode expansion.
832  {
833  IDoffset = 2;
834  }
835 
836  for (int i = 0; i < m_planes.size(); i++)
837  {
838  PlanesIDs.push_back(m_transposition->GetPlaneID(i) + IDoffset);
839  }
840 
841  m_planes[0]->GeneralGetFieldDefinitions(returnval, 1, HomoBasis, HomoLen,
842  strips, StripsIDs, PlanesIDs);
843  return returnval;
844 }
845 
847  std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
848 {
849  // Set up Homogeneous length details.
851 
852  std::vector<NekDouble> HomoLen;
853  HomoLen.push_back(m_lhom);
854 
855  std::vector<unsigned int> StripsIDs;
856 
857  bool strips;
858  m_session->MatchSolverInfo("HomoStrip", "True", strips, false);
859  if (strips)
860  {
861  StripsIDs.push_back(m_transposition->GetStripID());
862  }
863 
864  std::vector<unsigned int> PlanesIDs;
865  int IDoffset = 0;
866 
868  {
869  IDoffset = 2;
870  }
871 
872  for (int i = 0; i < m_planes.size(); i++)
873  {
874  PlanesIDs.push_back(m_transposition->GetPlaneID(i) + IDoffset);
875  }
876 
877  // enforce NumHomoDir == 1 by direct call
878  m_planes[0]->GeneralGetFieldDefinitions(fielddef, 1, HomoBasis, HomoLen,
879  strips, StripsIDs, PlanesIDs);
880 }
881 
882 /** This routine appends the data from the expansion list into
883  the output format where each element is given by looping
884  over its Fourier modes where as data in the expandion is
885  stored with all consecutive elements and then the Fourier
886  modes
887  */
890  std::vector<NekDouble> &fielddata, Array<OneD, NekDouble> &coeffs)
891 {
892  int i, n;
893  int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
894 
895  // Determine mapping from element ids to location in
896  // expansion list
897  if (m_elmtToExpId.size() == 0)
898  {
899  for (i = 0; i < m_planes[0]->GetExpSize(); ++i)
900  {
901  m_elmtToExpId[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
902  }
903  }
904 
905  for (i = 0; i < fielddef->m_elementIDs.size(); ++i)
906  {
907  int eid = m_elmtToExpId[fielddef->m_elementIDs[i]];
908  int datalen = (*m_exp)[eid]->GetNcoeffs();
909 
910  for (n = 0; n < m_planes.size(); ++n)
911  {
912  fielddata.insert(
913  fielddata.end(),
914  &coeffs[m_coeff_offset[eid] + n * ncoeffs_per_plane],
915  &coeffs[m_coeff_offset[eid] + n * ncoeffs_per_plane] + datalen);
916  }
917  }
918 }
919 
922  std::vector<NekDouble> &fielddata)
923 {
924  v_AppendFieldData(fielddef, fielddata, m_coeffs);
925 }
926 
927 // Extract the data in fielddata into the m_coeff list
930  std::vector<NekDouble> &fielddata, std::string &field,
931  Array<OneD, NekDouble> &coeffs, std::unordered_map<int, int> zIdToPlane)
932 {
933  int i, n;
934  int offset = 0;
935  int nzmodes = 1;
936  int datalen = fielddata.size() / fielddef->m_fields.size();
937  std::vector<unsigned int> fieldDefHomoZids;
938 
939  // Find data location according to field definition
940  for (i = 0; i < fielddef->m_fields.size(); ++i)
941  {
942  if (fielddef->m_fields[i] == field)
943  {
944  break;
945  }
946  offset += datalen;
947  }
948 
949  if (i == fielddef->m_fields.size())
950  {
951  cout << "Field " << field << "not found in data file. " << endl;
952  }
953  else
954  {
955 
956  int modes_offset = 0;
957  int planes_offset = 0;
958  Array<OneD, NekDouble> coeff_tmp;
959 
960  // Build map of plane IDs lying on this processor and determine
961  // mapping from element ids to location in expansion list.
962  if (zIdToPlane.size() == 0)
963  {
964  int IDoffset = m_homogeneousBasis->GetBasisType() ==
966  ? 2
967  : 0;
968  for (i = 0; i < m_planes.size(); ++i)
969  {
970  zIdToPlane[m_transposition->GetPlaneID(i) + IDoffset] = i;
971  }
972  }
973 
974  if (m_elmtToExpId.size() == 0)
975  {
976  for (i = 0; i < m_planes[0]->GetExpSize(); ++i)
977  {
978  m_elmtToExpId[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
979  }
980  }
981 
982  if (fielddef->m_numHomogeneousDir)
983  {
984  nzmodes = fielddef->m_homogeneousZIDs.size();
985  fieldDefHomoZids = fielddef->m_homogeneousZIDs;
986  }
987  else // input file is 2D and so set nzmodes to 1
988  {
989  nzmodes = 1;
990  fieldDefHomoZids.push_back(0);
991  }
992 
993  // calculate number of modes in the current partition
994  int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
995 
996  for (i = 0; i < fielddef->m_elementIDs.size(); ++i)
997  {
998  if (fielddef->m_uniOrder == true) // reset modes_offset to zero
999  {
1000  modes_offset = 0;
1001  }
1002 
1004  fielddef->m_shapeType, fielddef->m_numModes, modes_offset);
1005 
1006  auto it = m_elmtToExpId.find(fielddef->m_elementIDs[i]);
1007 
1008  // ensure element is on this partition for parallel case.
1009  if (it == m_elmtToExpId.end())
1010  {
1011  // increase offset for correct FieldData access
1012  offset += datalen * nzmodes;
1013  modes_offset +=
1014  (*m_exp)[0]->GetNumBases() + fielddef->m_numHomogeneousDir;
1015  continue;
1016  }
1017 
1018  int eid = it->second;
1019  bool sameBasis = true;
1020  for (int j = 0; j < fielddef->m_basis.size() - 1; ++j)
1021  {
1022  if (fielddef->m_basis[j] != (*m_exp)[eid]->GetBasisType(j))
1023  {
1024  sameBasis = false;
1025  break;
1026  }
1027  }
1028 
1029  for (n = 0; n < nzmodes; ++n, offset += datalen)
1030  {
1031 
1032  it = zIdToPlane.find(fieldDefHomoZids[n]);
1033 
1034  // Check to make sure this mode number lies in this field.
1035  if (it == zIdToPlane.end())
1036  {
1037  continue;
1038  }
1039 
1040  planes_offset = it->second;
1041  if (datalen == (*m_exp)[eid]->GetNcoeffs() && sameBasis)
1042  {
1043  Vmath::Vcopy(datalen, &fielddata[offset], 1,
1044  &coeffs[m_coeff_offset[eid] +
1045  planes_offset * ncoeffs_per_plane],
1046  1);
1047  }
1048  else // unpack data to new order
1049  {
1050  (*m_exp)[eid]->ExtractDataToCoeffs(
1051  &fielddata[offset], fielddef->m_numModes, modes_offset,
1052  &coeffs[m_coeff_offset[eid] +
1053  planes_offset * ncoeffs_per_plane],
1054  fielddef->m_basis);
1055  }
1056  }
1057  modes_offset +=
1058  (*m_exp)[0]->GetNumBases() + fielddef->m_numHomogeneousDir;
1059  }
1060  }
1061 }
1062 
1063 // Extract the data in fielddata into the m_coeff list
1065  const std::shared_ptr<ExpList> &fromExpList,
1066  const Array<OneD, const NekDouble> &fromCoeffs,
1067  Array<OneD, NekDouble> &toCoeffs)
1068 {
1069  int i;
1070  int fromNcoeffs_per_plane = fromExpList->GetPlane(0)->GetNcoeffs();
1071  int toNcoeffs_per_plane = m_planes[0]->GetNcoeffs();
1072  Array<OneD, NekDouble> tocoeffs_tmp, fromcoeffs_tmp;
1073 
1074  for (i = 0; i < m_planes.size(); ++i)
1075  {
1076  m_planes[i]->ExtractCoeffsToCoeffs(
1077  fromExpList->GetPlane(i),
1078  fromcoeffs_tmp = fromCoeffs + fromNcoeffs_per_plane * i,
1079  tocoeffs_tmp = toCoeffs + toNcoeffs_per_plane * i);
1080  }
1081 }
1082 
1084  int expansion, std::string var)
1085 {
1086  // If there is only one plane (e.g. HalfMode), we write a 2D plane.
1087  if (m_planes.size() == 1)
1088  {
1089  m_planes[0]->WriteVtkPieceData(outfile, expansion, var);
1090  return;
1091  }
1092 
1093  int i;
1094  int nq = (*m_exp)[expansion]->GetTotPoints();
1095  int npoints_per_plane = m_planes[0]->GetTotPoints();
1096 
1097  // If we are using Fourier points, output extra plane to fill domain
1098  int outputExtraPlane = 0;
1099  Array<OneD, NekDouble> extraPlane;
1100  if (m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier &&
1101  m_homogeneousBasis->GetPointsType() ==
1103  {
1104  outputExtraPlane = 1;
1105  // Get extra plane data
1106  if (m_StripZcomm->GetSize() == 1)
1107  {
1108  extraPlane = m_phys + m_phys_offset[expansion];
1109  }
1110  else
1111  {
1112  // Determine to and from rank for communication
1113  int size = m_StripZcomm->GetSize();
1114  int rank = m_StripZcomm->GetRank();
1115  int fromRank = (rank + 1) % size;
1116  int toRank = (rank == 0) ? size - 1 : rank - 1;
1117  // Communicate using SendRecv
1118  extraPlane = Array<OneD, NekDouble>(nq);
1119  Array<OneD, NekDouble> send(nq, m_phys + m_phys_offset[expansion]);
1120 
1121  m_StripZcomm->SendRecv(toRank, send, fromRank, extraPlane);
1122  }
1123  }
1124 
1125  // printing the fields of that zone
1126  outfile << " <DataArray type=\"Float64\" Name=\"" << var << "\">"
1127  << endl;
1128  outfile << " ";
1129  for (int n = 0; n < m_planes.size(); ++n)
1130  {
1131  const Array<OneD, NekDouble> phys =
1132  m_phys + m_phys_offset[expansion] + n * npoints_per_plane;
1133 
1134  for (i = 0; i < nq; ++i)
1135  {
1136  outfile << (fabs(phys[i]) < NekConstants::kNekZeroTol ? 0 : phys[i])
1137  << " ";
1138  }
1139  }
1140  if (outputExtraPlane)
1141  {
1142  for (i = 0; i < nq; ++i)
1143  {
1144  outfile << (fabs(extraPlane[i]) < NekConstants::kNekZeroTol
1145  ? 0
1146  : extraPlane[i])
1147  << " ";
1148  }
1149  }
1150  outfile << endl;
1151  outfile << " </DataArray>" << endl;
1152 }
1153 
1155  const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1156  Array<OneD, NekDouble> &outarray)
1157 {
1158  int cnt, cnt1;
1159  Array<OneD, NekDouble> tmparray;
1160  cnt = m_planes[0]->GetTotPoints();
1161  cnt1 = m_planes[0]->Get1DScaledTotPoints(scale);
1162 
1163  ASSERTL1(m_planes.size() * cnt1 <= outarray.size(),
1164  "size of outarray does not match internal estimage");
1165 
1166  for (int i = 0; i < m_planes.size(); i++)
1167  {
1168 
1169  m_planes[i]->PhysInterp1DScaled(scale, inarray + i * cnt,
1170  tmparray = outarray + i * cnt1);
1171  }
1172 }
1173 
1175  const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1176  Array<OneD, NekDouble> &outarray)
1177 {
1178  int cnt, cnt1;
1179  Array<OneD, NekDouble> tmparray;
1180  cnt = m_planes[0]->Get1DScaledTotPoints(scale);
1181  cnt1 = m_planes[0]->GetTotPoints();
1182  ASSERTL1(m_planes.size() * cnt <= inarray.size(),
1183  "size of outarray does not match internal estimage");
1184 
1185  for (int i = 0; i < m_planes.size(); i++)
1186  {
1187  m_planes[i]->PhysGalerkinProjection1DScaled(
1188  scale, inarray + i * cnt, tmparray = outarray + i * cnt1);
1189  }
1190 }
1192  const Array<OneD, const NekDouble> &inarray, Array<OneD, NekDouble> &out_d0,
1194 {
1195  int nT_pts = m_npoints; // number of total points = n. of
1196  // Fourier points * n. of points per plane (nT_pts)
1197  int nP_pts =
1198  nT_pts / m_planes.size(); // number of points per plane = n of
1199  // Fourier transform required (nP_pts)
1200 
1201  Array<OneD, NekDouble> temparray(nT_pts);
1202  Array<OneD, NekDouble> outarray(nT_pts);
1206 
1207  for (int i = 0; i < m_planes.size(); i++)
1208  {
1209  m_planes[i]->PhysDeriv(inarray + i * nP_pts, tmp2 = out_d0 + i * nP_pts,
1210  tmp3 = out_d1 + i * nP_pts);
1211  }
1212 
1213  if (out_d2 != NullNekDouble1DArray)
1214  {
1215  if (m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier ||
1216  m_homogeneousBasis->GetBasisType() ==
1218  m_homogeneousBasis->GetBasisType() ==
1220  m_homogeneousBasis->GetBasisType() ==
1222  {
1223  if (m_WaveSpace)
1224  {
1225  temparray = inarray;
1226  }
1227  else
1228  {
1229  HomogeneousFwdTrans(nT_pts, inarray, temparray);
1230  }
1231 
1232  NekDouble sign = -1.0;
1233  NekDouble beta;
1234 
1235  // Half Mode
1236  if (m_homogeneousBasis->GetBasisType() ==
1238  {
1239  beta = sign * 2 * M_PI * (m_transposition->GetK(0)) / m_lhom;
1240 
1241  Vmath::Smul(nP_pts, beta, temparray, 1, outarray, 1);
1242  }
1243  else if (m_homogeneousBasis->GetBasisType() ==
1245  {
1246  beta = -sign * 2 * M_PI * (m_transposition->GetK(0)) / m_lhom;
1247 
1248  Vmath::Smul(nP_pts, beta, temparray, 1, outarray, 1);
1249  }
1250 
1251  // Fully complex
1252  else
1253  {
1254  for (int i = 0; i < m_planes.size(); i++)
1255  {
1256  beta =
1257  -sign * 2 * M_PI * (m_transposition->GetK(i)) / m_lhom;
1258 
1259  Vmath::Smul(nP_pts, beta, tmp1 = temparray + i * nP_pts, 1,
1260  tmp2 = outarray + (i - int(sign)) * nP_pts, 1);
1261 
1262  sign = -1.0 * sign;
1263  }
1264  }
1265 
1266  if (m_WaveSpace)
1267  {
1268  out_d2 = outarray;
1269  }
1270  else
1271  {
1272  HomogeneousBwdTrans(nT_pts, outarray, out_d2);
1273  }
1274  }
1275  else
1276  {
1277  if (!m_session->DefinesSolverInfo("HomoStrip"))
1278  {
1279  ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,
1280  "Parallelisation in the homogeneous direction "
1281  "implemented just for Fourier basis");
1282  }
1283  else
1284  {
1285  ASSERTL0(m_StripZcomm->GetSize() == 1,
1286  "Parallelisation in the homogeneous direction "
1287  "implemented just for Fourier basis");
1288  }
1289 
1290  if (m_WaveSpace)
1291  {
1292  ASSERTL0(false, "Semi-phyisical time-stepping not "
1293  "implemented yet for non-Fourier "
1294  "basis");
1295  }
1296  else
1297  {
1298  StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
1299 
1300  m_transposition->Transpose(nT_pts, inarray, temparray, false,
1302 
1303  for (int i = 0; i < nP_pts; i++)
1304  {
1305  StdSeg.PhysDeriv(temparray + i * m_planes.size(),
1306  tmp2 = outarray + i * m_planes.size());
1307  }
1308 
1309  m_transposition->Transpose(nT_pts, outarray, out_d2, false,
1311 
1312  Vmath::Smul(nT_pts, 2.0 / m_lhom, out_d2, 1, out_d2, 1);
1313  }
1314  }
1315  }
1316 }
1317 
1319  Direction edir, const Array<OneD, const NekDouble> &inarray,
1320  Array<OneD, NekDouble> &out_d)
1321 
1322 {
1323  int nT_pts = m_npoints; // number of total points = n. of Fourier
1324  // points * n. of points per plane (nT_pts)
1325  int nP_pts =
1326  nT_pts / m_planes.size(); // number of points per plane = n of
1327  // Fourier transform required (nP_pts)
1328 
1329  int dir = (int)edir;
1330 
1331  Array<OneD, NekDouble> temparray(nT_pts);
1332  Array<OneD, NekDouble> outarray(nT_pts);
1335 
1336  if (dir < 2)
1337  {
1338  for (int i = 0; i < m_planes.size(); i++)
1339  {
1340  m_planes[i]->PhysDeriv(edir, inarray + i * nP_pts,
1341  tmp2 = out_d + i * nP_pts);
1342  }
1343  }
1344  else
1345  {
1346  if (m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier ||
1347  m_homogeneousBasis->GetBasisType() ==
1349  m_homogeneousBasis->GetBasisType() ==
1351  m_homogeneousBasis->GetBasisType() ==
1353  {
1354  if (m_WaveSpace)
1355  {
1356  temparray = inarray;
1357  }
1358  else
1359  {
1360  HomogeneousFwdTrans(nT_pts, inarray, temparray);
1361  }
1362 
1363  NekDouble sign = -1.0;
1364  NekDouble beta;
1365 
1366  // HalfMode
1367  if (m_homogeneousBasis->GetBasisType() ==
1369  {
1370  beta = 2 * sign * M_PI * (m_transposition->GetK(0)) / m_lhom;
1371 
1372  Vmath::Smul(nP_pts, beta, temparray, 1, outarray, 1);
1373  }
1374  else if (m_homogeneousBasis->GetBasisType() ==
1376  {
1377  beta = -2 * sign * M_PI * (m_transposition->GetK(0)) / m_lhom;
1378 
1379  Vmath::Smul(nP_pts, beta, temparray, 1, outarray, 1);
1380  }
1381  // Fully complex
1382  else
1383  {
1384  for (int i = 0; i < m_planes.size(); i++)
1385  {
1386  beta =
1387  -sign * 2 * M_PI * (m_transposition->GetK(i)) / m_lhom;
1388 
1389  Vmath::Smul(nP_pts, beta, tmp1 = temparray + i * nP_pts, 1,
1390  tmp2 = outarray + (i - int(sign)) * nP_pts, 1);
1391 
1392  sign = -1.0 * sign;
1393  }
1394  }
1395  if (m_WaveSpace)
1396  {
1397  out_d = outarray;
1398  }
1399  else
1400  {
1401  HomogeneousBwdTrans(nT_pts, outarray, out_d);
1402  }
1403  }
1404  else
1405  {
1406  if (!m_session->DefinesSolverInfo("HomoStrip"))
1407  {
1408  ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,
1409  "Parallelisation in the homogeneous direction "
1410  "implemented just for Fourier basis");
1411  }
1412  else
1413  {
1414  ASSERTL0(m_StripZcomm->GetSize() == 1,
1415  "Parallelisation in the homogeneous direction "
1416  "implemented just for Fourier basis");
1417  }
1418 
1419  if (m_WaveSpace)
1420  {
1421  ASSERTL0(false, "Semi-phyisical time-stepping not implemented "
1422  "yet for non-Fourier basis");
1423  }
1424  else
1425  {
1426  StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
1427 
1428  m_transposition->Transpose(nT_pts, inarray, temparray, false,
1430 
1431  for (int i = 0; i < nP_pts; i++)
1432  {
1433  StdSeg.PhysDeriv(temparray + i * m_planes.size(),
1434  tmp2 = outarray + i * m_planes.size());
1435  }
1436 
1437  m_transposition->Transpose(nT_pts, outarray, out_d, false,
1439 
1440  Vmath::Smul(nT_pts, 2.0 / m_lhom, out_d, 1, out_d, 1);
1441  }
1442  }
1443  }
1444 }
1445 
1447  const Array<OneD, const NekDouble> &inarray, Array<OneD, NekDouble> &out_d0,
1449 
1450 {
1451  v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
1452 }
1453 
1455  Direction edir, const Array<OneD, const NekDouble> &inarray,
1456  Array<OneD, NekDouble> &out_d)
1457 {
1458  v_PhysDeriv(edir, inarray, out_d);
1459 }
1460 
1462  void)
1463 {
1464  return m_transposition;
1465 }
1466 
1468 {
1469  return m_lhom;
1470 }
1471 
1473 {
1474  m_lhom = lhom;
1475 }
1476 
1478 {
1479  return m_transposition->GetPlanesIDs();
1480 }
1481 
1483  const Array<OneD, const NekDouble> &inarray)
1484 {
1485  NekDouble val = 0.0;
1486  int i = 0;
1487 
1488  for (i = 0; i < (*m_exp).size(); ++i)
1489  {
1490  val += (*m_exp)[i]->Integral(inarray + m_phys_offset[i]);
1491  }
1492  val *= m_lhom / m_homogeneousBasis->GetNumModes();
1493 
1494  m_comm->AllReduce(val, LibUtilities::ReduceSum);
1495 
1496  return val;
1497 }
1498 } // namespace MultiRegions
1499 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#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:50
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_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_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
ExpListHomogeneous1D(const ExpansionType type)
Default constructor.
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
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:103
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1076
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1492
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1120
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1049
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1111
void HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1807
std::unordered_map< int, int > m_elmtToExpId
Mapping from geometry ID of element to index inside m_exp.
Definition: ExpList.h:1132
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1051
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1122
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1092
void HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1798
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:848
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:369
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:209
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:574
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:248
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255