Nektar++
ExpListHomogeneous2D.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File ExpListHomogeneous2D.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 2-directions
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
38 #include <LibUtilities/Foundations/ManagerAccess.h> // for PointsManager, etc
39 #include <StdRegions/StdSegExp.h>
40 #include <StdRegions/StdQuadExp.h>
41 #include <LocalRegions/Expansion.h>
42 
43 using namespace std;
44 
45 namespace Nektar
46 {
47  namespace MultiRegions
48  {
49  // Forward declaration for typedefs
50  ExpListHomogeneous2D::ExpListHomogeneous2D():
51  ExpList(),
52  m_homogeneousBasis_y(LibUtilities::NullBasisSharedPtr),
53  m_homogeneousBasis_z(LibUtilities::NullBasisSharedPtr),
54  m_lhom_y(1),
55  m_lhom_z(1),
56  m_homogeneous2DBlockMat(MemoryManager<Homo2DBlockMatrixMap>::AllocateSharedPtr())
57  {
58  }
59 
61  const LibUtilities::BasisKey &HomoBasis_y,
62  const LibUtilities::BasisKey &HomoBasis_z,
63  const NekDouble lhom_y,
64  const NekDouble lhom_z,
65  const bool useFFT,
66  const bool dealiasing):
67  ExpList(pSession),
68  m_useFFT(useFFT),
69  m_lhom_y(lhom_y),
70  m_lhom_z(lhom_z),
72  m_dealiasing(dealiasing)
73  {
74  ASSERTL2(HomoBasis_y != LibUtilities::NullBasisKey,
75  "Homogeneous Basis in y direction is a null basis");
76  ASSERTL2(HomoBasis_z != LibUtilities::NullBasisKey,
77  "Homogeneous Basis in z direction is a null basis");
78 
81 
83 
84  m_Ycomm = m_comm->GetColumnComm()->GetRowComm();
85  m_Zcomm = m_comm->GetColumnComm()->GetRowComm();
86 
87  m_ny = m_homogeneousBasis_y->GetNumPoints()/m_Ycomm->GetSize();
88  m_nz = m_homogeneousBasis_z->GetNumPoints()/m_Zcomm->GetSize();
89 
91 
92  if(m_useFFT)
93  {
96  }
97 
98  if(m_dealiasing)
99  {
100  ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,"Remove dealiasing if you want to run in parallel");
101  SetPaddingBase();
102  }
103  }
104 
105 
106  /**
107  * @param In ExpListHomogeneous2D object to copy.
108  */
110  ExpList(In,false),
111  m_useFFT(In.m_useFFT),
112  m_FFT_y(In.m_FFT_y),
113  m_FFT_z(In.m_FFT_z),
115  m_Ycomm(In.m_Ycomm),
116  m_Zcomm(In.m_Ycomm),
119  m_lhom_y(In.m_lhom_y),
120  m_lhom_z(In.m_lhom_z),
122  m_ny(In.m_ny),
123  m_nz(In.m_nz),
127  MatFwdPAD(In.MatFwdPAD),
128  MatBwdPAD(In.MatBwdPAD)
129  {
130  m_lines = Array<OneD, ExpListSharedPtr>(In.m_lines.num_elements());
131  }
132 
134  const std::vector<unsigned int> &eIDs):
135  ExpList(In,eIDs,false),
136  m_useFFT(In.m_useFFT),
137  m_FFT_y(In.m_FFT_y),
138  m_FFT_z(In.m_FFT_z),
140  m_Ycomm(In.m_Ycomm),
141  m_Zcomm(In.m_Ycomm),
144  m_lhom_y(In.m_lhom_y),
145  m_lhom_z(In.m_lhom_z),
147  m_ny(In.m_ny),
148  m_nz(In.m_nz),
152  MatFwdPAD(In.MatFwdPAD),
153  MatBwdPAD(In.MatBwdPAD)
154  {
155  m_lines = Array<OneD, ExpListSharedPtr>(In.m_lines.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  Homogeneous2DTrans(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  Homogeneous2DTrans(inarray,outarray,false,coeffstate,Shuff,UnShuff);
183  }
184 
186  const Array<OneD, NekDouble> &inarray2,
187  Array<OneD, NekDouble> &outarray,
188  CoeffState coeffstate)
189  {
190  int npoints = outarray.num_elements(); // number of total physical points
191  int nlines = m_lines.num_elements(); // number of lines == number of Fourier modes = number of Fourier coeff = number of points per slab
192  int nslabs = npoints/nlines; // number of slabs = numebr of physical points per line
193 
194  Array<OneD, NekDouble> V1(npoints);
195  Array<OneD, NekDouble> V2(npoints);
196  Array<OneD, NekDouble> V1V2(npoints);
197  Array<OneD, NekDouble> ShufV1(npoints);
198  Array<OneD, NekDouble> ShufV2(npoints);
199  Array<OneD, NekDouble> ShufV1V2(npoints);
200 
201  if(m_WaveSpace)
202  {
203  V1 = inarray1;
204  V2 = inarray2;
205  }
206  else
207  {
208  HomogeneousFwdTrans(inarray1,V1,coeffstate);
209  HomogeneousFwdTrans(inarray2,V2,coeffstate);
210  }
211 
212  m_transposition->Transpose(V1,ShufV1,false,LibUtilities::eXtoYZ);
213  m_transposition->Transpose(V2,ShufV2,false,LibUtilities::eXtoYZ);
214 
215  Array<OneD, NekDouble> PadV1_slab_coeff(m_padsize_y*m_padsize_z,0.0);
216  Array<OneD, NekDouble> PadV2_slab_coeff(m_padsize_y*m_padsize_z,0.0);
217  Array<OneD, NekDouble> PadRe_slab_coeff(m_padsize_y*m_padsize_z,0.0);
218 
219  Array<OneD, NekDouble> PadV1_slab_phys(m_padsize_y*m_padsize_z,0.0);
220  Array<OneD, NekDouble> PadV2_slab_phys(m_padsize_y*m_padsize_z,0.0);
221  Array<OneD, NekDouble> PadRe_slab_phys(m_padsize_y*m_padsize_z,0.0);
222 
223  NekVector<NekDouble> PadIN_V1(m_padsize_y*m_padsize_z,PadV1_slab_coeff,eWrapper);
224  NekVector<NekDouble> PadOUT_V1(m_padsize_y*m_padsize_z,PadV1_slab_phys,eWrapper);
225 
226  NekVector<NekDouble> PadIN_V2(m_padsize_y*m_padsize_z,PadV2_slab_coeff,eWrapper);
227  NekVector<NekDouble> PadOUT_V2(m_padsize_y*m_padsize_z,PadV2_slab_phys,eWrapper);
228 
229  NekVector<NekDouble> PadIN_Re(m_padsize_y*m_padsize_z,PadRe_slab_phys,eWrapper);
230  NekVector<NekDouble> PadOUT_Re(m_padsize_y*m_padsize_z,PadRe_slab_coeff,eWrapper);
231 
232  //Looping on the slabs
233  for(int j = 0 ; j< nslabs ; j++)
234  {
235  //Copying the j-th slab of size N*M into a bigger slab of lenght 2*N*M
236  //We are in Fourier space
237  for(int i = 0 ; i< m_nz ; i++)
238  {
239  Vmath::Vcopy(m_ny,&(ShufV1[i*m_ny + j*nlines]),1,&(PadV1_slab_coeff[i*2*m_ny]),1);
240  Vmath::Vcopy(m_ny,&(ShufV2[i*m_ny + j*nlines]),1,&(PadV2_slab_coeff[i*2*m_ny]),1);
241  }
242 
243  //Moving to physical space using the padded system
244  PadOUT_V1 = (*MatBwdPAD)*PadIN_V1;
245  PadOUT_V2 = (*MatBwdPAD)*PadIN_V2;
246 
247  //Perfroming the vectors multiplication in physical
248  //space on the padded system
249  Vmath::Vmul(m_padsize_y*m_padsize_z,PadV1_slab_phys,1,PadV2_slab_phys,1,PadRe_slab_phys,1);
250 
251  //Moving back the result (V1*V2)_phys in Fourier
252  //space, padded system
253  PadOUT_Re = (*MatFwdPAD)*PadIN_Re;
254 
255  //Copying the first half of the padded pencil in the
256  //full vector (Fourier space)
257  for (int i = 0; i < m_nz; i++)
258  {
259  Vmath::Vcopy(m_ny,&(PadRe_slab_coeff[i*2*m_ny]),1,&(ShufV1V2[i*m_ny + j*nlines]),1);
260  }
261  }
262 
263  if(m_WaveSpace)
264  {
265  m_transposition->Transpose(ShufV1V2,outarray,false,LibUtilities::eYZtoX);
266  }
267  else
268  {
269  m_transposition->Transpose(ShufV1V2,V1V2,false,LibUtilities::eYZtoX);
270 
271  //Moving the results in physical space for the output
272  HomogeneousBwdTrans(V1V2,outarray,coeffstate);
273  }
274  }
275 
277  const Array<OneD, Array<OneD, NekDouble> > &inarray1,
278  const Array<OneD, Array<OneD, NekDouble> > &inarray2,
279  Array<OneD, Array<OneD, NekDouble> > &outarray,
280  CoeffState coeffstate)
281  {
282  boost::ignore_unused(coeffstate);
283 
284  // TODO Proper implementation of this
285  int ndim = inarray1.num_elements();
286  ASSERTL1( inarray2.num_elements() % ndim == 0,
287  "Wrong dimensions for DealiasedDotProd.");
288  int nvec = inarray2.num_elements() % ndim;
289  int npts = inarray1[0].num_elements();
290 
291  Array<OneD, NekDouble> out(npts);
292  for (int i = 0; i < nvec; i++)
293  {
294  Vmath::Zero(npts, outarray[i], 1);
295  for (int j = 0; j < ndim; j++)
296  {
297  DealiasedProd(inarray1[j], inarray2[i*ndim+j], out);
298  Vmath::Vadd(npts, outarray[i], 1, out, 1, outarray[i], 1);
299  }
300  }
301  }
302 
304  {
305  int cnt = 0, cnt1 = 0;
306  Array<OneD, NekDouble> tmparray;
307  int nlines = m_lines.num_elements();
308 
309  for(int n = 0; n < nlines; ++n)
310  {
311  m_lines[n]->FwdTrans(inarray+cnt, tmparray = outarray + cnt1,
312  coeffstate);
313  cnt += m_lines[n]->GetTotPoints();
314  cnt1 += m_lines[n]->GetNcoeffs();
315  }
316  if(!m_WaveSpace)
317  {
318  HomogeneousFwdTrans(outarray,outarray,coeffstate);
319  }
320  }
321 
323  {
324  int cnt = 0, cnt1 = 0;
325  Array<OneD, NekDouble> tmparray;
326  int nlines = m_lines.num_elements();
327 
328  for(int n = 0; n < nlines; ++n)
329  {
330  m_lines[n]->FwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
331 
332  cnt += m_lines[n]->GetTotPoints();
333  cnt1 += m_lines[n]->GetNcoeffs();
334  }
335  if(!m_WaveSpace)
336  {
337  HomogeneousFwdTrans(outarray,outarray);
338  }
339  }
340 
342  {
343  int cnt = 0, cnt1 = 0;
344  Array<OneD, NekDouble> tmparray;
345  int nlines = m_lines.num_elements();
346 
347  for(int n = 0; n < nlines; ++n)
348  {
349  m_lines[n]->BwdTrans(inarray+cnt, tmparray = outarray + cnt1,
350  coeffstate);
351  cnt += m_lines[n]->GetNcoeffs();
352  cnt1 += m_lines[n]->GetTotPoints();
353  }
354  if(!m_WaveSpace)
355  {
356  HomogeneousBwdTrans(outarray,outarray);
357  }
358  }
359 
361  {
362  int cnt = 0, cnt1 = 0;
363  Array<OneD, NekDouble> tmparray;
364  int nlines = m_lines.num_elements();
365 
366  for(int n = 0; n < nlines; ++n)
367  {
368  m_lines[n]->BwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
369 
370  cnt += m_lines[n]->GetNcoeffs();
371  cnt1 += m_lines[n]->GetTotPoints();
372  }
373  if(!m_WaveSpace)
374  {
375  HomogeneousBwdTrans(outarray,outarray);
376  }
377  }
378 
379 
381  {
382  int cnt = 0, cnt1 = 0;
383  Array<OneD, NekDouble> tmparray;
384  int nlines = m_lines.num_elements();
385 
386  for(int n = 0; n < nlines; ++n)
387  {
388  m_lines[n]->IProductWRTBase(inarray+cnt, tmparray = outarray + cnt1,coeffstate);
389 
390  cnt += m_lines[n]->GetNcoeffs();
391  cnt1 += m_lines[n]->GetTotPoints();
392  }
393  }
394 
396  {
397  int cnt = 0, cnt1 = 0;
398  Array<OneD, NekDouble> tmparray;
399  int nlines = m_lines.num_elements();
400 
401  for(int n = 0; n < nlines; ++n)
402  {
403  m_lines[n]->IProductWRTBase_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
404 
405  cnt += m_lines[n]->GetNcoeffs();
406  cnt1 += m_lines[n]->GetTotPoints();
407  }
408  }
409 
411  Array<OneD, NekDouble> &outarray,
412  bool IsForwards,
413  CoeffState coeffstate,
414  bool Shuff,
415  bool UnShuff)
416  {
417  boost::ignore_unused(Shuff, UnShuff);
418 
419  if(m_useFFT)
420  {
421 
422  int n = m_lines.num_elements(); //number of Fourier points in the Fourier directions (x-z grid)
423  int s = inarray.num_elements(); //number of total points = n. of Fourier points * n. of points per line
424  int p = s/n; //number of points per line = n of Fourier transform required
425 
426  Array<OneD, NekDouble> fft_in(s);
427  Array<OneD, NekDouble> fft_out(s);
428 
429  m_transposition->Transpose(inarray,fft_in,false,LibUtilities::eXtoYZ);
430 
431  if(IsForwards)
432  {
433  for(int i=0;i<(p*m_nz);i++)
434  {
435  m_FFT_y->FFTFwdTrans(m_tmpIN = fft_in + i*m_ny, m_tmpOUT = fft_out + i*m_ny);
436  }
437 
438  }
439  else
440  {
441  for(int i=0;i<(p*m_nz);i++)
442  {
443  m_FFT_y->FFTBwdTrans(m_tmpIN = fft_in + i*m_ny, m_tmpOUT = fft_out + i*m_ny);
444  }
445  }
446 
447  m_transposition->Transpose(fft_out,fft_in,false,LibUtilities::eYZtoZY);
448 
449  if(IsForwards)
450  {
451  for(int i=0;i<(p*m_ny);i++)
452  {
453  m_FFT_z->FFTFwdTrans(m_tmpIN = fft_in + i*m_nz, m_tmpOUT = fft_out + i*m_nz);
454  }
455 
456  }
457  else
458  {
459  for(int i=0;i<(p*m_ny);i++)
460  {
461  m_FFT_z->FFTBwdTrans(m_tmpIN = fft_in + i*m_nz, m_tmpOUT = fft_out + i*m_nz);
462  }
463  }
464 
465  //TODO: required ZYtoX routine
466  m_transposition->Transpose(fft_out,fft_in,false,LibUtilities::eZYtoYZ);
467 
468  m_transposition->Transpose(fft_in,outarray,false,LibUtilities::eYZtoX);
469 
470  }
471  else
472  {
473  DNekBlkMatSharedPtr blkmatY;
474  DNekBlkMatSharedPtr blkmatZ;
475 
476  if(inarray.num_elements() == m_npoints) //transform phys space
477  {
478  if(IsForwards)
479  {
482  }
483  else
484  {
487  }
488  }
489  else
490  {
491  if(IsForwards)
492  {
495  }
496  else
497  {
500  }
501  }
502 
503  int nrowsY = blkmatY->GetRows();
504  int ncolsY = blkmatY->GetColumns();
505 
506  Array<OneD, NekDouble> sortedinarrayY(ncolsY);
507  Array<OneD, NekDouble> sortedoutarrayY(nrowsY);
508 
509  int nrowsZ = blkmatZ->GetRows();
510  int ncolsZ = blkmatZ->GetColumns();
511 
512  Array<OneD, NekDouble> sortedinarrayZ(ncolsZ);
513  Array<OneD, NekDouble> sortedoutarrayZ(nrowsZ);
514 
515  NekVector<NekDouble> inY (ncolsY,sortedinarrayY,eWrapper);
516  NekVector<NekDouble> outY(nrowsY,sortedoutarrayY,eWrapper);
517 
518  NekVector<NekDouble> inZ (ncolsZ,sortedinarrayZ,eWrapper);
519  NekVector<NekDouble> outZ(nrowsZ,sortedoutarrayZ,eWrapper);
520 
521  m_transposition->Transpose(inarray,sortedinarrayY,!IsForwards,LibUtilities::eXtoYZ);
522 
523  outY = (*blkmatY)*inY;
524 
525  m_transposition->Transpose(sortedoutarrayY,sortedinarrayZ,false,LibUtilities::eYZtoZY);
526 
527  outZ = (*blkmatZ)*inZ;
528 
529  m_transposition->Transpose(sortedoutarrayZ,sortedoutarrayY,false,LibUtilities::eZYtoYZ);
530 
531  m_transposition->Transpose(sortedoutarrayY,outarray,false,LibUtilities::eYZtoX);
532  }
533  }
534 
536  {
537  auto matrixIter = m_homogeneous2DBlockMat->find(mattype);
538 
539  if(matrixIter == m_homogeneous2DBlockMat->end())
540  {
541  return ((*m_homogeneous2DBlockMat)[mattype] =
542  GenHomogeneous2DBlockMatrix(mattype,coeffstate));
543  }
544  else
545  {
546  return matrixIter->second;
547  }
548  }
549 
551  {
552  boost::ignore_unused(coeffstate);
553 
554  int i;
555  int n_exp = 0;
556 
557  DNekMatSharedPtr loc_mat;
558  DNekBlkMatSharedPtr BlkMatrix;
559 
561 
562  int NumPoints = 0;
563  int NumModes = 0;
564  int NumPencils = 0;
565 
566  if((mattype == eForwardsCoeffSpaceY1D) || (mattype == eBackwardsCoeffSpaceY1D)
567  ||(mattype == eForwardsPhysSpaceY1D) || (mattype == eBackwardsPhysSpaceY1D))
568  {
569  Basis = m_homogeneousBasis_y;
570  NumPoints = m_homogeneousBasis_y->GetNumModes();
571  NumModes = m_homogeneousBasis_y->GetNumPoints();
572  NumPencils = m_homogeneousBasis_z->GetNumPoints();
573  }
574  else
575  {
576  Basis = m_homogeneousBasis_z;
577  NumPoints = m_homogeneousBasis_z->GetNumModes();
578  NumModes = m_homogeneousBasis_z->GetNumPoints();
579  NumPencils = m_homogeneousBasis_y->GetNumPoints();
580  }
581 
582  if((mattype == eForwardsCoeffSpaceY1D) || (mattype == eForwardsCoeffSpaceZ1D)
583  ||(mattype == eBackwardsCoeffSpaceY1D)||(mattype == eBackwardsCoeffSpaceZ1D))
584  {
585  n_exp = m_lines[0]->GetNcoeffs();
586  }
587  else
588  {
589  n_exp = m_lines[0]->GetTotPoints(); // will operatore on m_phys
590  }
591 
592  Array<OneD,unsigned int> nrows(n_exp);
593  Array<OneD,unsigned int> ncols(n_exp);
594 
595  if((mattype == eForwardsCoeffSpaceY1D)||(mattype == eForwardsPhysSpaceY1D) ||
596  (mattype == eForwardsCoeffSpaceZ1D)||(mattype == eForwardsPhysSpaceZ1D))
597  {
598  nrows = Array<OneD, unsigned int>(n_exp*NumPencils,NumModes);
599  ncols = Array<OneD, unsigned int>(n_exp*NumPencils,NumPoints);
600  }
601  else
602  {
603  nrows = Array<OneD, unsigned int>(n_exp*NumPencils,NumPoints);
604  ncols = Array<OneD, unsigned int>(n_exp*NumPencils,NumModes);
605  }
606 
607  MatrixStorage blkmatStorage = eDIAGONAL;
608  BlkMatrix = MemoryManager<DNekBlkMat>::AllocateSharedPtr(nrows,ncols,blkmatStorage);
609 
610  StdRegions::StdSegExp StdSeg(Basis->GetBasisKey());
611 
612  if((mattype == eForwardsCoeffSpaceY1D)||(mattype == eForwardsPhysSpaceY1D) ||
613  (mattype == eForwardsCoeffSpaceZ1D)||(mattype == eForwardsPhysSpaceZ1D))
614  {
616  StdSeg.DetShapeType(),
617  StdSeg);
618 
619  loc_mat = StdSeg.GetStdMatrix(matkey);
620  }
621  else
622  {
624  StdSeg.DetShapeType(),
625  StdSeg);
626 
627  loc_mat = StdSeg.GetStdMatrix(matkey);
628  }
629 
630  // set up array of block matrices.
631  for(i = 0; i < (n_exp*NumPencils); ++i)
632  {
633  BlkMatrix->SetBlock(i,i,loc_mat);
634  }
635 
636  return BlkMatrix;
637  }
638 
639  std::vector<LibUtilities::FieldDefinitionsSharedPtr> ExpListHomogeneous2D::v_GetFieldDefinitions()
640  {
641  std::vector<LibUtilities::FieldDefinitionsSharedPtr> returnval;
642  // Set up Homogeneous length details.
644  HomoBasis[0] = m_homogeneousBasis_y;
645  HomoBasis[1] = m_homogeneousBasis_z;
646 
647  std::vector<NekDouble> HomoLen(2);
648  HomoLen[0] = m_lhom_y;
649  HomoLen[1] = m_lhom_z;
650 
651  int nhom_modes_y = m_homogeneousBasis_y->GetNumModes();
652  int nhom_modes_z = m_homogeneousBasis_z->GetNumModes();
653 
654  std::vector<unsigned int> sIDs
656 
657  std::vector<unsigned int> yIDs;
658  std::vector<unsigned int> zIDs;
659 
660  for(int n = 0; n < nhom_modes_z; ++n)
661  {
662  for(int m = 0; m < nhom_modes_y; ++m)
663  {
664  zIDs.push_back(n);
665  yIDs.push_back(m);
666  }
667  }
668 
669  m_lines[0]->GeneralGetFieldDefinitions(returnval, 2, HomoBasis,
670  HomoLen, false,
671  sIDs, zIDs, yIDs);
672  return returnval;
673  }
674 
675  void ExpListHomogeneous2D::v_GetFieldDefinitions(std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
676  {
677  // Set up Homogeneous length details.
679  HomoBasis[0] = m_homogeneousBasis_y;
680  HomoBasis[1] = m_homogeneousBasis_z;
681  std::vector<NekDouble> HomoLen(2);
682  HomoLen[0] = m_lhom_y;
683  HomoLen[1] = m_lhom_z;
684 
685  int nhom_modes_y = m_homogeneousBasis_y->GetNumModes();
686  int nhom_modes_z = m_homogeneousBasis_z->GetNumModes();
687 
688  std::vector<unsigned int> sIDs
690 
691  std::vector<unsigned int> yIDs;
692  std::vector<unsigned int> zIDs;
693 
694  for(int n = 0; n < nhom_modes_z; ++n)
695  {
696  for(int m = 0; m < nhom_modes_y; ++m)
697  {
698  zIDs.push_back(n);
699  yIDs.push_back(m);
700  }
701  }
702 
703  // enforce NumHomoDir == 1 by direct call
704  m_lines[0]->GeneralGetFieldDefinitions(fielddef, 2, HomoBasis,
705  HomoLen, false,
706  sIDs, zIDs, yIDs);
707  }
708 
710  {
711  int i,k;
712 
713  int NumMod_y = m_homogeneousBasis_y->GetNumModes();
714  int NumMod_z = m_homogeneousBasis_z->GetNumModes();
715 
716  int ncoeffs_per_line = m_lines[0]->GetNcoeffs();
717 
718  // Determine mapping from element ids to location in
719  // expansion list
720  map<int, int> ElmtID_to_ExpID;
721  for(i = 0; i < m_lines[0]->GetExpSize(); ++i)
722  {
723  ElmtID_to_ExpID[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
724  }
725 
726  for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
727  {
728  int eid = ElmtID_to_ExpID[fielddef->m_elementIDs[i]];
729  int datalen = (*m_exp)[eid]->GetNcoeffs();
730 
731  for(k = 0; k < (NumMod_y*NumMod_z); ++k)
732  {
733  fielddata.insert(fielddata.end(),&coeffs[m_coeff_offset[eid]+k*ncoeffs_per_line],&coeffs[m_coeff_offset[eid]+k*ncoeffs_per_line]+datalen);
734  }
735  }
736  }
737 
738  void ExpListHomogeneous2D::v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector<NekDouble> &fielddata)
739  {
740  v_AppendFieldData(fielddef,fielddata,m_coeffs);
741  }
742 
743  //Extract the data in fielddata into the m_coeff list
744  void ExpListHomogeneous2D::v_ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector<NekDouble> &fielddata, std::string &field, Array<OneD, NekDouble> &coeffs)
745  {
746  int i,k;
747  int offset = 0;
748  int datalen = fielddata.size()/fielddef->m_fields.size();
749  int ncoeffs_per_line = m_lines[0]->GetNcoeffs();
750  int NumMod_y = m_homogeneousBasis_y->GetNumModes();
751  int NumMod_z = m_homogeneousBasis_z->GetNumModes();
752 
753  // Find data location according to field definition
754  for(i = 0; i < fielddef->m_fields.size(); ++i)
755  {
756  if(fielddef->m_fields[i] == field)
757  {
758  break;
759  }
760  offset += datalen;
761  }
762 
763  ASSERTL0(i!= fielddef->m_fields.size(),"Field not found in data file");
764 
765  // Determine mapping from element ids to location in
766  // expansion list
767  map<int, int> ElmtID_to_ExpID;
768  for(i = 0; i < m_lines[0]->GetExpSize(); ++i)
769  {
770  ElmtID_to_ExpID[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
771  }
772 
773  for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
774  {
775  int eid = ElmtID_to_ExpID[fielddef->m_elementIDs[i]];
776  int datalen = (*m_exp)[eid]->GetNcoeffs();
777 
778  for(k = 0; k < (NumMod_y*NumMod_z); ++k)
779  {
780  Vmath::Vcopy(datalen,&fielddata[offset],1,&coeffs[m_coeff_offset[eid] + k*ncoeffs_per_line],1);
781  offset += datalen;
782  }
783  }
784  }
785 
786  void ExpListHomogeneous2D::v_WriteVtkPieceData(std::ostream &outfile, int expansion,
787  std::string var)
788  {
789  int i;
790  int nq = (*m_exp)[expansion]->GetTotPoints();
791  int npoints_per_line = m_lines[0]->GetTotPoints();
792 
793  // printing the fields of that zone
794  outfile << " <DataArray type=\"Float64\" Name=\""
795  << var << "\">" << endl;
796  outfile << " ";
797  for (int n = 0; n < m_lines.num_elements(); ++n)
798  {
799  const Array<OneD, NekDouble> phys = m_phys + m_phys_offset[expansion] + n*npoints_per_line;
800  for(i = 0; i < nq; ++i)
801  {
802  outfile << (fabs(phys[i]) < NekConstants::kNekZeroTol ? 0 : phys[i]) << " ";
803  }
804  }
805  outfile << endl;
806  outfile << " </DataArray>" << endl;
807  }
808 
810  Array<OneD, NekDouble> &out_d0,
811  Array<OneD, NekDouble> &out_d1,
812  Array<OneD, NekDouble> &out_d2)
813 
814  {
815  int nyzlines = m_lines.num_elements(); //number of Fourier points in the Fourier directions (nF_pts)
816  int npoints = inarray.num_elements(); //number of total points = n. of Fourier points * n. of points per line (nT_pts)
817  int n_points_line = npoints/nyzlines; //number of points per line
818 
819  Array<OneD, NekDouble> temparray(npoints);
820  Array<OneD, NekDouble> temparray1(npoints);
821  Array<OneD, NekDouble> temparray2(npoints);
825 
826  for( int i=0 ; i<nyzlines ; i++ )
827  {
828  m_lines[i]->PhysDeriv( tmp1 = inarray + i*n_points_line ,tmp2 = out_d0 + i*n_points_line);
829  }
830 
832  {
833  if(m_WaveSpace)
834  {
835  temparray = inarray;
836  }
837  else
838  {
839  HomogeneousFwdTrans(inarray,temparray);
840  }
841  NekDouble sign = -1.0;
842  NekDouble beta;
843 
844  //along y
845  for(int i = 0; i < m_ny; i++)
846  {
847  beta = -sign*2*M_PI*(i/2)/m_lhom_y;
848 
849  for(int j = 0; j < m_nz; j++)
850  {
851  Vmath::Smul(n_points_line,beta,tmp1 = temparray + n_points_line*(i+j*m_ny),1, tmp2 = temparray1 + n_points_line*((i-int(sign))+j*m_ny),1);
852  }
853 
854  sign = -1.0*sign;
855  }
856 
857  //along z
858  sign = -1.0;
859  for(int i = 0; i < m_nz; i++)
860  {
861  beta = -sign*2*M_PI*(i/2)/m_lhom_z;
862  Vmath::Smul(m_ny*n_points_line,beta,tmp1 = temparray + i*m_ny*n_points_line,1,tmp2 = temparray2 + (i-int(sign))*m_ny*n_points_line,1);
863  sign = -1.0*sign;
864  }
865  if(m_WaveSpace)
866  {
867  out_d1 = temparray1;
868  out_d2 = temparray2;
869  }
870  else
871  {
872  HomogeneousBwdTrans(temparray1,out_d1);
873  HomogeneousBwdTrans(temparray2,out_d2);
874  }
875  }
876  else
877  {
878  if(m_WaveSpace)
879  {
880  ASSERTL0(false,"Semi-phyisical time-stepping not implemented yet for non-Fourier basis")
881  }
882  else
883  {
884  StdRegions::StdQuadExp StdQuad(m_homogeneousBasis_y->GetBasisKey(),m_homogeneousBasis_z->GetBasisKey());
885 
886  m_transposition->Transpose(inarray,temparray,false,LibUtilities::eXtoYZ);
887 
888  for(int i = 0; i < n_points_line; i++)
889  {
890  StdQuad.PhysDeriv(tmp1 = temparray + i*nyzlines, tmp2 = temparray1 + i*nyzlines, tmp3 = temparray2 + i*nyzlines);
891  }
892 
893  m_transposition->Transpose(temparray1,out_d1,false,LibUtilities::eYZtoX);
894  m_transposition->Transpose(temparray2,out_d2,false,LibUtilities::eYZtoX);
895  Vmath::Smul(npoints,2.0/m_lhom_y,out_d1,1,out_d1,1);
896  Vmath::Smul(npoints,2.0/m_lhom_z,out_d2,1,out_d2,1);
897  }
898  }
899  }
900 
902  const Array<OneD, const NekDouble> &inarray,
903  Array<OneD, NekDouble> &out_d)
904 
905  {
906  int nyzlines = m_lines.num_elements(); //number of Fourier points in the Fourier directions (nF_pts)
907  int npoints = inarray.num_elements(); //number of total points = n. of Fourier points * n. of points per line (nT_pts)
908  int n_points_line = npoints/nyzlines; //number of points per line
909  //convert enum into int
910  int dir = (int)edir;
911 
912  Array<OneD, NekDouble> temparray(npoints);
913  Array<OneD, NekDouble> temparray1(npoints);
914  Array<OneD, NekDouble> temparray2(npoints);
918 
919  if (dir < 1)
920  {
921  for( int i=0 ; i<nyzlines ; i++)
922  {
923  m_lines[i]->PhysDeriv( tmp1 = inarray + i*n_points_line ,tmp2 = out_d + i*n_points_line);
924  }
925  }
926  else
927  {
929  {
930  if(m_WaveSpace)
931  {
932  temparray = inarray;
933  }
934  else
935  {
936  HomogeneousFwdTrans(inarray,temparray);
937  }
938  NekDouble sign = -1.0;
939  NekDouble beta;
940 
941  if (dir == 1)
942  {
943  //along y
944  for(int i = 0; i < m_ny; i++)
945  {
946  beta = -sign*2*M_PI*(i/2)/m_lhom_y;
947 
948  for(int j = 0; j < m_nz; j++)
949  {
950  Vmath::Smul(n_points_line,beta,tmp1 = temparray + n_points_line*(i+j*m_ny),1, tmp2 = temparray1 + n_points_line*((i-int(sign))+j*m_ny),1);
951  }
952  sign = -1.0*sign;
953  }
954  if(m_WaveSpace)
955  {
956  out_d = temparray1;
957  }
958  else
959  {
960  HomogeneousBwdTrans(temparray1,out_d);
961  }
962  }
963  else
964  {
965  //along z
966  for(int i = 0; i < m_nz; i++)
967  {
968  beta = -sign*2*M_PI*(i/2)/m_lhom_z;
969  Vmath::Smul(m_ny*n_points_line,beta,tmp1 = temparray + i*m_ny*n_points_line,1,tmp2 = temparray2 + (i-int(sign))*m_ny*n_points_line,1);
970  sign = -1.0*sign;
971  }
972  if(m_WaveSpace)
973  {
974  out_d = temparray2;
975  }
976  else
977  {
978  HomogeneousBwdTrans(temparray2,out_d);
979  }
980  }
981  }
982  else
983  {
984  if(m_WaveSpace)
985  {
986  ASSERTL0(false,"Semi-phyisical time-stepping not implemented yet for non-Fourier basis")
987  }
988  else
989  {
990  StdRegions::StdQuadExp StdQuad(m_homogeneousBasis_y->GetBasisKey(),m_homogeneousBasis_z->GetBasisKey());
991 
992  m_transposition->Transpose(inarray,temparray,false,LibUtilities::eXtoYZ);
993 
994  for(int i = 0; i < n_points_line; i++)
995  {
996  StdQuad.PhysDeriv(tmp1 = temparray + i*nyzlines, tmp2 = temparray1 + i*nyzlines, tmp3 = temparray2 + i*nyzlines);
997  }
998 
999  if (dir == 1)
1000  {
1001  m_transposition->Transpose(temparray1,out_d,false,LibUtilities::eYZtoX);
1002  Vmath::Smul(npoints,2.0/m_lhom_y,out_d,1,out_d,1);
1003  }
1004  else
1005  {
1006  m_transposition->Transpose(temparray2,out_d,false,LibUtilities::eYZtoX);
1007  Vmath::Smul(npoints,2.0/m_lhom_z,out_d,1,out_d,1);
1008  }
1009  }
1010  }
1011  }
1012  }
1013 
1015  Array<OneD, NekDouble> &out_d0,
1016  Array<OneD, NekDouble> &out_d1,
1017  Array<OneD, NekDouble> &out_d2)
1018 
1019  {
1020  v_PhysDeriv(inarray,out_d0,out_d1,out_d2);
1021  }
1022 
1024  const Array<OneD, const NekDouble> &inarray,
1025  Array<OneD, NekDouble> &out_d)
1026  {
1027  //convert int into enum
1028  v_PhysDeriv(edir,inarray,out_d);
1029  }
1030 
1032  {
1033  NekDouble size_y = 1.5*m_ny;
1034  NekDouble size_z = 1.5*m_nz;
1035  m_padsize_y = int(size_y);
1036  m_padsize_z = int(size_z);
1037 
1040 
1043 
1046 
1047  StdRegions::StdQuadExp StdQuad(m_paddingBasis_y->GetBasisKey(),m_paddingBasis_z->GetBasisKey());
1048 
1049  StdRegions::StdMatrixKey matkey1(StdRegions::eFwdTrans,StdQuad.DetShapeType(),StdQuad);
1050  StdRegions::StdMatrixKey matkey2(StdRegions::eBwdTrans,StdQuad.DetShapeType(),StdQuad);
1051 
1052  MatFwdPAD = StdQuad.GetStdMatrix(matkey1);
1053  MatBwdPAD = StdQuad.GetStdMatrix(matkey2);
1054  }
1055  } //end of namespace
1056 } //end of namespace
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
LibUtilities::TranspositionSharedPtr m_transposition
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
Represents a basis of a given type.
Definition: Basis.h:211
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:16
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void)
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
LibUtilities::BasisSharedPtr m_paddingBasis_z
Base expansion in z direction.
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:365
STL namespace.
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1069
Fourier Expansion .
Definition: BasisType.h:53
LibUtilities::NektarFFTSharedPtr m_FFT_z
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
NekDouble m_lhom_z
Width of homogeneous direction z.
virtual void v_DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
std::shared_ptr< Basis > BasisSharedPtr
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1052
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
NektarFFTFactory & GetNektarFFTFactory()
Definition: NektarFFT.cpp:69
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:71
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
BasisManagerT & BasisManager(void)
virtual void v_IProductWRTBase_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Class representing a segment element in reference space.
Definition: StdSegExp.h:53
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
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
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
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_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_HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
NekDouble m_lhom_y
Width of homogeneous direction y.
static std::vector< unsigned int > NullUnsignedIntVector
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Defines a specification for a set of points.
Definition: Points.h:59
double NekDouble
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
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)
int m_ny
Number of modes = number of poitns in y direction.
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
void Homogeneous2DTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool IsForwards, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
std::map< Homogeneous2DMatType, DNekBlkMatSharedPtr > Homo2DBlockMatrixMap
A map between homo matrix keys and their associated block matrices.
int m_nz
Number of modes = number of poitns in z direction.
LibUtilities::BasisSharedPtr m_paddingBasis_y
Base expansion in y direction.
virtual void v_BwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1020
void DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1958
#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.
LibUtilities::NektarFFTSharedPtr m_FFT_y
DNekBlkMatSharedPtr GenHomogeneous2DBlockMatrix(Homogeneous2DMatType mattype, CoeffState coeffstate=eLocal) const
virtual void v_FwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var)
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.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
virtual void v_HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
#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
DNekBlkMatSharedPtr GetHomogeneous2DBlockMatrix(Homogeneous2DMatType mattype, CoeffState coeffstate=eLocal) const
Describes the specification for a Basis.
Definition: Basis.h:49
std::shared_ptr< SessionReader > SessionReaderSharedPtr
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:302
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_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)