DASpalartAllmarasFv3FieldInversion.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 
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 defineTypeNameAndDebug(DASpalartAllmarasFv3FieldInversion, 0);
38 addToRunTimeSelectionTable(DATurbulenceModel, DASpalartAllmarasFv3FieldInversion, 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  Cv2_(dimensioned<scalar>::lookupOrAddToDict(
78  "Cv2",
79  this->coeffDict_,
80  5.0)),
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  nuTildaResPartDeriv_(
100  IOobject(
101  "nuTildaResPartDeriv",
102  mesh.time().timeName(),
103  mesh,
104  IOobject::NO_READ,
105  IOobject::NO_WRITE),
106  nuTildaRes_),
107  betaFieldInversion_(const_cast<volScalarField&>(
108  mesh.thisDb().lookupObject<volScalarField>("betaFieldInversion"))),
109  UData_(const_cast<volVectorField&>(
110  mesh.thisDb().lookupObject<volVectorField>("UData"))),
111  surfaceFriction_(const_cast<volScalarField&>(
112  mesh.thisDb().lookupObject<volScalarField>("surfaceFriction"))),
113  surfaceFrictionData_(const_cast<volScalarField&>(
114  mesh.thisDb().lookupObject<volScalarField>("surfaceFrictionData"))),
115  pData_(const_cast<volScalarField&>(
116  mesh.thisDb().lookupObject<volScalarField>("pData"))),
117  USingleComponentData_(const_cast<volScalarField&>(
118  mesh.thisDb().lookupObject<volScalarField>("USingleComponentData"))),
119  y_(mesh.thisDb().lookupObject<volScalarField>("yWall"))
120 {
121 
122  // initialize printInterval_ we need to check whether it is a steady state
123  // or unsteady primal solver
124  IOdictionary fvSchemes(
125  IOobject(
126  "fvSchemes",
127  mesh.time().system(),
128  mesh,
129  IOobject::MUST_READ,
130  IOobject::NO_WRITE,
131  false));
132  word ddtScheme = word(fvSchemes.subDict("ddtSchemes").lookup("default"));
133  if (ddtScheme == "steadyState")
134  {
136  daOption.getAllOptions().lookupOrDefault<label>("printInterval", 100);
137  }
138  else
139  {
141  daOption.getAllOptions().lookupOrDefault<label>("printIntervalUnsteady", 500);
142  }
143 }
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 // SA member functions. these functions are copied from
148 tmp<volScalarField> DASpalartAllmarasFv3FieldInversion::chi() const
149 {
150  return nuTilda_ / this->nu();
151 }
152 
154  const volScalarField& chi) const
155 {
156  const volScalarField chi3(pow3(chi));
157  return chi3 / (chi3 + pow3(Cv1_));
158 }
159 
161  const volScalarField& chi,
162  const volScalarField& fv1) const
163 {
164  return 1.0 / pow3(scalar(1) + chi / Cv2_);
165 }
166 
168  const volScalarField& chi,
169  const volScalarField& fv1) const
170 {
171 
172  const volScalarField chiByCv2((1 / Cv2_) * chi);
173 
174  return (scalar(1) + chi * fv1)
175  * (1 / Cv2_)
176  * (3 * (scalar(1) + chiByCv2) + sqr(chiByCv2))
177  / pow3(scalar(1) + chiByCv2);
178 }
179 
181  const volScalarField& Stilda) const
182 {
183  volScalarField r(
184  min(
185  nuTilda_
186  / (max(
187  Stilda,
188  dimensionedScalar("SMALL", Stilda.dimensions(), SMALL))
189  * sqr(kappa_ * y_)),
190  scalar(10.0)));
191  r.boundaryFieldRef() == 0.0;
192 
193  const volScalarField g(r + Cw2_ * (pow6(r) - r));
194 
195  return g * pow((1.0 + pow6(Cw3_)) / (pow6(g) + pow6(Cw3_)), 1.0 / 6.0);
196 }
197 
199 {
200  return tmp<volScalarField>(
201  new volScalarField("DnuTildaEff", (nuTilda_ + this->nu()) / sigmaNut_));
202 }
203 
204 // Augmented functions
206 {
207  /*
208  Description:
209  Update the name in modelStates based on the selected physical model at runtime
210 
211  Example:
212  In DAStateInfo, if the modelStates reads:
213 
214  modelStates = {"nut"}
215 
216  then for the SA model, calling correctModelStates(modelStates) will give:
217 
218  modelStates={"nuTilda"}
219 
220  while calling correctModelStates(modelStates) for the SST model will give
221 
222  modelStates={"k","omega"}
223 
224  We don't udpate the names for the radiation model becasue users are
225  supposed to set modelStates={"G"}
226  */
227 
228  // replace nut with nuTilda
229  forAll(modelStates, idxI)
230  {
231  word stateName = modelStates[idxI];
232  if (stateName == "nut")
233  {
234  modelStates[idxI] = "nuTilda";
235  }
236  }
237 }
238 
240 {
241  /*
242  Description:
243  Update nut based on other turbulence variables and update the BCs
244  Also update alphat if is present
245  */
246 
247  const volScalarField chi(this->chi());
248  const volScalarField fv1(this->fv1(chi));
249  nut_ = nuTilda_ * fv1;
250 
251  nut_.correctBoundaryConditions();
252 
253  // this is basically BasicTurbulenceModel::correctNut();
254  this->correctAlphat();
255 
256  return;
257 }
258 
260 {
261  /*
262  Description:
263  Update turbulence variable boundary values
264  */
265 
266  // correct the BCs for the perturbed fields
267  nuTilda_.correctBoundaryConditions();
268 }
269 
271 {
272  /*
273  Description:
274  Update nut based on nuTilda. Note: we need to update nut and its BC since we
275  may have perturbed other turbulence vars that affect the nut values
276  */
277 
278  this->correctNut();
279 }
280 
282 {
283  /*
284  Description:
285  Update the original variable connectivity for the adjoint state
286  residuals in stateCon. Basically, we modify/add state variables based on the
287  original model variables defined in stateCon.
288 
289  Input:
290 
291  stateResCon: the connectivity levels for a state residual, defined in Foam::DAJacCon
292 
293  Example:
294  If stateCon reads:
295  stateCon=
296  {
297  {"U", "p", "nut"},
298  {"p"}
299  }
300 
301  For the SA turbulence model, calling this function for will get a new stateCon
302  stateCon=
303  {
304  {"U", "p", "nuTilda"},
305  {"p"}
306  }
307 
308  For the SST turbulence model, calling this function will give
309  stateCon=
310  {
311  {"U", "p", "k", "omega"},
312  {"p", "U"}
313  }
314  ***NOTE***: we add a extra level of U connectivity because nut is
315  related to grad(U), k, and omega in SST!
316  */
317 
318  forAll(stateCon, idxI)
319  {
320  forAll(stateCon[idxI], idxJ)
321  {
322  word conStateName = stateCon[idxI][idxJ];
323  if (conStateName == "nut")
324  {
325  stateCon[idxI][idxJ] = "nuTilda";
326  }
327  }
328  }
329 }
330 
331 void DASpalartAllmarasFv3FieldInversion::addModelResidualCon(HashTable<List<List<word>>>& allCon) const
332 {
333  /*
334  Description:
335  Add the connectivity levels for all physical model residuals to allCon
336 
337  Input:
338  allCon: the connectivity levels for all state residual, defined in DAJacCon
339 
340  Example:
341  If stateCon reads:
342  allCon=
343  {
344  "URes":
345  {
346  {"U", "p", "nut"},
347  {"p"}
348  }
349  }
350 
351  For the SA turbulence model, calling this function for will get a new stateCon,
352  something like this:
353  allCon=
354  {
355  "URes":
356  {
357  {"U", "p", "nuTilda"},
358  {"p"}
359  },
360  "nuTildaRes":
361  {
362  {"U", "phi", "nuTilda"},
363  {"U"}
364  }
365  }
366 
367  */
368 
369  word pName;
370 
371  if (mesh_.thisDb().foundObject<volScalarField>("p"))
372  {
373  pName = "p";
374  }
375  else if (mesh_.thisDb().foundObject<volScalarField>("p_rgh"))
376  {
377  pName = "p_rgh";
378  }
379  else
380  {
381  FatalErrorIn(
382  "Neither p nor p_rgh was found in mesh.thisDb()!"
383  "addModelResidualCon failed to setup turbulence residuals!")
384  << exit(FatalError);
385  }
386 
387  // NOTE: for compressible flow, it depends on rho so we need to add T and p
388 #ifdef IncompressibleFlow
389  allCon.set(
390  "nuTildaRes",
391  {
392  {"U", "nuTilda", "phi"}, // lv0
393  {"U", "nuTilda"}, // lv1
394  {"nuTilda"} // lv2
395  });
396 #endif
397 
398 #ifdef CompressibleFlow
399  allCon.set(
400  "nuTildaRes",
401  {
402  {"U", "T", pName, "nuTilda", "phi"}, // lv0
403  {"U", "T", pName, "nuTilda"}, // lv1
404  {"T", pName, "nuTilda"} // lv2
405  });
406 #endif
407 }
408 
410 {
411  /*
412  Descroption:
413  Solve the residual equations and update the state. This function will be called
414  by the DASolver. It is needed because we want to control the output frequency
415  of the residual convergence every 100 steps. If using the correct from turbulence
416  it will output residual every step which will be too much of information.
417  */
418 
419  // We set the flag solveTurbState_ to 1 such that in the calcResiduals function
420  // we will solve and update nuTilda
421  solveTurbState_ = 1;
422  dictionary dummyOptions;
423  this->calcResiduals(dummyOptions);
424  // after it, we reset solveTurbState_ = 0 such that calcResiduals will not
425  // update nuTilda when calling from the adjoint class, i.e., solveAdjoint from DASolver.
426  solveTurbState_ = 0;
427 }
428 
430 {
431  /*
432  Descroption:
433  If solveTurbState_ == 1, this function solve and update nuTilda, and
434  is the same as calling turbulence.correct(). If solveTurbState_ == 0,
435  this function compute residuals for turbulence variables, e.g., nuTildaRes_
436 
437  Input:
438  options.isPC: 1 means computing residuals for preconditioner matrix.
439  This essentially use the first order scheme for div(phi,nuTilda)
440 
441  p_, U_, phi_, etc: State variables in OpenFOAM
442 
443  Output:
444  nuTildaRes_: If solveTurbState_ == 0, update the residual field variable
445 
446  nuTilda_: If solveTurbState_ == 1, update nuTilda
447  */
448 
449  // Copy and modify based on the "correct" function
450 
451  label printToScreen = this->isPrintTime(mesh_.time(), printInterval_);
452 
453  word divNuTildaScheme = "div(phi,nuTilda)";
454 
455  label isPC = 0;
456 
457  if (!solveTurbState_)
458  {
459  isPC = options.getLabel("isPC");
460 
461  if (isPC)
462  {
463  divNuTildaScheme = "div(pc)";
464  }
465  }
466 
467  const volScalarField chi(this->chi());
468  const volScalarField fv1(this->fv1(chi));
469 
470  const volScalarField Stilda(
471  this->fv3(chi, fv1) * ::sqrt(2.0) * mag(skew(fvc::grad(U_)))
472  + this->fv2(chi, fv1) * nuTilda_ / sqr(kappa_ * y_));
473 
474  tmp<fvScalarMatrix> nuTildaEqn(
475  fvm::ddt(phase_, rho_, nuTilda_)
476  + fvm::div(phaseRhoPhi_, nuTilda_, divNuTildaScheme)
477  - fvm::laplacian(phase_ * rho_ * DnuTildaEff(), nuTilda_)
478  - Cb2_ / sigmaNut_ * phase_ * rho_ * magSqr(fvc::grad(nuTilda_))
479  == Cb1_ * phase_ * rho_ * Stilda * nuTilda_ * betaFieldInversion_
480  - fvm::Sp(Cw1_ * phase_ * rho_ * fw(Stilda) * nuTilda_ / sqr(y_), nuTilda_));
481 
482  nuTildaEqn.ref().relax();
483 
484  if (solveTurbState_)
485  {
486 
487  // get the solver performance info such as initial
488  // and final residuals
489  SolverPerformance<scalar> solverNuTilda = solve(nuTildaEqn);
490  if (printToScreen)
491  {
492  Info << "nuTilda Initial residual: " << solverNuTilda.initialResidual() << endl
493  << " Final residual: " << solverNuTilda.finalResidual() << endl;
494  }
495 
496  DAUtility::boundVar(allOptions_, nuTilda_, printToScreen);
497  nuTilda_.correctBoundaryConditions();
498 
499  // ***************** NOTE*****************
500  // In the original SA, it is correctNut(fv1) and fv1 is not
501  // updated based on the latest nuTilda. We use correctNut which
502  // recompute fv1 with the latest nuTilda
503  this->correctNut();
504  }
505  else
506  {
507  // calculate residuals
508  nuTildaRes_ = nuTildaEqn.ref() & nuTilda_;
509  // need to normalize residuals
510  normalizeResiduals(nuTildaRes);
511  }
512 
513  return;
514 }
515 
517 {
518  /*
519  Description:
520  Return the value of the production term from the Spalart Allmaras model
521  */
522 
523  const volScalarField chi(this->chi());
524  const volScalarField fv1(this->fv1(chi));
525 
526  const volScalarField Stilda(
527  this->fv3(chi, fv1) * ::sqrt(2.0) * mag(skew(fvc::grad(U_)))
528  + this->fv2(chi, fv1) * nuTilda_ / sqr(kappa_ * y_));
529 
530  volScalarField SAProdTerm = Cb1_ * phase_ * rho_ * Stilda * nuTilda_;
531 
532  forAll(SAProdTerm, cellI)
533  {
534  prodTerm[cellI] = SAProdTerm[cellI];
535  }
536 }
537 
538 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
539 
540 } // End namespace Foam
541 
542 // ************************************************************************* //
Foam::DASpalartAllmarasFv3FieldInversion::Cw3_
dimensionedScalar Cw3_
Definition: DASpalartAllmarasFv3FieldInversion.H:63
Foam::DASpalartAllmarasFv3FieldInversion::Cb1_
dimensionedScalar Cb1_
Definition: DASpalartAllmarasFv3FieldInversion.H:59
Foam::DASpalartAllmarasFv3FieldInversion::fv3
tmp< volScalarField > fv3(const volScalarField &chi, const volScalarField &fv1) const
Definition: DASpalartAllmarasFv3FieldInversion.C:167
forAll
forAll(pseudoP.boundaryField(), patchI)
Definition: solvePseudoPEqn.H:10
DASpalartAllmarasFv3FieldInversion.H
Foam::DASpalartAllmarasFv3FieldInversion::DnuTildaEff
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
Definition: DASpalartAllmarasFv3FieldInversion.C:198
Foam::DAOption
Definition: DAOption.H:29
Foam::DASpalartAllmarasFv3FieldInversion::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute the turbulence residuals
Definition: DASpalartAllmarasFv3FieldInversion.C:429
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DASpalartAllmarasFv3FieldInversion::fv1
tmp< volScalarField > fv1(const volScalarField &chi) const
Definition: DASpalartAllmarasFv3FieldInversion.C:153
Foam::DASpalartAllmarasFv3FieldInversion::correctNut
virtual void correctNut()
update nut based on other turbulence variables and update the BCs
Definition: DASpalartAllmarasFv3FieldInversion.C:239
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
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::DASpalartAllmarasFv3FieldInversion::correctBoundaryConditions
virtual void correctBoundaryConditions()
update turbulence variable boundary values
Definition: DASpalartAllmarasFv3FieldInversion.C:259
Foam::DASpalartAllmarasFv3FieldInversion::correctModelStates
virtual void correctModelStates(wordList &modelStates) const
update the turbulence state for DAStateInfo::regStates_
Definition: DASpalartAllmarasFv3FieldInversion.C:205
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DASpalartAllmarasFv3FieldInversion::printInterval_
label printInterval_
time step interval to print residual
Definition: DASpalartAllmarasFv3FieldInversion.H:118
Foam::DASpalartAllmarasFv3FieldInversion::fw
tmp< volScalarField > fw(const volScalarField &Stilda) const
Definition: DASpalartAllmarasFv3FieldInversion.C:180
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::DASpalartAllmarasFv3FieldInversion::getTurbProdTerm
virtual void getTurbProdTerm(scalarList &prodTerm) const
return the value of the production term from the turbulence model
Definition: DASpalartAllmarasFv3FieldInversion.C:516
Foam::DASpalartAllmarasFv3FieldInversion::betaFieldInversion_
volScalarField & betaFieldInversion_
A beta field multiplying to the production term.
Definition: DASpalartAllmarasFv3FieldInversion.H:94
Foam::DASpalartAllmarasFv3FieldInversion::kappa_
dimensionedScalar kappa_
Definition: DASpalartAllmarasFv3FieldInversion.H:58
Foam::DASpalartAllmarasFv3FieldInversion::Cw1_
dimensionedScalar Cw1_
Definition: DASpalartAllmarasFv3FieldInversion.H:61
Foam::DASpalartAllmarasFv3FieldInversion::chi
tmp< volScalarField > chi() const
Definition: DASpalartAllmarasFv3FieldInversion.C:148
Foam::DATurbulenceModel::U_
volVectorField & U_
velocity
Definition: DATurbulenceModel.H:86
Foam::DASpalartAllmarasFv3FieldInversion::Cv1_
dimensionedScalar Cv1_
Definition: DASpalartAllmarasFv3FieldInversion.H:64
Foam::DATurbulenceModel::nut_
volScalarField & nut_
turbulence viscosity
Definition: DATurbulenceModel.H:83
Foam::DASpalartAllmarasFv3FieldInversion::nuTilda_
volScalarField & nuTilda_
Definition: DASpalartAllmarasFv3FieldInversion.H:88
Foam::DASpalartAllmarasFv3FieldInversion::Cw2_
dimensionedScalar Cw2_
Definition: DASpalartAllmarasFv3FieldInversion.H:62
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DASpalartAllmarasFv3FieldInversion::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DASpalartAllmarasFv3FieldInversion.C:270
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:262
Foam::DASpalartAllmarasFv3FieldInversion::DASpalartAllmarasFv3FieldInversion
DASpalartAllmarasFv3FieldInversion(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DASpalartAllmarasFv3FieldInversion.C:41
Foam::DASpalartAllmarasFv3FieldInversion::correct
virtual void correct()
solve the residual equations and update the state
Definition: DASpalartAllmarasFv3FieldInversion.C:409
Foam::DATurbulenceModel::correctAlphat
void correctAlphat()
update alphat
Definition: DATurbulenceModel.C:206
Foam::DASpalartAllmarasFv3FieldInversion::Cb2_
dimensionedScalar Cb2_
Definition: DASpalartAllmarasFv3FieldInversion.H:60
Foam::DASpalartAllmarasFv3FieldInversion::solveTurbState_
label solveTurbState_
whether to solve for turb states
Definition: DASpalartAllmarasFv3FieldInversion.H:115
solve
pseudoPEqn solve(solverDictP_)
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
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::DASpalartAllmarasFv3FieldInversion::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DASpalartAllmarasFv3FieldInversion.C:331
Foam::DASpalartAllmarasFv3FieldInversion::Cv2_
dimensionedScalar Cv2_
Definition: DASpalartAllmarasFv3FieldInversion.H:65
Foam::DASpalartAllmarasFv3FieldInversion::sigmaNut_
dimensionedScalar sigmaNut_
Definition: DASpalartAllmarasFv3FieldInversion.H:57
Foam::DASpalartAllmarasFv3FieldInversion::nuTildaRes_
volScalarField nuTildaRes_
Definition: DASpalartAllmarasFv3FieldInversion.H:89
Foam::DASpalartAllmarasFv3FieldInversion::y_
const volScalarField & y_
3D wall distance
Definition: DASpalartAllmarasFv3FieldInversion.H:112
Foam::DATurbulenceModel::phase_
volScalarField phase_
phase field
Definition: DATurbulenceModel.H:92
Foam::DASpalartAllmarasFv3FieldInversion::fv2
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
Definition: DASpalartAllmarasFv3FieldInversion.C:160
Foam::DASpalartAllmarasFv3FieldInversion::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DASpalartAllmarasFv3FieldInversion.C:281