Go to the documentation of this file.
48 #ifndef adjointTurbulenceModel_H
49 #define adjointTurbulenceModel_H
65 namespace incompressibleAdjoint
114 const word& adjointTurbulenceModelName
120 adjointTurbulenceModelName
133 const word& adjointTurbulenceModelName = typeName
145 const word& adjointTurbulenceModelName = typeName
181 lamTrans.
nu()() + turbVars().nutRef()
206 virtual bool read() = 0;
Manages the adjoint mean flow fields and their mean values.
class for managing incompressible objective functions.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A class for handling words, derived from Foam::string.
virtual bool read()=0
Read adjointLESProperties or adjointRASProperties dictionary.
A class for managing temporary objects.
static autoPtr< adjointTurbulenceModel > New(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName=typeName)
Return a reference to the selected turbulence model.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
TypeName("adjointTurbulenceModel")
Runtime type information.
A simple single-phase transport model based on viscosityModel.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
declareRunTimeNewSelectionTable(autoPtr, adjointTurbulenceModel, adjointTurbulenceModel,(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName),(primalVars, adjointVars, objManager, adjointTurbulenceModelName))
const singlePhaseTransportModel & laminarTransport() const
Return const reference to transport model.
Abstract base class for incompressible adjoint turbulence models (RAS, LES and laminar).
virtual bool writeData(Ostream &) const
Default dummy write function.
virtual void correct()=0
Solve the adjoint turbulence equations.
virtual void nullify()=0
Nullify all adjoint turbulence model fields and their old times.
Mesh data needed to do the Finite Volume discretisation.
virtual tmp< volVectorField > adjointMeanFlowSource()=0
Source term added to the adjoint mean flow due to the.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
virtual ~adjointTurbulenceModel()=default
Destructor.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
Macros to ease declaration of run-time selection tables.
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
incompressibleAdjointMeanFlowVars & adjointVars_
incompressibleVars & primalVars_
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const =0
Return the diffusion term for the momentum equation.
virtual tmp< volSymmTensorField > devReff() const =0
Return the effective stress tensor including the laminar stress.
tmp< volScalarField > nu() const
Return the laminar viscosity.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
virtual const volScalarField & nut() const
Return the turbulence viscosity.
Base class for solution control classes.