Protected Member Functions | |
SA functions | |
tmp< volScalarField > | chi () const |
tmp< volScalarField > | fv1 (const volScalarField &chi) const |
tmp< volScalarField > | fv2 (const volScalarField &chi, const volScalarField &fv1) const |
tmp< volScalarField > | fv3 (const volScalarField &chi, const volScalarField &fv1) const |
tmp< volScalarField > | fw (const volScalarField &Stilda) const |
Protected Attributes | |
SA parameters | |
dimensionedScalar | sigmaNut_ |
dimensionedScalar | kappa_ |
dimensionedScalar | Cb1_ |
dimensionedScalar | Cb2_ |
dimensionedScalar | Cw1_ |
dimensionedScalar | Cw2_ |
dimensionedScalar | Cw3_ |
dimensionedScalar | Cv1_ |
dimensionedScalar | Cv2_ |
Protected Attributes inherited from DATurbulenceModel | |
const fvMesh & | mesh_ |
fvMesh More... | |
const DAOption & | daOption_ |
DAOption object. More... | |
const dictionary & | allOptions_ |
all DAFoam options More... | |
volScalarField & | nut_ |
turbulence viscosity More... | |
volVectorField & | U_ |
velocity More... | |
surfaceScalarField & | phi_ |
face flux More... | |
volScalarField | phase_ |
phase field More... | |
surfaceScalarField & | phaseRhoPhi_ |
phase*phi*density field More... | |
IOdictionary | turbDict_ |
turbulence model property dict More... | |
dictionary | coeffDict_ |
turbulence model parameters dict More... | |
dimensionedScalar | kMin_ |
Lower limit of k. More... | |
dimensionedScalar | epsilonMin_ |
Lower limit of epsilon. More... | |
dimensionedScalar | omegaMin_ |
Lower limit for omega. More... | |
dimensionedScalar | nuTildaMin_ |
Lower limit for nuTilda. More... | |
scalar | Pr_ |
Prandtl number. More... | |
scalar | Prt_ = -9999.0 |
turbulent Prandtl number More... | |
Augmented variables for residual computation | |
volScalarField & | nuTilda_ |
volScalarField | nuTildaRes_ |
volScalarField | nuTildaResPartDeriv_ |
volScalarField & | betaFieldInversion_ |
A beta field multiplying to the production term. More... | |
volVectorField | UData_ |
reference velocity field More... | |
volScalarField | surfaceFriction_ |
a surface friction 'field' when using skin friction data for field inversion More... | |
volScalarField | surfaceFrictionData_ |
the reference field for surfaceFriction More... | |
volScalarField | pData_ |
the reference field for pressure More... | |
volScalarField | USingleComponentData_ |
reference field for profile data More... | |
const volScalarField & | y_ |
3D wall distance More... | |
label | solveTurbState_ = 0 |
whether to solve for turb states More... | |
label | printInterval_ |
time step interval to print residual More... | |
TypeName ("SpalartAllmarasFv3FieldInversion") | |
DASpalartAllmarasFv3FieldInversion (const word modelType, const fvMesh &mesh, const DAOption &daOption) | |
virtual | ~DASpalartAllmarasFv3FieldInversion () |
tmp< volScalarField > | DnuTildaEff () const |
Return the effective diffusivity for nuTilda. More... | |
virtual void | correctModelStates (wordList &modelStates) const |
update the turbulence state for DAStateInfo::regStates_ More... | |
virtual void | correctNut () |
update nut based on other turbulence variables and update the BCs More... | |
virtual void | correctBoundaryConditions () |
update turbulence variable boundary values More... | |
virtual void | updateIntermediateVariables () |
update any intermediate variables that are dependent on state variables and are used in calcResiduals More... | |
virtual void | correctStateResidualModelCon (List< List< word >> &stateCon) const |
update the original variable connectivity for the adjoint state residuals in stateCon More... | |
virtual void | addModelResidualCon (HashTable< List< List< word >>> &allCon) const |
add the model residual connectivity to stateCon More... | |
virtual void | calcResiduals (const dictionary &options) |
compute the turbulence residuals More... | |
virtual void | correct () |
solve the residual equations and update the state More... | |
virtual void | getTurbProdTerm (scalarList &prodTerm) const |
return the value of the production term from the turbulence model More... | |
Additional Inherited Members | |
Public Member Functions inherited from DATurbulenceModel | |
TypeName ("DATurbulenceModel") | |
declareRunTimeSelectionTable (autoPtr, DATurbulenceModel, dictionary,(const word modelType, const fvMesh &mesh, const DAOption &daOption),(modelType, mesh, daOption)) | |
DATurbulenceModel (const word modelType, const fvMesh &mesh, const DAOption &daOption) | |
virtual | ~DATurbulenceModel () |
void | correctWallDist () |
update wall distance for d_. Note: y_ will be automatically updated in mesh_ object More... | |
void | correctAlphat () |
update alphat More... | |
tmp< volSymmTensorField > | devRhoReff () const |
dev terms More... | |
tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) |
divDev terms More... | |
tmp< fvVectorMatrix > | divDevReff (volVectorField &U) |
divDev terms More... | |
tmp< volScalarField > | nuEff () const |
return effective viscosity More... | |
tmp< volScalarField > | getNut () |
get the nut field More... | |
tmp< volScalarField > | alphaEff () |
return effective thermal diffusivity More... | |
tmp< volScalarField > | nu () const |
get the nu field More... | |
tmp< volScalarField > | getAlpha () const |
get alpha field More... | |
tmp< volScalarField > | getRho () |
get the density field More... | |
tmp< volScalarField > | getPhase () |
get the phase field More... | |
scalar | getPrt () |
get the turbulent Prandtl number More... | |
tmp< Foam::volScalarField > | getMu () const |
get mu More... | |
bool | writeData (Ostream &os) const |
this is a virtual function for regIOobject More... | |
void | printYPlus () const |
print the min max and mean yPlus to screen More... | |
label | isPrintTime (const Time &runTime, const label printInterval) const |
virtual void | invTranProdNuTildaEqn (const volScalarField &mySource, volScalarField &pseudoNuTilda) |
Inverse transpose product, M_nuTilda^(-T) More... | |
virtual void | constructPseudoNuTildaEqn () |
virtual void | rhsSolvePseudoNuTildaEqn (const volScalarField &nuTildaSource) |
virtual void | calcLduResidualTurb (volScalarField &nuTildaRes) |
calculate the turbulence residual using LDU matrix More... | |
Static Public Member Functions inherited from DATurbulenceModel | |
static autoPtr< DATurbulenceModel > | New (const word modelType, const fvMesh &mesh, const DAOption &daOption) |
Definition at line 50 of file DASpalartAllmarasFv3FieldInversion.H.
DASpalartAllmarasFv3FieldInversion | ( | const word | modelType, |
const fvMesh & | mesh, | ||
const DAOption & | daOption | ||
) |
Definition at line 41 of file DASpalartAllmarasFv3FieldInversion.C.
References daOption(), mesh, and DASpalartAllmarasFv3FieldInversion::printInterval_.
|
inlinevirtual |
Definition at line 131 of file DASpalartAllmarasFv3FieldInversion.H.
|
protected |
Definition at line 148 of file DASpalartAllmarasFv3FieldInversion.C.
References DATurbulenceModel::nu(), and DASpalartAllmarasFv3FieldInversion::nuTilda_.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals(), DASpalartAllmarasFv3FieldInversion::correctNut(), DASpalartAllmarasFv3FieldInversion::fv1(), DASpalartAllmarasFv3FieldInversion::fv2(), DASpalartAllmarasFv3FieldInversion::fv3(), and DASpalartAllmarasFv3FieldInversion::getTurbProdTerm().
|
protected |
Definition at line 153 of file DASpalartAllmarasFv3FieldInversion.C.
References DASpalartAllmarasFv3FieldInversion::chi(), and DASpalartAllmarasFv3FieldInversion::Cv1_.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals(), DASpalartAllmarasFv3FieldInversion::correctNut(), DASpalartAllmarasFv3FieldInversion::fv3(), and DASpalartAllmarasFv3FieldInversion::getTurbProdTerm().
|
protected |
Definition at line 160 of file DASpalartAllmarasFv3FieldInversion.C.
References DASpalartAllmarasFv3FieldInversion::chi(), and DASpalartAllmarasFv3FieldInversion::Cv2_.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals(), and DASpalartAllmarasFv3FieldInversion::getTurbProdTerm().
|
protected |
Definition at line 167 of file DASpalartAllmarasFv3FieldInversion.C.
References DASpalartAllmarasFv3FieldInversion::chi(), DASpalartAllmarasFv3FieldInversion::Cv2_, and DASpalartAllmarasFv3FieldInversion::fv1().
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals(), and DASpalartAllmarasFv3FieldInversion::getTurbProdTerm().
|
protected |
Definition at line 180 of file DASpalartAllmarasFv3FieldInversion.C.
References DASpalartAllmarasFv3FieldInversion::Cw2_, DASpalartAllmarasFv3FieldInversion::Cw3_, DASpalartAllmarasFv3FieldInversion::kappa_, DASpalartAllmarasFv3FieldInversion::nuTilda_, and DASpalartAllmarasFv3FieldInversion::y_.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals().
TypeName | ( | "SpalartAllmarasFv3FieldInversion" | ) |
tmp< volScalarField > DnuTildaEff | ( | ) | const |
Return the effective diffusivity for nuTilda.
Definition at line 198 of file DASpalartAllmarasFv3FieldInversion.C.
References DATurbulenceModel::nu(), DASpalartAllmarasFv3FieldInversion::nuTilda_, and DASpalartAllmarasFv3FieldInversion::sigmaNut_.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals().
|
virtual |
update the turbulence state for DAStateInfo::regStates_
Implements DATurbulenceModel.
Definition at line 205 of file DASpalartAllmarasFv3FieldInversion.C.
References forAll().
|
virtual |
update nut based on other turbulence variables and update the BCs
Implements DATurbulenceModel.
Definition at line 239 of file DASpalartAllmarasFv3FieldInversion.C.
References DASpalartAllmarasFv3FieldInversion::chi(), DATurbulenceModel::correctAlphat(), DASpalartAllmarasFv3FieldInversion::fv1(), DATurbulenceModel::nut_, and DASpalartAllmarasFv3FieldInversion::nuTilda_.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals(), and DASpalartAllmarasFv3FieldInversion::updateIntermediateVariables().
|
virtual |
update turbulence variable boundary values
Implements DATurbulenceModel.
Definition at line 259 of file DASpalartAllmarasFv3FieldInversion.C.
References DASpalartAllmarasFv3FieldInversion::nuTilda_.
|
virtual |
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Implements DATurbulenceModel.
Definition at line 270 of file DASpalartAllmarasFv3FieldInversion.C.
References DASpalartAllmarasFv3FieldInversion::correctNut().
|
virtual |
update the original variable connectivity for the adjoint state residuals in stateCon
Implements DATurbulenceModel.
Definition at line 281 of file DASpalartAllmarasFv3FieldInversion.C.
References forAll().
|
virtual |
add the model residual connectivity to stateCon
Implements DATurbulenceModel.
Definition at line 331 of file DASpalartAllmarasFv3FieldInversion.C.
References DATurbulenceModel::mesh_.
|
virtual |
compute the turbulence residuals
Implements DATurbulenceModel.
Definition at line 429 of file DASpalartAllmarasFv3FieldInversion.C.
References DATurbulenceModel::allOptions_, DASpalartAllmarasFv3FieldInversion::betaFieldInversion_, DAUtility::boundVar(), DASpalartAllmarasFv3FieldInversion::Cb1_, DASpalartAllmarasFv3FieldInversion::Cb2_, DASpalartAllmarasFv3FieldInversion::chi(), DASpalartAllmarasFv3FieldInversion::correctNut(), DASpalartAllmarasFv3FieldInversion::Cw1_, DASpalartAllmarasFv3FieldInversion::DnuTildaEff(), DASpalartAllmarasFv3FieldInversion::fv1(), DASpalartAllmarasFv3FieldInversion::fv2(), DASpalartAllmarasFv3FieldInversion::fv3(), DASpalartAllmarasFv3FieldInversion::fw(), DATurbulenceModel::isPrintTime(), DASpalartAllmarasFv3FieldInversion::kappa_, DATurbulenceModel::mesh_, normalizeResiduals, DASpalartAllmarasFv3FieldInversion::nuTilda_, DASpalartAllmarasFv3FieldInversion::nuTildaRes_, DATurbulenceModel::phase_, DATurbulenceModel::phaseRhoPhi_, DASpalartAllmarasFv3FieldInversion::printInterval_, DASpalartAllmarasFv3FieldInversion::sigmaNut_, solve(), DASpalartAllmarasFv3FieldInversion::solveTurbState_, DATurbulenceModel::U_, and DASpalartAllmarasFv3FieldInversion::y_.
Referenced by DASpalartAllmarasFv3FieldInversion::correct().
|
virtual |
solve the residual equations and update the state
Implements DATurbulenceModel.
Definition at line 409 of file DASpalartAllmarasFv3FieldInversion.C.
References DASpalartAllmarasFv3FieldInversion::calcResiduals(), and DASpalartAllmarasFv3FieldInversion::solveTurbState_.
|
virtual |
return the value of the production term from the turbulence model
Reimplemented from DATurbulenceModel.
Definition at line 516 of file DASpalartAllmarasFv3FieldInversion.C.
References DASpalartAllmarasFv3FieldInversion::Cb1_, DASpalartAllmarasFv3FieldInversion::chi(), forAll(), DASpalartAllmarasFv3FieldInversion::fv1(), DASpalartAllmarasFv3FieldInversion::fv2(), DASpalartAllmarasFv3FieldInversion::fv3(), DASpalartAllmarasFv3FieldInversion::kappa_, DASpalartAllmarasFv3FieldInversion::nuTilda_, DATurbulenceModel::phase_, DATurbulenceModel::U_, and DASpalartAllmarasFv3FieldInversion::y_.
|
protected |
Definition at line 57 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals(), and DASpalartAllmarasFv3FieldInversion::DnuTildaEff().
|
protected |
Definition at line 58 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals(), DASpalartAllmarasFv3FieldInversion::fw(), and DASpalartAllmarasFv3FieldInversion::getTurbProdTerm().
|
protected |
Definition at line 59 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals(), and DASpalartAllmarasFv3FieldInversion::getTurbProdTerm().
|
protected |
Definition at line 60 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals().
|
protected |
Definition at line 61 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals().
|
protected |
Definition at line 62 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::fw().
|
protected |
Definition at line 63 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::fw().
|
protected |
Definition at line 64 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::fv1().
|
protected |
Definition at line 65 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::fv2(), and DASpalartAllmarasFv3FieldInversion::fv3().
|
protected |
Definition at line 88 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals(), DASpalartAllmarasFv3FieldInversion::chi(), DASpalartAllmarasFv3FieldInversion::correctBoundaryConditions(), DASpalartAllmarasFv3FieldInversion::correctNut(), DASpalartAllmarasFv3FieldInversion::DnuTildaEff(), DASpalartAllmarasFv3FieldInversion::fw(), and DASpalartAllmarasFv3FieldInversion::getTurbProdTerm().
|
protected |
Definition at line 89 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals().
|
protected |
Definition at line 90 of file DASpalartAllmarasFv3FieldInversion.H.
|
protected |
A beta field multiplying to the production term.
Definition at line 94 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals().
|
protected |
reference velocity field
Definition at line 97 of file DASpalartAllmarasFv3FieldInversion.H.
|
protected |
a surface friction 'field' when using skin friction data for field inversion
Definition at line 100 of file DASpalartAllmarasFv3FieldInversion.H.
|
protected |
the reference field for surfaceFriction
Definition at line 103 of file DASpalartAllmarasFv3FieldInversion.H.
|
protected |
the reference field for pressure
Definition at line 106 of file DASpalartAllmarasFv3FieldInversion.H.
|
protected |
reference field for profile data
Definition at line 109 of file DASpalartAllmarasFv3FieldInversion.H.
|
protected |
3D wall distance
Definition at line 112 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals(), DASpalartAllmarasFv3FieldInversion::fw(), and DASpalartAllmarasFv3FieldInversion::getTurbProdTerm().
|
protected |
whether to solve for turb states
Definition at line 115 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals(), and DASpalartAllmarasFv3FieldInversion::correct().
|
protected |
time step interval to print residual
Definition at line 118 of file DASpalartAllmarasFv3FieldInversion.H.
Referenced by DASpalartAllmarasFv3FieldInversion::calcResiduals(), and DASpalartAllmarasFv3FieldInversion::DASpalartAllmarasFv3FieldInversion().