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  */
56  const LibUtilities::SessionReaderSharedPtr& pSession):
57  IncNavierStokes(pSession)
58  {
59 
60  }
61 
63  {
64  int n;
65 
67 
68  // Set m_pressure to point to last field of m_fields;
69  if (boost::iequals(m_session->GetVariable(m_fields.num_elements()-1), "p"))
70  {
71  m_nConvectiveFields = m_fields.num_elements()-1;
73  }
74  else
75  {
76  ASSERTL0(false,"Need to set up pressure field definition");
77  }
78 
80 
81  // Integrate only the convective fields
82  for (n = 0; n < m_nConvectiveFields; ++n)
83  {
84  m_intVariables.push_back(n);
85  }
86 
87  // Load parameters for Spectral Vanishing Viscosity
88  m_session->MatchSolverInfo("SpectralVanishingViscosity","True",m_useSpecVanVisc,false);
89  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
90  m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
91 
92  // Needs to be set outside of next if so that it is turned off by default
93  m_session->MatchSolverInfo("SpectralVanishingViscosityHomo1D","True",m_useHomo1DSpecVanVisc,false);
94 
96  {
97  ASSERTL0(m_nConvectiveFields > 2,"Expect to have three velcoity fields with homogenous expansion");
98 
99  if(m_useHomo1DSpecVanVisc == false)
100  {
101  m_session->MatchSolverInfo("SpectralVanishingViscosity","True",m_useHomo1DSpecVanVisc,false);
102  }
103 
105  {
106  Array<OneD, unsigned int> planes;
107  planes = m_fields[0]->GetZIDs();
108 
109  int num_planes = planes.num_elements();
110  Array<OneD, NekDouble> SVV(num_planes,0.0);
111  NekDouble fac;
112  int kmodes = m_fields[0]->GetHomogeneousBasis()->GetNumModes();
113  int pstart;
114 
115  pstart = m_sVVCutoffRatio*kmodes;
116 
117  for(n = 0; n < num_planes; ++n)
118  {
119  if(planes[n] > pstart)
120  {
121  fac = (NekDouble)((planes[n] - kmodes)*(planes[n] - kmodes))/
122  ((NekDouble)((planes[n] - pstart)*(planes[n] - pstart)));
123  SVV[n] = m_sVVDiffCoeff*exp(-fac)/m_kinvis;
124  }
125  }
126 
127  for(int i = 0; i < m_velocity.num_elements(); ++i)
128  {
129  m_fields[m_velocity[i]]->SetHomo1DSpecVanVisc(SVV);
130  }
131  }
132 
133  }
134 
135  m_session->MatchSolverInfo("SmoothAdvection", "True", m_SmoothAdvection, false);
136 
137  if(m_subSteppingScheme) // Substepping
138  {
140  "Projection must be set to Mixed_CG_Discontinuous for "
141  "substepping");
142  }
143  else // Standard velocity correction scheme
144  {
145  // set explicit time-intregration class operators
147  }
148  m_extrapolation->SubSteppingTimeIntegration(m_intScheme->GetIntegrationMethod(), m_intScheme);
149  m_extrapolation->GenerateHOPBCMap();
150 
151  // set implicit time-intregration class operators
153  }
154 
155  /**
156  * Destructor
157  */
159  {
160  }
161 
162  /**
163  *
164  */
166  {
168 
170  {
172  s, "Substepping", LibUtilities::TimeIntegrationMethodMap[
173  m_subStepIntegrationScheme->GetIntegrationMethod()]);
174  }
175 
176  string dealias = m_homogen_dealiasing ? "Homogeneous1D" : "";
177  if (m_advObject->GetSpecHPDealiasing())
178  {
179  dealias += (dealias == "" ? "" : " + ") + string("spectral/hp");
180  }
181  if (dealias != "")
182  {
183  SolverUtils::AddSummaryItem(s, "Dealiasing", dealias);
184  }
185 
186  string smoothing = m_useSpecVanVisc ? "spectral/hp" : "";
188  {
189  smoothing += (smoothing == "" ? "" : " + ") + string("Homogeneous1D");
190  }
191  if (smoothing != "")
192  {
194  s, "Smoothing", "SVV (" + smoothing + " SVV (cut-off = "
195  + boost::lexical_cast<string>(m_sVVCutoffRatio)
196  + ", diff coeff = "
197  + boost::lexical_cast<string>(m_sVVDiffCoeff)+")");
198  }
199  }
200 
201  /**
202  *
203  */
205  {
206 
208 
209  // Set up Field Meta Data for output files
210  m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
211  m_fieldMetaDataMap["TimeStep"] = boost::lexical_cast<std::string>(m_timestep);
212 
213  for(int i = 0; i < m_nConvectiveFields; ++i)
214  {
215  m_fields[i]->LocalToGlobal();
216  m_fields[i]->ImposeDirichletConditions(m_fields[i]->UpdateCoeffs());
217  m_fields[i]->GlobalToLocal();
218  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
219  m_fields[i]->UpdatePhys());
220  }
221  }
222 
223 
224  /**
225  *
226  */
228  {
229  int nfields = m_fields.num_elements() - 1;
230  for (int k=0 ; k < nfields; ++k)
231  {
232  //Backward Transformation in physical space for time evolution
233  m_fields[k]->BwdTrans_IterPerExp(m_fields[k]->GetCoeffs(),
234  m_fields[k]->UpdatePhys());
235  }
236  }
237 
238  /**
239  *
240  */
242  {
243 
244  int nfields = m_fields.num_elements() - 1;
245  for (int k=0 ; k < nfields; ++k)
246  {
247  //Forward Transformation in physical space for time evolution
248  m_fields[k]->FwdTrans_IterPerExp(m_fields[k]->GetPhys(),m_fields[k]->UpdateCoeffs());
249  }
250  }
251 
252  /**
253  *
254  */
256  {
257  int vVar = m_session->GetVariables().size();
258  Array<OneD, bool> vChecks(vVar, false);
259  vChecks[vVar-1] = true;
260  return vChecks;
261  }
262 
263  /**
264  *
265  */
267  {
268  return m_session->GetVariables().size() - 1;
269  }
270 
271  /**
272  * Explicit part of the method - Advection, Forcing + HOPBCs
273  */
275  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
276  Array<OneD, Array<OneD, NekDouble> > &outarray,
277  const NekDouble time)
278  {
279  // Evaluate convection terms
281  inarray, outarray, m_time);
282 
283  // Smooth advection
285  {
286  for(int i = 0; i < m_nConvectiveFields; ++i)
287  {
288  m_pressure->SmoothField(outarray[i]);
289  }
290  }
291 
292  // Add forcing terms
293  std::vector<SolverUtils::ForcingSharedPtr>::const_iterator x;
294  for (x = m_forcing.begin(); x != m_forcing.end(); ++x)
295  {
296  (*x)->Apply(m_fields, inarray, outarray, time);
297  }
298 
299  // Calculate High-Order pressure boundary conditions
300  m_extrapolation->EvaluatePressureBCs(inarray,outarray,m_kinvis);
301  }
302 
303  /**
304  * Implicit part of the method - Poisson + nConv*Helmholtz
305  */
307  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
308  Array<OneD, Array<OneD, NekDouble> > &outarray,
309  const NekDouble time,
310  const NekDouble aii_Dt)
311  {
312  int i;
313  int phystot = m_fields[0]->GetTotPoints();
314 
316 
317  Array<OneD, Array< OneD, NekDouble> > F(m_nConvectiveFields);
318  for(i = 0; i < m_nConvectiveFields; ++i)
319  {
320  F[i] = Array<OneD, NekDouble> (phystot);
321  }
322 
323  // Enforcing boundary conditions on all fields
324  SetBoundaryConditions(time);
325 
326  // Substep the pressure boundary condition
327  m_extrapolation->SubStepSetPressureBCs(inarray,aii_Dt,m_kinvis);
328 
329  // Set up forcing term and coefficients for pressure Poisson equation
330  SetUpPressureForcing(inarray, F, aii_Dt);
331  factors[StdRegions::eFactorLambda] = 0.0;
332 
333  // Solver Pressure Poisson Equation
334  m_pressure->HelmSolve(F[0], m_pressure->UpdateCoeffs(), NullFlagList,
335  factors);
336 
337  // Set up forcing term and coefficients for Helmholtz problems
338  SetUpViscousForcing(inarray, F, aii_Dt);
339  factors[StdRegions::eFactorLambda] = 1.0/aii_Dt/m_kinvis;
340  if(m_useSpecVanVisc)
341  {
344  }
345 
346  // Solve Helmholtz system and put in Physical space
347  for(i = 0; i < m_nConvectiveFields; ++i)
348  {
349  m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(),
350  NullFlagList, factors);
351  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),outarray[i]);
352  }
353  }
354 
355  /**
356  * Forcing term for Poisson solver solver
357  */
359  const Array<OneD, const Array<OneD, NekDouble> > &fields,
360  Array<OneD, Array<OneD, NekDouble> > &Forcing,
361  const NekDouble aii_Dt)
362  {
363  int i;
364  int physTot = m_fields[0]->GetTotPoints();
365  int nvel = m_velocity.num_elements();
366  Array<OneD, NekDouble> wk(physTot, 0.0);
367 
368  Vmath::Zero(physTot,Forcing[0],1);
369 
370  for(i = 0; i < nvel; ++i)
371  {
372  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[i],fields[i], wk);
373  Vmath::Vadd(physTot,wk,1,Forcing[0],1,Forcing[0],1);
374  }
375  Vmath::Smul(physTot,1.0/aii_Dt,Forcing[0],1,Forcing[0],1);
376  }
377 
378  /**
379  * Forcing term for Helmholtz solver
380  */
382  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
383  Array<OneD, Array<OneD, NekDouble> > &Forcing,
384  const NekDouble aii_Dt)
385  {
386  NekDouble aii_dtinv = 1.0/aii_Dt;
387  int phystot = m_fields[0]->GetTotPoints();
388 
389  // Grad p
390  m_pressure->BwdTrans(m_pressure->GetCoeffs(),m_pressure->UpdatePhys());
391 
392  int nvel = m_velocity.num_elements();
393  if(nvel == 2)
394  {
395  m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1]);
396  }
397  else
398  {
399  m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1],
400  Forcing[2]);
401  }
402 
403  // Subtract inarray/(aii_dt) and divide by kinvis. Kinvis will
404  // need to be updated for the convected fields.
405  for(int i = 0; i < nvel; ++i)
406  {
407  Blas::Daxpy(phystot,-aii_dtinv,inarray[i],1,Forcing[i],1);
408  Blas::Dscal(phystot,1.0/m_kinvis,&(Forcing[i])[0],1);
409  }
410  }
411 
412 } //end of namespace