Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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 
76 
78  m_homogeneousBasis->GetNumPoints() /
79  m_StripZcomm->GetSize());
80 
81  if(m_useFFT)
82  {
84  "NekFFTW", m_homogeneousBasis->GetNumPoints());
85  }
86 
87  if(m_dealiasing)
88  {
89  if(m_useFFT)
90  {
91  NekDouble size = 1.5*m_homogeneousBasis->GetNumPoints();
92  m_padsize = int(size);
94  .CreateInstance("NekFFTW", m_padsize);
95  }
96  else
97  {
98  ASSERTL0(false, "Dealiasing available just in combination "
99  "with FFTW");
100  }
101  }
102  }
103 
104 
105  /**
106  * @param In ExpListHomogeneous1D object to copy.
107  */
109  ExpList(In,false),
110  m_transposition(In.m_transposition),
111  m_StripZcomm(In.m_StripZcomm),
112  m_useFFT(In.m_useFFT),
113  m_FFT(In.m_FFT),
114  m_tmpIN(In.m_tmpIN),
115  m_tmpOUT(In.m_tmpOUT),
116  m_homogeneousBasis(In.m_homogeneousBasis),
117  m_lhom(In.m_lhom),
118  m_homogeneous1DBlockMat(In.m_homogeneous1DBlockMat),
119  m_dealiasing(In.m_dealiasing),
120  m_padsize(In.m_padsize)
121  {
122  m_planes = Array<OneD, ExpListSharedPtr>(In.m_planes.num_elements());
123  }
124 
126  const std::vector<unsigned int> &eIDs):
127  ExpList(In,eIDs,false),
128  m_transposition(In.m_transposition),
129  m_useFFT(In.m_useFFT),
130  m_FFT(In.m_FFT),
131  m_tmpIN(In.m_tmpIN),
132  m_tmpOUT(In.m_tmpOUT),
133  m_homogeneousBasis(In.m_homogeneousBasis),
134  m_lhom(In.m_lhom),
135  m_homogeneous1DBlockMat(MemoryManager<Homo1DBlockMatrixMap>::AllocateSharedPtr()),
136  m_dealiasing(In.m_dealiasing),
137  m_padsize(In.m_padsize)
138  {
139  m_planes = Array<OneD, ExpListSharedPtr>(In.m_planes.num_elements());
140  }
141 
142  /**
143  * Destructor
144  */
146  {
147  }
148 
150  Array<OneD, NekDouble> &outarray,
151  CoeffState coeffstate,
152  bool Shuff,
153  bool UnShuff)
154  {
155  // Forwards trans
156  Homogeneous1DTrans(inarray,outarray,true,coeffstate,Shuff,UnShuff);
157  }
158 
160  Array<OneD, NekDouble> &outarray,
161  CoeffState coeffstate,
162  bool Shuff,
163  bool UnShuff)
164  {
165  // Backwards trans
166  Homogeneous1DTrans(inarray,outarray,false,coeffstate,Shuff,UnShuff);
167  }
168 
169  /**
170  * Dealiasing routine
171  */
173  const Array<OneD, NekDouble> &inarray2,
174  Array<OneD, NekDouble> &outarray,
175  CoeffState coeffstate)
176  {
177  // inarray1 = first term of the product in full physical space
178  // inarray2 = second term of the product in full physical space
179  // dealiased product stored in outarray
180 
181  int num_dofs = inarray1.num_elements();
182 
183  int N = m_homogeneousBasis->GetNumPoints();
184 
185  Array<OneD, NekDouble> V1(num_dofs);
186  Array<OneD, NekDouble> V2(num_dofs);
187  Array<OneD, NekDouble> V1V2(num_dofs);
188 
189  HomogeneousFwdTrans(inarray1,V1,coeffstate);
190  HomogeneousFwdTrans(inarray2,V2,coeffstate);
191 
192  int num_points_per_plane = num_dofs/m_planes.num_elements();
193  int num_proc;
194  if(!m_session->DefinesSolverInfo("HomoStrip"))
195  {
196  num_proc = m_comm->GetColumnComm()->GetSize();
197  }
198  else
199  {
200  num_proc = m_StripZcomm->GetSize();
201  }
202  int num_dfts_per_proc = num_points_per_plane / num_proc
203  + (num_points_per_plane % num_proc > 0);
204 
205  Array<OneD, NekDouble> ShufV1(num_dfts_per_proc*N,0.0);
206  Array<OneD, NekDouble> ShufV2(num_dfts_per_proc*N,0.0);
207  Array<OneD, NekDouble> ShufV1V2(num_dfts_per_proc*N,0.0);
208 
209  Array<OneD, NekDouble> ShufV1_PAD_coef(m_padsize,0.0);
210  Array<OneD, NekDouble> ShufV2_PAD_coef(m_padsize,0.0);
211  Array<OneD, NekDouble> ShufV1_PAD_phys(m_padsize,0.0);
212  Array<OneD, NekDouble> ShufV2_PAD_phys(m_padsize,0.0);
213 
214  Array<OneD, NekDouble> ShufV1V2_PAD_coef(m_padsize,0.0);
215  Array<OneD, NekDouble> ShufV1V2_PAD_phys(m_padsize,0.0);
216 
217  m_transposition->Transpose(V1, ShufV1, false, LibUtilities::eXYtoZ);
218  m_transposition->Transpose(V2, ShufV2, false, LibUtilities::eXYtoZ);
219 
220  // Looping on the pencils
221  for(int i = 0 ; i < num_dfts_per_proc ; i++)
222  {
223  // Copying the i-th pencil pf lenght N into a bigger
224  // pencil of lenght 2N We are in Fourier space
225  Vmath::Vcopy(N, &(ShufV1[i*N]), 1, &(ShufV1_PAD_coef[0]), 1);
226  Vmath::Vcopy(N, &(ShufV2[i*N]), 1, &(ShufV2_PAD_coef[0]), 1);
227 
228  // Moving to physical space using the padded system
229  m_FFT_deal->FFTBwdTrans(ShufV1_PAD_coef, ShufV1_PAD_phys);
230  m_FFT_deal->FFTBwdTrans(ShufV2_PAD_coef, ShufV2_PAD_phys);
231 
232  // Perfroming the vectors multiplication in physical space on
233  // the padded system
234  Vmath::Vmul(m_padsize, ShufV1_PAD_phys, 1,
235  ShufV2_PAD_phys, 1,
236  ShufV1V2_PAD_phys, 1);
237 
238  // Moving back the result (V1*V2)_phys in Fourier space, padded
239  // system
240  m_FFT_deal->FFTFwdTrans(ShufV1V2_PAD_phys, ShufV1V2_PAD_coef);
241 
242  // Copying the first half of the padded pencil in the full
243  // vector (Fourier space)
244  Vmath::Vcopy(N, &(ShufV1V2_PAD_coef[0]), 1,
245  &(ShufV1V2[i*N]), 1);
246  }
247 
248  m_transposition->Transpose(ShufV1V2, V1V2, false,
250 
251  // Moving the results in physical space for the output
252  HomogeneousBwdTrans(V1V2, outarray, coeffstate);
253  }
254 
255  /**
256  * Forward transform
257  */
259  {
260  int cnt = 0, cnt1 = 0;
261  Array<OneD, NekDouble> tmparray;
262 
263  for(int n = 0; n < m_planes.num_elements(); ++n)
264  {
265  m_planes[n]->FwdTrans(inarray+cnt, tmparray = outarray + cnt1,
266  coeffstate);
267  cnt += m_planes[n]->GetTotPoints();
268 
269  cnt1 += m_planes[n]->GetNcoeffs(); // need to skip ncoeffs
270  }
271  if(!m_WaveSpace)
272  {
273  HomogeneousFwdTrans(outarray,outarray,coeffstate);
274  }
275  }
276 
277  /**
278  * Forward transform element by element
279  */
281  {
282  int cnt = 0, cnt1 = 0;
283  Array<OneD, NekDouble> tmparray;
284 
285  //spectral element FwdTrans plane by plane
286  for(int n = 0; n < m_planes.num_elements(); ++n)
287  {
288  m_planes[n]->FwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
289  cnt += m_planes[n]->GetTotPoints();
290  cnt1 += m_planes[n]->GetNcoeffs();
291  }
292  if(!m_WaveSpace)
293  {
294  HomogeneousFwdTrans(outarray,outarray);
295  }
296  }
297 
298  /**
299  * Backward transform
300  */
302  {
303  int cnt = 0, cnt1 = 0;
304  Array<OneD, NekDouble> tmparray;
305 
306  for(int n = 0; n < m_planes.num_elements(); ++n)
307  {
308  m_planes[n]->BwdTrans(inarray+cnt, tmparray = outarray + cnt1,
309  coeffstate);
310  cnt += m_planes[n]->GetNcoeffs();
311  cnt1 += m_planes[n]->GetTotPoints();
312  }
313  if(!m_WaveSpace)
314  {
315  HomogeneousBwdTrans(outarray,outarray);
316  }
317  }
318 
319  /**
320  * Backward transform element by element
321  */
323  {
324  int cnt = 0, cnt1 = 0;
325  Array<OneD, NekDouble> tmparray;
326 
327  for(int n = 0; n < m_planes.num_elements(); ++n)
328  {
329  m_planes[n]->BwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
330 
331  cnt += m_planes[n]->GetNcoeffs();
332  cnt1 += m_planes[n]->GetTotPoints();
333  }
334  if(!m_WaveSpace)
335  {
336  HomogeneousBwdTrans(outarray,outarray);
337  }
338  }
339 
340  /**
341  * Inner product
342  */
344  {
345  int cnt = 0, cnt1 = 0;
346  Array<OneD, NekDouble> tmparray;
347 
348  for(int n = 0; n < m_planes.num_elements(); ++n)
349  {
350  m_planes[n]->IProductWRTBase(inarray+cnt, tmparray = outarray + cnt1,coeffstate);
351 
352  cnt1 += m_planes[n]->GetNcoeffs();
353  cnt += m_planes[n]->GetTotPoints();
354  }
355  }
356 
357  /**
358  * Inner product element by element
359  */
361  {
362  int cnt = 0, cnt1 = 0;
363  Array<OneD, NekDouble> tmparray;
364 
365  for(int n = 0; n < m_planes.num_elements(); ++n)
366  {
367  m_planes[n]->IProductWRTBase_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
368 
369  cnt1 += m_planes[n]->GetNcoeffs();
370  cnt += m_planes[n]->GetTotPoints();
371  }
372  }
373 
374  /**
375  * Homogeneous transform Bwd/Fwd (MVM and FFT)
376  */
378  bool IsForwards,
379  CoeffState coeffstate,
380  bool Shuff,
381  bool UnShuff)
382  {
383  int num_dofs;
384 
385  if(IsForwards)
386  {
387  num_dofs = inarray.num_elements();
388  }
389  else
390  {
391  num_dofs = outarray.num_elements();
392  }
393 
394  if(m_useFFT)
395  {
396  int num_points_per_plane = num_dofs/m_planes.num_elements();
397  int num_dfts_per_proc;
398  if(!m_session->DefinesSolverInfo("HomoStrip"))
399  {
400  int nP = m_comm->GetColumnComm()->GetSize();
401  num_dfts_per_proc = num_points_per_plane / nP
402  + (num_points_per_plane % nP > 0);
403  }
404  else
405  {
406  int nP = m_StripZcomm->GetSize();
407  num_dfts_per_proc = num_points_per_plane / nP
408  + (num_points_per_plane % nP > 0);
409  }
410 
411  Array<OneD, NekDouble> fft_in (num_dfts_per_proc*m_homogeneousBasis->GetNumPoints(),0.0);
412  Array<OneD, NekDouble> fft_out(num_dfts_per_proc*m_homogeneousBasis->GetNumPoints(),0.0);
413 
414  if(Shuff)
415  {
416  m_transposition->Transpose(inarray,fft_in,false,LibUtilities::eXYtoZ);
417  }
418  else
419  {
420  Vmath::Vcopy(num_dfts_per_proc*m_homogeneousBasis->GetNumPoints(),
421  inarray,1,fft_in,1);
422  }
423 
424  if(IsForwards)
425  {
426  for(int i = 0 ; i < num_dfts_per_proc ; i++)
427  {
428  m_FFT->FFTFwdTrans(m_tmpIN = fft_in + i*m_homogeneousBasis->GetNumPoints(), m_tmpOUT = fft_out + i*m_homogeneousBasis->GetNumPoints());
429  }
430  }
431  else
432  {
433  for(int i = 0 ; i < num_dfts_per_proc ; i++)
434  {
435  m_FFT->FFTBwdTrans(m_tmpIN = fft_in + i*m_homogeneousBasis->GetNumPoints(), m_tmpOUT = fft_out + i*m_homogeneousBasis->GetNumPoints());
436  }
437  }
438 
439  if(UnShuff)
440  {
441  m_transposition->Transpose(fft_out,outarray,false,LibUtilities::eZtoXY);
442  }
443  else
444  {
445  Vmath::Vcopy(num_dfts_per_proc*m_homogeneousBasis->GetNumPoints(),
446  fft_out,1,outarray,1);
447  }
448  }
449  else
450  {
451  DNekBlkMatSharedPtr blkmat;
452 
453  if(num_dofs == m_npoints) //transform phys space
454  {
455  if(IsForwards)
456  {
458  }
459  else
460  {
462  }
463  }
464  else
465  {
466  if(IsForwards)
467  {
469  }
470  else
471  {
473  }
474  }
475 
476  int nrows = blkmat->GetRows();
477  int ncols = blkmat->GetColumns();
478 
479  Array<OneD, NekDouble> sortedinarray(ncols,0.0);
480  Array<OneD, NekDouble> sortedoutarray(nrows,0.0);
481 
482  if(Shuff)
483  {
484  m_transposition->Transpose(inarray,sortedinarray,!IsForwards,LibUtilities::eXYtoZ);
485  }
486  else
487  {
488  Vmath::Vcopy(ncols,inarray,1,sortedinarray,1);
489  }
490 
491  // Create NekVectors from the given data arrays
492  NekVector<NekDouble> in (ncols,sortedinarray,eWrapper);
493  NekVector<NekDouble> out(nrows,sortedoutarray,eWrapper);
494 
495  // Perform matrix-vector multiply.
496  out = (*blkmat)*in;
497 
498  if(UnShuff)
499  {
500  m_transposition->Transpose(sortedoutarray,outarray,IsForwards,LibUtilities::eZtoXY);
501  }
502  else
503  {
504  Vmath::Vcopy(nrows,sortedoutarray,1,outarray,1);
505  }
506 
507  }
508  }
509 
511  {
512  Homo1DBlockMatrixMap::iterator matrixIter = m_homogeneous1DBlockMat->find(mattype);
513 
514  if(matrixIter == m_homogeneous1DBlockMat->end())
515  {
516  return ((*m_homogeneous1DBlockMat)[mattype] =
517  GenHomogeneous1DBlockMatrix(mattype,coeffstate));
518  }
519  else
520  {
521  return matrixIter->second;
522  }
523  }
524 
525 
527  {
528  DNekMatSharedPtr loc_mat;
529  DNekBlkMatSharedPtr BlkMatrix;
530  int n_exp = 0;
531  int num_trans_per_proc = 0;
532 
533  if((mattype == eForwardsCoeffSpace1D)
534  ||(mattype == eBackwardsCoeffSpace1D)) // will operate on m_coeffs
535  {
536  n_exp = m_planes[0]->GetNcoeffs();
537  }
538  else
539  {
540  n_exp = m_planes[0]->GetTotPoints(); // will operatore on m_phys
541  }
542 
543  num_trans_per_proc = n_exp/m_comm->GetColumnComm()->GetSize() + (n_exp%m_comm->GetColumnComm()->GetSize() > 0);
544 
545  Array<OneD,unsigned int> nrows(num_trans_per_proc);
546  Array<OneD,unsigned int> ncols(num_trans_per_proc);
547 
548  if((mattype == eForwardsCoeffSpace1D)||(mattype == eForwardsPhysSpace1D))
549  {
550  nrows = Array<OneD, unsigned int>(num_trans_per_proc,m_homogeneousBasis->GetNumModes());
551  ncols = Array<OneD, unsigned int>(num_trans_per_proc,m_homogeneousBasis->GetNumPoints());
552  }
553  else
554  {
555  nrows = Array<OneD, unsigned int>(num_trans_per_proc,m_homogeneousBasis->GetNumPoints());
556  ncols = Array<OneD, unsigned int>(num_trans_per_proc,m_homogeneousBasis->GetNumModes());
557  }
558 
559  MatrixStorage blkmatStorage = eDIAGONAL;
560  BlkMatrix = MemoryManager<DNekBlkMat>
561  ::AllocateSharedPtr(nrows,ncols,blkmatStorage);
562 
563  //Half Mode
565  {
566  StdRegions::StdPointExp StdPoint(m_homogeneousBasis->GetBasisKey());
567 
568  if((mattype == eForwardsCoeffSpace1D)||(mattype == eForwardsPhysSpace1D))
569  {
571  StdPoint.DetShapeType(),
572  StdPoint);
573 
574  loc_mat = StdPoint.GetStdMatrix(matkey);
575  }
576  else
577  {
579  StdPoint.DetShapeType(),
580  StdPoint);
581 
582  loc_mat = StdPoint.GetStdMatrix(matkey);
583  }
584  }
585  //other cases
586  else
587  {
588  StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
589 
590  if((mattype == eForwardsCoeffSpace1D)||(mattype == eForwardsPhysSpace1D))
591  {
593  StdSeg.DetShapeType(),
594  StdSeg);
595 
596  loc_mat = StdSeg.GetStdMatrix(matkey);
597  }
598  else
599  {
601  StdSeg.DetShapeType(),
602  StdSeg);
603 
604  loc_mat = StdSeg.GetStdMatrix(matkey);
605  }
606 
607  }
608 
609  // set up array of block matrices.
610  for(int i = 0; i < num_trans_per_proc; ++i)
611  {
612  BlkMatrix->SetBlock(i,i,loc_mat);
613  }
614 
615  return BlkMatrix;
616  }
617 
618  std::vector<LibUtilities::FieldDefinitionsSharedPtr> ExpListHomogeneous1D::v_GetFieldDefinitions()
619  {
620  std::vector<LibUtilities::FieldDefinitionsSharedPtr> returnval;
621 
622  // Set up Homogeneous length details.
624 
625  std::vector<NekDouble> HomoLen;
626  HomoLen.push_back(m_lhom);
627 
628  std::vector<unsigned int> StripsIDs;
629 
630  bool strips;
631  m_session->MatchSolverInfo("HomoStrip","True",strips,false);
632  if (strips)
633  {
634  StripsIDs.push_back(m_transposition->GetStripID());
635  }
636 
637  std::vector<unsigned int> PlanesIDs;
638  int IDoffset = 0;
639 
640  // introduce a 2 plane offset for single mode case so can
641  // be post-processed or used in MultiMode expansion.
643  {
644  IDoffset = 2;
645  }
646 
647  for(int i = 0; i < m_planes.num_elements(); i++)
648  {
649  PlanesIDs.push_back(m_transposition->GetPlaneID(i)+IDoffset);
650  }
651 
652  m_planes[0]->GeneralGetFieldDefinitions(returnval, 1, HomoBasis,
653  HomoLen, strips, StripsIDs, PlanesIDs);
654  return returnval;
655  }
656 
657  void ExpListHomogeneous1D::v_GetFieldDefinitions(std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
658  {
659  // Set up Homogeneous length details.
661 
662  std::vector<NekDouble> HomoLen;
663  HomoLen.push_back(m_lhom);
664 
665  std::vector<unsigned int> StripsIDs;
666 
667  bool strips;
668  m_session->MatchSolverInfo("HomoStrip","True",strips,false);
669  if (strips)
670  {
671  StripsIDs.push_back(m_transposition->GetStripID());
672  }
673 
674  std::vector<unsigned int> PlanesIDs;
675  int IDoffset = 0;
676 
678  {
679  IDoffset = 2;
680  }
681 
682  for(int i = 0; i < m_planes.num_elements(); i++)
683  {
684  PlanesIDs.push_back(m_transposition->GetPlaneID(i)+IDoffset);
685  }
686 
687  // enforce NumHomoDir == 1 by direct call
688  m_planes[0]->GeneralGetFieldDefinitions(fielddef, 1, HomoBasis,
689  HomoLen, strips, StripsIDs, PlanesIDs);
690  }
691 
692 
693  /** This routine appends the data from the expansion list into
694  the output format where each element is given by looping
695  over its Fourier modes where as data in the expandion is
696  stored with all consecutive elements and then the Fourier
697  modes
698  */
700  {
701  int i,n;
702  int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
703 
704  // Determine mapping from element ids to location in
705  // expansion list
706  map<int, int> ElmtID_to_ExpID;
707  for(i = 0; i < m_planes[0]->GetExpSize(); ++i)
708  {
709  ElmtID_to_ExpID[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
710  }
711 
712  for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
713  {
714  int eid = ElmtID_to_ExpID[fielddef->m_elementIDs[i]];
715  int datalen = (*m_exp)[eid]->GetNcoeffs();
716 
717  for(n = 0; n < m_planes.num_elements(); ++n)
718  {
719  fielddata.insert(fielddata.end(),&coeffs[m_coeff_offset[eid]+n*ncoeffs_per_plane],&coeffs[m_coeff_offset[eid]+n*ncoeffs_per_plane]+datalen);
720  }
721  }
722  }
723 
724  void ExpListHomogeneous1D::v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector<NekDouble> &fielddata)
725  {
726  v_AppendFieldData(fielddef,fielddata,m_coeffs);
727  }
728 
729  //Extract the data in fielddata into the m_coeff list
732  std::vector<NekDouble> &fielddata,
733  std::string &field,
734  Array<OneD, NekDouble> &coeffs)
735  {
736  int i,n;
737  int offset = 0;
738  int nzmodes = 1;
739  int datalen = fielddata.size()/fielddef->m_fields.size();
740  std::vector<unsigned int> fieldDefHomoZids;
741 
742 
743  // Find data location according to field definition
744  for(i = 0; i < fielddef->m_fields.size(); ++i)
745  {
746  if(fielddef->m_fields[i] == field)
747  {
748  break;
749  }
750  offset += datalen;
751  }
752 
753  if(i == fielddef->m_fields.size())
754  {
755  cout << "Field "<< field<< "not found in data file. " << endl;
756  }
757  else
758  {
759 
760  int modes_offset = 0;
761  int planes_offset = 0;
762  Array<OneD, NekDouble> coeff_tmp;
764 
765  // Build map of plane IDs lying on this processor.
766  std::map<int,int> homoZids;
767  for (i = 0; i < m_planes.num_elements(); ++i)
768  {
769  homoZids[m_transposition->GetPlaneID(i)] = i;
770  }
771 
772  if(fielddef->m_numHomogeneousDir)
773  {
774  nzmodes = fielddef->m_homogeneousZIDs.size();
775  fieldDefHomoZids = fielddef->m_homogeneousZIDs;
776  }
777  else // input file is 2D and so set nzmodes to 1
778  {
779  nzmodes = 1;
780  fieldDefHomoZids.push_back(0);
781  }
782 
783  // Determine mapping from element ids to location in expansion list.
784  map<int, int> ElmtID_to_ExpID;
785  for(i = 0; i < m_planes[0]->GetExpSize(); ++i)
786  {
787  ElmtID_to_ExpID[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
788  }
789 
790 
791  // calculate number of modes in the current partition
792  int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
793 
794  for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
795  {
796  if(fielddef->m_uniOrder == true) // reset modes_offset to zero
797  {
798  modes_offset = 0;
799  }
800 
801  int datalen = LibUtilities::GetNumberOfCoefficients(fielddef->m_shapeType,
802  fielddef->m_numModes,
803  modes_offset);
804 
805  it = ElmtID_to_ExpID.find(fielddef->m_elementIDs[i]);
806 
807  // ensure element is on this partition for parallel case.
808  if(it == ElmtID_to_ExpID.end())
809  {
810  // increase offset for correct FieldData access
811  offset += datalen*nzmodes;
812  modes_offset += (*m_exp)[0]->GetNumBases() +
813  fielddef->m_numHomogeneousDir;
814  continue;
815  }
816 
817  int eid = it->second;
818 
819 
820  for(n = 0; n < nzmodes; ++n, offset += datalen)
821  {
822 
823  it = homoZids.find(fieldDefHomoZids[n]);
824 
825  // Check to make sure this mode number lies in this field.
826  if (it == homoZids.end())
827  {
828  continue;
829  }
830 
831  planes_offset = it->second;
832  if(datalen == (*m_exp)[eid]->GetNcoeffs())
833  {
834  Vmath::Vcopy(datalen,&fielddata[offset],1,&coeffs[m_coeff_offset[eid]+planes_offset*ncoeffs_per_plane],1);
835  }
836  else // unpack data to new order
837  {
838  (*m_exp)[eid]->ExtractDataToCoeffs(&fielddata[offset], fielddef->m_numModes,modes_offset,&coeffs[m_coeff_offset[eid] + planes_offset*ncoeffs_per_plane]);
839  }
840  }
841  modes_offset += (*m_exp)[0]->GetNumBases() + fielddef->m_numHomogeneousDir;
842  }
843  }
844  }
845 
846  //Extract the data in fielddata into the m_coeff list
848  const boost::shared_ptr<ExpList> &fromExpList,const Array<OneD, const NekDouble> &fromCoeffs, Array<OneD, NekDouble> &toCoeffs)
849  {
850  int i;
851  int fromNcoeffs_per_plane = fromExpList->GetPlane(0)->GetNcoeffs();
852  int toNcoeffs_per_plane = m_planes[0]->GetNcoeffs();
853  Array<OneD, NekDouble> tocoeffs_tmp, fromcoeffs_tmp;
854 
855  for(i = 0; i < m_planes.num_elements(); ++i)
856  {
857  m_planes[i]->ExtractCoeffsToCoeffs(fromExpList->GetPlane(i),fromcoeffs_tmp = fromCoeffs + fromNcoeffs_per_plane*i, tocoeffs_tmp = toCoeffs + toNcoeffs_per_plane*i);
858  }
859  }
860 
861  void ExpListHomogeneous1D::v_WriteVtkPieceData(std::ostream &outfile, int expansion,
862  std::string var)
863  {
864  // If there is only one plane (e.g. HalfMode), we write a 2D plane.
865  if (m_planes.num_elements() == 1)
866  {
867  m_planes[0]->WriteVtkPieceData(outfile, expansion, var);
868  return;
869  }
870 
871  int i;
872  int nq = (*m_exp)[expansion]->GetTotPoints();
873  int npoints_per_plane = m_planes[0]->GetTotPoints();
874 
875  // If we are using Fourier points, output extra plane to fill domain
876  int outputExtraPlane = 0;
877  Array<OneD, NekDouble> extraPlane;
878  if ( m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier
879  && m_homogeneousBasis->GetPointsType() ==
881  {
882  outputExtraPlane = 1;
883  // Get extra plane data
884  if (m_StripZcomm->GetSize() == 1)
885  {
886  extraPlane = m_phys + m_phys_offset[expansion];
887  }
888  else
889  {
890  // Determine to and from rank for communication
891  int size = m_StripZcomm->GetSize();
892  int rank = m_StripZcomm->GetRank();
893  int fromRank = (rank+1) % size;
894  int toRank = (rank == 0) ? size-1 : rank-1;
895  // Communicate using SendRecv
896  extraPlane = Array<OneD, NekDouble>(nq);
897  Array<OneD, NekDouble> send (nq,
898  m_phys + m_phys_offset[expansion]);
899  m_StripZcomm->SendRecv(toRank, send,
900  fromRank, extraPlane);
901  }
902  }
903 
904  // printing the fields of that zone
905  outfile << " <DataArray type=\"Float64\" Name=\""
906  << var << "\">" << endl;
907  outfile << " ";
908  for (int n = 0; n < m_planes.num_elements(); ++n)
909  {
910  const Array<OneD, NekDouble> phys = m_phys + m_phys_offset[expansion] + n*npoints_per_plane;
911  for(i = 0; i < nq; ++i)
912  {
913  outfile << (fabs(phys[i]) < NekConstants::kNekZeroTol ? 0 : phys[i]) << " ";
914  }
915  }
916  if (outputExtraPlane)
917  {
918  for(i = 0; i < nq; ++i)
919  {
920  outfile << (fabs(extraPlane[i]) < NekConstants::kNekZeroTol ?
921  0 : extraPlane[i]) << " ";
922  }
923  }
924  outfile << endl;
925  outfile << " </DataArray>" << endl;
926  }
927 
929  {
930  int cnt,cnt1;
931  Array<OneD, NekDouble> tmparray;
932  cnt = m_planes[0]->GetTotPoints();
933  cnt1 = m_planes[0]->Get1DScaledTotPoints(scale);
934 
935  ASSERTL1(m_planes.num_elements()*cnt1 <= outarray.num_elements(),"size of outarray does not match internal estimage");
936 
937 
938  for(int i = 0; i < m_planes.num_elements(); i++)
939  {
940 
941  m_planes[i]->PhysInterp1DScaled(scale,inarray+i*cnt,
942  tmparray = outarray+i*cnt1);
943  }
944  }
945 
946 
948  {
949  int cnt,cnt1;
950  Array<OneD, NekDouble> tmparray;
951  cnt = m_planes[0]->Get1DScaledTotPoints(scale);
952  cnt1 = m_planes[0]->GetTotPoints();
953 
954  ASSERTL1(m_planes.num_elements()*cnt <= inarray.num_elements(),"size of outarray does not match internal estimage");
955 
956 
957  for(int i = 0; i < m_planes.num_elements(); i++)
958  {
959  m_planes[i]->PhysGalerkinProjection1DScaled(scale,inarray+i*cnt,
960  tmparray = outarray+i*cnt1);
961  }
962 
963  }
965  Array<OneD, NekDouble> &out_d0,
966  Array<OneD, NekDouble> &out_d1,
967  Array<OneD, NekDouble> &out_d2)
968  {
969  int nT_pts = inarray.num_elements(); //number of total points = n. of Fourier points * n. of points per plane (nT_pts)
970  int nP_pts = nT_pts/m_planes.num_elements(); //number of points per plane = n of Fourier transform required (nP_pts)
971 
972  Array<OneD, NekDouble> temparray(nT_pts);
973  Array<OneD, NekDouble> outarray(nT_pts);
977 
978  for(int i = 0; i < m_planes.num_elements(); i++)
979  {
980  m_planes[i]->PhysDeriv(inarray + i*nP_pts ,tmp2 = out_d0 + i*nP_pts , tmp3 = out_d1 + i*nP_pts );
981  }
982 
983  if(out_d2 != NullNekDouble1DArray)
984  {
987  {
988  if(m_WaveSpace)
989  {
990  temparray = inarray;
991  }
992  else
993  {
994  HomogeneousFwdTrans(inarray,temparray);
995  }
996 
997  NekDouble sign = -1.0;
998  NekDouble beta;
999 
1000  //Half Mode
1002  {
1003  beta = sign*2*M_PI*(m_transposition->GetK(0))/m_lhom;
1004 
1005  Vmath::Smul(nP_pts,beta,temparray,1,outarray,1);
1006  }
1007  else if(m_homogeneousBasis->GetBasisType() == LibUtilities::eFourierHalfModeIm)
1008  {
1009  beta = -sign*2*M_PI*(m_transposition->GetK(0))/m_lhom;
1010 
1011  Vmath::Smul(nP_pts,beta,temparray,1,outarray,1);
1012  }
1013 
1014  //Fully complex
1015  else
1016  {
1017  for(int i = 0; i < m_planes.num_elements(); i++)
1018  {
1019  beta = -sign*2*M_PI*(m_transposition->GetK(i))/m_lhom;
1020 
1021  Vmath::Smul(nP_pts,beta,tmp1 = temparray + i*nP_pts,1,tmp2 = outarray + (i-int(sign))*nP_pts,1);
1022 
1023  sign = -1.0*sign;
1024  }
1025  }
1026 
1027  if(m_WaveSpace)
1028  {
1029  out_d2 = outarray;
1030  }
1031  else
1032  {
1033  HomogeneousBwdTrans(outarray,out_d2);
1034  }
1035  }
1036  else
1037  {
1038  if(!m_session->DefinesSolverInfo("HomoStrip"))
1039  {
1040  ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,
1041  "Parallelisation in the homogeneous direction "
1042  "implemented just for Fourier basis");
1043  }
1044  else
1045  {
1046  ASSERTL0(m_StripZcomm->GetSize() == 1,
1047  "Parallelisation in the homogeneous direction "
1048  "implemented just for Fourier basis");
1049  }
1050 
1051  if(m_WaveSpace)
1052  {
1053  ASSERTL0(false, "Semi-phyisical time-stepping not "
1054  "implemented yet for non-Fourier "
1055  "basis");
1056  }
1057  else
1058  {
1059  StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
1060 
1061  m_transposition->Transpose(inarray,temparray,false,LibUtilities::eXYtoZ);
1062 
1063  for(int i = 0; i < nP_pts; i++)
1064  {
1065  StdSeg.PhysDeriv(temparray + i*m_planes.num_elements(), tmp2 = outarray + i*m_planes.num_elements());
1066  }
1067 
1068  m_transposition->Transpose(outarray,out_d2,false,LibUtilities::eZtoXY);
1069 
1070  Vmath::Smul(nT_pts,2.0/m_lhom,out_d2,1,out_d2,1);
1071  }
1072  }
1073  }
1074  }
1075 
1078 
1079  {
1080  int nT_pts = inarray.num_elements(); //number of total points = n. of Fourier points * n. of points per plane (nT_pts)
1081  int nP_pts = nT_pts/m_planes.num_elements(); //number of points per plane = n of Fourier transform required (nP_pts)
1082 
1083  int dir= (int)edir;
1084 
1085  Array<OneD, NekDouble> temparray(nT_pts);
1086  Array<OneD, NekDouble> outarray(nT_pts);
1089 
1090  if (dir < 2)
1091  {
1092  for(int i=0; i < m_planes.num_elements(); i++)
1093  {
1094  m_planes[i]->PhysDeriv(edir, inarray + i*nP_pts ,tmp2 = out_d + i*nP_pts);
1095  }
1096  }
1097  else
1098  {
1101  {
1102  if(m_WaveSpace)
1103  {
1104  temparray = inarray;
1105  }
1106  else
1107  {
1108  HomogeneousFwdTrans(inarray,temparray);
1109  }
1110 
1111  NekDouble sign = -1.0;
1112  NekDouble beta;
1113 
1114  //HalfMode
1116  {
1117  beta = 2*sign*M_PI*(m_transposition->GetK(0))/m_lhom;
1118 
1119  Vmath::Smul(nP_pts,beta,temparray,1,outarray,1);
1120  }
1121  else if(m_homogeneousBasis->GetBasisType() == LibUtilities::eFourierHalfModeIm)
1122  {
1123  beta = -2*sign*M_PI*(m_transposition->GetK(0))/m_lhom;
1124 
1125  Vmath::Smul(nP_pts,beta,temparray,1,outarray,1);
1126  }
1127  //Fully complex
1128  else
1129  {
1130  for(int i = 0; i < m_planes.num_elements(); i++)
1131  {
1132  beta = -sign*2*M_PI*(m_transposition->GetK(i))/m_lhom;
1133 
1134  Vmath::Smul(nP_pts,beta,tmp1 = temparray + i*nP_pts,1,tmp2 = outarray + (i-int(sign))*nP_pts,1);
1135 
1136  sign = -1.0*sign;
1137  }
1138  }
1139  if(m_WaveSpace)
1140  {
1141  out_d = outarray;
1142  }
1143  else
1144  {
1145  HomogeneousBwdTrans(outarray,out_d);
1146  }
1147  }
1148  else
1149  {
1150  if(!m_session->DefinesSolverInfo("HomoStrip"))
1151  {
1152  ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,
1153  "Parallelisation in the homogeneous direction "
1154  "implemented just for Fourier basis");
1155  }
1156  else
1157  {
1158  ASSERTL0(m_StripZcomm->GetSize() == 1,
1159  "Parallelisation in the homogeneous direction "
1160  "implemented just for Fourier basis");
1161  }
1162 
1163  if(m_WaveSpace)
1164  {
1165  ASSERTL0(false,"Semi-phyisical time-stepping not implemented yet for non-Fourier basis");
1166  }
1167  else
1168  {
1169  StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
1170 
1171  m_transposition->Transpose(inarray,temparray,false,LibUtilities::eXYtoZ);
1172 
1173  for(int i = 0; i < nP_pts; i++)
1174  {
1175  StdSeg.PhysDeriv(temparray + i*m_planes.num_elements(), tmp2 = outarray + i*m_planes.num_elements());
1176  }
1177 
1178  m_transposition->Transpose(outarray,out_d,false,LibUtilities::eZtoXY);
1179 
1180  Vmath::Smul(nT_pts,2.0/m_lhom,out_d,1,out_d,1);
1181  }
1182  }
1183  }
1184  }
1185 
1187  Array<OneD, NekDouble> &out_d0,
1188  Array<OneD, NekDouble> &out_d1,
1189  Array<OneD, NekDouble> &out_d2)
1190 
1191  {
1192  v_PhysDeriv(inarray,out_d0,out_d1,out_d2);
1193  }
1194 
1196  const Array<OneD, const NekDouble> &inarray,
1197  Array<OneD, NekDouble> &out_d)
1198  {
1199  v_PhysDeriv(edir,inarray,out_d);
1200  }
1201 
1203  {
1204  return m_transposition;
1205  }
1206 
1208  {
1209  return m_lhom;
1210  }
1211 
1213  {
1214  return m_transposition->GetPlanesIDs();
1215  }
1216  } //end of namespace
1217 } //end of namespace
1218 
1219 
1220 /**
1221 * $Log: v $
1222 *
1223 **/
1224 
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:188
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
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:22
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:118
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:956
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:939
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:988
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
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:64
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:199
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:977
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:991
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:910
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_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:907
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:240
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.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:169
LibUtilities::NektarFFTSharedPtr m_FFT