37#ifndef adjointSolver_H
38#define adjointSolver_H
69 adjointSolver(
const adjointSolver&) =
delete;
72 void operator=(
const adjointSolver&) =
delete;
277 autoPtr<boundaryVectorField>& dSfdbMult,
278 autoPtr<boundaryVectorField>& dnfdbMult,
279 autoPtr<boundaryVectorField>& dxdbDirectMult,
280 autoPtr<pointBoundaryVectorField>& pointDxDirectDbMult,
291 autoPtr<boundaryVectorField>& bcDxDbMult,
314 const word& designVariablesName,
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
virtual void clearSensitivities()
Clears the sensitivity field known by the adjoint solver.
virtual void accumulateDivDxDbMultiplier(autoPtr< scalarField > &divDxDbMult, const scalar dt)
Compute the multiplier for div(dxdb).
objectiveManager objectiveManager_
Object to manage objective functions.
virtual tmp< volScalarField > yWall() const
Return the distance field, to be used in the solution of the adjoint eikonal PDE.
tmp< scalarField > sensitivities_
Sensitivities field.
void allocateSensitivities()
Allocate the sensitivity derivatives.
virtual void accumulateGradDxDbMultiplier(volTensorField &gradDxDbMult, const scalar dt)
Compute the multiplier for grad(dxdb).
virtual tmp< volScalarField > adjointEikonalSource()
Return the source the adjoint eikonal equation.
virtual const scalarField & getObjectiveSensitivities(autoPtr< designVariables > &designVars)
Grab a reference to the computed sensitivities.
virtual dimensionSet daDimensions() const
Return the dimensions of the adjoint distance field.
virtual bool writeData(Ostream &os) const
Write the sensitivity derivatives.
static autoPtr< adjointSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName, const word &solverName)
Return a reference to the selected turbulence model.
bool isDoubleSidedConstraint()
Is the solving referring to a double-sided constraint.
virtual void topOSensMultiplier(scalarField &betaMult, const word &designVariablesName, const scalar dt)
Compute the multiplier of beta.
autoPtr< adjointSensitivity > adjointSensitivity_
Sensitivity Derivatives engine.
virtual void accumulateGeometryVariationsMultipliers(autoPtr< boundaryVectorField > &dSfdbMult, autoPtr< boundaryVectorField > &dnfdbMult, autoPtr< boundaryVectorField > &dxdbDirectMult, autoPtr< pointBoundaryVectorField > &pointDxDirectDbMult, const labelHashSet &sensitivityPatchIDs, const scalar dt)
Accumulate the multipliers of geometric quantities defined at the boundary, usually through an object...
virtual ~adjointSolver()=default
Destructor.
virtual bool readDict(const dictionary &dict)
const primalSolver & getPrimalSolver() const
Return a const-reference to the primal solver corresponding to this adjoint solver.
dictionary designVarsDict() const
Return the dictionary corresponding to the design variables.
virtual void accumulateOptionsDxDbMultiplier(vectorField &optionsDxDbMult, const scalar dt)
Contributions from fvOptions that inlcude geometric aspects in them and change when the geometry is d...
bool computeSensitivities_
Are sensitivities computed.
virtual bool includeDistance() const
Does the adjoint to an equation computing distances need to taken into consideration.
virtual void updatePrimalBasedQuantities()
Update primal based quantities, e.g. the primal fields in adjoint turbulence models.
const word & primalSolverName() const
Return the primal solver name.
TypeName("adjointSolver")
Run-time type information.
virtual void preCalculateSensitivities()
Actions to be performed before calculating sensitivities.
bool isConstraint()
Is the solving referring to a constraint.
virtual void accumulateBCSensitivityIntegrand(autoPtr< boundaryVectorField > &bcDxDbMult, const labelHashSet &sensitivityPatchIDs, const scalar dt)
Contributions from boundary functions that inlcude geometric aspects in them and change when the geom...
const word primalSolverName_
Name of primal solver.
declareRunTimeNewSelectionTable(autoPtr, adjointSolver, adjointSolver,(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName, const word &solverName),(mesh, managerType, dict, primalSolverName, solverName))
bool isDoubleSidedConstraint_
Is the adjoint solver used to tackle a double-sided constraint.
bool isConstraint_
Is the adjoint solver used to tackle a constraint.
virtual dimensionSet maDimensions() const
Return the dimensions of the adjoint grid displacement variable.
const objectiveManager & getObjectiveManager() const
Return a const reference to the objective manager.
virtual void computeObjectiveSensitivities(autoPtr< designVariables > &designVars)
Compute sensitivities of the underlaying objectives.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Abstract base class for defining design variables.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
Class for managing objective functions.
Base class for primal solvers.
const word & managerType() const
Return the manager type.
const dictionary & dict() const
Return the solver dictionary.
const fvMesh & mesh() const
Return the solver mesh.
const word & solverName() const
Return the solver name.
solver(const solver &)=delete
No copy construct.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
OBJstream os(runTime.globalPath()/outputName)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Macros to ease declaration of run-time selection tables.
#define declareRunTimeNewSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection for derived classes.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.