Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // 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: Evaluation of the linearised advective term
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 #include <StdRegions/StdSegExp.h>
38 
48 
49 
50 namespace Nektar
51 {
54  "Linearised",
56 
57 /**
58  * Constructor. Creates ...
59  *
60  * \param
61  * \param
62  */
63 
65  Advection()
66 {
67 }
68 
69 
73 {
74  Advection::v_InitObject(pSession, pFields);
75 
76  m_session = pSession;
78  ::AllocateSharedPtr(m_session, pFields[0]->GetGraph());
79  m_spacedim = pFields[0]->GetGraph()->GetSpaceDimension();
80  m_expdim = pFields[0]->GetGraph()->GetMeshDimension();
81 
82  //Setting parameters for homogeneous problems
83  m_HomoDirec = 0;
84  m_useFFT = false;
86  m_singleMode = false;
87  m_halfMode = false;
88  m_multipleModes = false;
89 
90  if(m_session->DefinesSolverInfo("HOMOGENEOUS"))
91  {
92  std::string HomoStr = m_session->GetSolverInfo("HOMOGENEOUS");
93  m_spacedim = 3;
94 
95  if((HomoStr == "HOMOGENEOUS1D")||(HomoStr == "Homogeneous1D")||
96  (HomoStr == "1D")||(HomoStr == "Homo1D"))
97  {
99  m_LhomZ = m_session->GetParameter("LZ");
100  m_HomoDirec = 1;
101 
102  ASSERTL0(m_session->DefinesSolverInfo("ModeType"),
103  "Need to specify ModeType as HalfMode,SingleMode or "
104  "MultipleModes");
105 
106  m_session->MatchSolverInfo("ModeType", "SingleMode",
107  m_singleMode, false);
108  m_session->MatchSolverInfo("ModeType", "HalfMode",
109  m_halfMode, false);
110  m_session->MatchSolverInfo("ModeType", "MultipleModes",
111  m_multipleModes, false);
112 
113  if(m_singleMode)
114  {
115  m_npointsZ = 2;
116  }
117  else if(m_halfMode)
118  {
119  m_npointsZ = 1;
120  }
121  else if(m_multipleModes)
122  {
123  m_npointsZ = m_session->GetParameter("HomModesZ");
124  }
125  }
126 
127  if((HomoStr == "HOMOGENEOUS2D")||(HomoStr == "Homogeneous2D")||
128  (HomoStr == "2D")||(HomoStr == "Homo2D"))
129  {
131  m_session->LoadParameter("HomModesY", m_npointsY);
132  m_session->LoadParameter("LY", m_LhomY);
133  m_session->LoadParameter("HomModesZ", m_npointsZ);
134  m_session->LoadParameter("LZ", m_LhomZ);
135  m_HomoDirec = 2;
136  }
137 
138  if((HomoStr == "HOMOGENEOUS3D")||(HomoStr == "Homogeneous3D")||
139  (HomoStr == "3D")||(HomoStr == "Homo3D"))
140  {
142  m_session->LoadParameter("HomModesX",m_npointsX);
143  m_session->LoadParameter("LX", m_LhomX );
144  m_session->LoadParameter("HomModesY",m_npointsY);
145  m_session->LoadParameter("LY", m_LhomY );
146  m_session->LoadParameter("HomModesZ",m_npointsZ);
147  m_session->LoadParameter("LZ", m_LhomZ );
148  m_HomoDirec = 3;
149  }
150 
151  if(m_session->DefinesSolverInfo("USEFFT"))
152  {
153  m_useFFT = true;
154  }
155  }
156  else
157  {
158  m_npointsZ = 1; // set to default value so can use to identify 2d or 3D (homogeneous) expansions
159  }
160 
161  int nvar = m_session->GetVariables().size();
163  for (int i = 0; i < nvar; ++i)
164  {
165  m_baseflow[i] = Array<OneD, NekDouble>(pFields[i]->GetTotPoints(), 0.0);
166  }
167 
168  ASSERTL0(m_session->DefinesFunction("BaseFlow"),
169  "Base flow must be defined for linearised forms.");
170  string file = m_session->GetFunctionFilename("BaseFlow", 0);
171 
172 
173  //Periodic base flows
174  if(m_session->DefinesParameter("N_slices"))
175  {
176  m_session->LoadParameter("N_slices",m_slices);
177  if(m_slices>1)
178  {
179  DFT(file,pFields,m_slices);
180  }
181  else
182  {
183  ASSERTL0(false, "Number of slices must be a positive number "
184  "greater than 1");
185  }
186  }
187  //Steady base-flow
188  else
189  {
190  m_slices=1;
191 
192  //BaseFlow from file
193  if (m_session->GetFunctionType("BaseFlow", m_session->GetVariable(0))
195  {
196  ImportFldBase(file,pFields,1);
197 
198  }
199  //analytic base flow
200  else
201  {
202  int nq = pFields[0]->GetNpoints();
203  Array<OneD,NekDouble> x0(nq);
204  Array<OneD,NekDouble> x1(nq);
205  Array<OneD,NekDouble> x2(nq);
206 
207  // get the coordinates (assuming all fields have the same
208  // discretisation)
209  pFields[0]->GetCoords(x0,x1,x2);
210 
211  for(unsigned int i = 0 ; i < pFields.num_elements(); i++)
212  {
214  = m_session->GetFunction("BaseFlow", i);
215 
216  ifunc->Evaluate(x0,x1,x2,m_baseflow[i]);
217  }
218  }
219  }
220 
221  if(m_session->DefinesParameter("period"))
222  {
223  m_period=m_session->GetParameter("period");
224  }
225  else
226  {
227  m_period=(m_session->GetParameter("TimeStep")*m_slices)/(m_slices-1.);
228  }
229 
230 }
231 
233 {
234 }
235 
236 
237 //Advection function
238 
240  const int nConvectiveFields,
242  const Array<OneD, Array<OneD, NekDouble> > &advVel,
243  const Array<OneD, Array<OneD, NekDouble> > &inarray,
244  Array<OneD, Array<OneD, NekDouble> > &outarray,
245  const NekDouble &time)
246 {
247  ASSERTL1(nConvectiveFields == inarray.num_elements(),
248  "Number of convective fields and Inarray are not compatible");
249 
250  int nPointsTot = fields[0]->GetNpoints();
251  int ndim = advVel.num_elements();
252 
254  = Array<OneD, NekDouble> (nPointsTot*nConvectiveFields);
255 
256  Array<OneD, NekDouble> grad0,grad1,grad2;
257 
258  // Evaluation of the gradient of each component of the base flow
259  // \nabla U
260  Array<OneD, NekDouble> grad_base_u0,grad_base_u1,grad_base_u2;
261 
262  // \nabla V
263  Array<OneD, NekDouble> grad_base_v0,grad_base_v1,grad_base_v2;
264 
265  // \nabla W
266  Array<OneD, NekDouble> grad_base_w0,grad_base_w1,grad_base_w2;
267 
268  grad0 = Array<OneD, NekDouble> (nPointsTot);
269  grad_base_u0 = Array<OneD, NekDouble> (nPointsTot);
270  grad_base_v0 = Array<OneD, NekDouble> (nPointsTot);
271  grad_base_w0 = Array<OneD, NekDouble> (nPointsTot);
272 
273  // Evaluation of the base flow for periodic cases
274  if (m_slices > 1)
275  {
276  ASSERTL0(m_session->GetFunctionType("BaseFlow", 0)
278  "Base flow should be a sequence of files.");
279 
280  for (int i = 0; i < ndim; ++i)
281  {
283  time, m_period);
284  }
285  }
286 
287 
288  //Evaluate the linearised advection term
289  switch(ndim)
290  {
291  case 1: // 1D
292  ASSERTL0(false,"Not set up for 1D");
293  break;
294  case 2: //2D
295  {
296  grad1 = Array<OneD, NekDouble> (nPointsTot);
297  grad_base_u1 = Array<OneD, NekDouble> (nPointsTot);
298  grad_base_v1 = Array<OneD, NekDouble> (nPointsTot);
299 
300  //Derivates of the base flow
301  fields[0]-> PhysDeriv(m_baseflow[0], grad_base_u0, grad_base_u1);
302  fields[0]-> PhysDeriv(m_baseflow[1], grad_base_v0, grad_base_v1);
303 
304  //x-equation
305  fields[0]->PhysDeriv(inarray[0],grad0,grad1);
306  // Evaluate U du'/dx
307  Vmath::Vmul (nPointsTot,grad0,1,m_baseflow[0],1,outarray[0],1);
308  //Evaluate U du'/dx+ V du'/dy
309  Vmath::Vvtvp(nPointsTot,grad1,1,m_baseflow[1],1,outarray[0],1,
310  outarray[0],1);
311  //Evaluate (U du'/dx+ V du'/dy)+u' dU/dx
312  Vmath::Vvtvp(nPointsTot,grad_base_u0,1,advVel[0],1,outarray[0],1,
313  outarray[0],1);
314  //Evaluate (U du'/dx+ V du'/dy +u' dU/dx)+v' dU/dy
315  Vmath::Vvtvp(nPointsTot,grad_base_u1,1,advVel[1],1,outarray[0],1,
316  outarray[0],1);
317  Vmath::Neg(nPointsTot,outarray[0],1);
318 
319  //y-equation
320  fields[0]->PhysDeriv(inarray[1],grad0,grad1);
321  // Evaluate U dv'/dx
322  Vmath::Vmul (nPointsTot,grad0,1,m_baseflow[0],1,outarray[1],1);
323  //Evaluate U dv'/dx+ V dv'/dy
324  Vmath::Vvtvp(nPointsTot,grad1,1,m_baseflow[1],1,outarray[1],1,
325  outarray[1],1);
326  //Evaluate (U dv'/dx+ V dv'/dy)+u' dV/dx
327  Vmath::Vvtvp(nPointsTot,grad_base_v0,1,advVel[0],1,outarray[1],1,
328  outarray[1],1);
329  //Evaluate (U dv'/dx+ V dv'/dy +u' dv/dx)+v' dV/dy
330  Vmath::Vvtvp(nPointsTot,grad_base_v1,1,advVel[1],1,outarray[1],1,
331  outarray[1],1);
332  Vmath::Neg(nPointsTot,outarray[1],1);
333  }
334  break;
335  case 3: //3D
336  {
337  grad1 = Array<OneD, NekDouble> (nPointsTot);
338  grad2 = Array<OneD, NekDouble> (nPointsTot);
339 
340  grad_base_u1 = Array<OneD, NekDouble> (nPointsTot);
341  grad_base_v1 = Array<OneD, NekDouble> (nPointsTot);
342  grad_base_w1 = Array<OneD, NekDouble> (nPointsTot,0.0);
343 
344  grad_base_u2 = Array<OneD, NekDouble> (nPointsTot,0.0);
345  grad_base_v2 = Array<OneD, NekDouble> (nPointsTot,0.0);
346  grad_base_w2 = Array<OneD, NekDouble> (nPointsTot,0.0);
347 
348  // note base flow should be moved to GetBaseFlow method
349  if(m_halfMode) // can assume W = 0 and d/dz == 0
350  {
351  // note base flow should be moved to GetBaseFlow method
352  fields[0]->PhysDeriv(m_baseflow[0],
353  grad_base_u0, grad_base_u1);
354  fields[0]->PhysDeriv(m_baseflow[1],
355  grad_base_v0, grad_base_v1);
356  }
357  else
358  {
359  if(m_multipleModes)
360  {
361  // Differentiate base flow in physical space
362  bool oldwavespace = fields[0]->GetWaveSpace();
363  fields[0]->SetWaveSpace(false);
364  fields[0]->PhysDeriv(m_baseflow[0], grad_base_u0,
365  grad_base_u1, grad_base_u2);
366  fields[0]->PhysDeriv(m_baseflow[1], grad_base_v0,
367  grad_base_v1, grad_base_v2);
368  fields[0]->PhysDeriv(m_baseflow[2], grad_base_w0,
369  grad_base_w1, grad_base_w2);
370  fields[0]->SetWaveSpace(oldwavespace);
371 
372 
373  }
374  else // has to be single mode where d/dz = 0
375  {
376  fields[0]->PhysDeriv(m_baseflow[0], grad_base_u0,
377  grad_base_u1);
378  fields[0]->PhysDeriv(m_baseflow[1], grad_base_v0,
379  grad_base_v1);
380  fields[0]->PhysDeriv(m_baseflow[2], grad_base_w0,
381  grad_base_w1);
382  }
383  }
384 
385  //x-equation
386  //
387  // Could probably clean up independent looping if clean up
388  // base flow derivative evaluation
389  if(m_multipleModes)
390  {
391  fields[0]->PhysDeriv(inarray[0], grad0, grad1, grad2);
392  // transform gradients into physical fouier space
393  fields[0]->HomogeneousBwdTrans(grad0, grad0);
394  fields[0]->HomogeneousBwdTrans(grad1, grad1);
395  fields[0]->HomogeneousBwdTrans(grad2, grad2);
396  }
397  else
398  {
399  if(m_halfMode) //W = 0 so no need for d/dz
400  {
401  fields[0]->PhysDeriv(inarray[0], grad0, grad1);
402  }
403  else // Single mode
404  {
405  fields[0]->PhysDeriv(inarray[0], grad0, grad1, grad2);
406  }
407  }
408  //Evaluate: U du'/dx
409  Vmath::Vmul (nPointsTot, grad0, 1, m_baseflow[0], 1,
410  outarray[0], 1);
411  //Evaluate and add: V du'/dy
412  Vmath::Vvtvp(nPointsTot, grad1, 1, m_baseflow[1], 1,
413  outarray[0], 1, outarray[0], 1);
414  //Evaluate and add: u' dU/dx
415  Vmath::Vvtvp(nPointsTot,grad_base_u0,1,advVel[0],1,
416  outarray[0],1,outarray[0],1);
417  //Evaluate and add: v' dU/dy
418  Vmath::Vvtvp(nPointsTot,grad_base_u1,1,advVel[1],1,
419  outarray[0],1,outarray[0],1);
420 
421  if(!m_halfMode)
422  {
423  //Evaluate an add W du'/dz
424  Vmath::Vvtvp(nPointsTot, grad2, 1, m_baseflow[2],
425  1,outarray[0], 1, outarray[0], 1);
426  }
427 
428  if(m_multipleModes)
429  {
430  //Evaluate and add w' dU/dz
431  Vmath::Vvtvp(nPointsTot,grad_base_u2,1,
432  advVel[2],1,outarray[0],1,outarray[0],1);
433  fields[0]->HomogeneousFwdTrans(outarray[0],outarray[0]);
434  }
435  Vmath::Neg(nPointsTot,outarray[0],1);
436 
437  //y-equation------------
438  if(m_multipleModes)
439  {
440  fields[0]->PhysDeriv(inarray[1], grad0, grad1, grad2);
441  // transform gradients into physical fouier space
442  fields[0]->HomogeneousBwdTrans(grad0, grad0);
443  fields[0]->HomogeneousBwdTrans(grad1, grad1);
444  fields[0]->HomogeneousBwdTrans(grad2, grad2);
445  }
446  else
447  {
448  if(m_halfMode) //W = 0 so no need for d/dz
449  {
450  fields[0]->PhysDeriv(inarray[1], grad0, grad1);
451  }
452  else // Single mode
453  {
454  fields[0]->PhysDeriv(inarray[1], grad0, grad1, grad2);
455  }
456  }
457 
458  //Evaluate U dv'/dx
459  Vmath::Vmul (nPointsTot, grad0, 1, m_baseflow[0], 1,
460  outarray[1], 1);
461  //Evaluate V dv'/dy
462  Vmath::Vvtvp(nPointsTot, grad1, 1, m_baseflow[1], 1,
463  outarray[1], 1, outarray[1], 1);
464  //Evaluate u' dV/dx
465  Vmath::Vvtvp(nPointsTot,grad_base_v0,1,advVel[0],1
466  ,outarray[1],1,outarray[1],1);
467  //Evaluate v' dV/dy
468  Vmath::Vvtvp(nPointsTot,grad_base_v1,1,advVel[1],1,
469  outarray[1],1,outarray[1],1);
470 
471  if(!m_halfMode)
472  {
473  //Evaluate W du'/dz
474  Vmath::Vvtvp(nPointsTot, grad2, 1, m_baseflow[2], 1,
475  outarray[1], 1, outarray[1], 1);
476  }
477 
478  if(m_multipleModes)
479  {
480  //Evaluate w' dV/dz
481  Vmath::Vvtvp(nPointsTot,grad_base_v2,1,advVel[2],1,
482  outarray[1],1,outarray[1],1);
483  fields[0]->HomogeneousFwdTrans(outarray[1],outarray[1]);
484  }
485 
486  Vmath::Neg(nPointsTot,outarray[1],1);
487 
488  //z-equation ------------
489  if(m_multipleModes)
490  {
491  fields[0]->PhysDeriv(inarray[2], grad0, grad1, grad2);
492  // transform gradients into physical fouier space
493  fields[0]->HomogeneousBwdTrans(grad0, grad0);
494  fields[0]->HomogeneousBwdTrans(grad1, grad1);
495  fields[0]->HomogeneousBwdTrans(grad2, grad2);
496  }
497  else
498  {
499  if(m_halfMode) //W = 0 so no need for d/dz
500  {
501  fields[0]->PhysDeriv(inarray[2], grad0, grad1);
502  }
503  else // Single mode
504  {
505  fields[0]->PhysDeriv(inarray[2], grad0, grad1, grad2);
506  }
507  }
508 
509  //Evaluate U dw'/dx
510  Vmath::Vmul (nPointsTot, grad0, 1, m_baseflow[0], 1,
511  outarray[2], 1);
512  //Evaluate V dw'/dx
513  Vmath::Vvtvp(nPointsTot, grad1, 1, m_baseflow[1], 1,
514  outarray[2], 1, outarray[2], 1);
515 
516  if(!m_halfMode)// if halfmode W = 0
517  {
518  //Evaluate u' dW/dx
519  Vmath::Vvtvp(nPointsTot,grad_base_w0,1,advVel[0],1,
520  outarray[2],1,outarray[2],1);
521  //Evaluate v' dW/dy
522  Vmath::Vvtvp(nPointsTot,grad_base_w1,1,advVel[1],1,
523  outarray[2],1,outarray[2],1);
524 
525  //Evaluate W dw'/dz
526  Vmath::Vvtvp(nPointsTot, grad2, 1, m_baseflow[2], 1,
527  outarray[2], 1, outarray[2], 1);
528  }
529 
530  if(m_multipleModes)
531  {
532  //Evaluate w' dW/dz
533  Vmath::Vvtvp(nPointsTot,grad_base_w2,1,advVel[2],1,
534  outarray[2],1,outarray[2],1);
535  fields[0]->HomogeneousFwdTrans(outarray[2],outarray[2]);
536  }
537  Vmath::Neg(nPointsTot,outarray[2],1);
538  }
539  }
540 }
541 
543  const Array<OneD, Array<OneD, NekDouble> > &inarray)
544 {
545  if (m_session->GetSolverInfo("EqType") == "UnsteadyNavierStokes")
546  {
547  // The SFD method is only applied to the velocity variables in
548  // incompressible
549  ASSERTL1(inarray.num_elements() == (m_baseflow.num_elements() - 1),
550  "Number of base flow variables does not match what is "
551  "expected.");
552  }
553  else
554  {
555  ASSERTL1(inarray.num_elements() == (m_baseflow.num_elements()),
556  "Number of base flow variables does not match what is expected.");
557  }
558 
559  int npts = inarray[0].num_elements();
560 
561  for (int i = 0; i < inarray.num_elements(); ++i)
562  {
563  ASSERTL1(npts == m_baseflow[i].num_elements(),
564  "Size of base flow array does not match expected.");
565  Vmath::Vcopy(npts, inarray[i], 1, m_baseflow[i], 1);
566  }
567 }
568 
569 
570 /**
571  * Import field from infile and load into \a m_fields. This routine will
572  * also perform a \a BwdTrans to ensure data is in both the physical and
573  * coefficient storage.
574  * @param infile Filename to read.
575  */
576 void LinearisedAdvection::ImportFldBase(std::string pInfile,
578 {
579  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
580  std::vector<std::vector<NekDouble> > FieldData;
581 
582  int nqtot = m_baseflow[0].num_elements();
583  int nvar = m_session->GetVariables().size();
584  int s;
585  Array<OneD, NekDouble> tmp_coeff(pFields[0]->GetNcoeffs(), 0.0);
586 
587  int numexp = pFields[0]->GetExpSize();
588  Array<OneD,int> ElementGIDs(numexp);
589 
590  // Define list of global element ids
591  for(int i = 0; i < numexp; ++i)
592  {
593  ElementGIDs[i] = pFields[0]->GetExp(i)->GetGeom()->GetGlobalID();
594  }
595 
596  //Get Homogeneous
599  m_session->GetComm());
600  fld->Import(pInfile, FieldDef, FieldData,
602  ElementGIDs);
603 
604 
605  if(m_session->DefinesSolverInfo("HOMOGENEOUS"))
606  {
607  std::string HomoStr = m_session->GetSolverInfo("HOMOGENEOUS");
608  }
609 
610  // copy FieldData into m_fields
611  for(int j = 0; j < nvar; ++j)
612  {
613  for(int i = 0; i < FieldDef.size(); ++i)
614  {
615  if((m_session->DefinesSolverInfo("HOMOGENEOUS") &&
616  (m_session->GetSolverInfo("HOMOGENEOUS")=="HOMOGENEOUS1D" ||
617  m_session->GetSolverInfo("HOMOGENEOUS")=="1D" ||
618  m_session->GetSolverInfo("HOMOGENEOUS")=="Homo1D")) &&
619  m_multipleModes==false)
620  {
621  // w-component must be ignored and set to zero.
622  // SJS I do not believe this is always the case
623  if (j != nvar - 2)
624  {
625  // p component (it is 4th variable of the 3D and
626  // corresponds 3rd variable of 2D)
627  s = (j == nvar - 1) ? 2 : j;
628 
629  //extraction of the 2D
630  pFields[j]->ExtractDataToCoeffs(
631  FieldDef[i],
632  FieldData[i],
633  FieldDef[i]->m_fields[s],
634  tmp_coeff);
635 
636  }
637 
638  //Zero higher modes than being considered in
639  //multi-mode expansion
640  if (m_npointsZ > 2)
641  {
642  // @Chris not clear why this is being done. Seems
643  // to be zeroing all but the first two complex
644  // modes?
645  int ncplane = (pFields[0]->GetNcoeffs()) / m_npointsZ;
646  Vmath::Zero(ncplane*(m_npointsZ-2),
647  &tmp_coeff[2*ncplane], 1);
648  }
649  }
650  // 2D cases
651  else
652  {
653  bool flag = FieldDef[i]->m_fields[j] ==
654  m_session->GetVariable(j);
655 
656  ASSERTL0(flag, (std::string("Order of ") + pInfile
657  + std::string(" data and that defined in "
658  "m_boundaryconditions differs")).c_str());
659 
660  pFields[j]->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
661  FieldDef[i]->m_fields[j],
662  tmp_coeff);
663  }
664  }
665 
666  if(m_singleMode || m_halfMode)
667  {
668  //pFields[j]->SetWaveSpace(true);
669 
670  pFields[j]->GetPlane(0)->BwdTrans(tmp_coeff, m_baseflow[j]);
671 
672  if(m_singleMode)
673  {
674  //copy the bwd trans into the second plane for single
675  //Mode Analysis
676  int ncplane=(pFields[0]->GetNpoints())/m_npointsZ;
677  Vmath::Vcopy(ncplane,&m_baseflow[j][0],1,&m_baseflow[j][ncplane],1);
678  }
679  }
680  else // fully 3D base flow - put in physical space.
681  {
682  bool oldwavespace = pFields[j]->GetWaveSpace();
683  pFields[j]->SetWaveSpace(false);
684  pFields[j]->BwdTrans(tmp_coeff, m_baseflow[j]);
685  pFields[j]->SetWaveSpace(oldwavespace);
686  }
687  }
688 
689  if(m_session->DefinesParameter("N_slices"))
690  {
691  int nConvectiveFields = pFields.num_elements()-1;
692 
693  for(int i=0; i<nConvectiveFields;++i)
694  {
695 
696  Vmath::Vcopy(nqtot, &m_baseflow[i][0], 1, &m_interp[i][slice*nqtot], 1);
697  }
698 
699  }
700 }
701 
702 
704  const NekDouble m_slices,
705  const Array<OneD, const NekDouble> &inarray,
706  Array<OneD, NekDouble> &outarray,
707  const NekDouble m_time,
708  const NekDouble m_period)
709 {
710  int npoints = m_baseflow[0].num_elements();
711  NekDouble BetaT = 2*M_PI*fmod (m_time, m_period) / m_period;
712  NekDouble phase;
713  Array<OneD, NekDouble> auxiliary(npoints);
714 
715  Vmath::Vcopy(npoints,&inarray[0],1,&outarray[0],1);
716  Vmath::Svtvp(npoints, cos(0.5*m_slices*BetaT),&inarray[npoints],1,&outarray[0],1,&outarray[0],1);
717 
718  for (int i = 2; i < m_slices; i += 2)
719  {
720  phase = (i>>1) * BetaT;
721 
722  Vmath::Svtvp(npoints, cos(phase),&inarray[i*npoints],1,&outarray[0],1,&outarray[0],1);
723  Vmath::Svtvp(npoints, -sin(phase), &inarray[(i+1)*npoints], 1, &outarray[0], 1,&outarray[0],1);
724  }
725 
726 }
727 
729 {
730  DNekMatSharedPtr loc_mat;
731  DNekBlkMatSharedPtr BlkMatrix;
732  int n_exp = 0;
733 
734  n_exp = m_baseflow[0].num_elements(); // will operatore on m_phys
735 
736  Array<OneD,unsigned int> nrows(n_exp);
737  Array<OneD,unsigned int> ncols(n_exp);
738 
739  nrows = Array<OneD, unsigned int>(n_exp,m_slices);
740  ncols = Array<OneD, unsigned int>(n_exp,m_slices);
741 
742  MatrixStorage blkmatStorage = eDIAGONAL;
743  BlkMatrix = MemoryManager<DNekBlkMat>
744  ::AllocateSharedPtr(nrows,ncols,blkmatStorage);
745 
746 
749  StdRegions::StdSegExp StdSeg(BK);
750 
752  StdSeg.DetShapeType(),
753  StdSeg);
754 
755  loc_mat = StdSeg.GetStdMatrix(matkey);
756 
757  // set up array of block matrices.
758  for(int i = 0; i < n_exp; ++i)
759  {
760  BlkMatrix->SetBlock(i,i,loc_mat);
761  }
762 
763  return BlkMatrix;
764 }
765 
766 //Discrete Fourier Transform for Floquet analysis
767 void LinearisedAdvection::DFT(const string file,
769  const NekDouble m_slices)
770 {
771  int ConvectedFields = m_baseflow.num_elements()-1;
772  int npoints = m_baseflow[0].num_elements();
773  m_interp = Array<OneD, Array<OneD, NekDouble> > (ConvectedFields);
774 
775  for (int i = 0; i < ConvectedFields; ++i)
776  {
777  m_interp[i] = Array<OneD,NekDouble>(npoints*m_slices, 0.0);
778  }
779 
780  // Import the slides into auxiliary vector
781  // The base flow should be stored in the form "filename_%d.ext"
782  // A subdirectory can also be included, such as "dir/filename_%d.ext"
783  size_t found = file.find("%d");
784  ASSERTL0(found != string::npos && file.find("%d", found+1) == string::npos,
785  "Since N_slices is specified, the filename provided for function "
786  "'BaseFlow' must include exactly one instance of the format "
787  "specifier '%d', to index the time-slices.");
788  char* buffer = new char[file.length() + 8];
789  for (int i = 0; i < m_slices; ++i)
790  {
791  sprintf(buffer, file.c_str(), i);
792  ImportFldBase(buffer,pFields,i);
793  }
794  delete[] buffer;
795 
796 
797  // Discrete Fourier Transform of the fields
798  for(int k=0; k< ConvectedFields;++k)
799  {
800 #ifdef NEKTAR_USING_FFTW
801 
802  //Discrete Fourier Transform using FFTW
803  Array<OneD, NekDouble> fft_in(npoints*m_slices);
804  Array<OneD, NekDouble> fft_out(npoints*m_slices);
805 
808 
809  //Shuffle the data
810  for(int j= 0; j < m_slices; ++j)
811  {
812  Vmath::Vcopy(npoints,&m_interp[k][j*npoints],1,&(fft_in[j]),m_slices);
813  }
814 
815  m_FFT = LibUtilities::GetNektarFFTFactory().CreateInstance("NekFFTW", m_slices);
816 
817  //FFT Transform
818  for(int i=0; i<npoints; i++)
819  {
820  m_FFT->FFTFwdTrans(m_tmpIN =fft_in + i*m_slices, m_tmpOUT =fft_out + i*m_slices);
821 
822  }
823 
824  //Reshuffle data
825  for(int s = 0; s < m_slices; ++s)
826  {
827  Vmath::Vcopy(npoints,&fft_out[s],m_slices,&m_interp[k][s*npoints],1);
828 
829  }
830 
831  Vmath::Zero(fft_in.num_elements(),&fft_in[0],1);
832  Vmath::Zero(fft_out.num_elements(),&fft_out[0],1);
833 #else
834  //Discrete Fourier Transform using MVM
835  DNekBlkMatSharedPtr blkmat;
837 
838  int nrows = blkmat->GetRows();
839  int ncols = blkmat->GetColumns();
840 
841  Array<OneD, NekDouble> sortedinarray(ncols);
842  Array<OneD, NekDouble> sortedoutarray(nrows);
843 
844  //Shuffle the data
845  for(int j= 0; j < m_slices; ++j)
846  {
847  Vmath::Vcopy(npoints,&m_interp[k][j*npoints],1,&(sortedinarray[j]),m_slices);
848  }
849 
850  // Create NekVectors from the given data arrays
851  NekVector<NekDouble> in (ncols,sortedinarray,eWrapper);
852  NekVector<NekDouble> out(nrows,sortedoutarray,eWrapper);
853 
854  // Perform matrix-vector multiply.
855  out = (*blkmat)*in;
856 
857  //Reshuffle data
858  for(int s = 0; s < m_slices; ++s)
859  {
860  Vmath::Vcopy(npoints,&sortedoutarray[s],m_slices,&m_interp[k][s*npoints],1);
861  }
862 
863  for(int r=0; r<sortedinarray.num_elements();++r)
864  {
865  sortedinarray[0]=0;
866  sortedoutarray[0]=0;
867  }
868 
869 #endif
870 
871  //scaling of the Fourier coefficients
872  NekDouble j=-1;
873  for (int i = 2; i < m_slices; i += 2)
874  {
875  Vmath::Smul(2*npoints,j,&m_interp[k][i*npoints],1,&m_interp[k][i*npoints],1);
876  j=-j;
877 
878  }
879 
880  }
881 
882 }
883 
884 } //end of namespace
885 
enum HomogeneousType m_HomogeneousType
void DFT(const string file, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const NekDouble m_slices)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Array< OneD, NekDouble > m_tmpIN
int m_HomoDirec
number of homogenous directions
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 boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
bool m_singleMode
flag to determine if use single mode or not
Array< OneD, Array< OneD, NekDouble > > m_baseflow
Storage for base flow.
LibUtilities::NektarFFTSharedPtr m_FFT
auxiliary variables
void ImportFldBase(std::string pInfile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, int slice)
Import Base flow.
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:471
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:428
array buffer
Definition: GsLib.hpp:56
int m_slices
number of slices
Array< OneD, NekDouble > m_tmpOUT
Fourier Expansion .
Definition: BasisType.h:52
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 void v_SetBaseFlow(const Array< OneD, Array< OneD, NekDouble > > &inarray)
Overrides the base flow used during linearised advection.
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:700
Array< OneD, Array< OneD, NekDouble > > m_interp
interpolation vector
Class representing a segment element in reference space.
Definition: StdSegExp.h:54
void UpdateBase(const NekDouble m_slices, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble m_time, const NekDouble m_period)
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:64
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
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)
Advects a vector field.
int m_npointsZ
number of points in Z direction (if homogeneous)
static std::string npts
Definition: InputFld.cpp:43
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:46
Defines a specification for a set of points.
Definition: Points.h:58
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
boost::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:225
double NekDouble
NekDouble m_period
period length
DNekBlkMatSharedPtr GetFloquetBlockMatrix(FloquetMatType mattype, bool UseContCoeffs=false) const
int m_npointsY
number of points in Y direction (if homogeneous)
boost::shared_ptr< Equation > EquationSharedPtr
bool m_multipleModes
flag to determine if use multiple mode or not
static SolverUtils::AdvectionSharedPtr create(std::string)
Creates an instance of this class.
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
NekDouble m_LhomY
physical length in Y direction (if homogeneous)
LibUtilities::SessionReaderSharedPtr m_session
boost::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:72
int m_npointsX
number of points in X direction (if homogeneous)
NekDouble m_LhomX
physical length in X direction (if homogeneous)
static std::string className
Name of class.
bool m_halfMode
flag to determine if use half mode or not
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
NekDouble m_LhomZ
physical length in Z direction (if homogeneous)
#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
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:54
bool m_useFFT
flag to determine if use or not the FFT for transformations
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
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215