40#ifndef objectiveIncompressible_H
41#define objectiveIncompressible_H
56class objectiveIncompressible
113 objectiveIncompressible(
const objectiveIncompressible&) =
delete;
116 void operator=(
const objectiveIncompressible&) =
delete;
130 objectiveIncompressible,
135 const word& adjointSolverName,
136 const word& primalSolverName
138 (
mesh,
dict, adjointSolverName, primalSolverName)
145 objectiveIncompressible
149 const word& adjointSolverName,
150 const word& primalSolverName
161 const word& adjointSolverName,
162 const word& primalSolverName
173 virtual scalar
J() = 0;
354 virtual bool write(
const bool valid =
true)
const;
380#include "objectiveIncompressibleI.H"
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,...
Mesh data needed to do the Finite Volume discretisation.
Abstract base class for incompressible turbulence models.
Base class for solution control classes.
bool hasdJdp() const noexcept
const volScalarField & dJdTMvar1()
Contribution to field adjoint turbulence model variable 1.
autoPtr< volScalarField > dJdTMvar2Ptr_
Second turbulence model variable.
const volScalarField & dJdT()
Contribution to field adjoint energy eq.
virtual void addSource(fvScalarMatrix &matrix)
Scalar sources are more ambigious (adjoint pressure, turbulence model, energy, etc equations),...
virtual void update_boundarydJdTMvar2()
virtual void update_gradDxDbMultiplier()
Update grad( dx/db multiplier). Volume-based sensitivity term.
static autoPtr< objectiveIncompressible > New(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Return a reference to the selected turbulence model.
virtual void update_dJdTMvar1()
bool hasdJdv() const noexcept
Inline functions for checking whether pointers are set or not.
bool hasdJdT() const noexcept
autoPtr< boundaryTensorField > bdJdGradUPtr_
Term multiplying gradU variations.
TypeName("incompressible")
Runtime type information.
autoPtr< boundaryScalarField > bdJdTMvar2Ptr_
Adjoint outlet turbulence model var 2.
virtual void update_boundarydJdvt()
autoPtr< volScalarField > dJdpPtr_
bool hasBoundarydJdv() const noexcept
virtual void update_dJdbField()
void update_dJdTMvar(autoPtr< volScalarField > &dJdTMvarPtr, tmp< volScalarField >(incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const, const volScalarField &JacobianMultiplier, const labelList &zones)
Compute dJdTMVar{1,2}.
bool hasBoundarydJdvn() const noexcept
virtual void update_divDxDbMultiplier()
Update div( dx/db multiplier). Volume-based sensitivity term.
virtual void addSource(fvVectorMatrix &matrix)
Vector sources can be given only to the adjoint momentum equations. Implemented in base objectiveInco...
const volScalarField & dJdTMvar2()
Contribution to field adjoint turbulence model variable 2.
void allocatedJdTurbulence()
Allocate fields related to the differentiation of turbulence models, if necessary.
virtual void update_dJdTMvar2()
virtual void update_dJdv()
Update vol and boundary fields and derivatives.
virtual void nullify()
Update objective function derivatives.
const boundaryVectorField & boundarydJdv()
Objective partial deriv wrt velocity for all patches.
virtual ~objectiveIncompressible()=default
Destructor.
const volVectorField & dJdv()
Contribution to field adjoint momentum eqs.
bool hasBoundarydJdT() const noexcept
virtual void update_boundarydJdp()
declareRunTimeSelectionTable(autoPtr, objectiveIncompressible, dictionary,(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName),(mesh, dict, adjointSolverName, primalSolverName))
virtual void update_meanValues()
Some objectives need to store some auxiliary values. If averaging is enabled, update these mean value...
autoPtr< volScalarField > dJdTMvar1Ptr_
First turbulence model variable.
const volScalarField & dJdp()
Contribution to field adjoint continuity eq.
autoPtr< boundaryScalarField > bdJdTMvar1Ptr_
Adjoint outlet turbulence model var 1.
bool hasdJdTMVar1() const noexcept
const boundaryScalarField & boundarydJdT()
Objective partial deriv wrt temperature for all patches.
const boundaryVectorField & boundarydJdp()
Objective partial deriv wrt pressure (times normal) for all patches.
const boundaryVectorField & boundarydJdvt()
Objective partial deriv wrt tangent velocity for all patches.
virtual scalar J()=0
Return the objective function value.
virtual void update_boundarydJdTMvar1()
autoPtr< boundaryScalarField > bdJdnutPtr_
Jacobian wrt to nut.
autoPtr< boundaryScalarField > bdJdTPtr_
Adjoint outlet temperature.
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
virtual void update_boundarydJdGradU()
bool hasdJdTMVar2() const noexcept
autoPtr< volScalarField > dJdTPtr_
const boundaryTensorField & boundarydJdGradU()
Objective partial deriv wrt gradU.
autoPtr< boundaryVectorField > bdJdpPtr_
Adjoint (intlet,wall) velocity.
virtual void doNormalization()
Normalize all fields allocated by the objective.
virtual void update_boundarydJdnut()
const boundaryScalarField & boundarydJdnut()
Objective partial deriv wrt nut for all patches.
virtual void update()
Update objective function derivatives.
virtual void update_boundarydJdv()
const boundaryScalarField & boundarydJdTMvar1()
Objective partial deriv wrt turbulence model var 1 for all patches.
virtual void update_dJdT()
const boundaryScalarField & boundarydJdTMvar2()
Objective partial deriv wrt turbulence model var 2 for all patches.
bool hasBoundarydJdTMVar2() const noexcept
autoPtr< boundaryVectorField > bdJdvtPtr_
Adjoint outlet velocity.
bool hasBoundarydJdTMVar1() const noexcept
const incompressibleVars & vars_
virtual void update_dJdb()
bool hasBoundarydJdGradU() const noexcept
virtual void update_boundarydJdvn()
autoPtr< boundaryVectorField > bdJdvPtr_
bool hasBoundarydJdnut() const noexcept
void checkCellZonesSize(const labelList &zoneIDs) const
Check if cellZones provided include at least one cell.
const boundaryScalarField & boundarydJdvn()
Objective partial deriv wrt normal velocity for all patches.
bool hasBoundarydJdp() const noexcept
virtual void update_boundarydJdT()
virtual void update_dJdp()
autoPtr< volVectorField > dJdvPtr_
bool hasBoundarydJdvt() const noexcept
const dictionary & dict() const
Return objective dictionary.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
const labelIOList & zoneIDs
volScalarField::Boundary boundaryScalarField
volFields
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
fvPatchField< tensor > fvPatchTensorField
volTensorField::Boundary boundaryTensorField
fvMatrix< vector > fvVectorMatrix
volVectorField::Boundary boundaryVectorField
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
#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.