46#ifndef adjointRASModel_H
47#define adjointRASModel_H
80 public adjointTurbulenceModel,
88 adjointRASModel(
const adjointRASModel&) =
delete;
91 void operator=(
const adjointRASModel&) =
delete;
172 const word& adjointTurbulenceModelName
178 adjointTurbulenceModelName
192 const word& adjointTurbulenceModelName =
193 adjointTurbulenceModel::typeName
205 const word& adjointTurbulenceModelName =
206 adjointTurbulenceModel::typeName
340 const word& designVarsName
Bound the given scalar field if it has gone unbounded.
Useful typenames for fields defined only at the boundaries.
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
dictionary()
Default construct, a top-level empty dictionary.
Manages the adjoint mean flow fields and their mean values.
virtual tmp< scalarField > diffusionCoeffVar1(label patchI) const
Diffusion coefficient of the first primal and adjoint turbulence model equation. Needed for some adjo...
virtual tmp< volVectorField > adjointMeanFlowSource()=0
Source terms to the adjoint momentum equation due to the differentiation of the turbulence model.
virtual tmp< volScalarField > nutJacobianTMVar2() const
Jacobian of nut wrt the second turbulence model variable.
volScalarField & getAdjointTMVariable1Inst()
Return non-constant reference to adjoint turbulence model variable 1.
virtual const boundaryVectorField & wallShapeSensitivities()=0
Sensitivity terms for shape optimisation, emerging from the turbulence model differentiation.
dictionary coeffDict_
Model coefficients dictionary.
virtual const boundaryVectorField & adjointMomentumBCSource() const =0
Source for the outlet adjoint momentum BC coming from differentiating the turbulence model.
virtual tmp< volScalarField > distanceSensitivities()=0
Sensitivity terms resulting from the differentiation of the distance field. Misses dxdb,...
bool includeDistance_
Does the turbulence model include distances and should the adjoint to the distance field be computed.
virtual tmp< volTensorField > FISensitivityTerm()=0
Term contributing to the computation of FI-based sensitivities.
autoPtr< boundaryVectorField > wallShapeSensitivitiesPtr_
Wall sensitivity term for shape optimisation.
autoPtr< volScalarField > adjointTMVariable1MeanPtr_
Adjoint turbulence model variable 1, mean value.
autoPtr< volScalarField > adjointTMVariable2MeanPtr_
Adjoint turbulence model variable 2, mean value.
volScalarField & getAdjointTMVariable2Inst()
Return non-constant reference to adjoint turbulence model variable 2.
virtual void correct()
Solve the adjoint turbulence equations.
virtual tmp< volSymmTensorField > devReff() const =0
Return the effective stress tensor including the laminar stress.
const word & adjointSolverName() const
Const access to the adjoint solver name.
nearWallDist y_
Near wall distance boundary field.
Switch adjointTurbulence_
Turbulence on/off flag.
virtual tmp< volSymmTensorField > devReff(const volVectorField &U) const =0
Return the effective stress tensor based on a given velocity field.
void restoreInitValues()
Restore field values to the initial ones.
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const =0
Return the diffusion term for the momentum equation.
virtual tmp< scalarField > topologySensitivities(const word &designVarsName) const =0
Term contributing to the computation of topology optimisation sensitivities.
bool changedPrimalSolution_
Has the primal solution changed?
autoPtr< boundaryVectorField > adjMomentumBCSourcePtr_
Source to the adjoint momentum BC emerging from differentiating the turbulence model.
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt the first turbulence model variable.
const nearWallDist & y() const
Return the near wall distances.
virtual void printCoeffs()
Print model coefficients.
bool includeDistance() const
Should the adjoint to the eikonal equation be computed.
autoPtr< volScalarField > & getAdjointTMVariable2InstPtr()
Return non-constant autoPtr to adjoint turbulence model variable 2.
autoPtr< volScalarField > adjointTMVariable1Ptr_
Adjoint turbulence model variable 1.
autoPtr< volScalarField > adjointTMVariable2Ptr_
Adjoint turbulence model variable 2.
void computeMeanFields()
Average adjoint fields on the fly.
declareRunTimeSelectionTable(autoPtr, adjointRASModel, dictionary,(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName),(primalVars, adjointVars, objManager, adjointTurbulenceModelName))
void resetMeanFields()
Reset mean fields to zero.
void setChangedPrimalSolution()
Set flag of changed primal solution to true.
objectiveManager & objectiveManager_
Reference to the objectiveManager.
volScalarField & getAdjointTMVariable2()
Return non-constant reference to adjoint turbulence model variable 2.
volScalarField & getAdjointTMVariable1()
Return non-constant reference to adjoint turbulence model variable 1.
TypeName("adjointRASModel")
Runtime type information.
const wordList & getAdjointTMVariablesBaseNames() const
Return reference to the adjoint turbulence model variables base names.
wordList adjointTMVariablesBaseNames_
Base names of the adjoint fields.
virtual tmp< volVectorField > nutJacobianU(tmp< volScalarField > &dNutdUMult) const
Jacobian of nut wrt the flow velocity.
virtual tmp< scalarField > diffusionCoeffVar2(label patchI) const
Diffusion coefficient of the second primal and adjoint turbulence model equation. Needed for some adj...
autoPtr< boundaryVectorField > wallFloCoSensitivitiesPtr_
Wall sensitivity term for flow control optimisation.
autoPtr< volScalarField > & getAdjointTMVariable1InstPtr()
Return non-constant autoPtr to adjoint turbulence model variable 1.
virtual void nullify()=0
Nullify all adjoint turbulence model fields and their old times.
virtual ~adjointRASModel()=default
Destructor.
virtual const boundaryVectorField & wallFloCoSensitivities()=0
Sensitivity terms for flow control, emerging from the turbulence model differentiation.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Switch printCoeffs_
Flag to print the model coeffs at run-time.
void setMeanFields()
Set mean fields.
virtual bool read()
Read adjointRASProperties dictionary.
static autoPtr< adjointRASModel > New(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName=adjointTurbulenceModel::typeName)
Return a reference to the selected adjointRAS model.
const word & primalSolverName() const
Const access to the primal solver name.
incompressibleAdjointMeanFlowVars & adjointVars_
incompressibleVars & primalVars_
Base class for solution control classes.
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest.
Class for managing objective functions.
A class for managing temporary objects.
const word & solverName() const
Return solver name.
A class for handling words, derived from Foam::string.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Namespace for incompressible adjoint turbulence models.
List< word > wordList
List of word.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
volVectorField::Boundary boundaryVectorField
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.