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