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