Go to the documentation of this file.
69 #ifndef adjointSpalartAllmaras_H
70 #define adjointSpalartAllmaras_H
79 namespace incompressibleAdjoint
81 namespace adjointRASModels
306 const word& adjointTurbulenceModelName
307 = adjointTurbulenceModel::typeName,
308 const word& modelName = typeName
Manages the adjoint mean flow fields and their mean values.
class for managing incompressible objective functions.
tmp< volScalarField > fv1(const volScalarField &chi) const
volVectorField gradNuTilda_
tmp< volScalarField > r(const volScalarField &Stilda) const
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
Continuous adjoint to the Spalart-Allmaras one-eqn mixing-length model for incompressible flows.
tmp< volScalarField > dStilda_dNuTilda(const volScalarField &Omega, const volScalarField &fv2, const volScalarField &dFv2dChi) const
dimensionedScalar sigmaNut_
TypeName("adjointSpalartAllmaras")
Runtime type information.
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
virtual tmp< volVectorField > adjointMeanFlowSource()
tmp< volScalarField > dStilda_dOmega(const volScalarField &Omega, const volScalarField &fv2) const
tmp< volScalarField > dnut_dNuTilda(const volScalarField &fv1, const volScalarField &dFv1dChi) const
tmp< volScalarField > dfw_dOmega(const volScalarField &Stilda, const volScalarField &dfwdr, const volScalarField &dStildadOmega) const
virtual const boundaryVectorField & adjointMomentumBCSource() const
tmp< volScalarField > dfw_dDelta(const volScalarField &Stilda, const volScalarField &dfwdr, const volScalarField &dStildadDelta) const
tmp< volScalarField > dr_dNuTilda(const volScalarField &Stilda) const
dimensionedScalar minStilda_
tmp< volScalarField > dFv2_dChi(const volScalarField &chi, const volScalarField &fv1, const volScalarField &dFv1dChi) const
virtual tmp< volScalarField > distanceSensitivities()
tmp< volScalarField > dr_dDelta(const volScalarField &Stilda) const
volScalarField & nuaTilda()
Access to the adjoint Spalart - Allmaras field.
tmp< volScalarField > dfw_dNuTilda(const volScalarField &Stilda, const volScalarField &dfwdr, const volScalarField &dStildadNuTilda) const
tmp< volScalarField > dfw_dr(const volScalarField &Stilda) const
virtual bool read()
Read adjointRASProperties dictionary.
volSymmTensorField symmAdjointProductionU_
tmp< volVectorField > conservativeMomentumSource()
Conservative source term for the adjoint momentum equations.
tmp< volScalarField > dD_dNuTilda(const volScalarField &fw, const volScalarField &dfwdNuTilda) const
virtual tmp< volTensorField > FISensitivityTerm()
Term contributing to the computation of FI-based sensitivities.
const volScalarField & nuTilda() const
References to the primal turbulence model variables.
virtual void nullify()
Nullify all adjoint turbulence model fields and their old times.
virtual tmp< scalarField > diffusionCoeffVar1(label patchI) const
const volScalarField & nut() const
Return the turbulence viscosity.
const volScalarField & y_
Wall distance.
Abstract base class for incompressible turbulence models.
tmp< volScalarField > dFv1_dChi(const volScalarField &chi) const
tmp< volScalarField > allocateMask()
Allocate the mask field.
tmp< volScalarField > dr_dStilda(const volScalarField &Stilda) const
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1) const
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt the first turbulence model variable.
tmp< volScalarField > fw(const volScalarField &Stilda) const
volTensorField momentumSourceMult_
void updatePrimalRelatedFields()
Update the constant primal-related fields.
autoPtr< volScalarField > adjointTMVariable1Ptr_
Adjoint turbulence model variable 1.
tmp< volScalarField > chi() const
bool limitAdjointProduction_
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor including the laminar stress.
tmp< volScalarField > dP_dNuTilda(const volScalarField &dStildadNuTilda) const
volScalarField mask_
Field for masking (convergence enhancement)
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const
Return the diffusion term for the momentum equation.
tmp< volScalarField > dStilda_dDelta(const volScalarField &Omega, const volScalarField &fv2) const
virtual const boundaryVectorField & wallShapeSensitivities()
Sensitivity terms for shape optimisation, emerging from.
virtual void correct()
Solve the adjoint turbulence equations.
volScalarField productionDestructionSource_
tmp< volScalarField > DnuTildaEff() const
Generic GeometricField class.
virtual ~adjointSpalartAllmaras()=default
Destructor.
Base class for solution control classes.
virtual const boundaryVectorField & wallFloCoSensitivities()
Sensitivity terms for flow control, emerging from the.