Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
IncNavierStokes.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File IncNavierStokes.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: Incompressible Navier Stokes class definition built on
33 // ADRBase class
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
37 #include <iomanip>
38 #include <boost/algorithm/string.hpp>
39 
46 
47 namespace Nektar
48 {
49 
50  /**
51  * Constructor. Creates ...
52  *
53  * \param
54  * \param
55  */
57  UnsteadySystem(pSession),
58  AdvectionSystem(pSession),
59  m_SmoothAdvection(false),
60  m_steadyStateSteps(0)
61  {
62  }
63 
65  {
66  AdvectionSystem::v_InitObject();
67 
68  int i,j;
69  int numfields = m_fields.num_elements();
70  std::string velids[] = {"u","v","w"};
71 
72  // Set up Velocity field to point to the first m_expdim of m_fields;
74 
75  for(i = 0; i < m_spacedim; ++i)
76  {
77  for(j = 0; j < numfields; ++j)
78  {
79  std::string var = m_boundaryConditions->GetVariable(j);
80  if(boost::iequals(velids[i], var))
81  {
82  m_velocity[i] = j;
83  break;
84  }
85 
86  ASSERTL0(j != numfields, "Failed to find field: " + var);
87  }
88  }
89 
90  // Set up equation type enum using kEquationTypeStr
91  for(i = 0; i < (int) eEquationTypeSize; ++i)
92  {
93  bool match;
94  m_session->MatchSolverInfo("EQTYPE",kEquationTypeStr[i],match,false);
95  if(match)
96  {
98  break;
99  }
100  }
101  ASSERTL0(i != eEquationTypeSize,"EQTYPE not found in SOLVERINFO section");
102 
103  // This probably should to into specific implementations
104  // Equation specific Setups
105  switch(m_equationType)
106  {
107  case eSteadyStokes:
108  case eSteadyOseen:
109  case eSteadyNavierStokes:
110  case eSteadyLinearisedNS:
111  break;
113  case eUnsteadyStokes:
114  {
115  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
116  m_session->LoadParameter("IO_CFLSteps", m_cflsteps, 0);
117  m_session->LoadParameter("SteadyStateSteps", m_steadyStateSteps, 0);
118  m_session->LoadParameter("SteadyStateTol", m_steadyStateTol, 1e-6);
119  }
120  break;
121  case eNoEquationType:
122  default:
123  ASSERTL0(false,"Unknown or undefined equation type");
124  }
125 
126  m_session->LoadParameter("Kinvis", m_kinvis);
127 
128  // Default advection type per solver
129  std::string vConvectiveType;
130  switch(m_equationType)
131  {
132  case eUnsteadyStokes:
133  case eSteadyLinearisedNS:
134  vConvectiveType = "NoAdvection";
135  break;
137  case eSteadyNavierStokes:
138  vConvectiveType = "Convective";
139  break;
141  vConvectiveType = "Linearised";
142  break;
143  default:
144  break;
145  }
146 
147  // Check if advection type overridden
148  if (m_session->DefinesTag("AdvectiveType") && m_equationType != eUnsteadyStokes &&
150  {
151  vConvectiveType = m_session->GetTag("AdvectiveType");
152  }
153 
154  // Initialise advection
155  m_advObject = SolverUtils::GetAdvectionFactory().CreateInstance(vConvectiveType, vConvectiveType);
156  m_advObject->InitObject( m_session, m_fields);
157 
158  // Forcing terms
161 
162  // check to see if any Robin boundary conditions and if so set
163  // up m_field to boundary condition maps;
167 
168  for (i = 0; i < m_fields.num_elements(); ++i)
169  {
170  bool Set = false;
171 
174  int radpts = 0;
175 
176  BndConds = m_fields[i]->GetBndConditions();
177  BndExp = m_fields[i]->GetBndCondExpansions();
178  for(int n = 0; n < BndConds.num_elements(); ++n)
179  {
180  if(boost::iequals(BndConds[n]->GetUserDefined(),"Radiation"))
181  {
182  ASSERTL0(BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eRobin,
183  "Radiation boundary condition must be of type Robin <R>");
184 
185  if(Set == false)
186  {
187  m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],m_fieldsBCToTraceID[i]);
188  Set = true;
189  }
190  radpts += BndExp[n]->GetTotPoints();
191  }
192  if(boost::iequals(BndConds[n]->GetUserDefined(),"ZeroNormalComponent"))
193  {
194  ASSERTL0(BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eDirichlet,
195  "Zero Normal Component boundary condition option must be of type Dirichlet <D>");
196 
197  if(Set == false)
198  {
199  m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],m_fieldsBCToTraceID[i]);
200  Set = true;
201  }
202  }
203  }
204 
205  m_fieldsRadiationFactor[i] = Array<OneD, NekDouble>(radpts);
206 
207  radpts = 0; // reset to use as a counter
208 
209  for(int n = 0; n < BndConds.num_elements(); ++n)
210  {
211  if(boost::iequals(BndConds[n]->GetUserDefined(),"Radiation"))
212  {
213 
214  int npoints = BndExp[n]->GetNpoints();
215  Array<OneD, NekDouble> x0(npoints,0.0);
216  Array<OneD, NekDouble> x1(npoints,0.0);
217  Array<OneD, NekDouble> x2(npoints,0.0);
218  Array<OneD, NekDouble> tmpArray;
219 
220  BndExp[n]->GetCoords(x0,x1,x2);
221 
222  LibUtilities::Equation coeff =
223  boost::static_pointer_cast<
225  >(BndConds[n])->m_robinPrimitiveCoeff;
226 
227  coeff.Evaluate(x0,x1,x2,m_time,
228  tmpArray = m_fieldsRadiationFactor[i]+ radpts);
229  //Vmath::Neg(npoints,tmpArray = m_fieldsRadiationFactor[i]+ radpts,1);
230  radpts += npoints;
231  }
232  }
233  }
234 
235  // Set up Field Meta Data for output files
236  m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
237  m_fieldMetaDataMap["TimeStep"] = boost::lexical_cast<std::string>(m_timestep);
238  }
239 
240  /**
241  * Destructor
242  */
244  {
245  }
246 
247 
248  /**
249  *
250  */
252  Array<OneD, Array<OneD, NekDouble> > &physfield,
254  {
255  ASSERTL1(flux.num_elements() == m_velocity.num_elements(),"Dimension of flux array and velocity array do not match");
256 
257  for(int j = 0; j < flux.num_elements(); ++j)
258  {
259  Vmath::Vmul(GetNpoints(), physfield[i], 1, m_fields[m_velocity[j]]->GetPhys(), 1, flux[j], 1);
260  }
261  }
262 
263  /**
264  * Calcualate numerical fluxes
265  */
267  Array<OneD, Array<OneD, NekDouble> > &numflux)
268  {
269  /// Counter variable
270  int i;
271 
272  /// Number of trace points
273  int nTracePts = GetTraceNpoints();
274 
275  /// Number of spatial dimensions
276  int nDimensions = m_spacedim;
277 
278  /// Forward state array
279  Array<OneD, NekDouble> Fwd(2*nTracePts);
280 
281  /// Backward state array
282  Array<OneD, NekDouble> Bwd = Fwd + nTracePts;
283 
284  /// Normal velocity array
285  Array<OneD, NekDouble> Vn (nTracePts, 0.0);
286 
287  // Extract velocity field along the trace space and multiply by trace normals
288  for(i = 0; i < nDimensions; ++i)
289  {
290  m_fields[0]->ExtractTracePhys(m_fields[m_velocity[i]]->GetPhys(), Fwd);
291  Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Fwd, 1, Vn, 1, Vn, 1);
292  }
293 
294  /// Compute the numerical fluxes at the trace points
295  for(i = 0; i < numflux.num_elements(); ++i)
296  {
297  /// Extract forwards/backwards trace spaces
298  m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
299 
300  /// Upwind between elements
301  m_fields[i]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux[i]);
302 
303  /// Calculate the numerical fluxes multipling Fwd or Bwd
304  /// by the normal advection velocity
305  Vmath::Vmul(nTracePts, numflux[i], 1, Vn, 1, numflux[i], 1);
306  }
307  }
308 
309  /**
310  * Evaluation -N(V) for all fields except pressure using m_velocity
311  */
313  Array<OneD, Array<OneD, NekDouble> > &outarray,
315  {
316  int i;
317  int nqtot = m_fields[0]->GetTotPoints();
318  int VelDim = m_velocity.num_elements();
319  Array<OneD, Array<OneD, NekDouble> > velocity(VelDim);
321 
322  for(i = 0; i < VelDim; ++i)
323  {
324  if(m_fields[i]->GetWaveSpace() && !m_singleMode && !m_halfMode)
325  {
326  velocity[i] = Array<OneD, NekDouble>(nqtot,0.0);
327  m_fields[i]->HomogeneousBwdTrans(inarray[m_velocity[i]],velocity[i]);
328  }
329  else
330  {
331  velocity[i] = inarray[m_velocity[i]];
332  }
333  }
334 
335  // Set up Derivative work space;
336  if(wk.num_elements())
337  {
338  ASSERTL0(wk.num_elements() >= nqtot*VelDim,
339  "Workspace is not sufficient");
340  Deriv = wk;
341  }
342  else
343  {
344  Deriv = Array<OneD, NekDouble> (nqtot*VelDim);
345  }
346 
348  velocity, inarray, outarray, m_time);
349  }
350 
351  /**
352  * Time dependent boundary conditions updating
353  */
355  {
356  int i, n;
357  std::string varName;
358  int nvariables = m_fields.num_elements();
359 
360  for (i = 0; i < nvariables; ++i)
361  {
362  for(n = 0; n < m_fields[i]->GetBndConditions().num_elements(); ++n)
363  {
364  if(m_fields[i]->GetBndConditions()[n]->IsTimeDependent())
365  {
366  varName = m_session->GetVariable(i);
367  m_fields[i]->EvaluateBoundaryConditions(time, varName);
368  }
369 
370  }
371 
372  // Set Radiation conditions if required
374  }
375 
377  }
378 
379  /**
380  * Probably should be pushed back into ContField?
381  */
383  {
384  int i,n;
385 
388 
389  BndConds = m_fields[fieldid]->GetBndConditions();
390  BndExp = m_fields[fieldid]->GetBndCondExpansions();
391 
394 
395  int cnt;
396  int elmtid,nq,offset, boundary;
397  Array<OneD, NekDouble> Bvals, U;
398  int cnt1 = 0;
399 
400  for(cnt = n = 0; n < BndConds.num_elements(); ++n)
401  {
402  std::string type = BndConds[n]->GetUserDefined();
403 
404  if((BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eRobin)&&(boost::iequals(type,"Radiation")))
405  {
406  for(i = 0; i < BndExp[n]->GetExpSize(); ++i,cnt++)
407  {
408  elmtid = m_fieldsBCToElmtID[m_velocity[fieldid]][cnt];
409  elmt = m_fields[fieldid]->GetExp(elmtid);
410  offset = m_fields[fieldid]->GetPhys_Offset(elmtid);
411 
412  U = m_fields[fieldid]->UpdatePhys() + offset;
413  Bc = BndExp[n]->GetExp(i);
414 
415  boundary = m_fieldsBCToTraceID[fieldid][cnt];
416 
417  // Get edge values and put into ubc
418  nq = Bc->GetTotPoints();
419  Array<OneD, NekDouble> ubc(nq);
420  elmt->GetTracePhysVals(boundary,Bc,U,ubc);
421 
422  Vmath::Vmul(nq,&m_fieldsRadiationFactor[fieldid][cnt1 +
423  BndExp[n]->GetPhys_Offset(i)],1,&ubc[0],1,&ubc[0],1);
424 
425  Bvals = BndExp[n]->UpdateCoeffs()+BndExp[n]->GetCoeff_Offset(i);
426 
427  Bc->IProductWRTBase(ubc,Bvals);
428  }
429  cnt1 += BndExp[n]->GetTotPoints();
430  }
431  else
432  {
433  cnt += BndExp[n]->GetExpSize();
434  }
435  }
436  }
437 
438 
440  {
441  // use static trip since cannot use UserDefinedTag for zero
442  // velocity and have time dependent conditions
443  static bool Setup = false;
444 
445  if(Setup == true)
446  {
447  return;
448  }
449  Setup = true;
450 
451  int i,n;
452 
454  BndConds(m_spacedim);
456  BndExp(m_spacedim);
457 
458 
459  for(i = 0; i < m_spacedim; ++i)
460  {
461  BndConds[i] = m_fields[m_velocity[i]]->GetBndConditions();
462  BndExp[i] = m_fields[m_velocity[i]]->GetBndCondExpansions();
463  }
464 
466 
467  int cnt;
468  int elmtid,nq, boundary;
469 
471  Array<OneD, NekDouble> Bphys,Bcoeffs;
472 
473  int fieldid = m_velocity[0];
474 
475  for(cnt = n = 0; n < BndConds[0].num_elements(); ++n)
476  {
477  if((BndConds[0][n]->GetBoundaryConditionType() == SpatialDomains::eDirichlet)&& (boost::iequals(BndConds[0][n]->GetUserDefined(),"ZeroNormalComponent")))
478  {
479  for(i = 0; i < BndExp[0][n]->GetExpSize(); ++i,cnt++)
480  {
481  elmtid = m_fieldsBCToElmtID[fieldid][cnt];
482  elmt = m_fields[0]->GetExp(elmtid);
483  boundary = m_fieldsBCToTraceID[fieldid][cnt];
484 
485  normals = elmt->GetSurfaceNormal(boundary);
486 
487  nq = BndExp[0][n]->GetExp(i)->GetTotPoints();
488  Array<OneD, NekDouble> normvel(nq,0.0);
489 
490  for(int k = 0; k < m_spacedim; ++k)
491  {
492  Bphys = BndExp[k][n]->UpdatePhys()+
493  BndExp[k][n]->GetPhys_Offset(i);
494  Bc = BndExp[k][n]->GetExp(i);
495  Vmath::Vvtvp(nq,normals[k],1,Bphys,1,normvel,1,
496  normvel,1);
497  }
498 
499  // negate normvel for next step
500  Vmath::Neg(nq,normvel,1);
501 
502  for(int k = 0; k < m_spacedim; ++k)
503  {
504  Bphys = BndExp[k][n]->UpdatePhys()+
505  BndExp[k][n]->GetPhys_Offset(i);
506  Bcoeffs = BndExp[k][n]->UpdateCoeffs()+
507  BndExp[k][n]->GetCoeff_Offset(i);
508  Bc = BndExp[k][n]->GetExp(i);
509  Vmath::Vvtvp(nq,normvel,1,normals[k],1,Bphys,1,
510  Bphys,1);
511  Bc->FwdTrans_BndConstrained(Bphys,Bcoeffs);
512  }
513  }
514  }
515  else
516  {
517  cnt += BndExp[0][n]->GetExpSize();
518  }
519  }
520  }
521 
522 
523  /**
524  * Add an additional forcing term programmatically.
525  */
527  {
528  m_forcing.push_back(pForce);
529  }
530 
531 
532  /**
533  * Decide if at a steady state if the discrerte L2 sum of the
534  * coefficients is the same as the previous step to within the
535  * tolerance m_steadyStateTol;
536  */
538  {
539  static NekDouble previousL2 = 0.0;
540  bool returnval = false;
541 
542  NekDouble L2 = 0.0;
543 
544  // calculate L2 discrete summation
545  int ncoeffs = m_fields[0]->GetNcoeffs();
546 
547  for(int i = 0; i < m_fields.num_elements(); ++i)
548  {
549  L2 += Vmath::Dot(ncoeffs,m_fields[i]->GetCoeffs(),1,m_fields[i]->GetCoeffs(),1);
550  }
551 
552  if(fabs(L2-previousL2) < ncoeffs*m_steadyStateTol)
553  {
554  returnval = true;
555  }
556 
557  previousL2 = L2;
558 
559  return returnval;
560  }
561 
562  /**
563  *
564  */
566  {
567  int n_vel = m_velocity.num_elements();
568  int n_element = m_fields[0]->GetExpSize();
569 
570  const Array<OneD, int> ExpOrder = GetNumExpModesPerExp();
571  Array<OneD, int> ExpOrderList (n_element, ExpOrder);
572 
573  const NekDouble cLambda = 0.2; // Spencer book pag. 317
574 
575  Array<OneD, NekDouble> cfl (n_element, 0.0);
576  Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
578 
579  if(m_HomogeneousType == eHomogeneous1D) // just do check on 2D info
580  {
581  velfields = Array<OneD, Array<OneD, NekDouble> >(2);
582 
583  for(int i = 0; i < 2; ++i)
584  {
585  velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
586  }
587  }
588  else
589  {
590  velfields = Array<OneD, Array<OneD, NekDouble> >(n_vel);
591 
592  for(int i = 0; i < n_vel; ++i)
593  {
594  velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
595  }
596  }
597 
598  stdVelocity = m_extrapolation->GetMaxStdVelocity(velfields);
599 
600  for(int el = 0; el < n_element; ++el)
601  {
602  cfl[el] = m_timestep*(stdVelocity[el] * cLambda *
603  (ExpOrder[el]-1) * (ExpOrder[el]-1));
604  }
605 
606  return cfl;
607  }
608 
609  /**
610  *
611  */
613  {
614  int n_element = m_fields[0]->GetExpSize();
615 
617 
618  elmtid = Vmath::Imax(n_element,cfl,1);
619  NekDouble CFL,CFL_loc;
620 
621  CFL = CFL_loc = cfl[elmtid];
622  m_comm->AllReduce(CFL,LibUtilities::ReduceMax);
623 
624  // unshuffle elmt id if data is not stored in consecutive order.
625  elmtid = m_fields[0]->GetExp(elmtid)->GetGeom()->GetGlobalID();
626  if(CFL != CFL_loc)
627  {
628  elmtid = -1;
629  }
630 
631  m_comm->AllReduce(elmtid,LibUtilities::ReduceMax);
632 
633  // express element id with respect to plane
635  {
636  elmtid = elmtid%m_fields[0]->GetPlane(0)->GetExpSize();
637  }
638  return CFL;
639  }
640 
641 
642  /**
643  * Perform the extrapolation.
644  */
646  {
647  m_extrapolation->SubStepSaveFields(step);
648  m_extrapolation->SubStepAdvance(m_intSoln,step,m_time);
649  return false;
650  }
651 
652 
653  /**
654  * Estimate CFL and perform steady-state check
655  */
657  {
658  if(m_cflsteps && !((step+1)%m_cflsteps))
659  {
660  int elmtid;
661  NekDouble cfl = GetCFLEstimate(elmtid);
662 
663  if(m_comm->GetRank() == 0)
664  {
665  cout << "CFL (zero plane): "<< cfl << " (in elmt "
666  << elmtid << ")" << endl;
667  }
668  }
669 
670  if(m_steadyStateSteps && step && (!((step+1)%m_steadyStateSteps)))
671  {
672  if(CalcSteadyState() == true)
673  {
674  cout << "Reached Steady State to tolerance "
675  << m_steadyStateTol << endl;
676  return true;
677  }
678  }
679 
680  return false;
681  }
682 } //end of namespace
683 
EquationType m_equationType
equation type;
bool m_singleMode
Flag to determine if single homogeneous mode is used.
void SetZeroNormalVelocity()
Set Normal Velocity Component to Zero.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void SetBoundaryConditions(NekDouble time)
time dependent boundary conditions updating
Array< OneD, Array< OneD, int > > m_fieldsBCToTraceID
Mapping from BCs to Elmt Edge IDs.
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.
virtual bool v_PreIntegrate(int step)
NekDouble m_time
Current time of simulation.
void SetRadiationBoundaryForcing(int fieldid)
Set Radiation forcing term.
NekDouble m_kinvis
Kinematic viscosity.
bool CalcSteadyState(void)
evaluate steady state
NekDouble m_timestep
Time step size.
bool m_halfMode
Flag to determine if half homogeneous mode is used.
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp()
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
ExtrapolateSharedPtr m_extrapolation
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
SOLVER_UTILS_EXPORT typedef boost::shared_ptr< Forcing > ForcingSharedPtr
A shared pointer to an EquationSystem object.
Definition: Forcing.h:51
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtr > Load(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
Definition: Forcing.cpp:84
virtual bool v_PostIntegrate(int step)
int m_nConvectiveFields
Number of fields to be convected;.
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
int Imax(int n, const T *x, const int incx)
Return the index of the maximum element in x.
Definition: Vmath.cpp:741
LibUtilities::CommSharedPtr m_comm
Communicator.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
Array< OneD, Array< OneD, NekDouble > > m_fieldsRadiationFactor
RHS Factor for Radiation Condition.
int m_cflsteps
dump cfl estimate
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
NekDouble Evaluate() const
Definition: Equation.h:102
int m_steadyStateSteps
Check for steady state at step interval.
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
Pointer to boundary conditions object.
int m_spacedim
Spatial dimension (>= expansion dim).
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
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:382
double NekDouble
virtual void v_NumericalFlux(Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
T Dot(int n, const T *w, const T *x)
vvtvp (vector times vector times vector): z = w*x*y
Definition: Vmath.cpp:900
SOLVER_UTILS_EXPORT int GetPhys_Offset(int n)
SOLVER_UTILS_EXPORT int GetNpoints()
void EvaluateAdvectionTerms(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
void AddForcing(const SolverUtils::ForcingSharedPtr &pForce)
Array< OneD, Array< OneD, int > > m_fieldsBCToElmtID
Mapping from BCs to Elmt IDs.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
const std::string kEquationTypeStr[]
int m_infosteps
Number of time steps between outputting status information.
NekDouble m_steadyStateTol
Tolerance to which steady state should be evaluated at.
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
virtual void v_GetFluxVector(const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
virtual void v_InitObject()
Init object for UnsteadySystem class.
virtual int v_GetForceDimension()=0
enum HomogeneousType m_HomogeneousType
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
NekDouble GetCFLEstimate(int &elmtid)
Array< OneD, NekDouble > GetElmtCFLVals(void)
IncNavierStokes(const LibUtilities::SessionReaderSharedPtr &pSession)
Constructor.