37#ifndef incompressibleAdjointSolver_H
38#define incompressibleAdjointSolver_H
54class incompressibleAdjointSolver
63 incompressibleAdjointSolver
65 const incompressibleAdjointSolver&
69 void operator=(
const incompressibleAdjointSolver&) =
delete;
109 incompressibleAdjointSolver,
125 incompressibleAdjointSolver
261 const word& designVariablesName,
271 virtual bool write(
const bool valid =
true)
const
273 if (
mesh_.time().writeTime())
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
const word & primalSolverName() const
Return the primal solver name.
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,...
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.
autoPtr< ATCModel > ATCModel_
Adjoint Transpose Convection options.
bool hasBCdxdbMult(const labelHashSet &sensitivityPatchIDs)
Compute, if necessary, and return the hasBCdxdbMult_ bool.
virtual tmp< volScalarField > yWall() const
Return the distance field, to be used in the solution of the adjoint eikonal PDE.
static autoPtr< incompressibleAdjointSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName, const word &solverName)
Return a reference to the selected incompressible adjoint solver.
TypeName("incompressible")
Run-time type information.
declareRunTimeSelectionTable(autoPtr, incompressibleAdjointSolver, dictionary,(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName, const word &solverName),(mesh, managerType, dict, primalSolverName, solverName))
const incompressibleVars & getPrimalVars() const
Access to the incompressible primal variables set.
virtual tmp< volScalarField > adjointEikonalSource()
Return the source the adjoint eikonal equation.
const autoPtr< ATCModel > & getATCModel() const
Access to the ATC model.
virtual dimensionSet daDimensions() const
Return the dimensions of the adjoint distance field.
virtual ~incompressibleAdjointSolver()=default
Destructor.
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 void accumulateDivDxDbMultiplier(autoPtr< scalarField > &divDxDbMult, const scalar dt)
Compute the multiplier for div(dxdb).
virtual bool write(const bool valid=true) const
In case of multi-point runs with turbulent flows, output dummy turbulence fields with the base names,...
virtual const incompressibleAdjointVars & getAdjointVars() const
Access to the incompressible adjoint variables set.
virtual bool includeDistance() const
Should the adjoint to the eikonal equation be solved.
virtual void updatePrimalBasedQuantities()
Update primal based quantities, e.g. the primal fields in adjoint turbulence models.
virtual void accumulateGradDxDbMultiplier(volTensorField &gradDxDbMult, const scalar dt)
Compute the multiplier for grad(dxdb).
virtual bool writeNow() const
In case of multi-point runs with turbulent flows, output dummy turbulence fields with the base names,...
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...
virtual void topOSensMultiplier(scalarField &betaMult, const word &designVariablesName, const scalar dt)
Compute the multiplier of beta.
Switch hasBCdxdbMult_
Auxiliary bool to avoid a potentially expensive part of the sensitivity computation.
virtual void accumulateOptionsDxDbMultiplier(vectorField &optionsDxDbMult, const scalar dt)
Contributions from fvOptions that inlcude geometric aspects in them and change when the geometry is d...
virtual dimensionSet maDimensions() const
Return the dimensions of the adjoint grid displacement variable.
incompressibleVars & primalVars_
Primal variable set.
Class including all adjoint fields for incompressible flows.
Base class for solution control classes.
bool write() const
Write dummy turbulent fields to allow for continuation in multi-point, turbulent runs.
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.
fvMesh & mesh_
Reference to the mesh database.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
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.
#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.