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 
36 #include <LibUtilities/Foundations/ManagerAccess.h> // for PointsManager, etc
37 #include <StdRegions/StdSegExp.h>
38 #include <StdRegions/StdPointExp.h>
39 #include <LocalRegions/Expansion.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),
51  m_homogeneousBasis(LibUtilities::NullBasisSharedPtr),
52  m_lhom(1),
53  m_homogeneous1DBlockMat(MemoryManager<Homo1DBlockMatrixMap>::AllocateSharedPtr())
54  {
55  }
56 
58  const ExpansionType type,
60  const LibUtilities::BasisKey &HomoBasis,
61  const NekDouble lhom,
62  const bool useFFT,
63  const bool dealiasing):
64  ExpList(type),
65  m_useFFT(useFFT),
66  m_lhom(lhom),
67  m_homogeneous1DBlockMat(MemoryManager<Homo1DBlockMatrixMap>::AllocateSharedPtr()),
68  m_dealiasing(dealiasing)
69  {
70  m_session = pSession;
71  m_comm = pSession->GetComm();
72 
74  "Homogeneous Basis is a null basis");
75 
77 
78  m_StripZcomm = m_session->DefinesSolverInfo("HomoStrip") ?
79  m_comm->GetColumnComm()->GetColumnComm() :
80  m_comm->GetColumnComm();
81 
82  ASSERTL0( m_homogeneousBasis->GetNumPoints() %
83  m_StripZcomm->GetSize() == 0,
84  "HomModesZ should be a multiple of npz.");
85 
86  if ( (m_homogeneousBasis->GetBasisType() !=
88  && (m_homogeneousBasis->GetBasisType() !=
90  {
91  ASSERTL0(
92  (m_homogeneousBasis->GetNumPoints() /
93  m_StripZcomm->GetSize()) % 2 == 0,
94  "HomModesZ/npz should be an even integer.");
95  }
96 
98  ::AllocateSharedPtr(HomoBasis, m_comm,
99  m_StripZcomm);
100 
102  m_homogeneousBasis->GetNumPoints() /
103  m_StripZcomm->GetSize());
104 
105  if(m_useFFT)
106  {
108  "NekFFTW", m_homogeneousBasis->GetNumPoints());
109  }
110 
111  if(m_dealiasing)
112  {
113  if(m_useFFT)
114  {
115  int size = m_homogeneousBasis->GetNumPoints() +
116  m_homogeneousBasis->GetNumPoints() / 2;
117  m_padsize = size + (size % 2);
119  .CreateInstance("NekFFTW", m_padsize);
120  }
121  else
122  {
123  ASSERTL0(false, "Dealiasing available just in combination "
124  "with FFTW");
125  }
126  }
127  }
128 
129 
130  /**
131  * @param In ExpListHomogeneous1D object to copy.
132  */
134  ExpList(In,false),
135  m_transposition(In.m_transposition),
136  m_StripZcomm(In.m_StripZcomm),
137  m_useFFT(In.m_useFFT),
138  m_FFT(In.m_FFT),
139  m_tmpIN(In.m_tmpIN),
140  m_tmpOUT(In.m_tmpOUT),
141  m_homogeneousBasis(In.m_homogeneousBasis),
142  m_lhom(In.m_lhom),
143  m_homogeneous1DBlockMat(In.m_homogeneous1DBlockMat),
144  m_dealiasing(In.m_dealiasing),
145  m_padsize(In.m_padsize)
146  {
148  }
149 
151  const std::vector<unsigned int> &eIDs):
152  ExpList(In,eIDs,false),
153  m_transposition(In.m_transposition),
154  m_useFFT(In.m_useFFT),
155  m_FFT(In.m_FFT),
156  m_tmpIN(In.m_tmpIN),
157  m_tmpOUT(In.m_tmpOUT),
158  m_homogeneousBasis(In.m_homogeneousBasis),
159  m_lhom(In.m_lhom),
160  m_homogeneous1DBlockMat(MemoryManager<Homo1DBlockMatrixMap>::AllocateSharedPtr()),
161  m_dealiasing(In.m_dealiasing),
162  m_padsize(In.m_padsize)
163  {
165  }
166 
167  /**
168  * Destructor
169  */
171  {
172  }
173 
175  Array<OneD, NekDouble> &outarray,
176  bool Shuff,
177  bool UnShuff)
178  {
179  // Forwards trans
180  Homogeneous1DTrans(inarray,outarray,true,Shuff,UnShuff);
181  }
182 
184  Array<OneD, NekDouble> &outarray,
185  bool Shuff,
186  bool UnShuff)
187  {
188  // Backwards trans
189  Homogeneous1DTrans(inarray,outarray,false,Shuff,UnShuff);
190  }
191 
192  /**
193  * Dealiasing routine
194  *
195  * @param inarray1 First term of the product
196  * @param inarray2 Second term of the product
197  * @param outarray Dealiased product
198  */
200  const Array<OneD, NekDouble> &inarray2,
201  Array<OneD, NekDouble> &outarray)
202  {
203  int num_dofs = inarray1.size();
204  int N = m_homogeneousBasis->GetNumPoints();
205 
206  Array<OneD, NekDouble> V1(num_dofs);
207  Array<OneD, NekDouble> V2(num_dofs);
208  Array<OneD, NekDouble> V1V2(num_dofs);
209 
210  if(m_WaveSpace)
211  {
212  V1 = inarray1;
213  V2 = inarray2;
214  }
215  else
216  {
217  HomogeneousFwdTrans(inarray1,V1);
218  HomogeneousFwdTrans(inarray2,V2);
219  }
220 
221  int num_points_per_plane = num_dofs/m_planes.size();
222  int num_proc;
223  if(!m_session->DefinesSolverInfo("HomoStrip"))
224  {
225  num_proc = m_comm->GetColumnComm()->GetSize();
226  }
227  else
228  {
229  num_proc = m_StripZcomm->GetSize();
230  }
231  int num_dfts_per_proc = num_points_per_plane / num_proc
232  + (num_points_per_plane % num_proc > 0);
233 
234  Array<OneD, NekDouble> ShufV1(num_dfts_per_proc*N,0.0);
235  Array<OneD, NekDouble> ShufV2(num_dfts_per_proc*N,0.0);
236  Array<OneD, NekDouble> ShufV1V2(num_dfts_per_proc*N,0.0);
237 
238  Array<OneD, NekDouble> ShufV1_PAD_coef(m_padsize,0.0);
239  Array<OneD, NekDouble> ShufV2_PAD_coef(m_padsize,0.0);
240  Array<OneD, NekDouble> ShufV1_PAD_phys(m_padsize,0.0);
241  Array<OneD, NekDouble> ShufV2_PAD_phys(m_padsize,0.0);
242 
243  Array<OneD, NekDouble> ShufV1V2_PAD_coef(m_padsize,0.0);
244  Array<OneD, NekDouble> ShufV1V2_PAD_phys(m_padsize,0.0);
245 
246  m_transposition->Transpose(V1, ShufV1, false, LibUtilities::eXYtoZ);
247  m_transposition->Transpose(V2, ShufV2, false, LibUtilities::eXYtoZ);
248 
249  // Looping on the pencils
250  for(int i = 0 ; i < num_dfts_per_proc ; i++)
251  {
252  // Copying the i-th pencil pf lenght N into a bigger
253  // pencil of lenght 2N We are in Fourier space
254  Vmath::Vcopy(N, &(ShufV1[i*N]), 1, &(ShufV1_PAD_coef[0]), 1);
255  Vmath::Vcopy(N, &(ShufV2[i*N]), 1, &(ShufV2_PAD_coef[0]), 1);
256 
257  // Moving to physical space using the padded system
258  m_FFT_deal->FFTBwdTrans(ShufV1_PAD_coef, ShufV1_PAD_phys);
259  m_FFT_deal->FFTBwdTrans(ShufV2_PAD_coef, ShufV2_PAD_phys);
260 
261  // Perfroming the vectors multiplication in physical space on
262  // the padded system
263  Vmath::Vmul(m_padsize, ShufV1_PAD_phys, 1,
264  ShufV2_PAD_phys, 1,
265  ShufV1V2_PAD_phys, 1);
266 
267  // Moving back the result (V1*V2)_phys in Fourier space, padded
268  // system
269  m_FFT_deal->FFTFwdTrans(ShufV1V2_PAD_phys, ShufV1V2_PAD_coef);
270 
271  // Copying the first half of the padded pencil in the full
272  // vector (Fourier space)
273  Vmath::Vcopy(N, &(ShufV1V2_PAD_coef[0]), 1,
274  &(ShufV1V2[i*N]), 1);
275  }
276 
277  // Moving the results to the output
278  if (m_WaveSpace)
279  {
280  m_transposition->Transpose(ShufV1V2, outarray, false,
282  }
283  else
284  {
285  m_transposition->Transpose(ShufV1V2, V1V2, false,
287  HomogeneousBwdTrans(V1V2, outarray);
288  }
289  }
290 
291  /**
292  * Dealiasing routine for dot product
293  *
294  * @param inarray1 First term of the product with dimension ndim
295  * (e.g. U)
296  * @param inarray2 Second term of the product with dimension ndim*nvec
297  * (e.g. grad U)
298  * @param outarray Dealiased product with dimension nvec
299  */
301  const Array<OneD, Array<OneD, NekDouble> > &inarray1,
302  const Array<OneD, Array<OneD, NekDouble> > &inarray2,
303  Array<OneD, Array<OneD, NekDouble> > &outarray)
304  {
305  int ndim = inarray1.size();
306  ASSERTL1( inarray2.size() % ndim == 0,
307  "Wrong dimensions for DealiasedDotProd.");
308  int nvec = inarray2.size() / ndim;
309 
310  int num_dofs = inarray1[0].size();
311  int N = m_homogeneousBasis->GetNumPoints();
312 
313  int num_points_per_plane = num_dofs/m_planes.size();
314  int num_proc;
315  if(!m_session->DefinesSolverInfo("HomoStrip"))
316  {
317  num_proc = m_comm->GetColumnComm()->GetSize();
318  }
319  else
320  {
321  num_proc = m_StripZcomm->GetSize();
322  }
323  int num_dfts_per_proc = num_points_per_plane / num_proc
324  + (num_points_per_plane % num_proc > 0);
325 
326  // Get inputs in Fourier space
328  Array<OneD, Array<OneD, NekDouble> > V2(ndim*nvec);
329  if(m_WaveSpace)
330  {
331  for (int i = 0; i < ndim; i++)
332  {
333  V1[i] = inarray1[i];
334  }
335  for (int i = 0; i < ndim*nvec; i++)
336  {
337  V2[i] = inarray2[i];
338  }
339  }
340  else
341  {
342  for (int i = 0; i < ndim; i++)
343  {
344  V1[i] = Array<OneD, NekDouble> (num_dofs);
345  HomogeneousFwdTrans(inarray1[i],V1[i]);
346  }
347  for (int i = 0; i < ndim*nvec; i++)
348  {
349  V2[i] = Array<OneD, NekDouble> (num_dofs);
350  HomogeneousFwdTrans(inarray2[i],V2[i]);
351  }
352  }
353 
354  // Allocate variables for ffts
355  Array<OneD, Array<OneD, NekDouble> > ShufV1(ndim);
356  Array<OneD, NekDouble> ShufV1_PAD_coef(m_padsize,0.0);
357  Array<OneD, Array<OneD, NekDouble> > ShufV1_PAD_phys(ndim);
358  for (int i = 0; i < ndim; i++)
359  {
360  ShufV1[i] = Array<OneD, NekDouble>
361  (num_dfts_per_proc*N,0.0);
362  ShufV1_PAD_phys[i] = Array<OneD, NekDouble>
363  (m_padsize,0.0);
364  }
365 
366  Array<OneD, Array<OneD, NekDouble> > ShufV2(ndim*nvec);
367  Array<OneD, NekDouble> ShufV2_PAD_coef(m_padsize,0.0);
368  Array<OneD, Array<OneD, NekDouble> > ShufV2_PAD_phys(ndim*nvec);
369  for (int i = 0; i < ndim*nvec; i++)
370  {
371  ShufV2[i] = Array<OneD, NekDouble>
372  (num_dfts_per_proc*N,0.0);
373  ShufV2_PAD_phys[i] = Array<OneD, NekDouble>
374  (m_padsize,0.0);
375  }
376 
377  Array<OneD, Array<OneD, NekDouble> > ShufV1V2(nvec);
378  Array<OneD, NekDouble> ShufV1V2_PAD_coef(m_padsize,0.0);
379  Array<OneD, NekDouble> ShufV1V2_PAD_phys(m_padsize,0.0);
380  for (int i = 0; i < nvec; i++)
381  {
382  ShufV1V2[i] = Array<OneD, NekDouble>
383  (num_dfts_per_proc*N,0.0);
384  }
385 
386  // Transpose inputs
387  for (int i = 0; i < ndim; i++)
388  {
389  m_transposition->Transpose(V1[i], ShufV1[i], false,
391  }
392  for (int i = 0; i < ndim*nvec; i++)
393  {
394  m_transposition->Transpose(V2[i], ShufV2[i], false,
396  }
397 
398  // Looping on the pencils
399  for(int i = 0 ; i < num_dfts_per_proc ; i++)
400  {
401  for (int j = 0; j < ndim; j++)
402  {
403  // Copying the i-th pencil pf lenght N into a bigger
404  // pencil of lenght 1.5N We are in Fourier space
405  Vmath::Vcopy(N, &(ShufV1[j][i*N]), 1,
406  &(ShufV1_PAD_coef[0]), 1);
407  // Moving to physical space using the padded system
408  m_FFT_deal->FFTBwdTrans(ShufV1_PAD_coef, ShufV1_PAD_phys[j]);
409  }
410  for (int j = 0; j < ndim*nvec; j++)
411  {
412  Vmath::Vcopy(N, &(ShufV2[j][i*N]), 1,
413  &(ShufV2_PAD_coef[0]), 1);
414  m_FFT_deal->FFTBwdTrans(ShufV2_PAD_coef, ShufV2_PAD_phys[j]);
415  }
416 
417  // Performing the vectors multiplication in physical space on
418  // the padded system
419  for (int j = 0; j < nvec; j++)
420  {
421  Vmath::Zero(m_padsize, ShufV1V2_PAD_phys, 1);
422  for (int k = 0; k < ndim; k++)
423  {
424  Vmath::Vvtvp(m_padsize, ShufV1_PAD_phys[k], 1,
425  ShufV2_PAD_phys[j*ndim+k], 1,
426  ShufV1V2_PAD_phys, 1,
427  ShufV1V2_PAD_phys, 1);
428  }
429  // Moving back the result (V1*V2)_phys in Fourier space,
430  // padded system
431  m_FFT_deal->FFTFwdTrans(ShufV1V2_PAD_phys, ShufV1V2_PAD_coef);
432  // Copying the first half of the padded pencil in the full
433  // vector (Fourier space)
434  Vmath::Vcopy(N, &(ShufV1V2_PAD_coef[0]), 1,
435  &(ShufV1V2[j][i*N]), 1);
436  }
437  }
438 
439  // Moving the results to the output
440  if (m_WaveSpace)
441  {
442  for (int j = 0; j < nvec; j++)
443  {
444  m_transposition->Transpose(ShufV1V2[j], outarray[j],
445  false,
447  }
448  }
449  else
450  {
451  Array<OneD, NekDouble> V1V2(num_dofs);
452  for (int j = 0; j < nvec; j++)
453  {
454  m_transposition->Transpose(ShufV1V2[j], V1V2, false,
456  HomogeneousBwdTrans(V1V2, outarray[j]);
457  }
458  }
459  }
460 
461  /**
462  * Forward transform
463  */
465  {
466  int cnt = 0, cnt1 = 0;
467  Array<OneD, NekDouble> tmparray;
468 
469  for(int n = 0; n < m_planes.size(); ++n)
470  {
471  m_planes[n]->FwdTrans(inarray+cnt, tmparray = outarray + cnt1);
472  cnt += m_planes[n]->GetTotPoints();
473 
474  cnt1 += m_planes[n]->GetNcoeffs(); // need to skip ncoeffs
475  }
476  if(!m_WaveSpace)
477  {
478  HomogeneousFwdTrans(outarray,outarray);
479  }
480  }
481 
482  /**
483  * Forward transform element by element
484  */
486  const NekDouble> &inarray,
487  Array<OneD, NekDouble> &outarray)
488  {
489  int cnt = 0, cnt1 = 0;
490  Array<OneD, NekDouble> tmparray;
491 
492  //spectral element FwdTrans plane by plane
493  for(int n = 0; n < m_planes.size(); ++n)
494  {
495  m_planes[n]->FwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
496  cnt += m_planes[n]->GetTotPoints();
497  cnt1 += m_planes[n]->GetNcoeffs();
498  }
499  if(!m_WaveSpace)
500  {
501  HomogeneousFwdTrans(outarray,outarray);
502  }
503  }
504 
505  /**
506  * Forward transform element by element with boundaries constrained
507  */
509  const NekDouble> &inarray,
510  Array<OneD, NekDouble> &outarray)
511  {
512  int cnt = 0, cnt1 = 0;
513  Array<OneD, NekDouble> tmparray;
514 
515  //spectral element FwdTrans plane by plane
516  for(int n = 0; n < m_planes.size(); ++n)
517  {
518  m_planes[n]->FwdTrans_BndConstrained(inarray+cnt, tmparray = outarray + cnt1);
519  cnt += m_planes[n]->GetTotPoints();
520  cnt1 += m_planes[n]->GetNcoeffs();
521  }
522  if(!m_WaveSpace)
523  {
524  HomogeneousFwdTrans(outarray,outarray);
525  }
526  }
527 
528  /**
529  * Backward transform
530  */
532  {
533  int cnt = 0, cnt1 = 0;
534  Array<OneD, NekDouble> tmparray;
535 
536  for(int n = 0; n < m_planes.size(); ++n)
537  {
538  m_planes[n]->BwdTrans(inarray+cnt, tmparray = outarray + cnt1);
539  cnt += m_planes[n]->GetNcoeffs();
540  cnt1 += m_planes[n]->GetTotPoints();
541  }
542  if(!m_WaveSpace)
543  {
544  HomogeneousBwdTrans(outarray,outarray);
545  }
546  }
547 
548  /**
549  * Backward transform element by element
550  */
552  {
553  int cnt = 0, cnt1 = 0;
554  Array<OneD, NekDouble> tmparray;
555 
556  for(int n = 0; n < m_planes.size(); ++n)
557  {
558  m_planes[n]->BwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
559 
560  cnt += m_planes[n]->GetNcoeffs();
561  cnt1 += m_planes[n]->GetTotPoints();
562  }
563  if(!m_WaveSpace)
564  {
565  HomogeneousBwdTrans(outarray,outarray);
566  }
567  }
568 
569  /**
570  * Inner product
571  */
573  {
574  int cnt = 0, cnt1 = 0;
575  Array<OneD, NekDouble> tmparray, tmpIn;
576 
577  if(m_WaveSpace)
578  {
579  tmpIn = inarray;
580  }
581  else
582  {
583  tmpIn = Array<OneD, NekDouble> (inarray.size(), 0.0);
584  HomogeneousFwdTrans(inarray,tmpIn);
585  }
586 
587  for(int n = 0; n < m_planes.size(); ++n)
588  {
589  m_planes[n]->IProductWRTBase(tmpIn+cnt, tmparray = outarray + cnt1);
590 
591  cnt1 += m_planes[n]->GetNcoeffs();
592  cnt += m_planes[n]->GetTotPoints();
593  }
594  }
595 
596  /**
597  * Inner product element by element
598  */
600  {
601  int cnt = 0, cnt1 = 0;
602  Array<OneD, NekDouble> tmparray, tmpIn;
603 
604  if(m_WaveSpace)
605  {
606  tmpIn = inarray;
607  }
608  else
609  {
610  tmpIn = Array<OneD, NekDouble> (inarray.size(), 0.0);
611  HomogeneousFwdTrans(inarray,tmpIn);
612  }
613 
614  for(int n = 0; n < m_planes.size(); ++n)
615  {
616  m_planes[n]->IProductWRTBase_IterPerExp(tmpIn+cnt, tmparray = outarray + cnt1);
617 
618  cnt1 += m_planes[n]->GetNcoeffs();
619  cnt += m_planes[n]->GetTotPoints();
620  }
621  }
622 
623  /**
624  * Homogeneous transform Bwd/Fwd (MVM and FFT)
625  */
627  bool IsForwards,
628  bool Shuff,
629  bool UnShuff)
630  {
631  int num_dofs;
632 
633  if(IsForwards)
634  {
635  num_dofs = inarray.size();
636  }
637  else
638  {
639  num_dofs = outarray.size();
640  }
641 
642  if(m_useFFT)
643  {
644  int num_points_per_plane = num_dofs/m_planes.size();
645  int num_dfts_per_proc;
646  if(!m_session->DefinesSolverInfo("HomoStrip"))
647  {
648  int nP = m_comm->GetColumnComm()->GetSize();
649  num_dfts_per_proc = num_points_per_plane / nP
650  + (num_points_per_plane % nP > 0);
651  }
652  else
653  {
654  int nP = m_StripZcomm->GetSize();
655  num_dfts_per_proc = num_points_per_plane / nP
656  + (num_points_per_plane % nP > 0);
657  }
658 
659  Array<OneD, NekDouble> fft_in (num_dfts_per_proc*m_homogeneousBasis->GetNumPoints(),0.0);
660  Array<OneD, NekDouble> fft_out(num_dfts_per_proc*m_homogeneousBasis->GetNumPoints(),0.0);
661 
662  if(Shuff)
663  {
664  m_transposition->Transpose(inarray,fft_in,false,LibUtilities::eXYtoZ);
665  }
666  else
667  {
668  Vmath::Vcopy(num_dfts_per_proc*m_homogeneousBasis->GetNumPoints(),
669  inarray,1,fft_in,1);
670  }
671 
672  if(IsForwards)
673  {
674  for(int i = 0 ; i < num_dfts_per_proc ; i++)
675  {
676  m_FFT->FFTFwdTrans(m_tmpIN = fft_in + i*m_homogeneousBasis->GetNumPoints(), m_tmpOUT = fft_out + i*m_homogeneousBasis->GetNumPoints());
677  }
678  }
679  else
680  {
681  for(int i = 0 ; i < num_dfts_per_proc ; i++)
682  {
683  m_FFT->FFTBwdTrans(m_tmpIN = fft_in + i*m_homogeneousBasis->GetNumPoints(), m_tmpOUT = fft_out + i*m_homogeneousBasis->GetNumPoints());
684  }
685  }
686 
687  if(UnShuff)
688  {
689  m_transposition->Transpose(fft_out,outarray,false,LibUtilities::eZtoXY);
690  }
691  else
692  {
693  Vmath::Vcopy(num_dfts_per_proc*m_homogeneousBasis->GetNumPoints(),
694  fft_out,1,outarray,1);
695  }
696  }
697  else
698  {
699  DNekBlkMatSharedPtr blkmat;
700 
701  if(num_dofs == m_npoints) //transform phys space
702  {
703  if(IsForwards)
704  {
706  }
707  else
708  {
710  }
711  }
712  else
713  {
714  if(IsForwards)
715  {
717  }
718  else
719  {
721  }
722  }
723 
724  int nrows = blkmat->GetRows();
725  int ncols = blkmat->GetColumns();
726 
727  Array<OneD, NekDouble> sortedinarray(ncols,0.0);
728  Array<OneD, NekDouble> sortedoutarray(nrows,0.0);
729 
730  if(Shuff)
731  {
732  m_transposition->Transpose(inarray,sortedinarray,!IsForwards,LibUtilities::eXYtoZ);
733  }
734  else
735  {
736  Vmath::Vcopy(ncols,inarray,1,sortedinarray,1);
737  }
738 
739  // Create NekVectors from the given data arrays
740  NekVector<NekDouble> in (ncols,sortedinarray,eWrapper);
741  NekVector<NekDouble> out(nrows,sortedoutarray,eWrapper);
742 
743  // Perform matrix-vector multiply.
744  out = (*blkmat)*in;
745 
746  if(UnShuff)
747  {
748  m_transposition->Transpose(sortedoutarray,outarray,IsForwards,LibUtilities::eZtoXY);
749  }
750  else
751  {
752  Vmath::Vcopy(nrows,sortedoutarray,1,outarray,1);
753  }
754 
755  }
756  }
757 
759  {
760  auto matrixIter = m_homogeneous1DBlockMat->find(mattype);
761 
762  if(matrixIter == m_homogeneous1DBlockMat->end())
763  {
764  return ((*m_homogeneous1DBlockMat)[mattype] =
765  GenHomogeneous1DBlockMatrix(mattype));
766  }
767  else
768  {
769  return matrixIter->second;
770  }
771  }
772 
773 
775  {
776  DNekMatSharedPtr loc_mat;
777  DNekBlkMatSharedPtr BlkMatrix;
778  int n_exp = 0;
779  int num_trans_per_proc = 0;
780 
781  if((mattype == eForwardsCoeffSpace1D)
782  ||(mattype == eBackwardsCoeffSpace1D)) // will operate on m_coeffs
783  {
784  n_exp = m_planes[0]->GetNcoeffs();
785  }
786  else
787  {
788  n_exp = m_planes[0]->GetTotPoints(); // will operatore on m_phys
789  }
790 
791  num_trans_per_proc = n_exp/m_comm->GetColumnComm()->GetSize() + (n_exp%m_comm->GetColumnComm()->GetSize() > 0);
792 
793  Array<OneD,unsigned int> nrows(num_trans_per_proc);
794  Array<OneD,unsigned int> ncols(num_trans_per_proc);
795 
796  if((mattype == eForwardsCoeffSpace1D)||(mattype == eForwardsPhysSpace1D))
797  {
798  nrows = Array<OneD, unsigned int>(num_trans_per_proc,m_homogeneousBasis->GetNumModes());
799  ncols = Array<OneD, unsigned int>(num_trans_per_proc,m_homogeneousBasis->GetNumPoints());
800  }
801  else
802  {
803  nrows = Array<OneD, unsigned int>(num_trans_per_proc,m_homogeneousBasis->GetNumPoints());
804  ncols = Array<OneD, unsigned int>(num_trans_per_proc,m_homogeneousBasis->GetNumModes());
805  }
806 
807  MatrixStorage blkmatStorage = eDIAGONAL;
808  BlkMatrix = MemoryManager<DNekBlkMat>
809  ::AllocateSharedPtr(nrows,ncols,blkmatStorage);
810 
811  //Half Mode
813  {
814  StdRegions::StdPointExp StdPoint(m_homogeneousBasis->GetBasisKey());
815 
816  if((mattype == eForwardsCoeffSpace1D)||(mattype == eForwardsPhysSpace1D))
817  {
819  StdPoint.DetShapeType(),
820  StdPoint);
821 
822  loc_mat = StdPoint.GetStdMatrix(matkey);
823  }
824  else
825  {
827  StdPoint.DetShapeType(),
828  StdPoint);
829 
830  loc_mat = StdPoint.GetStdMatrix(matkey);
831  }
832  }
833  //other cases
834  else
835  {
836  StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
837 
838  if((mattype == eForwardsCoeffSpace1D)||(mattype == eForwardsPhysSpace1D))
839  {
841  StdSeg.DetShapeType(),
842  StdSeg);
843 
844  loc_mat = StdSeg.GetStdMatrix(matkey);
845  }
846  else
847  {
849  StdSeg.DetShapeType(),
850  StdSeg);
851 
852  loc_mat = StdSeg.GetStdMatrix(matkey);
853  }
854 
855  }
856 
857  // set up array of block matrices.
858  for(int i = 0; i < num_trans_per_proc; ++i)
859  {
860  BlkMatrix->SetBlock(i,i,loc_mat);
861  }
862 
863  return BlkMatrix;
864  }
865 
866  std::vector<LibUtilities::FieldDefinitionsSharedPtr> ExpListHomogeneous1D::v_GetFieldDefinitions()
867  {
868  std::vector<LibUtilities::FieldDefinitionsSharedPtr> returnval;
869 
870  // Set up Homogeneous length details.
872 
873  std::vector<NekDouble> HomoLen;
874  HomoLen.push_back(m_lhom);
875 
876  std::vector<unsigned int> StripsIDs;
877 
878  bool strips;
879  m_session->MatchSolverInfo("HomoStrip","True",strips,false);
880  if (strips)
881  {
882  StripsIDs.push_back(m_transposition->GetStripID());
883  }
884 
885  std::vector<unsigned int> PlanesIDs;
886  int IDoffset = 0;
887 
888  // introduce a 2 plane offset for single mode case so can
889  // be post-processed or used in MultiMode expansion.
891  {
892  IDoffset = 2;
893  }
894 
895  for(int i = 0; i < m_planes.size(); i++)
896  {
897  PlanesIDs.push_back(m_transposition->GetPlaneID(i)+IDoffset);
898  }
899 
900  m_planes[0]->GeneralGetFieldDefinitions(returnval, 1, HomoBasis,
901  HomoLen, strips, StripsIDs, PlanesIDs);
902  return returnval;
903  }
904 
905  void ExpListHomogeneous1D::v_GetFieldDefinitions(std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
906  {
907  // Set up Homogeneous length details.
909 
910  std::vector<NekDouble> HomoLen;
911  HomoLen.push_back(m_lhom);
912 
913  std::vector<unsigned int> StripsIDs;
914 
915  bool strips;
916  m_session->MatchSolverInfo("HomoStrip","True",strips,false);
917  if (strips)
918  {
919  StripsIDs.push_back(m_transposition->GetStripID());
920  }
921 
922  std::vector<unsigned int> PlanesIDs;
923  int IDoffset = 0;
924 
926  {
927  IDoffset = 2;
928  }
929 
930  for(int i = 0; i < m_planes.size(); i++)
931  {
932  PlanesIDs.push_back(m_transposition->GetPlaneID(i)+IDoffset);
933  }
934 
935  // enforce NumHomoDir == 1 by direct call
936  m_planes[0]->GeneralGetFieldDefinitions(fielddef, 1, HomoBasis,
937  HomoLen, strips, StripsIDs, PlanesIDs);
938  }
939 
940 
941  /** This routine appends the data from the expansion list into
942  the output format where each element is given by looping
943  over its Fourier modes where as data in the expandion is
944  stored with all consecutive elements and then the Fourier
945  modes
946  */
948  {
949  int i,n;
950  int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
951 
952  // Determine mapping from element ids to location in
953  // expansion list
954  if (m_elmtToExpId.size() == 0)
955  {
956  for(i = 0; i < m_planes[0]->GetExpSize(); ++i)
957  {
958  m_elmtToExpId[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
959  }
960  }
961 
962  for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
963  {
964  int eid = m_elmtToExpId[fielddef->m_elementIDs[i]];
965  int datalen = (*m_exp)[eid]->GetNcoeffs();
966 
967  for(n = 0; n < m_planes.size(); ++n)
968  {
969  fielddata.insert(fielddata.end(),&coeffs[m_coeff_offset[eid]+n*ncoeffs_per_plane],&coeffs[m_coeff_offset[eid]+n*ncoeffs_per_plane]+datalen);
970  }
971  }
972  }
973 
974  void ExpListHomogeneous1D::v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector<NekDouble> &fielddata)
975  {
976  v_AppendFieldData(fielddef,fielddata,m_coeffs);
977  }
978 
979  //Extract the data in fielddata into the m_coeff list
982  std::vector<NekDouble> &fielddata,
983  std::string &field,
984  Array<OneD, NekDouble> &coeffs)
985  {
986  int i,n;
987  int offset = 0;
988  int nzmodes = 1;
989  int datalen = fielddata.size()/fielddef->m_fields.size();
990  std::vector<unsigned int> fieldDefHomoZids;
991 
992 
993  // Find data location according to field definition
994  for(i = 0; i < fielddef->m_fields.size(); ++i)
995  {
996  if(fielddef->m_fields[i] == field)
997  {
998  break;
999  }
1000  offset += datalen;
1001  }
1002 
1003  if(i == fielddef->m_fields.size())
1004  {
1005  cout << "Field "<< field<< "not found in data file. " << endl;
1006  }
1007  else
1008  {
1009 
1010  int modes_offset = 0;
1011  int planes_offset = 0;
1012  Array<OneD, NekDouble> coeff_tmp;
1013 
1014  // Build map of plane IDs lying on this processor and determine
1015  // mapping from element ids to location in expansion list.
1016  if (m_zIdToPlane.size() == 0)
1017  {
1018  for (i = 0; i < m_planes.size(); ++i)
1019  {
1020  m_zIdToPlane[m_transposition->GetPlaneID(i)] = i;
1021  }
1022 
1023  for (i = 0; i < m_planes[0]->GetExpSize(); ++i)
1024  {
1025  m_elmtToExpId[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
1026  }
1027  }
1028 
1029  if(fielddef->m_numHomogeneousDir)
1030  {
1031  nzmodes = fielddef->m_homogeneousZIDs.size();
1032  fieldDefHomoZids = fielddef->m_homogeneousZIDs;
1033  }
1034  else // input file is 2D and so set nzmodes to 1
1035  {
1036  nzmodes = 1;
1037  fieldDefHomoZids.push_back(0);
1038  }
1039 
1040  // calculate number of modes in the current partition
1041  int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
1042 
1043  for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
1044  {
1045  if(fielddef->m_uniOrder == true) // reset modes_offset to zero
1046  {
1047  modes_offset = 0;
1048  }
1049 
1050  int datalen = LibUtilities::GetNumberOfCoefficients(fielddef->m_shapeType,
1051  fielddef->m_numModes,
1052  modes_offset);
1053 
1054  auto it = m_elmtToExpId.find(fielddef->m_elementIDs[i]);
1055 
1056  // ensure element is on this partition for parallel case.
1057  if(it == m_elmtToExpId.end())
1058  {
1059  // increase offset for correct FieldData access
1060  offset += datalen*nzmodes;
1061  modes_offset += (*m_exp)[0]->GetNumBases() +
1062  fielddef->m_numHomogeneousDir;
1063  continue;
1064  }
1065 
1066  int eid = it->second;
1067  bool sameBasis = true;
1068  for (int j = 0; j < fielddef->m_basis.size()-1; ++j)
1069  {
1070  if (fielddef->m_basis[j] != (*m_exp)[eid]->GetBasisType(j))
1071  {
1072  sameBasis = false;
1073  break;
1074  }
1075  }
1076 
1077  for(n = 0; n < nzmodes; ++n, offset += datalen)
1078  {
1079 
1080  it = m_zIdToPlane.find(fieldDefHomoZids[n]);
1081 
1082  // Check to make sure this mode number lies in this field.
1083  if (it == m_zIdToPlane.end())
1084  {
1085  continue;
1086  }
1087 
1088  planes_offset = it->second;
1089  if(datalen == (*m_exp)[eid]->GetNcoeffs() && sameBasis)
1090  {
1091  Vmath::Vcopy(datalen,&fielddata[offset],1,&coeffs[m_coeff_offset[eid]+planes_offset*ncoeffs_per_plane],1);
1092  }
1093  else // unpack data to new order
1094  {
1095  (*m_exp)[eid]->ExtractDataToCoeffs(&fielddata[offset], fielddef->m_numModes,modes_offset,&coeffs[m_coeff_offset[eid] + planes_offset*ncoeffs_per_plane], fielddef->m_basis);
1096  }
1097  }
1098  modes_offset += (*m_exp)[0]->GetNumBases() + fielddef->m_numHomogeneousDir;
1099  }
1100  }
1101  }
1102 
1103  //Extract the data in fielddata into the m_coeff list
1105  const std::shared_ptr<ExpList> &fromExpList,const Array<OneD, const NekDouble> &fromCoeffs, Array<OneD, NekDouble> &toCoeffs)
1106  {
1107  int i;
1108  int fromNcoeffs_per_plane = fromExpList->GetPlane(0)->GetNcoeffs();
1109  int toNcoeffs_per_plane = m_planes[0]->GetNcoeffs();
1110  Array<OneD, NekDouble> tocoeffs_tmp, fromcoeffs_tmp;
1111 
1112  for(i = 0; i < m_planes.size(); ++i)
1113  {
1114  m_planes[i]->ExtractCoeffsToCoeffs(fromExpList->GetPlane(i),fromcoeffs_tmp = fromCoeffs + fromNcoeffs_per_plane*i, tocoeffs_tmp = toCoeffs + toNcoeffs_per_plane*i);
1115  }
1116  }
1117 
1118  void ExpListHomogeneous1D::v_WriteVtkPieceData(std::ostream &outfile, int expansion,
1119  std::string var)
1120  {
1121  // If there is only one plane (e.g. HalfMode), we write a 2D plane.
1122  if (m_planes.size() == 1)
1123  {
1124  m_planes[0]->WriteVtkPieceData(outfile, expansion, var);
1125  return;
1126  }
1127 
1128  int i;
1129  int nq = (*m_exp)[expansion]->GetTotPoints();
1130  int npoints_per_plane = m_planes[0]->GetTotPoints();
1131 
1132  // If we are using Fourier points, output extra plane to fill domain
1133  int outputExtraPlane = 0;
1134  Array<OneD, NekDouble> extraPlane;
1135  if ( m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier
1136  && m_homogeneousBasis->GetPointsType() ==
1138  {
1139  outputExtraPlane = 1;
1140  // Get extra plane data
1141  if (m_StripZcomm->GetSize() == 1)
1142  {
1143  extraPlane = m_phys + m_phys_offset[expansion];
1144  }
1145  else
1146  {
1147  // Determine to and from rank for communication
1148  int size = m_StripZcomm->GetSize();
1149  int rank = m_StripZcomm->GetRank();
1150  int fromRank = (rank+1) % size;
1151  int toRank = (rank == 0) ? size-1 : rank-1;
1152  // Communicate using SendRecv
1153  extraPlane = Array<OneD, NekDouble>(nq);
1154  Array<OneD, NekDouble> send (nq,
1155  m_phys + m_phys_offset[expansion]);
1156  m_StripZcomm->SendRecv(toRank, send,
1157  fromRank, extraPlane);
1158  }
1159  }
1160 
1161  // printing the fields of that zone
1162  outfile << " <DataArray type=\"Float64\" Name=\""
1163  << var << "\">" << endl;
1164  outfile << " ";
1165  for (int n = 0; n < m_planes.size(); ++n)
1166  {
1167  const Array<OneD, NekDouble> phys = m_phys + m_phys_offset[expansion] + n*npoints_per_plane;
1168  for(i = 0; i < nq; ++i)
1169  {
1170  outfile << (fabs(phys[i]) < NekConstants::kNekZeroTol ? 0 : phys[i]) << " ";
1171  }
1172  }
1173  if (outputExtraPlane)
1174  {
1175  for(i = 0; i < nq; ++i)
1176  {
1177  outfile << (fabs(extraPlane[i]) < NekConstants::kNekZeroTol ?
1178  0 : extraPlane[i]) << " ";
1179  }
1180  }
1181  outfile << endl;
1182  outfile << " </DataArray>" << endl;
1183  }
1184 
1186  {
1187  int cnt,cnt1;
1188  Array<OneD, NekDouble> tmparray;
1189  cnt = m_planes[0]->GetTotPoints();
1190  cnt1 = m_planes[0]->Get1DScaledTotPoints(scale);
1191 
1192  ASSERTL1(m_planes.size()*cnt1 <= outarray.size(),"size of outarray does not match internal estimage");
1193 
1194 
1195  for(int i = 0; i < m_planes.size(); i++)
1196  {
1197 
1198  m_planes[i]->PhysInterp1DScaled(scale,inarray+i*cnt,
1199  tmparray = outarray+i*cnt1);
1200  }
1201  }
1202 
1203 
1205  {
1206  int cnt,cnt1;
1207  Array<OneD, NekDouble> tmparray;
1208  cnt = m_planes[0]->Get1DScaledTotPoints(scale);
1209  cnt1 = m_planes[0]->GetTotPoints();
1210 
1211  ASSERTL1(m_planes.size()*cnt <= inarray.size(),"size of outarray does not match internal estimage");
1212 
1213 
1214  for(int i = 0; i < m_planes.size(); i++)
1215  {
1216  m_planes[i]->PhysGalerkinProjection1DScaled(scale,inarray+i*cnt,
1217  tmparray = outarray+i*cnt1);
1218  }
1219 
1220  }
1222  Array<OneD, NekDouble> &out_d0,
1223  Array<OneD, NekDouble> &out_d1,
1224  Array<OneD, NekDouble> &out_d2)
1225  {
1226  int nT_pts = inarray.size(); //number of total points = n. of Fourier points * n. of points per plane (nT_pts)
1227  int nP_pts = nT_pts/m_planes.size(); //number of points per plane = n of Fourier transform required (nP_pts)
1228 
1229  Array<OneD, NekDouble> temparray(nT_pts);
1230  Array<OneD, NekDouble> outarray(nT_pts);
1234 
1235  for(int i = 0; i < m_planes.size(); i++)
1236  {
1237  m_planes[i]->PhysDeriv(inarray + i*nP_pts ,tmp2 = out_d0 + i*nP_pts , tmp3 = out_d1 + i*nP_pts );
1238  }
1239 
1240  if(out_d2 != NullNekDouble1DArray)
1241  {
1244  {
1245  if(m_WaveSpace)
1246  {
1247  temparray = inarray;
1248  }
1249  else
1250  {
1251  HomogeneousFwdTrans(inarray,temparray);
1252  }
1253 
1254  NekDouble sign = -1.0;
1255  NekDouble beta;
1256 
1257  //Half Mode
1259  {
1260  beta = sign*2*M_PI*(m_transposition->GetK(0))/m_lhom;
1261 
1262  Vmath::Smul(nP_pts,beta,temparray,1,outarray,1);
1263  }
1264  else if(m_homogeneousBasis->GetBasisType() == LibUtilities::eFourierHalfModeIm)
1265  {
1266  beta = -sign*2*M_PI*(m_transposition->GetK(0))/m_lhom;
1267 
1268  Vmath::Smul(nP_pts,beta,temparray,1,outarray,1);
1269  }
1270 
1271  //Fully complex
1272  else
1273  {
1274  for(int i = 0; i < m_planes.size(); i++)
1275  {
1276  beta = -sign*2*M_PI*(m_transposition->GetK(i))/m_lhom;
1277 
1278  Vmath::Smul(nP_pts,beta,tmp1 = temparray + i*nP_pts,1,tmp2 = outarray + (i-int(sign))*nP_pts,1);
1279 
1280  sign = -1.0*sign;
1281  }
1282  }
1283 
1284  if(m_WaveSpace)
1285  {
1286  out_d2 = outarray;
1287  }
1288  else
1289  {
1290  HomogeneousBwdTrans(outarray,out_d2);
1291  }
1292  }
1293  else
1294  {
1295  if(!m_session->DefinesSolverInfo("HomoStrip"))
1296  {
1297  ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,
1298  "Parallelisation in the homogeneous direction "
1299  "implemented just for Fourier basis");
1300  }
1301  else
1302  {
1303  ASSERTL0(m_StripZcomm->GetSize() == 1,
1304  "Parallelisation in the homogeneous direction "
1305  "implemented just for Fourier basis");
1306  }
1307 
1308  if(m_WaveSpace)
1309  {
1310  ASSERTL0(false, "Semi-phyisical time-stepping not "
1311  "implemented yet for non-Fourier "
1312  "basis");
1313  }
1314  else
1315  {
1316  StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
1317 
1318  m_transposition->Transpose(inarray,temparray,false,LibUtilities::eXYtoZ);
1319 
1320  for(int i = 0; i < nP_pts; i++)
1321  {
1322  StdSeg.PhysDeriv(temparray + i*m_planes.size(), tmp2 = outarray + i*m_planes.size());
1323  }
1324 
1325  m_transposition->Transpose(outarray,out_d2,false,LibUtilities::eZtoXY);
1326 
1327  Vmath::Smul(nT_pts,2.0/m_lhom,out_d2,1,out_d2,1);
1328  }
1329  }
1330  }
1331  }
1332 
1335 
1336  {
1337  int nT_pts = inarray.size(); //number of total points = n. of Fourier points * n. of points per plane (nT_pts)
1338  int nP_pts = nT_pts/m_planes.size(); //number of points per plane = n of Fourier transform required (nP_pts)
1339 
1340  int dir= (int)edir;
1341 
1342  Array<OneD, NekDouble> temparray(nT_pts);
1343  Array<OneD, NekDouble> outarray(nT_pts);
1346 
1347  if (dir < 2)
1348  {
1349  for(int i=0; i < m_planes.size(); i++)
1350  {
1351  m_planes[i]->PhysDeriv(edir, inarray + i*nP_pts ,tmp2 = out_d + i*nP_pts);
1352  }
1353  }
1354  else
1355  {
1358  {
1359  if(m_WaveSpace)
1360  {
1361  temparray = inarray;
1362  }
1363  else
1364  {
1365  HomogeneousFwdTrans(inarray,temparray);
1366  }
1367 
1368  NekDouble sign = -1.0;
1369  NekDouble beta;
1370 
1371  //HalfMode
1373  {
1374  beta = 2*sign*M_PI*(m_transposition->GetK(0))/m_lhom;
1375 
1376  Vmath::Smul(nP_pts,beta,temparray,1,outarray,1);
1377  }
1378  else if(m_homogeneousBasis->GetBasisType() == LibUtilities::eFourierHalfModeIm)
1379  {
1380  beta = -2*sign*M_PI*(m_transposition->GetK(0))/m_lhom;
1381 
1382  Vmath::Smul(nP_pts,beta,temparray,1,outarray,1);
1383  }
1384  //Fully complex
1385  else
1386  {
1387  for(int i = 0; i < m_planes.size(); i++)
1388  {
1389  beta = -sign*2*M_PI*(m_transposition->GetK(i))/m_lhom;
1390 
1391  Vmath::Smul(nP_pts,beta,tmp1 = temparray + i*nP_pts,1,tmp2 = outarray + (i-int(sign))*nP_pts,1);
1392 
1393  sign = -1.0*sign;
1394  }
1395  }
1396  if(m_WaveSpace)
1397  {
1398  out_d = outarray;
1399  }
1400  else
1401  {
1402  HomogeneousBwdTrans(outarray,out_d);
1403  }
1404  }
1405  else
1406  {
1407  if(!m_session->DefinesSolverInfo("HomoStrip"))
1408  {
1409  ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,
1410  "Parallelisation in the homogeneous direction "
1411  "implemented just for Fourier basis");
1412  }
1413  else
1414  {
1415  ASSERTL0(m_StripZcomm->GetSize() == 1,
1416  "Parallelisation in the homogeneous direction "
1417  "implemented just for Fourier basis");
1418  }
1419 
1420  if(m_WaveSpace)
1421  {
1422  ASSERTL0(false,"Semi-phyisical time-stepping not implemented yet for non-Fourier basis");
1423  }
1424  else
1425  {
1426  StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
1427 
1428  m_transposition->Transpose(inarray,temparray,false,LibUtilities::eXYtoZ);
1429 
1430  for(int i = 0; i < nP_pts; i++)
1431  {
1432  StdSeg.PhysDeriv(temparray + i*m_planes.size(), tmp2 = outarray + i*m_planes.size());
1433  }
1434 
1435  m_transposition->Transpose(outarray,out_d,false,LibUtilities::eZtoXY);
1436 
1437  Vmath::Smul(nT_pts,2.0/m_lhom,out_d,1,out_d,1);
1438  }
1439  }
1440  }
1441  }
1442 
1444  Array<OneD, NekDouble> &out_d0,
1445  Array<OneD, NekDouble> &out_d1,
1446  Array<OneD, NekDouble> &out_d2)
1447 
1448  {
1449  v_PhysDeriv(inarray,out_d0,out_d1,out_d2);
1450  }
1451 
1453  const Array<OneD, const NekDouble> &inarray,
1454  Array<OneD, NekDouble> &out_d)
1455  {
1456  v_PhysDeriv(edir,inarray,out_d);
1457  }
1458 
1460  {
1461  return m_transposition;
1462  }
1463 
1465  {
1466  return m_lhom;
1467  }
1468 
1470  {
1471  m_lhom = lhom;
1472  }
1473 
1475  {
1476  return m_transposition->GetPlanesIDs();
1477  }
1478 
1480  const NekDouble> &inarray)
1481  {
1482  NekDouble val = 0.0;
1483  int i = 0;
1484 
1485  for (i = 0; i < (*m_exp).size(); ++i)
1486  {
1487  val += (*m_exp)[i]->Integral(inarray + m_phys_offset[i]);
1488  }
1489  val *= m_lhom/m_homogeneousBasis->GetNumModes();
1490 
1491  m_comm->AllReduce(val, LibUtilities::ReduceSum);
1492 
1493  return val;
1494  }
1495  } //end of namespace
1496 } //end of namespace
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:274
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:15
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:145
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 void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
NekDouble m_lhom
Width of homogeneous direction.
LibUtilities::TranspositionSharedPtr m_transposition
virtual void v_ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
virtual void v_DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_DealiasedDotProd(const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray)
Array< OneD, ExpListSharedPtr > m_planes
std::unordered_map< int, int > m_zIdToPlane
NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
virtual void v_FwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
LibUtilities::NektarFFTSharedPtr m_FFT_deal
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual Array< OneD, const unsigned int > v_GetZIDs(void)
void Homogeneous1DTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool IsForwards, bool Shuff=true, bool UnShuff=true)
virtual void v_SetHomoLen(const NekDouble lhom)
LibUtilities::NektarFFTSharedPtr m_FFT
virtual void v_IProductWRTBase_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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)
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void)
DNekBlkMatSharedPtr GenHomogeneous1DBlockMatrix(Homogeneous1DMatType mattype) const
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var)
virtual void v_HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
virtual void v_ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
Extract data from raw field data into expansion list.
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
DNekBlkMatSharedPtr GetHomogeneous1DBlockMatrix(Homogeneous1DMatType mattype) const
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_BwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void)
ExpListHomogeneous1D(const ExpansionType type)
Default constructor.
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:107
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1252
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1800
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1304
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1220
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1290
std::unordered_map< int, int > m_elmtToExpId
Mapping from geometry ID of element to index inside m_exp.
Definition: ExpList.h:1321
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1223
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1307
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1269
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:617
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:376
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.
Definition: StdSegExp.h:54
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:370
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:179
std::shared_ptr< Transposition > TranspositionSharedPtr
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:65
int GetNumberOfCoefficients(ShapeType shape, std::vector< unsigned int > &modes, int offset)
Definition: ShapeType.hpp:313
@ eFourierSingleMode
Fourier ModifiedExpansion with just the first mode .
Definition: BasisType.h:59
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode
Definition: BasisType.h:61
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode
Definition: BasisType.h:60
@ eFourier
Fourier Expansion .
Definition: BasisType.h:53
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:1
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:71
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
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:192
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:513
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:225
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199