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