DASpalartAllmaras.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6  This file is modified from OpenFOAM's source code
7  src/TurbulenceModels/turbulenceModels/RAS/SpalartAllmaras/SpalartAllmaras.C
8 
9  OpenFOAM: The Open Source CFD Toolbox
10 
11  Copyright (C): 2011-2016 OpenFOAM Foundation
12 
13  OpenFOAM License:
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "DASpalartAllmaras.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 defineTypeNameAndDebug(DASpalartAllmaras, 0);
38 addToRunTimeSelectionTable(DATurbulenceModel, DASpalartAllmaras, dictionary);
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42  const word modelType,
43  const fvMesh& mesh,
44  const DAOption& daOption)
45  : DATurbulenceModel(modelType, mesh, daOption),
46  // SA parameters
47  sigmaNut_(dimensioned<scalar>::lookupOrAddToDict(
48  "sigmaNut",
49  this->coeffDict_,
50  0.66666)),
51  kappa_(dimensioned<scalar>::lookupOrAddToDict(
52  "kappa",
53  this->coeffDict_,
54  0.41)),
55 
56  Cb1_(dimensioned<scalar>::lookupOrAddToDict(
57  "Cb1",
58  this->coeffDict_,
59  0.1355)),
60  Cb2_(dimensioned<scalar>::lookupOrAddToDict(
61  "Cb2",
62  this->coeffDict_,
63  0.622)),
64  Cw1_(Cb1_ / sqr(kappa_) + (1.0 + Cb2_) / sigmaNut_),
65  Cw2_(dimensioned<scalar>::lookupOrAddToDict(
66  "Cw2",
67  this->coeffDict_,
68  0.3)),
69  Cw3_(dimensioned<scalar>::lookupOrAddToDict(
70  "Cw3",
71  this->coeffDict_,
72  2.0)),
73  Cv1_(dimensioned<scalar>::lookupOrAddToDict(
74  "Cv1",
75  this->coeffDict_,
76  7.1)),
77  Cs_(dimensioned<scalar>::lookupOrAddToDict(
78  "Cs",
79  this->coeffDict_,
80  0.3)),
81  // Augmented variables
82  nuTilda_(const_cast<volScalarField&>(
83  mesh.thisDb().lookupObject<volScalarField>("nuTilda"))),
84  nuTildaRes_(
85  IOobject(
86  "nuTildaRes",
87  mesh.time().timeName(),
88  mesh,
89  IOobject::NO_READ,
90  IOobject::NO_WRITE),
91  mesh,
92 #ifdef CompressibleFlow
93  dimensionedScalar("nuTildaRes", dimensionSet(1, -1, -2, 0, 0, 0, 0), 0.0),
94 #endif
95 #ifdef IncompressibleFlow
96  dimensionedScalar("nuTildaRes", dimensionSet(0, 2, -2, 0, 0, 0, 0), 0.0),
97 #endif
98  zeroGradientFvPatchField<scalar>::typeName),
99  y_(mesh.thisDb().lookupObject<volScalarField>("yWall"))
100 {
101 
102  // initialize printInterval_ we need to check whether it is a steady state
103  // or unsteady primal solver
104  IOdictionary fvSchemes(
105  IOobject(
106  "fvSchemes",
107  mesh.time().system(),
108  mesh,
109  IOobject::MUST_READ,
110  IOobject::NO_WRITE,
111  false));
112  word ddtScheme = word(fvSchemes.subDict("ddtSchemes").lookup("default"));
113  if (ddtScheme == "steadyState")
114  {
116  daOption.getAllOptions().lookupOrDefault<label>("printInterval", 100);
117  }
118  else
119  {
121  daOption.getAllOptions().lookupOrDefault<label>("printIntervalUnsteady", 500);
122  }
123 }
124 
125 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
126 
127 // SA member functions. these functions are copied from
128 // src/TurbulenceModels/turbulenceModels/RAS/SpalartAllmaras/SpalartAllmaras.C
129 tmp<volScalarField> DASpalartAllmaras::chi() const
130 {
131  return nuTilda_ / this->nu();
132 }
133 
134 tmp<volScalarField> DASpalartAllmaras::fv1(
135  const volScalarField& chi) const
136 {
137  const volScalarField chi3(pow3(chi));
138  return chi3 / (chi3 + pow3(Cv1_));
139 }
140 
141 tmp<volScalarField> DASpalartAllmaras::fv2(
142  const volScalarField& chi,
143  const volScalarField& fv1) const
144 {
145  return 1.0 - chi / (1.0 + chi * fv1);
146 }
147 
148 tmp<volScalarField> DASpalartAllmaras::Stilda(
149  const volScalarField& chi,
150  const volScalarField& fv1) const
151 {
152  volScalarField Omega(::sqrt(2.0) * mag(skew(fvc::grad(U_))));
153 
154  return (
155  max(
156  Omega
157  + fv2(chi, fv1) * nuTilda_ / sqr(kappa_ * y_),
158  Cs_ * Omega));
159 }
160 
161 tmp<volScalarField> DASpalartAllmaras::fw(
162  const volScalarField& Stilda) const
163 {
164  volScalarField r(
165  min(
166  nuTilda_
167  / (max(
168  Stilda,
169  dimensionedScalar("SMALL", Stilda.dimensions(), SMALL))
170  * sqr(kappa_ * y_)),
171  scalar(10.0)));
172  r.boundaryFieldRef() == 0.0;
173 
174  const volScalarField g(r + Cw2_ * (pow6(r) - r));
175 
176  return g * pow((1.0 + pow6(Cw3_)) / (pow6(g) + pow6(Cw3_)), 1.0 / 6.0);
177 }
178 
179 tmp<volScalarField> DASpalartAllmaras::DnuTildaEff() const
180 {
181  return tmp<volScalarField>(
182  new volScalarField("DnuTildaEff", (nuTilda_ + this->nu()) / sigmaNut_));
183 }
184 
185 // Augmented functions
186 void DASpalartAllmaras::correctModelStates(wordList& modelStates) const
187 {
188  /*
189  Description:
190  Update the name in modelStates based on the selected physical model at runtime
191 
192  Example:
193  In DAStateInfo, if the modelStates reads:
194 
195  modelStates = {"nut"}
196 
197  then for the SA model, calling correctModelStates(modelStates) will give:
198 
199  modelStates={"nuTilda"}
200 
201  while calling correctModelStates(modelStates) for the SST model will give
202 
203  modelStates={"k","omega"}
204 
205  We don't udpate the names for the radiation model becasue users are
206  supposed to set modelStates={"G"}
207  */
208 
209  // replace nut with nuTilda
210  forAll(modelStates, idxI)
211  {
212  word stateName = modelStates[idxI];
213  if (stateName == "nut")
214  {
215  modelStates[idxI] = "nuTilda";
216  }
217  }
218 }
219 
221 {
222  /*
223  Description:
224  Update nut based on other turbulence variables and update the BCs
225  Also update alphat if is present
226  */
227 
228  const volScalarField chi(this->chi());
229  const volScalarField fv1(this->fv1(chi));
230  nut_ = nuTilda_ * fv1;
231 
232  nut_.correctBoundaryConditions();
233 
234  // this is basically BasicTurbulenceModel::correctNut();
235  this->correctAlphat();
236 
237  return;
238 }
239 
241 {
242  /*
243  Description:
244  Update turbulence variable boundary values
245  */
246 
247  // correct the BCs for the perturbed fields
248  nuTilda_.correctBoundaryConditions();
249 }
250 
252 {
253  /*
254  Description:
255  Update nut based on nuTilda. Note: we need to update nut and its BC since we
256  may have perturbed other turbulence vars that affect the nut values
257  */
258 
259  this->correctNut();
260 }
261 
262 void DASpalartAllmaras::correctStateResidualModelCon(List<List<word>>& stateCon) const
263 {
264  /*
265  Description:
266  Update the original variable connectivity for the adjoint state
267  residuals in stateCon. Basically, we modify/add state variables based on the
268  original model variables defined in stateCon.
269 
270  Input:
271 
272  stateResCon: the connectivity levels for a state residual, defined in Foam::DAJacCon
273 
274  Example:
275  If stateCon reads:
276  stateCon=
277  {
278  {"U", "p", "nut"},
279  {"p"}
280  }
281 
282  For the SA turbulence model, calling this function for will get a new stateCon
283  stateCon=
284  {
285  {"U", "p", "nuTilda"},
286  {"p"}
287  }
288 
289  For the SST turbulence model, calling this function will give
290  stateCon=
291  {
292  {"U", "p", "k", "omega"},
293  {"p", "U"}
294  }
295  ***NOTE***: we add a extra level of U connectivity because nut is
296  related to grad(U), k, and omega in SST!
297  */
298 
299  forAll(stateCon, idxI)
300  {
301  forAll(stateCon[idxI], idxJ)
302  {
303  word conStateName = stateCon[idxI][idxJ];
304  if (conStateName == "nut")
305  {
306  stateCon[idxI][idxJ] = "nuTilda";
307  }
308  }
309  }
310 }
311 
312 void DASpalartAllmaras::addModelResidualCon(HashTable<List<List<word>>>& allCon) const
313 {
314  /*
315  Description:
316  Add the connectivity levels for all physical model residuals to allCon
317 
318  Input:
319  allCon: the connectivity levels for all state residual, defined in DAJacCon
320 
321  Example:
322  If stateCon reads:
323  allCon=
324  {
325  "URes":
326  {
327  {"U", "p", "nut"},
328  {"p"}
329  }
330  }
331 
332  For the SA turbulence model, calling this function for will get a new stateCon,
333  something like this:
334  allCon=
335  {
336  "URes":
337  {
338  {"U", "p", "nuTilda"},
339  {"p"}
340  },
341  "nuTildaRes":
342  {
343  {"U", "phi", "nuTilda"},
344  {"U"}
345  }
346  }
347 
348  */
349 
350  word pName;
351 
352  if (mesh_.thisDb().foundObject<volScalarField>("p"))
353  {
354  pName = "p";
355  }
356  else if (mesh_.thisDb().foundObject<volScalarField>("p_rgh"))
357  {
358  pName = "p_rgh";
359  }
360  else
361  {
362  FatalErrorIn(
363  "Neither p nor p_rgh was found in mesh.thisDb()!"
364  "addModelResidualCon failed to setup turbulence residuals!")
365  << exit(FatalError);
366  }
367 
368  // NOTE: for compressible flow, it depends on rho so we need to add T and p
369 #ifdef IncompressibleFlow
370  allCon.set(
371  "nuTildaRes",
372  {
373  {"U", "nuTilda", "phi"}, // lv0
374  {"U", "nuTilda"}, // lv1
375  {"nuTilda"} // lv2
376  });
377 #endif
378 
379 #ifdef CompressibleFlow
380  allCon.set(
381  "nuTildaRes",
382  {
383  {"U", "T", pName, "nuTilda", "phi"}, // lv0
384  {"U", "T", pName, "nuTilda"}, // lv1
385  {"T", pName, "nuTilda"} // lv2
386  });
387 #endif
388 }
389 
391 {
392  /*
393  Descroption:
394  Solve the residual equations and update the state. This function will be called
395  by the DASolver. It is needed because we want to control the output frequency
396  of the residual convergence every 100 steps. If using the correct from turbulence
397  it will output residual every step which will be too much of information.
398  */
399 
400  // We set the flag solveTurbState_ to 1 such that in the calcResiduals function
401  // we will solve and update nuTilda
402  solveTurbState_ = 1;
403  dictionary dummyOptions;
404  this->calcResiduals(dummyOptions);
405  // after it, we reset solveTurbState_ = 0 such that calcResiduals will not
406  // update nuTilda when calling from the adjoint class, i.e., solveAdjoint from DASolver.
407  solveTurbState_ = 0;
408 }
409 
410 void DASpalartAllmaras::calcResiduals(const dictionary& options)
411 {
412  /*
413  Descroption:
414  If solveTurbState_ == 1, this function solve and update nuTilda, and
415  is the same as calling turbulence.correct(). If solveTurbState_ == 0,
416  this function compute residuals for turbulence variables, e.g., nuTildaRes_
417 
418  Input:
419  options.isPC: 1 means computing residuals for preconditioner matrix.
420  This essentially use the first order scheme for div(phi,nuTilda)
421 
422  p_, U_, phi_, etc: State variables in OpenFOAM
423 
424  Output:
425  nuTildaRes_: If solveTurbState_ == 0, update the residual field variable
426 
427  nuTilda_: If solveTurbState_ == 1, update nuTilda
428  */
429 
430  // Copy and modify based on the "correct" function
431 
432  label printToScreen = this->isPrintTime(mesh_.time(), printInterval_);
433 
434  word divNuTildaScheme = "div(phi,nuTilda)";
435 
436  label isPC = 0;
437 
438  if (!solveTurbState_)
439  {
440  isPC = options.getLabel("isPC");
441 
442  if (isPC)
443  {
444  divNuTildaScheme = "div(pc)";
445  }
446  }
447 
448  const volScalarField chi(this->chi());
449  const volScalarField fv1(this->fv1(chi));
450 
451  const volScalarField Stilda(this->Stilda(chi, fv1));
452 
453  tmp<fvScalarMatrix> nuTildaEqn(
454  fvm::ddt(phase_, rho_, nuTilda_)
455  + fvm::div(phaseRhoPhi_, nuTilda_, divNuTildaScheme)
456  - fvm::laplacian(phase_ * rho_ * DnuTildaEff(), nuTilda_)
457  - Cb2_ / sigmaNut_ * phase_ * rho_ * magSqr(fvc::grad(nuTilda_))
458  == Cb1_ * phase_ * rho_ * Stilda * nuTilda_
459  - fvm::Sp(Cw1_ * phase_ * rho_ * fw(Stilda) * nuTilda_ / sqr(y_), nuTilda_));
460 
461  nuTildaEqn.ref().relax();
462 
463  if (solveTurbState_)
464  {
465 
466  // get the solver performance info such as initial
467  // and final residuals
468  SolverPerformance<scalar> solverNuTilda = solve(nuTildaEqn);
469  if (printToScreen)
470  {
471  Info << "nuTilda Initial residual: " << solverNuTilda.initialResidual() << endl
472  << " Final residual: " << solverNuTilda.finalResidual() << endl;
473  }
474 
475  DAUtility::boundVar(allOptions_, nuTilda_, printToScreen);
476  nuTilda_.correctBoundaryConditions();
477 
478  // ***************** NOTE*****************
479  // In the original SA, it is correctNut(fv1) and fv1 is not
480  // updated based on the latest nuTilda. We use correctNut which
481  // recompute fv1 with the latest nuTilda
482  this->correctNut();
483  }
484  else
485  {
486  // calculate residuals
487  nuTildaRes_ = nuTildaEqn.ref() & nuTilda_;
488  // need to normalize residuals
489  normalizeResiduals(nuTildaRes);
490  }
491 
492  return;
493 }
494 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
495 
496 } // End namespace Foam
497 
498 // ************************************************************************* //
Foam::DASpalartAllmaras::sigmaNut_
dimensionedScalar sigmaNut_
Definition: DASpalartAllmaras.H:55
forAll
forAll(pseudoP.boundaryField(), patchI)
Definition: solvePseudoPEqn.H:10
DASpalartAllmaras.H
Foam::DASpalartAllmaras::Cw1_
dimensionedScalar Cw1_
Definition: DASpalartAllmaras.H:59
Foam::DASpalartAllmaras::DASpalartAllmaras
DASpalartAllmaras(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DASpalartAllmaras.C:41
Foam::DASpalartAllmaras::solveTurbState_
label solveTurbState_
whether to solve for turb states
Definition: DASpalartAllmaras.H:94
Foam::DASpalartAllmaras::Cv1_
dimensionedScalar Cv1_
Definition: DASpalartAllmaras.H:62
Foam::DAOption
Definition: DAOption.H:29
Foam::DASpalartAllmaras::correct
virtual void correct()
solve the residual equations and update the state
Definition: DASpalartAllmaras.C:390
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DASpalartAllmaras::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DASpalartAllmaras.C:262
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
Foam::DASpalartAllmaras::Cb1_
dimensionedScalar Cb1_
Definition: DASpalartAllmaras.H:57
Foam::DATurbulenceModel::phaseRhoPhi_
surfaceScalarField & phaseRhoPhi_
phase*phi*density field
Definition: DATurbulenceModel.H:95
Foam::DATurbulenceModel::allOptions_
const dictionary & allOptions_
all DAFoam options
Definition: DATurbulenceModel.H:80
Foam::DASpalartAllmaras::Cw2_
dimensionedScalar Cw2_
Definition: DASpalartAllmaras.H:60
Foam::DASpalartAllmaras::kappa_
dimensionedScalar kappa_
Definition: DASpalartAllmaras.H:56
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAUtility::boundVar
static void boundVar(const dictionary &allOptions, volScalarField &var, const label printToScreen)
bound a volScalar variable based on parametes defined in DAOption::allOptions_
Definition: DAUtility.C:419
Foam::DASpalartAllmaras::printInterval_
label printInterval_
time step interval to print residual
Definition: DASpalartAllmaras.H:97
Foam::DATurbulenceModel::U_
volVectorField & U_
velocity
Definition: DATurbulenceModel.H:86
Foam::DATurbulenceModel::nut_
volScalarField & nut_
turbulence viscosity
Definition: DATurbulenceModel.H:83
Foam::DASpalartAllmaras::Cs_
dimensionedScalar Cs_
Definition: DASpalartAllmaras.H:63
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DASpalartAllmaras::chi
tmp< volScalarField > chi() const
Definition: DASpalartAllmaras.C:129
Foam::DASpalartAllmaras::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute the turbulence residuals
Definition: DASpalartAllmaras.C:410
Foam::DASpalartAllmaras::Stilda
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1) const
Definition: DASpalartAllmaras.C:148
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:262
Foam::DASpalartAllmaras::nuTildaRes_
volScalarField nuTildaRes_
Definition: DASpalartAllmaras.H:87
Foam::DASpalartAllmaras::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DASpalartAllmaras.C:312
Foam::DATurbulenceModel::correctAlphat
void correctAlphat()
update alphat
Definition: DATurbulenceModel.C:206
Foam::DASpalartAllmaras::nuTilda_
volScalarField & nuTilda_
Definition: DASpalartAllmaras.H:86
solve
pseudoPEqn solve(solverDictP_)
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::DASpalartAllmaras::Cw3_
dimensionedScalar Cw3_
Definition: DASpalartAllmaras.H:61
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::DASpalartAllmaras::correctBoundaryConditions
virtual void correctBoundaryConditions()
update turbulence variable boundary values
Definition: DASpalartAllmaras.C:240
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
Foam::DATurbulenceModel::mesh_
const fvMesh & mesh_
fvMesh
Definition: DATurbulenceModel.H:74
Foam::DATurbulenceModel::isPrintTime
label isPrintTime(const Time &runTime, const label printInterval) const
Definition: DATurbulenceModel.C:464
Foam::DASpalartAllmaras::correctNut
virtual void correctNut()
update nut based on other turbulence variables and update the BCs
Definition: DASpalartAllmaras.C:220
Foam::DASpalartAllmaras::fv1
tmp< volScalarField > fv1(const volScalarField &chi) const
Definition: DASpalartAllmaras.C:134
Foam::DASpalartAllmaras::DnuTildaEff
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
Definition: DASpalartAllmaras.C:179
Foam::DASpalartAllmaras::correctModelStates
virtual void correctModelStates(wordList &modelStates) const
update the turbulence state for DAStateInfo::regStates_
Definition: DASpalartAllmaras.C:186
Foam::DASpalartAllmaras::Cb2_
dimensionedScalar Cb2_
Definition: DASpalartAllmaras.H:58
Foam::DASpalartAllmaras::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DASpalartAllmaras.C:251
Foam::DASpalartAllmaras::fv2
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
Definition: DASpalartAllmaras.C:141
Foam::DASpalartAllmaras::fw
tmp< volScalarField > fw(const volScalarField &Stilda) const
Definition: DASpalartAllmaras.C:161
Foam::DATurbulenceModel::phase_
volScalarField phase_
phase field
Definition: DATurbulenceModel.H:92
Foam::DASpalartAllmaras::y_
const volScalarField & y_
3D wall distance
Definition: DASpalartAllmaras.H:91