Nektar++
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_subSteppingScheme(false),
60  m_SmoothAdvection(false),
61  m_steadyStateSteps(0)
62  {
63  }
64 
66  {
67  AdvectionSystem::v_InitObject();
68 
69  int i,j;
70  int numfields = m_fields.num_elements();
71  std::string velids[] = {"u","v","w"};
72 
73  // Set up Velocity field to point to the first m_expdim of m_fields;
75 
76  for(i = 0; i < m_spacedim; ++i)
77  {
78  for(j = 0; j < numfields; ++j)
79  {
80  std::string var = m_boundaryConditions->GetVariable(j);
81  if(boost::iequals(velids[i], var))
82  {
83  m_velocity[i] = j;
84  break;
85  }
86 
87  ASSERTL0(j != numfields, "Failed to find field: " + var);
88  }
89  }
90 
91  // Set up equation type enum using kEquationTypeStr
92  for(i = 0; i < (int) eEquationTypeSize; ++i)
93  {
94  bool match;
95  m_session->MatchSolverInfo("EQTYPE",kEquationTypeStr[i],match,false);
96  if(match)
97  {
99  break;
100  }
101  }
102  ASSERTL0(i != eEquationTypeSize,"EQTYPE not found in SOLVERINFO section");
103 
104  // This probably should to into specific implementations
105  // Equation specific Setups
106  switch(m_equationType)
107  {
108  case eSteadyStokes:
109  case eSteadyOseen:
110  case eSteadyNavierStokes:
111  case eSteadyLinearisedNS:
112  break;
114  case eUnsteadyStokes:
115  {
116  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
117  m_session->LoadParameter("IO_CFLSteps", m_cflsteps, 0);
118  m_session->LoadParameter("SteadyStateSteps", m_steadyStateSteps, 0);
119  m_session->LoadParameter("SteadyStateTol", m_steadyStateTol, 1e-6);
120 
121  // check to see if any user defined boundary condition is
122  // indeed implemented
123 
124  for(int n = 0; n < m_fields[0]->GetBndConditions().num_elements(); ++n)
125  {
126  std::string type =m_fields[0]->GetBndConditions()[n]->GetUserDefined();
127  if(!type.empty())
128  // Time Dependent Boundary Condition (if no user
129  // defined then this is empty)
130  ASSERTL0 (boost::iequals(type,"Wall_Forces") ||
131  boost::iequals(type,"TimeDependent") ||
132  boost::iequals(type,"MovingBody") ||
133  boost::iequals(type,"Radiation") ||
134  boost::iequals(type,"I") ||
135  boost::iequals(type,"HOutflow"),
136  "Unknown USERDEFINEDTYPE boundary condition");
137  }
138  }
139  break;
140  case eNoEquationType:
141  default:
142  ASSERTL0(false,"Unknown or undefined equation type");
143  }
144 
145  m_session->LoadParameter("Kinvis", m_kinvis);
146 
147  // Default advection type per solver
148  std::string vConvectiveType;
149  switch(m_equationType)
150  {
151  case eUnsteadyStokes:
152  case eSteadyLinearisedNS:
153  vConvectiveType = "NoAdvection";
154  break;
156  case eSteadyNavierStokes:
157  vConvectiveType = "Convective";
158  break;
160  vConvectiveType = "Linearised";
161  break;
162  default:
163  break;
164  }
165 
166  // Check if advection type overridden
167  if (m_session->DefinesTag("AdvectiveType") && m_equationType != eUnsteadyStokes &&
169  {
170  vConvectiveType = m_session->GetTag("AdvectiveType");
171  }
172 
173  // Initialise advection
174  m_advObject = SolverUtils::GetAdvectionFactory().CreateInstance(vConvectiveType, vConvectiveType);
175  m_advObject->InitObject( m_session, m_fields);
176 
177  // Forcing terms
180 
181  // check to see if any Robin boundary conditions and if so set
182  // up m_field to boundary condition maps;
186 
187  for (i = 0; i < m_fields.num_elements(); ++i)
188  {
189  bool Set = false;
190 
193  int radpts = 0;
194 
195  BndConds = m_fields[i]->GetBndConditions();
196  BndExp = m_fields[i]->GetBndCondExpansions();
197  for(int n = 0; n < BndConds.num_elements(); ++n)
198  {
199  if(boost::iequals(BndConds[n]->GetUserDefined(),"Radiation"))
200  {
201  ASSERTL0(BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eRobin,
202  "Radiation boundary condition must be of type Robin <R>");
203 
204  if(Set == false)
205  {
206  m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],m_fieldsBCToTraceID[i]);
207  Set = true;
208  }
209  radpts += BndExp[n]->GetTotPoints();
210  }
211  }
212 
213  m_fieldsRadiationFactor[i] = Array<OneD, NekDouble>(radpts);
214 
215  radpts = 0; // reset to use as a counter
216 
217  for(int n = 0; n < BndConds.num_elements(); ++n)
218  {
219  if(boost::iequals(BndConds[n]->GetUserDefined(),"Radiation"))
220  {
221 
222  int npoints = BndExp[n]->GetNpoints();
223  Array<OneD, NekDouble> x0(npoints,0.0);
224  Array<OneD, NekDouble> x1(npoints,0.0);
225  Array<OneD, NekDouble> x2(npoints,0.0);
226  Array<OneD, NekDouble> tmpArray;
227 
228  BndExp[n]->GetCoords(x0,x1,x2);
229 
230  LibUtilities::Equation coeff =
231  boost::static_pointer_cast<
233  >(BndConds[n])->m_robinPrimitiveCoeff;
234 
235  coeff.Evaluate(x0,x1,x2,m_time,
236  tmpArray = m_fieldsRadiationFactor[i]+ radpts);
237  //Vmath::Neg(npoints,tmpArray = m_fieldsRadiationFactor[i]+ radpts,1);
238  radpts += npoints;
239  }
240  }
241  }
242 
243  // Set up Field Meta Data for output files
244  m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
245  m_fieldMetaDataMap["TimeStep"] = boost::lexical_cast<std::string>(m_timestep);
246  }
247 
248  /**
249  * Destructor
250  */
252  {
253  }
254 
255 
256  /**
257  *
258  */
260  Array<OneD, Array<OneD, NekDouble> > &physfield,
262  {
263  ASSERTL1(flux.num_elements() == m_velocity.num_elements(),"Dimension of flux array and velocity array do not match");
264 
265  for(int j = 0; j < flux.num_elements(); ++j)
266  {
267  Vmath::Vmul(GetNpoints(), physfield[i], 1, m_fields[m_velocity[j]]->GetPhys(), 1, flux[j], 1);
268  }
269  }
270 
271  /**
272  * Calcualate numerical fluxes
273  */
275  Array<OneD, Array<OneD, NekDouble> > &numflux)
276  {
277  /// Counter variable
278  int i;
279 
280  /// Number of trace points
281  int nTracePts = GetTraceNpoints();
282 
283  /// Number of spatial dimensions
284  int nDimensions = m_spacedim;
285 
286  /// Forward state array
287  Array<OneD, NekDouble> Fwd(2*nTracePts);
288 
289  /// Backward state array
290  Array<OneD, NekDouble> Bwd = Fwd + nTracePts;
291 
292  /// Normal velocity array
293  Array<OneD, NekDouble> Vn (nTracePts, 0.0);
294 
295  // Extract velocity field along the trace space and multiply by trace normals
296  for(i = 0; i < nDimensions; ++i)
297  {
298  m_fields[0]->ExtractTracePhys(m_fields[m_velocity[i]]->GetPhys(), Fwd);
299  Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Fwd, 1, Vn, 1, Vn, 1);
300  }
301 
302  /// Compute the numerical fluxes at the trace points
303  for(i = 0; i < numflux.num_elements(); ++i)
304  {
305  /// Extract forwards/backwards trace spaces
306  m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
307 
308  /// Upwind between elements
309  m_fields[i]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux[i]);
310 
311  /// Calculate the numerical fluxes multipling Fwd or Bwd
312  /// by the normal advection velocity
313  Vmath::Vmul(nTracePts, numflux[i], 1, Vn, 1, numflux[i], 1);
314  }
315  }
316 
317  /**
318  * Evaluation -N(V) for all fields except pressure using m_velocity
319  */
321  Array<OneD, Array<OneD, NekDouble> > &outarray,
323  {
324  int i;
325  int nqtot = m_fields[0]->GetTotPoints();
326  int VelDim = m_velocity.num_elements();
327  Array<OneD, Array<OneD, NekDouble> > velocity(VelDim);
329 
330  for(i = 0; i < VelDim; ++i)
331  {
332  if(m_fields[i]->GetWaveSpace() && !m_singleMode && !m_halfMode)
333  {
334  velocity[i] = Array<OneD, NekDouble>(nqtot,0.0);
335  m_fields[i]->HomogeneousBwdTrans(inarray[m_velocity[i]],velocity[i]);
336  }
337  else
338  {
339  velocity[i] = inarray[m_velocity[i]];
340  }
341  }
342 
343  // Set up Derivative work space;
344  if(wk.num_elements())
345  {
346  ASSERTL0(wk.num_elements() >= nqtot*VelDim,
347  "Workspace is not sufficient");
348  Deriv = wk;
349  }
350  else
351  {
352  Deriv = Array<OneD, NekDouble> (nqtot*VelDim);
353  }
354 
356  velocity, inarray, outarray, m_time);
357  }
358 
359  /**
360  * Time dependent boundary conditions updating
361  */
363  {
364  int i, n;
365  std::string varName;
366  int nvariables = m_fields.num_elements();
367 
368  for (i = 0; i < nvariables; ++i)
369  {
370  for(n = 0; n < m_fields[i]->GetBndConditions().num_elements(); ++n)
371  {
372  if(m_fields[i]->GetBndConditions()[n]->IsTimeDependent() ||
373  m_fields[i]->GetBndConditions()[n]->GetUserDefined() ==
374  "MovingBody")
375  {
376  varName = m_session->GetVariable(i);
377  m_fields[i]->EvaluateBoundaryConditions(time, varName);
378  }
379 
380  }
381 
382  // Set Radiation conditions if required
384  }
385  }
386 
387  /**
388  * Probably should be pushed back into ContField?
389  */
391  {
392  int i,n;
393 
396 
397  BndConds = m_fields[fieldid]->GetBndConditions();
398  BndExp = m_fields[fieldid]->GetBndCondExpansions();
399 
402 
403  int cnt;
404  int elmtid,nq,offset, boundary;
405  Array<OneD, NekDouble> Bvals, U;
406  int cnt1 = 0;
407 
408  for(cnt = n = 0; n < BndConds.num_elements(); ++n)
409  {
410  std::string type = BndConds[n]->GetUserDefined();
411 
412  if((BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eRobin)&&(boost::iequals(type,"Radiation")))
413  {
414  for(i = 0; i < BndExp[n]->GetExpSize(); ++i,cnt++)
415  {
416  elmtid = m_fieldsBCToElmtID[fieldid][cnt];
417  elmt = m_fields[fieldid]->GetExp(elmtid);
418  offset = m_fields[fieldid]->GetPhys_Offset(elmtid);
419 
420  U = m_fields[fieldid]->UpdatePhys() + offset;
421  Bc = BndExp[n]->GetExp(i);
422 
423  boundary = m_fieldsBCToTraceID[fieldid][cnt];
424 
425  // Get edge values and put into ubc
426  nq = Bc->GetTotPoints();
427  Array<OneD, NekDouble> ubc(nq);
428  elmt->GetTracePhysVals(boundary,Bc,U,ubc);
429 
430  Vmath::Vmul(nq,&m_fieldsRadiationFactor[fieldid][cnt1 +
431  BndExp[n]->GetPhys_Offset(i)],1,&ubc[0],1,&ubc[0],1);
432 
433  Bvals = BndExp[n]->UpdateCoeffs()+BndExp[n]->GetCoeff_Offset(i);
434 
435  Bc->IProductWRTBase(ubc,Bvals);
436  }
437  cnt1 += BndExp[n]->GetTotPoints();
438  }
439  else
440  {
441  cnt += BndExp[n]->GetExpSize();
442  }
443  }
444  }
445 
446 
447  /**
448  * Add an additional forcing term programmatically.
449  */
451  {
452  m_forcing.push_back(pForce);
453  }
454 
455 
456  /**
457  * Decide if at a steady state if the discrerte L2 sum of the
458  * coefficients is the same as the previous step to within the
459  * tolerance m_steadyStateTol;
460  */
462  {
463  static NekDouble previousL2 = 0.0;
464  bool returnval = false;
465 
466  NekDouble L2 = 0.0;
467 
468  // calculate L2 discrete summation
469  int ncoeffs = m_fields[0]->GetNcoeffs();
470 
471  for(int i = 0; i < m_fields.num_elements(); ++i)
472  {
473  L2 += Vmath::Dot(ncoeffs,m_fields[i]->GetCoeffs(),1,m_fields[i]->GetCoeffs(),1);
474  }
475 
476  if(fabs(L2-previousL2) < ncoeffs*m_steadyStateTol)
477  {
478  returnval = true;
479  }
480 
481  previousL2 = L2;
482 
483  return returnval;
484  }
485 
486  /**
487  *
488  */
490  {
491  int n_vel = m_velocity.num_elements();
492  int n_element = m_fields[0]->GetExpSize();
493 
494  const Array<OneD, int> ExpOrder = GetNumExpModesPerExp();
495  Array<OneD, int> ExpOrderList (n_element, ExpOrder);
496 
497  const NekDouble cLambda = 0.2; // Spencer book pag. 317
498 
499  Array<OneD, NekDouble> cfl (n_element, 0.0);
500  Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
502 
503  if(m_HomogeneousType == eHomogeneous1D) // just do check on 2D info
504  {
505  velfields = Array<OneD, Array<OneD, NekDouble> >(2);
506 
507  for(int i = 0; i < 2; ++i)
508  {
509  velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
510  }
511  }
512  else
513  {
514  velfields = Array<OneD, Array<OneD, NekDouble> >(n_vel);
515 
516  for(int i = 0; i < n_vel; ++i)
517  {
518  velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
519  }
520  }
521 
522  stdVelocity = m_extrapolation->GetMaxStdVelocity(velfields);
523 
524  for(int el = 0; el < n_element; ++el)
525  {
526  cfl[el] = m_timestep*(stdVelocity[el] * cLambda *
527  (ExpOrder[el]-1) * (ExpOrder[el]-1));
528  }
529 
530  return cfl;
531  }
532 
533  /**
534  *
535  */
537  {
538  int n_element = m_fields[0]->GetExpSize();
539 
541 
542  elmtid = Vmath::Imax(n_element,cfl,1);
543  NekDouble CFL,CFL_loc;
544 
545  CFL = CFL_loc = cfl[elmtid];
546  m_comm->AllReduce(CFL,LibUtilities::ReduceMax);
547 
548  // unshuffle elmt id if data is not stored in consecutive order.
549  elmtid = m_fields[0]->GetExp(elmtid)->GetGeom()->GetGlobalID();
550  if(CFL != CFL_loc)
551  {
552  elmtid = -1;
553  }
554 
555  m_comm->AllReduce(elmtid,LibUtilities::ReduceMax);
556 
557  // express element id with respect to plane
559  {
560  elmtid = elmtid%m_fields[0]->GetPlane(0)->GetExpSize();
561  }
562  return CFL;
563  }
564 
565 
566  /**
567  * Perform the extrapolation.
568  */
570  {
571  m_extrapolation->SubStepSaveFields(step);
572  m_extrapolation->SubStepAdvance(m_intSoln,step,m_time);
573  return false;
574  }
575 
576 
577  /**
578  * Estimate CFL and perform steady-state check
579  */
581  {
582  if(m_cflsteps && !((step+1)%m_cflsteps))
583  {
584  int elmtid;
585  NekDouble cfl = GetCFLEstimate(elmtid);
586 
587  if(m_comm->GetRank() == 0)
588  {
589  cout << "CFL (zero plane): "<< cfl << " (in elmt "
590  << elmtid << ")" << endl;
591  }
592  }
593 
594  if(m_steadyStateSteps && step && (!((step+1)%m_steadyStateSteps)))
595  {
596  if(CalcSteadyState() == true)
597  {
598  cout << "Reached Steady State to tolerance "
599  << m_steadyStateTol << endl;
600  return true;
601  }
602  }
603 
604  return false;
605  }
606 } //end of namespace
607 
EquationType m_equationType
equation type;
bool m_singleMode
Flag to determine if single homogeneous mode is used.
#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:50
int Imax(int n, const T *x, const int incx)
Return the index of the maximum element in x.
Definition: Vmath.cpp:732
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
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:891
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.