Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
UnsteadyAdvectionDiffusion.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File UnsteadyAdvectionDiffusion.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: Unsteady advection-diffusion solve routines
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <iostream>
37 
39 
40 using namespace std;
41 
42 namespace Nektar
43 {
44  string UnsteadyAdvectionDiffusion::className
46  "UnsteadyAdvectionDiffusion",
47  UnsteadyAdvectionDiffusion::create);
48 
49  UnsteadyAdvectionDiffusion::UnsteadyAdvectionDiffusion(
51  : UnsteadySystem(pSession),
52  AdvectionSystem(pSession)
53  {
54  m_planeNumber = 0;
55  }
56 
57  /**
58  * @brief Initialisation object for the unsteady linear advection
59  * diffusion equation.
60  */
62  {
64 
65  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
66  m_session->LoadParameter("epsilon", m_epsilon, 0.0);
67 
68  // turn on substepping
69  m_session->MatchSolverInfo("Extrapolation", "SubStepping",
70  m_subSteppingScheme, false);
71 
72  // Define Velocity fields
74  std::vector<std::string> vel;
75  vel.push_back("Vx");
76  vel.push_back("Vy");
77  vel.push_back("Vz");
78  vel.resize(m_spacedim);
79 
80  EvaluateFunction(vel, m_velocity, "AdvectionVelocity");
81 
82  m_session->MatchSolverInfo(
83  "SpectralVanishingViscosity", "True", m_useSpecVanVisc, false);
84 
86  {
87  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
88  m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
89  }
90 
91  // Type of advection and diffusion classes to be used
92  switch(m_projectionType)
93  {
94  // Discontinuous field
96  {
97  // Do not forwards transform initial condition
98  m_homoInitialFwd = false;
99 
100  // Advection term
101  string advName;
102  string riemName;
103  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
105  CreateInstance(advName, advName);
106  m_advObject->SetFluxVector(&UnsteadyAdvectionDiffusion::
107  GetFluxVectorAdv, this);
108  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
110  CreateInstance(riemName);
111  m_riemannSolver->SetScalar("Vn", &UnsteadyAdvectionDiffusion::
112  GetNormalVelocity, this);
113  m_advObject->SetRiemannSolver(m_riemannSolver);
114  m_advObject->InitObject (m_session, m_fields);
115 
116  // Diffusion term
118  {
119  std::string diffName;
120  m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
122  CreateInstance(diffName, diffName);
123  m_diffusion->SetFluxVector(&UnsteadyAdvectionDiffusion::
124  GetFluxVectorDiff, this);
125  m_diffusion->InitObject(m_session, m_fields);
126  }
127 
128  ASSERTL0(m_subSteppingScheme == false,"SubSteppingScheme is not set up for DG projection");
129  break;
130  }
131  // Continuous field
134  {
135  // Advection term
136  std::string advName;
137  m_session->LoadSolverInfo("AdvectionType", advName,
138  "NonConservative");
140  CreateInstance(advName, advName);
141  m_advObject->SetFluxVector(&UnsteadyAdvectionDiffusion::
142  GetFluxVectorAdv, this);
143 
144  if(advName.compare("WeakDG") == 0)
145  {
146  string riemName;
147  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
149  CreateInstance(riemName);
150  m_riemannSolver->SetScalar("Vn",
152  GetNormalVelocity, this);
153  m_advObject->SetRiemannSolver(m_riemannSolver);
154  m_advObject->InitObject (m_session, m_fields);
155  }
156 
157  // In case of Galerkin explicit diffusion gives an error
159  {
160  ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
161  }
162  // In case of Galerkin implicit diffusion: do nothing
163  break;
164  }
165  default:
166  {
167  ASSERTL0(false, "Unsupported projection type.");
168  break;
169  }
170  }
171 
177 
178  if(m_subSteppingScheme) // Substepping
179  {
181  "Projection must be set to Mixed_CG_Discontinuous for "
182  "substepping");
184  m_intScheme->GetIntegrationMethod(), m_intScheme);
185  }
186  }
187 
188  /**
189  * @brief Unsteady linear advection diffusion equation destructor.
190  */
192  {
193  }
194 
195  /**
196  * @brief Get the normal velocity for the unsteady linear advection
197  * diffusion equation.
198  */
200  {
201  return GetNormalVel(m_velocity);
202  }
203 
204 
206  const Array<OneD, const Array<OneD, NekDouble> > &velfield)
207  {
208  // Number of trace (interface) points
209  int i;
210  int nTracePts = GetTraceNpoints();
211 
212  // Auxiliary variable to compute the normal velocity
213  Array<OneD, NekDouble> tmp(nTracePts);
214  m_traceVn = Array<OneD, NekDouble>(nTracePts, 0.0);
215 
216  // Reset the normal velocity
217  Vmath::Zero(nTracePts, m_traceVn, 1);
218 
219  for (i = 0; i < velfield.num_elements(); ++i)
220  {
221  m_fields[0]->ExtractTracePhys(velfield[i], tmp);
222 
223  Vmath::Vvtvp(nTracePts,
224  m_traceNormals[i], 1,
225  tmp, 1,
226  m_traceVn, 1,
227  m_traceVn, 1);
228  }
229 
230  return m_traceVn;
231  }
232 
233  /**
234  * @brief Compute the right-hand side for the unsteady linear advection
235  * diffusion problem.
236  *
237  * @param inarray Given fields.
238  * @param outarray Calculated solution.
239  * @param time Time.
240  */
242  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
243  Array<OneD, Array<OneD, NekDouble> >&outarray,
244  const NekDouble time)
245  {
246  // Number of fields (variables of the problem)
247  int nVariables = inarray.num_elements();
248 
249  // Number of solution points
250  int nSolutionPts = GetNpoints();
251 
252  // RHS computation using the new advection base class
253  m_advObject->Advect(nVariables, m_fields, m_velocity,
254  inarray, outarray, time);
255 
256  // Negate the RHS
257  for (int i = 0; i < nVariables; ++i)
258  {
259  Vmath::Neg(nSolutionPts, outarray[i], 1);
260  }
261 
263  {
264  Array<OneD, Array<OneD, NekDouble> > outarrayDiff(nVariables);
265  for (int i = 0; i < nVariables; ++i)
266  {
267  outarrayDiff[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
268  }
269 
270  m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
271 
272  for (int i = 0; i < nVariables; ++i)
273  {
274  Vmath::Vadd(nSolutionPts, &outarray[i][0], 1,
275  &outarrayDiff[i][0], 1, &outarray[i][0], 1);
276  }
277  }
278 
279  }
280 
281  /**
282  * @brief Compute the projection for the unsteady advection
283  * diffusion problem.
284  *
285  * @param inarray Given fields.
286  * @param outarray Calculated solution.
287  * @param time Time.
288  */
290  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
291  Array<OneD, Array<OneD, NekDouble> > &outarray,
292  const NekDouble time)
293  {
294  int i;
295  int nvariables = inarray.num_elements();
296  SetBoundaryConditions(time);
297  switch(m_projectionType)
298  {
300  {
301  // Just copy over array
302  int npoints = GetNpoints();
303 
304  for(i = 0; i < nvariables; ++i)
305  {
306  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
307  }
308  break;
309  }
312  {
314 
315  for(i = 0; i < nvariables; ++i)
316  {
317  m_fields[i]->FwdTrans(inarray[i], coeffs);
318  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
319  }
320  break;
321  }
322  default:
323  {
324  ASSERTL0(false, "Unknown projection scheme");
325  break;
326  }
327  }
328  }
329 
330 
331  /* @brief Compute the diffusion term implicitly.
332  *
333  * @param inarray Given fields.
334  * @param outarray Calculated solution.
335  * @param time Time.
336  * @param lambda Diffusion coefficient.
337  */
339  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
340  Array<OneD, Array<OneD, NekDouble> >&outarray,
341  const NekDouble time,
342  const NekDouble lambda)
343  {
344  int nvariables = inarray.num_elements();
345  int nq = m_fields[0]->GetNpoints();
346 
348  factors[StdRegions::eFactorLambda] = 1.0/lambda/m_epsilon;
349 
350  if(m_useSpecVanVisc)
351  {
354  }
356  {
357  factors[StdRegions::eFactorTau] = 1.0;
358  }
359 
360  Array<OneD, Array< OneD, NekDouble> > F(nvariables);
361  for (int n = 0; n < nvariables; ++n)
362  {
363  F[n] = Array<OneD, NekDouble> (nq);
364  }
365 
366  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
367  // inarray = input: \hat{rhs} -> output: \hat{Y}
368  // outarray = output: nabla^2 \hat{Y}
369  // where \hat = modal coeffs
370  for (int i = 0; i < nvariables; ++i)
371  {
372  // Multiply 1.0/timestep/lambda
373  Vmath::Smul(nq, -factors[StdRegions::eFactorLambda],
374  inarray[i], 1, F[i], 1);
375  }
376 
377  //Setting boundary conditions
378  SetBoundaryConditions(time);
379 
380  for (int i = 0; i < nvariables; ++i)
381  {
382  // Solve a system of equations with Helmholtz solver
383  m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(),
384  NullFlagList, factors);
385 
386  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
387  }
388  }
389 
390  /**
391  * @brief Return the flux vector for the advection part.
392  *
393  * @param physfield Fields.
394  * @param flux Resulting flux.
395  */
397  const Array<OneD, Array<OneD, NekDouble> > &physfield,
399  {
400  ASSERTL1(flux[0].num_elements() == m_velocity.num_elements(),
401  "Dimension of flux array and velocity array do not match");
402 
403  const int nq = m_fields[0]->GetNpoints();
404 
405  for (int i = 0; i < flux.num_elements(); ++i)
406  {
407  for (int j = 0; j < flux[0].num_elements(); ++j)
408  {
409  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1,
410  flux[i][j], 1);
411  }
412  }
413  }
414 
415  /**
416  * @brief Return the flux vector for the diffusion part.
417  *
418  * @param i Equation number.
419  * @param j Spatial direction.
420  * @param physfield Fields.
421  * @param derivatives First order derivatives.
422  * @param flux Resulting flux.
423  */
425  const int i,
426  const int j,
427  const Array<OneD, Array<OneD, NekDouble> > &physfield,
428  Array<OneD, Array<OneD, NekDouble> > &derivatives,
430  {
431  for (int k = 0; k < flux.num_elements(); ++k)
432  {
433  Vmath::Zero(GetNpoints(), flux[k], 1);
434  }
435  Vmath::Vcopy(GetNpoints(), physfield[i], 1, flux[j], 1);
436  }
437 
440  {
442  }
443 
444 
445  /**
446  * Perform the extrapolation.
447  */
449  {
451  {
453  }
454 
455  return false;
456  }
457 
458 
459  /**
460  *
461  */
463  const LibUtilities::TimeIntegrationSolutionSharedPtr &integrationSoln,
464  int nstep,
465  NekDouble time)
466  {
467  int n;
468  int nsubsteps;
469 
470  NekDouble dt;
471 
472  Array<OneD, Array<OneD, NekDouble> > fields, velfields;
473 
474  static int ncalls = 1;
475  int nint = min(ncalls++, m_intSteps);
476 
479 
480  LibUtilities::CommSharedPtr comm = m_session->GetComm();
481 
482  // Get the proper time step with CFL control
483  dt = GetSubstepTimeStep();
484 
485  nsubsteps = (m_timestep > dt)? ((int)(m_timestep/dt)+1):1;
486  nsubsteps = max(m_minsubsteps, nsubsteps);
487 
488  dt = m_timestep/nsubsteps;
489 
490  if (m_infosteps && !((nstep+1)%m_infosteps) && comm->GetRank() == 0)
491  {
492  cout << "Sub-integrating using "<< nsubsteps
493  << " steps over Dt = " << m_timestep
494  << " (SubStep CFL=" << m_cflSafetyFactor << ")"<< endl;
495  }
496 
497  for (int m = 0; m < nint; ++m)
498  {
499  // We need to update the fields held by the m_integrationSoln
500  fields = integrationSoln->UpdateSolutionVector()[m];
501 
502  // Initialise NS solver which is set up to use a GLM method
503  // with calls to EvaluateAdvection_SetPressureBCs and
504  // SolveUnsteadyStokesSystem
506  SubIntegrationSoln = m_subStepIntegrationScheme->
507  InitializeScheme(dt, fields, time, m_subStepIntegrationOps);
508 
509  for(n = 0; n < nsubsteps; ++n)
510  {
511  fields = m_subStepIntegrationScheme->TimeIntegrate(n, dt, SubIntegrationSoln,
513  }
514 
515  // Reset time integrated solution in m_integrationSoln
516  integrationSoln->SetSolVector(m,fields);
517  }
518  }
519 
520 
521  /**
522  *
523  */
525  {
526  int n_element = m_fields[0]->GetExpSize();
527 
528  const Array<OneD, int> ExpOrder=m_fields[0]->EvalBasisNumModesMaxPerExp();
529  Array<OneD, int> ExpOrderList (n_element, ExpOrder);
530 
531  const NekDouble cLambda = 0.2; // Spencer book pag. 317
532 
533  Array<OneD, NekDouble> tstep (n_element, 0.0);
534  Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
535 
536  stdVelocity = GetMaxStdVelocity(m_velocity);
537 
538  for(int el = 0; el < n_element; ++el)
539  {
540  tstep[el] = m_cflSafetyFactor /
541  (stdVelocity[el] * cLambda *
542  (ExpOrder[el]-1) * (ExpOrder[el]-1));
543  }
544 
545  NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
546  m_session->GetComm()->AllReduce(TimeStep,LibUtilities::ReduceMin);
547 
548  return TimeStep;
549  }
550 
552  int intMethod,
553  const LibUtilities::TimeIntegrationWrapperSharedPtr &IntegrationScheme)
554  {
555  // Set to 1 for first step and it will then be increased in
556  // time advance routines
557  switch(intMethod)
558  {
561  {
563 
564  }
565  break;
567  {
569  }
570  break;
571  default:
572  ASSERTL0(0,"Integration method not suitable: Options include BackwardEuler or BDFImplicitOrder1");
573  break;
574  }
575  m_intSteps = IntegrationScheme->GetIntegrationSteps();
576 
577  // set explicit time-integration class operators
580  }
581 
582  /**
583  * Explicit Advection terms used by SubStepAdvance time integration
584  */
586  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
587  Array<OneD, Array<OneD, NekDouble> > &outarray,
588  const NekDouble time)
589  {
590  int i;
591  int nVariables = inarray.num_elements();
592 
593  /// Get the number of coefficients
594  int ncoeffs = m_fields[0]->GetNcoeffs();
595 
596  /// Define an auxiliary variable to compute the RHS
597  Array<OneD, Array<OneD, NekDouble> > WeakAdv(nVariables);
598  WeakAdv[0] = Array<OneD, NekDouble> (ncoeffs*nVariables);
599  for(i = 1; i < nVariables; ++i)
600  {
601  WeakAdv[i] = WeakAdv[i-1] + ncoeffs;
602  }
603 
604  // Currently assume velocity field is time independent and does not therefore
605  // need extrapolating.
606  // RHS computation using the advection base class
607  m_advObject->Advect(nVariables, m_fields, m_velocity,
608  inarray, outarray, time);
609 
610  for(i = 0; i < nVariables; ++i)
611  {
612  m_fields[i]->IProductWRTBase(outarray[i],WeakAdv[i]);
613  // negation requried due to sign of DoAdvection term to be consistent
614  Vmath::Neg(ncoeffs, WeakAdv[i], 1);
615  }
616 
617  AddAdvectionPenaltyFlux(m_velocity, inarray, WeakAdv);
618 
619 
620  /// Operations to compute the RHS
621  for(i = 0; i < nVariables; ++i)
622  {
623  // Negate the RHS
624  Vmath::Neg(ncoeffs, WeakAdv[i], 1);
625 
626  /// Multiply the flux by the inverse of the mass matrix
627  m_fields[i]->MultiplyByElmtInvMass(WeakAdv[i], WeakAdv[i]);
628 
629  /// Store in outarray the physical values of the RHS
630  m_fields[i]->BwdTrans(WeakAdv[i], outarray[i]);
631  }
632  }
633 
634  /**
635  * Projection used by SubStepAdvance time integration
636  */
638  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
639  Array<OneD, Array<OneD, NekDouble> > &outarray,
640  const NekDouble time)
641  {
642  ASSERTL1(inarray.num_elements() == outarray.num_elements(),"Inarray and outarray of different sizes ");
643 
644  for(int i = 0; i < inarray.num_elements(); ++i)
645  {
646  Vmath::Vcopy(inarray[i].num_elements(),inarray[i],1,outarray[i],1);
647  }
648  }
649 
651  const Array<OneD, const Array<OneD, NekDouble> > &velfield,
652  const Array<OneD, const Array<OneD, NekDouble> > &physfield,
653  Array<OneD, Array<OneD, NekDouble> > &Outarray)
654  {
655  ASSERTL1(physfield.num_elements() == Outarray.num_elements(),
656  "Physfield and outarray are of different dimensions");
657 
658  int i;
659 
660  /// Number of trace points
661  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
662 
663  /// Forward state array
664  Array<OneD, NekDouble> Fwd(3*nTracePts);
665 
666  /// Backward state array
667  Array<OneD, NekDouble> Bwd = Fwd + nTracePts;
668 
669  /// upwind numerical flux state array
670  Array<OneD, NekDouble> numflux = Bwd + nTracePts;
671 
672  /// Normal velocity array
673  Array<OneD, NekDouble> Vn = GetNormalVel(velfield);
674 
675  for(i = 0; i < physfield.num_elements(); ++i)
676  {
677  /// Extract forwards/backwards trace spaces
678  /// Note: Needs to have correct i value to get boundary conditions
679  m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
680 
681  /// Upwind between elements
682  m_fields[0]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux);
683 
684  /// Construct difference between numflux and Fwd,Bwd
685  Vmath::Vsub(nTracePts, numflux, 1, Fwd, 1, Fwd, 1);
686  Vmath::Vsub(nTracePts, numflux, 1, Bwd, 1, Bwd, 1);
687 
688  /// Calculate the numerical fluxes multipling Fwd, Bwd and
689  /// numflux by the normal advection velocity
690  Vmath::Vmul(nTracePts, Fwd, 1, Vn, 1, Fwd, 1);
691  Vmath::Vmul(nTracePts, Bwd, 1, Vn, 1, Bwd, 1);
692 
693  m_fields[0]->AddFwdBwdTraceIntegral(Fwd,Bwd,Outarray[i]);
694  }
695  }
696 
697 
699  const Array<OneD, Array<OneD,NekDouble> > inarray)
700  {
701 
702  int n_points_0 = m_fields[0]->GetExp(0)->GetTotPoints();
703  int n_element = m_fields[0]->GetExpSize();
704  int nvel = inarray.num_elements();
705  int cnt;
706 
707  ASSERTL0(nvel >= 2, "Method not implemented for 1D");
708 
709  NekDouble pntVelocity;
710 
711  // Getting the standard velocity vector on the 2D normal space
712  Array<OneD, Array<OneD, NekDouble> > stdVelocity(nvel);
713  Array<OneD, NekDouble> maxV(n_element, 0.0);
715 
716  for (int i = 0; i < nvel; ++i)
717  {
718  stdVelocity[i] = Array<OneD, NekDouble>(n_points_0);
719  }
720 
721  if (nvel == 2)
722  {
723  cnt = 0.0;
724  for (int el = 0; el < n_element; ++el)
725  {
726  int n_points = m_fields[0]->GetExp(el)->GetTotPoints();
727  ptsKeys = m_fields[0]->GetExp(el)->GetPointsKeys();
728 
729  // reset local space if necessary
730  if(n_points != n_points_0)
731  {
732  for (int j = 0; j < nvel; ++j)
733  {
734  stdVelocity[j] = Array<OneD, NekDouble>(n_points);
735  }
736  n_points_0 = n_points;
737  }
738 
740  m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetDerivFactors(ptsKeys);
741 
742  if (m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetGtype()
744  {
745  for (int i = 0; i < n_points; i++)
746  {
747  stdVelocity[0][i] = gmat[0][i]*inarray[0][i+cnt]
748  + gmat[2][i]*inarray[1][i+cnt];
749 
750  stdVelocity[1][i] = gmat[1][i]*inarray[0][i+cnt]
751  + gmat[3][i]*inarray[1][i+cnt];
752  }
753  }
754  else
755  {
756  for (int i = 0; i < n_points; i++)
757  {
758  stdVelocity[0][i] = gmat[0][0]*inarray[0][i+cnt]
759  + gmat[2][0]*inarray[1][i+cnt];
760 
761  stdVelocity[1][i] = gmat[1][0]*inarray[0][i+cnt]
762  + gmat[3][0]*inarray[1][i+cnt];
763  }
764  }
765 
766  cnt += n_points;
767 
768 
769  for (int i = 0; i < n_points; i++)
770  {
771  pntVelocity = stdVelocity[0][i]*stdVelocity[0][i]
772  + stdVelocity[1][i]*stdVelocity[1][i];
773 
774  if (pntVelocity>maxV[el])
775  {
776  maxV[el] = pntVelocity;
777  }
778  }
779  maxV[el] = sqrt(maxV[el]);
780  }
781  }
782  else
783  {
784  cnt = 0;
785  for (int el = 0; el < n_element; ++el)
786  {
787 
788  int n_points = m_fields[0]->GetExp(el)->GetTotPoints();
789  ptsKeys = m_fields[0]->GetExp(el)->GetPointsKeys();
790 
791  // reset local space if necessary
792  if(n_points != n_points_0)
793  {
794  for (int j = 0; j < nvel; ++j)
795  {
796  stdVelocity[j] = Array<OneD, NekDouble>(n_points);
797  }
798  n_points_0 = n_points;
799  }
800 
802  m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetDerivFactors(ptsKeys);
803 
804  if (m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetGtype()
806  {
807  for (int i = 0; i < n_points; i++)
808  {
809  stdVelocity[0][i] = gmat[0][i]*inarray[0][i+cnt]
810  + gmat[3][i]*inarray[1][i+cnt]
811  + gmat[6][i]*inarray[2][i+cnt];
812 
813  stdVelocity[1][i] = gmat[1][i]*inarray[0][i+cnt]
814  + gmat[4][i]*inarray[1][i+cnt]
815  + gmat[7][i]*inarray[2][i+cnt];
816 
817  stdVelocity[2][i] = gmat[2][i]*inarray[0][i+cnt]
818  + gmat[5][i]*inarray[1][i+cnt]
819  + gmat[8][i]*inarray[2][i+cnt];
820  }
821  }
822  else
823  {
824  for (int i = 0; i < n_points; i++)
825  {
826  stdVelocity[0][i] = gmat[0][0]*inarray[0][i+cnt]
827  + gmat[3][0]*inarray[1][i+cnt]
828  + gmat[6][0]*inarray[2][i+cnt];
829 
830  stdVelocity[1][i] = gmat[1][0]*inarray[0][i+cnt]
831  + gmat[4][0]*inarray[1][i+cnt]
832  + gmat[7][0]*inarray[2][i+cnt];
833 
834  stdVelocity[2][i] = gmat[2][0]*inarray[0][i+cnt]
835  + gmat[5][0]*inarray[1][i+cnt]
836  + gmat[8][0]*inarray[2][i+cnt];
837  }
838  }
839 
840  cnt += n_points;
841 
842  for (int i = 0; i < n_points; i++)
843  {
844  pntVelocity = stdVelocity[0][i]*stdVelocity[0][i]
845  + stdVelocity[1][i]*stdVelocity[1][i]
846  + stdVelocity[2][i]*stdVelocity[2][i];
847 
848  if (pntVelocity > maxV[el])
849  {
850  maxV[el] = pntVelocity;
851  }
852  }
853 
854  maxV[el] = sqrt(maxV[el]);
855  //cout << maxV[el]*maxV[el] << endl;
856  }
857  }
858 
859  return maxV;
860  }
861 }
void SubStepProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
void GetFluxVectorAdv(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Evaluate the flux at each solution point for the advection part.
Array< OneD, NekDouble > GetMaxStdVelocity(const Array< OneD, Array< OneD, NekDouble > > inarray)
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:242
BDF multi-step scheme of order 1 (implicit)
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
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
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
NekDouble m_time
Current time of simulation.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
LibUtilities::TimeIntegrationWrapperSharedPtr m_subStepIntegrationScheme
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
NekDouble m_timestep
Time step size.
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:47
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
Definition: Vmath.cpp:871
Array< OneD, Array< OneD, NekDouble > > m_velocity
boost::shared_ptr< TimeIntegrationWrapper > TimeIntegrationWrapperSharedPtr
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:42
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:442
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
STL namespace.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
void SubStepAdvection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Array< OneD, NekDouble > & GetNormalVel(const Array< OneD, const Array< OneD, NekDouble > > &velfield)
Get the normal velocity based on input velfield.
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity based on m_velocity.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Perform the projection.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
virtual void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
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:213
virtual bool v_PreIntegrate(int step)
PreIntegration step for substepping.
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
Base class for unsteady solvers.
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
Print Summary.
RiemannSolverFactory & GetRiemannSolverFactory()
int m_spacedim
Spatial dimension (>= expansion dim).
void SetUpSubSteppingTimeIntegration(int intMethod, const LibUtilities::TimeIntegrationWrapperSharedPtr &IntegrationScheme)
SolverUtils::DiffusionSharedPtr m_diffusion
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:46
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:396
double NekDouble
TimeIntegrationWrapperFactory & GetTimeIntegrationWrapperFactory()
SOLVER_UTILS_EXPORT void EvaluateFunction(Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
Evaluates a function as specified in the session file.
EquationSystemFactory & GetEquationSystemFactory()
void SubStepAdvance(const LibUtilities::TimeIntegrationSolutionSharedPtr &integrationSoln, int nstep, NekDouble time)
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:343
BDF multi-step scheme of order 2 (implicit)
virtual void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
Solve implicitly the diffusion term.
SOLVER_UTILS_EXPORT int GetNpoints()
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetExpSize()
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
void AddAdvectionPenaltyFlux(const Array< OneD, const Array< OneD, NekDouble > > &velfield, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &Outarray)
SOLVER_UTILS_EXPORT int GetNcoeffs()
boost::shared_ptr< TimeIntegrationSolution > TimeIntegrationSolutionSharedPtr
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
void GetFluxVectorDiff(const int i, const int j, const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &derivatives, Array< OneD, Array< OneD, NekDouble > > &flux)
Evaluate the flux at each solution point for the diffusion part.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:373
A base class for PDEs which include an advection component.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
Geometry is curved or has non-constant factors.
virtual void v_InitObject()
Initialise the object.
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:299
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:183
static FlagList NullFlagList
An empty flag list.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215