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