Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
VelocityCorrectionScheme.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File VelocityCorrection.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: Velocity Correction Scheme for the Incompressible
33 // Navier Stokes equations
34 ///////////////////////////////////////////////////////////////////////////////
35 
38 #include <SolverUtils/Core/Misc.h>
39 
40 #include <boost/algorithm/string.hpp>
41 
42 namespace Nektar
43 {
46  "VelocityCorrectionScheme",
48 
49  /**
50  * Constructor. Creates ...
51  *
52  * \param
53  * \param
54  */
57  : UnsteadySystem(pSession),
58  IncNavierStokes(pSession),
59  m_varCoeffLap(StdRegions::NullVarCoeffMap)
60  {
61 
62  }
63 
65  {
66  int n;
67 
69  m_explicitDiffusion = false;
70 
71  // Set m_pressure to point to last field of m_fields;
72  if (boost::iequals(m_session->GetVariable(m_fields.num_elements()-1), "p"))
73  {
74  m_nConvectiveFields = m_fields.num_elements()-1;
76  }
77  else
78  {
79  ASSERTL0(false,"Need to set up pressure field definition");
80  }
81 
82  // creation of the extrapolation object
84  {
85  std::string vExtrapolation = "Standard";
86 
87  if (m_session->DefinesSolverInfo("Extrapolation"))
88  {
89  vExtrapolation = m_session->GetSolverInfo("Extrapolation");
90  }
91 
93  vExtrapolation,
94  m_session,
95  m_fields,
96  m_pressure,
97  m_velocity,
98  m_advObject);
99  }
100 
101  // Integrate only the convective fields
102  for (n = 0; n < m_nConvectiveFields; ++n)
103  {
104  m_intVariables.push_back(n);
105  }
106 
109 
110  // Load parameters for Spectral Vanishing Viscosity
111  m_session->MatchSolverInfo("SpectralVanishingViscosity","True",
112  m_useSpecVanVisc,false);
113  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
114  m_session->LoadParameter("SVVDiffCoeff", m_sVVDiffCoeff, 0.1);
115  // Needs to be set outside of next if so that it is turned off by default
116  m_session->MatchSolverInfo("SpectralVanishingViscosityHomo1D","True",
117  m_useHomo1DSpecVanVisc,false);
118 
119 
120  m_session->MatchSolverInfo("SPECTRALHPDEALIASING","True",
121  m_specHP_dealiasing,false);
122 
123 
125  {
126  ASSERTL0(m_nConvectiveFields > 2,"Expect to have three velocity fields with homogenous expansion");
127 
128  if(m_useHomo1DSpecVanVisc == false)
129  {
130  m_session->MatchSolverInfo("SpectralVanishingViscosity","True",m_useHomo1DSpecVanVisc,false);
131  }
132 
134  {
136  planes = m_fields[0]->GetZIDs();
137 
138  int num_planes = planes.num_elements();
139  Array<OneD, NekDouble> SVV(num_planes,0.0);
140  NekDouble fac;
141  int kmodes = m_fields[0]->GetHomogeneousBasis()->GetNumModes();
142  int pstart;
143 
144  pstart = m_sVVCutoffRatio*kmodes;
145 
146  for(n = 0; n < num_planes; ++n)
147  {
148  if(planes[n] > pstart)
149  {
150  fac = (NekDouble)((planes[n] - kmodes)*(planes[n] - kmodes))/
151  ((NekDouble)((planes[n] - pstart)*(planes[n] - pstart)));
152  SVV[n] = m_sVVDiffCoeff*exp(-fac)/m_kinvis;
153  }
154  }
155 
156  for(int i = 0; i < m_velocity.num_elements(); ++i)
157  {
158  m_fields[m_velocity[i]]->SetHomo1DSpecVanVisc(SVV);
159  }
160  }
161 
162  }
163 
164  m_session->MatchSolverInfo("SmoothAdvection", "True", m_SmoothAdvection, false);
165 
166  // set explicit time-intregration class operators
168 
169  m_extrapolation->SubSteppingTimeIntegration(m_intScheme->GetIntegrationMethod(), m_intScheme);
170  m_extrapolation->GenerateHOPBCMap();
171 
172  // set implicit time-intregration class operators
174  }
175 
176  /**
177  * Destructor
178  */
180  {
181  }
182 
183  /**
184  *
185  */
187  {
188  UnsteadySystem::v_GenerateSummary(s);
189 
190  if (m_extrapolation->GetSubStepIntegrationMethod() !=
192  {
193  SolverUtils::AddSummaryItem(s, "Substepping",
195  m_extrapolation->GetSubStepIntegrationMethod()]);
196  }
197 
198  string dealias = m_homogen_dealiasing ? "Homogeneous1D" : "";
200  {
201  dealias += (dealias == "" ? "" : " + ") + string("spectral/hp");
202  }
203  if (dealias != "")
204  {
205  SolverUtils::AddSummaryItem(s, "Dealiasing", dealias);
206  }
207 
208  string smoothing = m_useSpecVanVisc ? "spectral/hp" : "";
210  {
211  smoothing += (smoothing == "" ? "" : " + ") + string("Homogeneous1D");
212  }
213  if (smoothing != "")
214  {
216  s, "Smoothing", "SVV (" + smoothing + " SVV (cut-off = "
217  + boost::lexical_cast<string>(m_sVVCutoffRatio)
218  + ", diff coeff = "
219  + boost::lexical_cast<string>(m_sVVDiffCoeff)+")");
220  }
221  }
222 
223  /**
224  *
225  */
227  {
228 
229  UnsteadySystem::v_DoInitialise();
230 
231  // Set up Field Meta Data for output files
232  m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
233  m_fieldMetaDataMap["TimeStep"] = boost::lexical_cast<std::string>(m_timestep);
234 
235  // set boundary conditions here so that any normal component
236  // correction are imposed before they are imposed on intiial
237  // field below
239 
240  for(int i = 0; i < m_nConvectiveFields; ++i)
241  {
242  m_fields[i]->LocalToGlobal();
243  m_fields[i]->ImposeDirichletConditions(m_fields[i]->UpdateCoeffs());
244  m_fields[i]->GlobalToLocal();
245  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
246  m_fields[i]->UpdatePhys());
247  }
248  }
249 
250 
251  /**
252  *
253  */
255  {
256  int nfields = m_fields.num_elements() - 1;
257  for (int k=0 ; k < nfields; ++k)
258  {
259  //Backward Transformation in physical space for time evolution
260  m_fields[k]->BwdTrans_IterPerExp(m_fields[k]->GetCoeffs(),
261  m_fields[k]->UpdatePhys());
262  }
263  }
264 
265  /**
266  *
267  */
269  {
270 
271  int nfields = m_fields.num_elements() - 1;
272  for (int k=0 ; k < nfields; ++k)
273  {
274  //Forward Transformation in physical space for time evolution
275  m_fields[k]->FwdTrans_IterPerExp(m_fields[k]->GetPhys(),m_fields[k]->UpdateCoeffs());
276  }
277  }
278 
279  /**
280  *
281  */
283  {
284  int vVar = m_session->GetVariables().size();
285  Array<OneD, bool> vChecks(vVar, false);
286  vChecks[vVar-1] = true;
287  return vChecks;
288  }
289 
290  /**
291  *
292  */
294  {
295  return m_session->GetVariables().size() - 1;
296  }
297 
298  /**
299  * Explicit part of the method - Advection, Forcing + HOPBCs
300  */
302  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
303  Array<OneD, Array<OneD, NekDouble> > &outarray,
304  const NekDouble time)
305  {
306  EvaluateAdvectionTerms(inarray, outarray);
307 
308  // Smooth advection
310  {
311  for(int i = 0; i < m_nConvectiveFields; ++i)
312  {
313  m_pressure->SmoothField(outarray[i]);
314  }
315  }
316 
317  // Add forcing terms
318  std::vector<SolverUtils::ForcingSharedPtr>::const_iterator x;
319  for (x = m_forcing.begin(); x != m_forcing.end(); ++x)
320  {
321  (*x)->Apply(m_fields, inarray, outarray, time);
322  }
323 
324  // Calculate High-Order pressure boundary conditions
325  m_extrapolation->EvaluatePressureBCs(inarray,outarray,m_kinvis);
326  }
327 
328  /**
329  * Implicit part of the method - Poisson + nConv*Helmholtz
330  */
332  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
333  Array<OneD, Array<OneD, NekDouble> > &outarray,
334  const NekDouble time,
335  const NekDouble aii_Dt)
336  {
337  int i;
338  int phystot = m_fields[0]->GetTotPoints();
339 
341  for(i = 0; i < m_nConvectiveFields; ++i)
342  {
343  F[i] = Array<OneD, NekDouble> (phystot);
344  }
345 
346  // Enforcing boundary conditions on all fields
347  SetBoundaryConditions(time);
348 
349  // Substep the pressure boundary condition if using substepping
350  m_extrapolation->SubStepSetPressureBCs(inarray,aii_Dt,m_kinvis);
351 
352  // Set up forcing term for pressure Poisson equation
353  SetUpPressureForcing(inarray, F, aii_Dt);
354 
355  // Solve Pressure System
356  SolvePressure (F[0]);
357 
358  // Set up forcing term for Helmholtz problems
359  SetUpViscousForcing(inarray, F, aii_Dt);
360 
361  // Solve velocity system
362  SolveViscous( F, outarray, aii_Dt);
363  }
364 
365  /**
366  * Forcing term for Poisson solver solver
367  */
369  const Array<OneD, const Array<OneD, NekDouble> > &fields,
370  Array<OneD, Array<OneD, NekDouble> > &Forcing,
371  const NekDouble aii_Dt)
372  {
373  int i;
374  int physTot = m_fields[0]->GetTotPoints();
375  int nvel = m_velocity.num_elements();
376  Array<OneD, NekDouble> wk(physTot, 0.0);
377 
378  Vmath::Zero(physTot,Forcing[0],1);
379 
380  for(i = 0; i < nvel; ++i)
381  {
382  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[i],fields[i], wk);
383  Vmath::Vadd(physTot,wk,1,Forcing[0],1,Forcing[0],1);
384  }
385 
386  Vmath::Smul(physTot,1.0/aii_Dt,Forcing[0],1,Forcing[0],1);
387  }
388 
389  /**
390  * Forcing term for Helmholtz solver
391  */
393  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
394  Array<OneD, Array<OneD, NekDouble> > &Forcing,
395  const NekDouble aii_Dt)
396  {
397  NekDouble aii_dtinv = 1.0/aii_Dt;
398  int phystot = m_fields[0]->GetTotPoints();
399 
400  // Grad p
401  m_pressure->BwdTrans(m_pressure->GetCoeffs(),m_pressure->UpdatePhys());
402 
403  int nvel = m_velocity.num_elements();
404  if(nvel == 2)
405  {
406  m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1]);
407  }
408  else
409  {
410  m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0],
411  Forcing[1], Forcing[2]);
412  }
413 
414  // Subtract inarray/(aii_dt) and divide by kinvis. Kinvis will
415  // need to be updated for the convected fields.
416  for(int i = 0; i < nvel; ++i)
417  {
418  Blas::Daxpy(phystot,-aii_dtinv,inarray[i],1,Forcing[i],1);
419  Blas::Dscal(phystot,1.0/m_kinvis,&(Forcing[i])[0],1);
420  }
421  }
422 
423 
424  /**
425  * Solve pressure system
426  */
428  const Array<OneD, NekDouble> &Forcing)
429  {
431  // Setup coefficient for equation
432  factors[StdRegions::eFactorLambda] = 0.0;
433 
434  // Solver Pressure Poisson Equation
435  m_pressure->HelmSolve(Forcing, m_pressure->UpdateCoeffs(), NullFlagList,
436  factors);
437  }
438 
439  /**
440  * Solve velocity system
441  */
443  const Array<OneD, const Array<OneD, NekDouble> > &Forcing,
444  Array<OneD, Array<OneD, NekDouble> > &outarray,
445  const NekDouble aii_Dt)
446  {
448  // Setup coefficients for equation
449  factors[StdRegions::eFactorLambda] = 1.0/aii_Dt/m_kinvis;
450  if(m_useSpecVanVisc)
451  {
454  }
455 
456  // Solve Helmholtz system and put in Physical space
457  for(int i = 0; i < m_nConvectiveFields; ++i)
458  {
459  m_fields[i]->HelmSolve(Forcing[i], m_fields[i]->UpdateCoeffs(),
460  NullFlagList, factors);
461  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),outarray[i]);
462  }
463  }
464 
465 } //end of namespace
EquationType m_equationType
equation type;
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Creates an instance of this class.
virtual Array< OneD, bool > v_GetSystemSingularChecks()
virtual void v_DoInitialise(void)
Sets up initial conditions.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void SetBoundaryConditions(NekDouble time)
time dependent boundary conditions updating
void EvaluateAdvection_SetPressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
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
void SolveUnsteadyStokesSystem(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
NekDouble m_time
Current time of simulation.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
NekDouble m_kinvis
Kinematic viscosity.
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
void SolvePressure(const Array< OneD, NekDouble > &Forcing)
ExtrapolateFactory & GetExtrapolateFactory()
Definition: Extrapolate.cpp:48
VelocityCorrectionScheme(const LibUtilities::SessionReaderSharedPtr &pSession)
Constructor.
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
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
Print a summary of time stepping parameters.
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
bool m_useHomo1DSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
NekDouble m_sVVDiffCoeff
Diffusion coefficient of SVV modes.
ExtrapolateSharedPtr m_extrapolation
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
static std::string className
Name of class.
const char *const TimeIntegrationMethodMap[]
int m_nConvectiveFields
Number of fields to be convected;.
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
void SolveViscous(const Array< OneD, const Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
void SetUpViscousForcing(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
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:199
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
virtual void v_SetUpViscousForcing(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
void SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
bool m_SmoothAdvection
bool to identify if advection term smoothing is requested
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:50
virtual void v_TransCoeffToPhys(void)
Virtual function for transformation to physical space.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
double NekDouble
virtual void v_SolvePressure(const Array< OneD, NekDouble > &Forcing)
static const NekDouble kNekUnsetDouble
bool m_homogen_dealiasing
Flag to determine if dealiasing is used for homogeneous simulations.
EquationSystemFactory & GetEquationSystemFactory()
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
virtual void v_SolveViscous(const Array< OneD, const Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
void EvaluateAdvectionTerms(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
NekDouble m_sVVCutoffRatio
cutt off ratio from which to start decayhing modes
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
Array< OneD, NekDouble > m_saved_aii_Dt
Save aiiDt value to use as a trip to reset global matrix setup.
bool m_useSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
virtual void v_TransPhysToCoeff(void)
Virtual function for transformation to coefficient space.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
This class is the base class for Navier Stokes problems.
virtual void v_InitObject()
Init object for UnsteadySystem class.
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
virtual void v_SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
virtual void v_InitObject()
Init object for UnsteadySystem class.
virtual void v_EvaluateAdvection_SetPressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
enum HomogeneousType m_HomogeneousType
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:285
static FlagList NullFlagList
An empty flag list.
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:227
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215