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