Nektar++
LinearisedAdvection.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: LinearisedAdvection.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: Evaluation of the linearised advective term
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
36 #include <StdRegions/StdSegExp.h>
37 #include <boost/format.hpp>
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 string LinearisedAdvection::className =
45  "Linearised", LinearisedAdvection::create,
46  "Linearised Non-Conservative");
47 
48 /**
49  * Constructor. Creates ...
50  *
51  * \param
52  * \param
53  */
54 
55 LinearisedAdvection::LinearisedAdvection() : Advection()
56 {
57 }
58 
62 {
63  Advection::v_InitObject(pSession, pFields);
64 
65  m_session = pSession;
68  m_session, pFields[0]->GetGraph());
69  m_spacedim = pFields[0]->GetGraph()->GetSpaceDimension();
70  m_expdim = pFields[0]->GetGraph()->GetMeshDimension();
71 
72  // Setting parameters for homogeneous problems
73  m_HomoDirec = 0;
74  m_useFFT = false;
76  m_singleMode = false;
77  m_halfMode = false;
78  m_multipleModes = false;
79 
80  if (m_session->DefinesSolverInfo("HOMOGENEOUS"))
81  {
82  std::string HomoStr = m_session->GetSolverInfo("HOMOGENEOUS");
83  m_spacedim = 3;
84 
85  if ((HomoStr == "HOMOGENEOUS1D") || (HomoStr == "Homogeneous1D") ||
86  (HomoStr == "1D") || (HomoStr == "Homo1D"))
87  {
89  m_LhomZ = m_session->GetParameter("LZ");
90  m_HomoDirec = 1;
91 
92  ASSERTL0(m_session->DefinesSolverInfo("ModeType"),
93  "Need to specify ModeType as HalfMode,SingleMode or "
94  "MultipleModes");
95 
96  m_session->MatchSolverInfo("ModeType", "SingleMode", m_singleMode,
97  false);
98  m_session->MatchSolverInfo("ModeType", "HalfMode", m_halfMode,
99  false);
100  m_session->MatchSolverInfo("ModeType", "MultipleModes",
101  m_multipleModes, false);
102 
103  if (m_singleMode)
104  {
105  m_npointsZ = 2;
106  }
107  else if (m_halfMode)
108  {
109  m_npointsZ = 1;
110  }
111  else if (m_multipleModes)
112  {
113  m_npointsZ = m_session->GetParameter("HomModesZ");
114  }
115  }
116 
117  if ((HomoStr == "HOMOGENEOUS2D") || (HomoStr == "Homogeneous2D") ||
118  (HomoStr == "2D") || (HomoStr == "Homo2D"))
119  {
121  m_session->LoadParameter("HomModesY", m_npointsY);
122  m_session->LoadParameter("LY", m_LhomY);
123  m_session->LoadParameter("HomModesZ", m_npointsZ);
124  m_session->LoadParameter("LZ", m_LhomZ);
125  m_HomoDirec = 2;
126  }
127 
128  if ((HomoStr == "HOMOGENEOUS3D") || (HomoStr == "Homogeneous3D") ||
129  (HomoStr == "3D") || (HomoStr == "Homo3D"))
130  {
132  m_session->LoadParameter("HomModesX", m_npointsX);
133  m_session->LoadParameter("LX", m_LhomX);
134  m_session->LoadParameter("HomModesY", m_npointsY);
135  m_session->LoadParameter("LY", m_LhomY);
136  m_session->LoadParameter("HomModesZ", m_npointsZ);
137  m_session->LoadParameter("LZ", m_LhomZ);
138  m_HomoDirec = 3;
139  }
140 
141  if (m_session->DefinesSolverInfo("USEFFT"))
142  {
143  m_useFFT = true;
144  }
145  }
146  else
147  {
148  m_npointsZ = 1; // set to default value so can use to identify 2d or 3D
149  // (homogeneous) expansions
150  }
151 
152  size_t nvar = m_session->GetVariables().size();
154  for (size_t i = 0; i < nvar; ++i)
155  {
156  m_baseflow[i] = Array<OneD, NekDouble>(pFields[i]->GetTotPoints(), 0.0);
157  }
158 
159  size_t nBaseDerivs = (m_halfMode || m_singleMode) ? 2 : m_spacedim;
160  m_gradBase = Array<OneD, Array<OneD, NekDouble>>(nvar * nBaseDerivs);
161  for (size_t i = 0; i < nvar; ++i)
162  {
163  for (size_t j = 0; j < nBaseDerivs; ++j)
164  {
165  m_gradBase[i * nBaseDerivs + j] =
166  Array<OneD, NekDouble>(pFields[i]->GetTotPoints(), 0.0);
167  }
168  }
169 
170  ASSERTL0(m_session->DefinesFunction("BaseFlow"),
171  "Base flow must be defined for linearised forms.");
172  string file = m_session->GetFunctionFilename("BaseFlow", 0);
173 
174  // Periodic base flows
175  if (m_session->DefinesParameter("N_slices"))
176  {
177  m_session->LoadParameter("N_slices", m_slices);
178  if (m_slices > 1)
179  {
180  ASSERTL0(m_session->GetFunctionType("BaseFlow", 0) ==
182  "Base flow should be a sequence of files.");
183  m_session->LoadParameter("BaseFlow_interporder", m_interporder, 0);
185  "BaseFlow_interporder should be smaller than or equal to "
186  "N_slices.");
188  m_session->LoadParameter("N_start", m_start, 0);
189  m_session->LoadParameter("N_skip", m_skip, 1);
190  DFT(file, pFields);
191  }
192  else
193  {
194  ASSERTL0(false, "Number of slices must be a positive number "
195  "greater than 1");
196  }
197  }
198  // Steady base-flow
199  else
200  {
201  m_slices = 1;
202 
203  // BaseFlow from file
204  if (m_session->GetFunctionType("BaseFlow", m_session->GetVariable(0)) ==
206  {
207  ImportFldBase(file, pFields, 1);
208  }
209  // analytic base flow
210  else
211  {
212  size_t nq = pFields[0]->GetNpoints();
213  Array<OneD, NekDouble> x0(nq);
214  Array<OneD, NekDouble> x1(nq);
215  Array<OneD, NekDouble> x2(nq);
216 
217  // get the coordinates (assuming all fields have the same
218  // discretisation)
219  pFields[0]->GetCoords(x0, x1, x2);
220 
221  for (size_t i = 0; i < pFields.size(); i++)
222  {
224  m_session->GetFunction("BaseFlow", i);
225 
226  ifunc->Evaluate(x0, x1, x2, m_baseflow[i]);
227  }
228  }
229  }
230 
231  for (size_t i = 0; i < nvar; ++i)
232  {
233  UpdateGradBase(i, pFields[i]);
234  }
235 
236  if (m_session->DefinesParameter("period"))
237  {
238  m_period = m_session->GetParameter("period");
239  }
240  else
241  {
242  m_period =
243  (m_session->GetParameter("TimeStep") * m_slices) / (m_slices - 1.);
244  }
245  if (m_session->GetComm()->GetRank() == 0)
246  {
247  cout << "baseflow info : interpolation order " << m_interporder
248  << ", period " << m_period << ", periodicity ";
249  if (m_isperiodic)
250  {
251  cout << "yes\n";
252  }
253  else
254  {
255  cout << "no\n";
256  }
257  cout << "baseflow info : files from " << m_start << " to "
258  << (m_start + (m_slices - 1) * m_skip) << " (skip " << m_skip
259  << ") with " << (m_slices - (m_interporder > 1))
260  << " time intervals" << endl;
261  }
262 }
263 
265 {
266 }
267 
268 // Advection function
269 
271  const int nConvectiveFields,
273  const Array<OneD, Array<OneD, NekDouble>> &advVel,
274  const Array<OneD, Array<OneD, NekDouble>> &inarray,
275  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time,
276  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
277  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
278 {
279  boost::ignore_unused(pFwd, pBwd);
280  ASSERTL1(nConvectiveFields == inarray.size(),
281  "Number of convective fields and Inarray are not compatible");
282 
283  size_t nPointsTot = fields[0]->GetNpoints();
284  size_t ndim = advVel.size();
285  size_t nBaseDerivs = (m_halfMode || m_singleMode) ? 2 : m_spacedim;
286  size_t nDerivs = (m_halfMode) ? 2 : m_spacedim;
287 
288  Array<OneD, Array<OneD, NekDouble>> velocity(ndim);
289  for (size_t i = 0; i < ndim; ++i)
290  {
291  if (fields[i]->GetWaveSpace() && !m_singleMode && !m_halfMode)
292  {
293  velocity[i] = Array<OneD, NekDouble>(nPointsTot, 0.0);
294  fields[i]->HomogeneousBwdTrans(nPointsTot, advVel[i], velocity[i]);
295  }
296  else
297  {
298  velocity[i] = advVel[i];
299  }
300  }
301 
302  Array<OneD, Array<OneD, NekDouble>> grad(nDerivs);
303  for (size_t i = 0; i < nDerivs; ++i)
304  {
305  grad[i] = Array<OneD, NekDouble>(nPointsTot);
306  }
307 
308  // Evaluation of the base flow for periodic cases
309  if (m_slices > 1)
310  {
311  for (size_t i = 0; i < ndim; ++i)
312  {
313  UpdateBase(m_interp[i], m_baseflow[i], time);
314  UpdateGradBase(i, fields[i]);
315  }
316  }
317 
318  // Evaluate the linearised advection term
319  for (size_t i = 0; i < (size_t)nConvectiveFields; ++i)
320  {
321  // Calculate gradient
322  switch (nDerivs)
323  {
324  case 1:
325  {
326  fields[i]->PhysDeriv(inarray[i], grad[0]);
327  }
328  break;
329  case 2:
330  {
331  fields[i]->PhysDeriv(inarray[i], grad[0], grad[1]);
332  }
333  break;
334  case 3:
335  {
336  fields[i]->PhysDeriv(inarray[i], grad[0], grad[1], grad[2]);
337  if (m_multipleModes)
338  {
339  // transform gradients into physical Fourier space
340  fields[i]->HomogeneousBwdTrans(nPointsTot, grad[0],
341  grad[0]);
342  fields[i]->HomogeneousBwdTrans(nPointsTot, grad[1],
343  grad[1]);
344  fields[i]->HomogeneousBwdTrans(nPointsTot, grad[2],
345  grad[2]);
346  }
347  }
348  break;
349  }
350 
351  // Calculate U_j du'_i/dx_j
352  Vmath::Vmul(nPointsTot, grad[0], 1, m_baseflow[0], 1, outarray[i], 1);
353  for (size_t j = 1; j < nDerivs; ++j)
354  {
355  Vmath::Vvtvp(nPointsTot, grad[j], 1, m_baseflow[j], 1, outarray[i],
356  1, outarray[i], 1);
357  }
358 
359  // Add u'_j dU_i/dx_j
360  size_t lim = (m_halfMode || m_singleMode) ? 2 : ndim;
361  if (m_halfMode && i == 2)
362  {
363  lim = 0;
364  }
365  for (size_t j = 0; j < lim; ++j)
366  {
367  Vmath::Vvtvp(nPointsTot, m_gradBase[i * nBaseDerivs + j], 1,
368  velocity[j], 1, outarray[i], 1, outarray[i], 1);
369  }
370 
371  if (m_multipleModes)
372  {
373  fields[i]->HomogeneousFwdTrans(nPointsTot, outarray[i],
374  outarray[i]);
375  }
376  Vmath::Neg(nPointsTot, outarray[i], 1);
377  }
378 }
379 
381  const Array<OneD, Array<OneD, NekDouble>> &inarray,
383 {
384  if (m_session->GetSolverInfo("EqType") == "UnsteadyNavierStokes")
385  {
386  // The SFD method is only applied to the velocity variables in
387  // incompressible
388  ASSERTL1(inarray.size() == (m_baseflow.size() - 1),
389  "Number of base flow variables does not match what is "
390  "expected.");
391  }
392  else
393  {
394  ASSERTL1(
395  inarray.size() == (m_baseflow.size()),
396  "Number of base flow variables does not match what is expected.");
397  }
398 
399  size_t npts = inarray[0].size();
400 
401  for (size_t i = 0; i < inarray.size(); ++i)
402  {
403  ASSERTL1(npts == m_baseflow[i].size(),
404  "Size of base flow array does not match expected.");
405  Vmath::Vcopy(npts, inarray[i], 1, m_baseflow[i], 1);
406  UpdateGradBase(i, fields[i]);
407  }
408 }
409 
410 /**
411  * Import field from infile and load into \a m_fields. This routine will
412  * also perform a \a BwdTrans to ensure data is in both the physical and
413  * coefficient storage.
414  * @param pInFile Filename to read.
415  * @param pFields Array of expansion lists
416  */
418  std::string pInfile, Array<OneD, MultiRegions::ExpListSharedPtr> &pFields,
419  int pSlice)
420 {
421  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
422  std::vector<std::vector<NekDouble>> FieldData;
423 
424  size_t nqtot = m_baseflow[0].size();
425  Array<OneD, NekDouble> tmp_coeff(pFields[0]->GetNcoeffs(), 0.0);
426 
427  size_t numexp = pFields[0]->GetExpSize();
428  Array<OneD, int> ElementGIDs(numexp);
429 
430  // Define list of global element ids
431  for (size_t i = 0; i < numexp; ++i)
432  {
433  ElementGIDs[i] = pFields[0]->GetExp(i)->GetGeom()->GetGlobalID();
434  }
435 
436  // Get Homogeneous
439  fld->Import(pInfile, FieldDef, FieldData,
441 
442  size_t nSessionVar = m_session->GetVariables().size();
443  size_t nSessionConvVar = nSessionVar - 1;
444  size_t nFileVar = FieldDef[0]->m_fields.size();
445 
446  std::unordered_map<int, int> zIdToPlane;
447  if (m_singleMode || m_halfMode)
448  {
449  zIdToPlane[0] = 0;
450  }
451 
452  for (size_t j = 0; j < nFileVar; ++j)
453  {
454  size_t k = 0;
455  for (; k < nSessionConvVar; ++k)
456  {
457  if (m_session->GetVariable(k) == FieldDef[0]->m_fields[j])
458  {
459  break;
460  }
461  }
462  if (k == nSessionConvVar)
463  {
464  continue;
465  }
466  for (size_t i = 0; i < FieldDef.size(); ++i)
467  {
468  pFields[j]->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
469  FieldDef[i]->m_fields[j], tmp_coeff,
470  zIdToPlane);
471  }
472 
473  if (m_singleMode || m_halfMode)
474  {
475  pFields[j]->GetPlane(0)->BwdTrans(tmp_coeff, m_baseflow[k]);
476 
477  if (m_singleMode)
478  {
479  // copy the bwd trans into the second plane for single
480  // Mode Analysis
481  int ncplane = (pFields[0]->GetNpoints()) / m_npointsZ;
482  Vmath::Vcopy(ncplane, &m_baseflow[k][0], 1,
483  &m_baseflow[k][ncplane], 1);
484  }
485  }
486  else // fully 3D base flow - put in physical space.
487  {
488  bool oldwavespace = pFields[j]->GetWaveSpace();
489  pFields[j]->SetWaveSpace(false);
490  pFields[j]->BwdTrans(tmp_coeff, m_baseflow[k]);
491  pFields[j]->SetWaveSpace(oldwavespace);
492  }
493  }
494 
495  // If time-periodic, put loaded data into the slice storage.
496  if (m_slices > 1)
497  {
498  for (size_t i = 0; i < nSessionConvVar; ++i)
499  {
500  Vmath::Vcopy(nqtot, &m_baseflow[i][0], 1,
501  &m_interp[i][pSlice * nqtot], 1);
502  }
503  }
504 }
505 
507  const Array<OneD, const NekDouble> &inarray,
508  Array<OneD, NekDouble> &outarray, const NekDouble time)
509 {
510  int npoints = m_baseflow[0].size();
511  if (m_isperiodic)
512  {
513  NekDouble BetaT = 2 * M_PI * fmod(time, m_period) / m_period;
514  NekDouble phase;
515  Array<OneD, NekDouble> auxiliary(npoints);
516 
517  Vmath::Vcopy(npoints, &inarray[0], 1, &outarray[0], 1);
518  Vmath::Svtvp(npoints, cos(0.5 * m_slices * BetaT), &inarray[npoints], 1,
519  &outarray[0], 1, &outarray[0], 1);
520 
521  for (int i = 2; i < m_slices; i += 2)
522  {
523  phase = (i >> 1) * BetaT;
524 
525  Vmath::Svtvp(npoints, cos(phase), &inarray[i * npoints], 1,
526  &outarray[0], 1, &outarray[0], 1);
527  Vmath::Svtvp(npoints, -sin(phase), &inarray[(i + 1) * npoints], 1,
528  &outarray[0], 1, &outarray[0], 1);
529  }
530  }
531  else
532  {
533  NekDouble x = time;
534  x = x / m_period * (m_slices - 1);
535  int ix = x;
536  if (ix < 0)
537  {
538  ix = 0;
539  }
540  if (ix > m_slices - 2)
541  {
542  ix = m_slices - 2;
543  }
544  int padleft = m_interporder / 2 - 1;
545  if (padleft > ix)
546  {
547  padleft = ix;
548  }
549  int padright = m_interporder - 1 - padleft;
550  if (padright > m_slices - 1 - ix)
551  {
552  padright = m_slices - 1 - ix;
553  }
554  padleft = m_interporder - 1 - padright;
556  for (int i = 0; i < m_interporder; ++i)
557  {
558  for (int j = 0; j < m_interporder; ++j)
559  {
560  if (i != j)
561  {
562  coeff[i] *= (x - ix + padleft - (NekDouble)j) /
563  ((NekDouble)i - (NekDouble)j);
564  }
565  }
566  }
567  Vmath::Zero(npoints, &outarray[0], 1);
568  for (int i = ix - padleft; i < ix + padright + 1; ++i)
569  {
570  Vmath::Svtvp(npoints, coeff[i - ix + padleft],
571  &inarray[i * npoints], 1, &outarray[0], 1,
572  &outarray[0], 1);
573  }
574  }
575 }
576 
578  const int var, const MultiRegions::ExpListSharedPtr &field)
579 {
580  int nBaseDerivs = (m_halfMode || m_singleMode) ? 2 : m_spacedim;
581  switch (m_spacedim)
582  {
583  case 1: // 1D
584  {
585  field->PhysDeriv(m_baseflow[var],
586  m_gradBase[var * nBaseDerivs + 0]);
587  }
588  break;
589  case 2: // 2D
590  {
591  field->PhysDeriv(m_baseflow[var], m_gradBase[var * nBaseDerivs + 0],
592  m_gradBase[var * nBaseDerivs + 1]);
593  }
594  break;
595  case 3:
596  {
597  if (m_halfMode) // can assume W = 0 and d/dz == 0
598  {
599  if (var < 2)
600  {
601  field->PhysDeriv(m_baseflow[var],
602  m_gradBase[var * nBaseDerivs + 0],
603  m_gradBase[var * nBaseDerivs + 1]);
604  }
605  }
606  else if (m_singleMode) // single mode where d/dz = 0
607  {
608  field->PhysDeriv(m_baseflow[var],
609  m_gradBase[var * nBaseDerivs + 0],
610  m_gradBase[var * nBaseDerivs + 1]);
611  }
612  else
613  {
614  // Differentiate base flow in physical space
615  bool oldwavespace = field->GetWaveSpace();
616  field->SetWaveSpace(false);
617  field->PhysDeriv(m_baseflow[var],
618  m_gradBase[var * nBaseDerivs + 0],
619  m_gradBase[var * nBaseDerivs + 1],
620  m_gradBase[var * nBaseDerivs + 2]);
621  field->SetWaveSpace(oldwavespace);
622  }
623  }
624  break;
625  }
626 }
627 
629  FloquetMatType mattype, bool UseContCoeffs) const
630 {
631  boost::ignore_unused(mattype, UseContCoeffs);
632 
633  DNekMatSharedPtr loc_mat;
634  DNekBlkMatSharedPtr BlkMatrix;
635  size_t n_exp = 0;
636 
637  n_exp = m_baseflow[0].size(); // will operatore on m_phys
638 
639  Array<OneD, unsigned int> nrows(n_exp);
640  Array<OneD, unsigned int> ncols(n_exp);
641 
642  nrows = Array<OneD, unsigned int>(n_exp, m_slices);
643  ncols = Array<OneD, unsigned int>(n_exp, m_slices);
644 
645  MatrixStorage blkmatStorage = eDIAGONAL;
646  BlkMatrix = MemoryManager<DNekBlkMat>::AllocateSharedPtr(nrows, ncols,
647  blkmatStorage);
648 
652  StdRegions::StdSegExp StdSeg(BK);
653 
655  StdSeg.DetShapeType(), StdSeg);
656 
657  loc_mat = StdSeg.GetStdMatrix(matkey);
658 
659  // set up array of block matrices.
660  for (size_t i = 0; i < n_exp; ++i)
661  {
662  BlkMatrix->SetBlock(i, i, loc_mat);
663  }
664 
665  return BlkMatrix;
666 }
667 
668 // Discrete Fourier Transform for Floquet analysis
670  const string file, Array<OneD, MultiRegions::ExpListSharedPtr> &pFields)
671 {
672  size_t ConvectedFields = m_baseflow.size() - 1;
673  size_t npoints = m_baseflow[0].size();
674  m_interp = Array<OneD, Array<OneD, NekDouble>>(ConvectedFields);
675 
676  for (size_t i = 0; i < ConvectedFields; ++i)
677  {
678  m_interp[i] = Array<OneD, NekDouble>(npoints * m_slices, 0.0);
679  }
680 
681  // Import the slides into auxiliary vector
682  // The base flow should be stored in the form "filename_%d.ext"
683  // A subdirectory can also be included, such as "dir/filename_%d.ext"
684  size_t found = file.find("%d");
685  ASSERTL0(found != string::npos &&
686  file.find("%d", found + 1) == string::npos,
687  "Since N_slices is specified, the filename provided for function "
688  "'BaseFlow' must include exactly one instance of the format "
689  "specifier '%d', to index the time-slices.");
690  size_t nstart = m_start;
691  for (size_t i = nstart; i < nstart + m_slices * m_skip; i += m_skip)
692  {
693  boost::format filename(file);
694  filename % i;
695  ImportFldBase(filename.str(), pFields, (i - nstart) / m_skip);
696  if (m_session->GetComm()->GetRank() == 0)
697  {
698  cout << "read base flow file " << filename.str() << endl;
699  }
700  }
701  if (!m_isperiodic)
702  {
703  return;
704  }
705 
706  // Discrete Fourier Transform of the fields
707  for (size_t k = 0; k < ConvectedFields; ++k)
708  {
709 #ifdef NEKTAR_USING_FFTW
710 
711  // Discrete Fourier Transform using FFTW
712  Array<OneD, NekDouble> fft_in(npoints * m_slices);
713  Array<OneD, NekDouble> fft_out(npoints * m_slices);
714 
717 
718  // Shuffle the data
719  for (int j = 0; j < m_slices; ++j)
720  {
721  Vmath::Vcopy(npoints, &m_interp[k][j * npoints], 1, &(fft_in[j]),
722  m_slices);
723  }
724 
726  m_slices);
727 
728  // FFT Transform
729  for (size_t i = 0; i < npoints; i++)
730  {
731  m_FFT->FFTFwdTrans(m_tmpIN = fft_in + i * m_slices,
732  m_tmpOUT = fft_out + i * m_slices);
733  }
734 
735  // Reshuffle data
736  for (int s = 0; s < m_slices; ++s)
737  {
738  Vmath::Vcopy(npoints, &fft_out[s], m_slices,
739  &m_interp[k][s * npoints], 1);
740  }
741 
742  Vmath::Zero(fft_in.size(), &fft_in[0], 1);
743  Vmath::Zero(fft_out.size(), &fft_out[0], 1);
744 #else
745  // Discrete Fourier Transform using MVM
746  DNekBlkMatSharedPtr blkmat;
748 
749  int nrows = blkmat->GetRows();
750  int ncols = blkmat->GetColumns();
751 
752  Array<OneD, NekDouble> sortedinarray(ncols);
753  Array<OneD, NekDouble> sortedoutarray(nrows);
754 
755  // Shuffle the data
756  for (int j = 0; j < m_slices; ++j)
757  {
758  Vmath::Vcopy(npoints, &m_interp[k][j * npoints], 1,
759  &(sortedinarray[j]), m_slices);
760  }
761 
762  // Create NekVectors from the given data arrays
763  NekVector<NekDouble> in(ncols, sortedinarray, eWrapper);
764  NekVector<NekDouble> out(nrows, sortedoutarray, eWrapper);
765 
766  // Perform matrix-vector multiply.
767  out = (*blkmat) * in;
768 
769  // Reshuffle data
770  for (int s = 0; s < m_slices; ++s)
771  {
772  Vmath::Vcopy(npoints, &sortedoutarray[s], m_slices,
773  &m_interp[k][s * npoints], 1);
774  }
775 
776  for (size_t r = 0; r < sortedinarray.size(); ++r)
777  {
778  sortedinarray[0] = 0;
779  sortedoutarray[0] = 0;
780  }
781 
782 #endif
783  }
784 }
785 
786 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
Describes the specification for a Basis.
Definition: Basis.h:50
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition: FieldIO.cpp:226
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
Defines a specification for a set of points.
Definition: Points.h:59
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
int m_npointsZ
number of points in Z direction (if homogeneous)
NekDouble m_LhomZ
physical length in Z direction (if homogeneous)
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
Initialises the advection object.
enum HomogeneousType m_HomogeneousType
void UpdateBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble time)
NekDouble m_period
period length
int m_npointsX
number of points in X direction (if homogeneous)
int m_npointsY
number of points in Y direction (if homogeneous)
int m_HomoDirec
number of homogenous directions
bool m_useFFT
flag to determine if use or not the FFT for transformations
LibUtilities::SessionReaderSharedPtr m_session
Array< OneD, NekDouble > m_tmpIN
virtual void v_Advect(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &advVel, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray) override
Advects a vector field.
void UpdateGradBase(const int var, const MultiRegions::ExpListSharedPtr &field)
DNekBlkMatSharedPtr GetFloquetBlockMatrix(FloquetMatType mattype, bool UseContCoeffs=false) const
bool m_singleMode
flag to determine if use single mode or not
Array< OneD, Array< OneD, NekDouble > > m_baseflow
Storage for base flow.
Array< OneD, NekDouble > m_tmpOUT
Array< OneD, Array< OneD, NekDouble > > m_gradBase
bool m_multipleModes
flag to determine if use multiple mode or not
void DFT(const std::string file, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
LibUtilities::NektarFFTSharedPtr m_FFT
auxiliary variables
bool m_halfMode
flag to determine if use half mode or not
NekDouble m_LhomY
physical length in Y direction (if homogeneous)
void ImportFldBase(std::string pInfile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, int slice)
Import Base flow.
NekDouble m_LhomX
physical length in X direction (if homogeneous)
virtual void v_SetBaseFlow(const Array< OneD, Array< OneD, NekDouble >> &inarray, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields) override
Overrides the base flow used during linearised advection.
Array< OneD, Array< OneD, NekDouble > > m_interp
interpolation vector
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
An abstract base class encapsulating the concept of advection of a vector field.
Definition: Advection.h:70
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
Definition: Advection.cpp:353
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:609
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
Class representing a segment element in reference space All interface of this class sits in StdExpans...
Definition: StdSegExp.h:48
NektarFFTFactory & GetNektarFFTFactory()
Definition: NektarFFT.cpp:69
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:327
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:129
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:53
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:76
@ eFourier
Fourier Expansion .
Definition: BasisType.h:57
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:77
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
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:209
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:622
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:574
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255