Go to the documentation of this file.
37 #ifndef adjointSolver_H
38 #define adjointSolver_H
110 const word& managerType,
124 const word& managerType,
136 const word& managerType,
Base class for adjoint solvers.
class for managing incompressible objective functions.
bool computeSensitivities_
Are sensitivities computed.
bool isConstraint_
Is the adjoint solver used to tackle a constraint.
const word primalSolverName_
Name of primal solver.
A class for handling words, derived from Foam::string.
TypeName("adjointSolver")
Run-time type information.
const word & primalSolverName() const
Return the primal solver name.
virtual bool isConstraint()
Is the solving referring to a constraint.
virtual void updatePrimalBasedQuantities()
virtual bool readDict(const dictionary &dict)
declareRunTimeNewSelectionTable(autoPtr, adjointSolver, adjointSolver,(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName),(mesh, managerType, dict, primalSolverName))
autoPtr< objectiveManager > objectiveManagerPtr_
Object to manage objective functions.
virtual const scalarField & getObjectiveSensitivities()=0
Grab a reference to the computed sensitivities.
virtual void computeObjectiveSensitivities()=0
Compute sensitivities of the underlaying objectives.
virtual sensitivity & getSensitivityBase()=0
Return the base sensitivity object.
Base class for solution control classes.
virtual ~adjointSolver()=default
Destructor.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual void clearSensitivities()
Clears the sensitivity field known by the adjoint solver.
Mesh data needed to do the Finite Volume discretisation.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
const objectiveManager & getObjectiveManager() const
Return a const reference to the objective manager.
const fvMesh & mesh() const
Return the solver mesh.
Macros to ease declaration of run-time selection tables.
Base class for primal solvers.
autoPtr< scalarField > sensitivities_
Sensitivities field.
Abstract base class for adjoint sensitivities.
const primalSolver & getPrimalSolver() const
virtual const dictionary & dict() const
Return the solver dictionary.
static autoPtr< adjointSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName)
Return a reference to the selected turbulence model.