Go to the documentation of this file.
164 void operator=(
const objective&) =
delete;
186 const word& adjointSolverName,
187 const word& primalSolverName
189 (
mesh,
dict, adjointSolverName, primalSolverName)
200 const word& adjointSolverName,
201 const word& primalSolverName
212 const word& objectiveType,
213 const word& adjointSolverName,
214 const word& primalSolverName
227 virtual scalar
J() = 0;
306 virtual void update() = 0;
353 virtual void write()
const;
virtual void nullify()
Nullify adjoint contributions.
autoPtr< boundaryVectorField > bdxdbDirectMultPtr_
autoPtr< boundaryVectorField > bdJdbPtr_
Term from material derivative.
const vectorField3 & boundaryEdgeMultiplier()
Multiplier located at patch boundary edges.
bool hasBoundaryEdgeContribution() const
autoPtr< scalar > integrationStartTimePtr_
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
A class for handling words, derived from Foam::string.
A class for handling file names.
declareRunTimeNewSelectionTable(autoPtr, objective, objective,(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName),(mesh, dict, adjointSolverName, primalSolverName))
autoPtr< boundaryTensorField > bdJdStressPtr_
For use with discrete-like sensitivities.
const word objectiveName_
const boundaryVectorField & dndbMultiplier()
Multiplier of delta(n dS)/delta b for all patches.
const boundaryVectorField & dxdbDirectMultiplier()
Multiplier of delta(x)/delta b for all patches.
static autoPtr< objective > New(const fvMesh &mesh, const dictionary &dict, const word &objectiveType, const word &adjointSolverName, const word &primalSolverName)
Return a reference to the selected turbulence model.
void setInstantValueFilePtr() const
Set the output file ptr for the instantaneous value.
autoPtr< boundaryVectorField > bdxdbMultPtr_
Term multiplying delta(x)/delta b at the boundary.
void setObjectiveFilePtr() const
Set the output file ptr.
virtual void update_divDxDbMultiplier()
Update div( dx/db multiplier). Volume-based sensitivity term.
Base class for solver control classes.
virtual void update_boundarydJdb()
Update objective function derivative term.
virtual void updateNormalizationFactor()
virtual void update_dSdbMultiplier()
Update d (normal dS) / db multiplier. Surface-based sensitivity term.
virtual void update_dxdbMultiplier()
Update d (x) / db multiplier. Surface-based sensitivity term.
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
void setMeanValueFilePtr() const
Set the output file ptr for the mean value.
const word adjointSolverName_
virtual void writeMeanValue() const
Write mean objective function history.
virtual ~objective()=default
Destructor.
virtual void writeInstantaneousValue() const
Write objective function history at each primal solver iteration.
autoPtr< volScalarField > dJdbPtr_
bool hasIntegrationStartTime() const
bool hasIntegrationEndTime() const
bool hasDivDxDbMult() const
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
autoPtr< vectorField3 > bEdgeContribution_
void accumulateJMean()
Accumulate contribution for the mean objective value.
const volScalarField & divDxDbMultiplier()
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
const volTensorField & gradDxDbMultiplier()
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
autoPtr< boundaryVectorField > bdndbMultPtr_
Term multiplying delta(n)/delta b.
virtual void update_boundaryEdgeContribution()
Update boundary edge contributions.
const boundaryVectorField & dxdbMultiplier()
Multiplier of delta(x)/delta b for all patches.
const word primalSolverName_
virtual void update_gradDxDbMultiplier()
Update grad( dx/db multiplier). Volume-based sensitivity term.
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
autoPtr< OFstream > meanValueFilePtr_
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
autoPtr< scalar > integrationEndTimePtr_
const boundaryTensorField & boundarydJdStress()
Objective partial deriv wrt stress tensor.
autoPtr< OFstream > instantValueFilePtr_
Mesh data needed to do the Finite Volume discretisation.
bool hasBoundarydJdb() const
bool hasBoundarydJdStress() const
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
virtual void update_dndbMultiplier()
Update d (normal) / db multiplier. Surface-based sensitivity term.
const word & objectiveName() const
virtual void update_dxdbDirectMultiplier()
Useful typenames for fields defined only at the boundaries.
Macros to ease declaration of run-time selection tables.
const boundaryVectorField & dSdbMultiplier()
Multiplier of delta(n dS)/delta b for all patches.
fileName objFunctionFolder_
Output file variables.
bool hasGradDxDbMult() const
const boundaryVectorField & boundarydJdb()
Contribution to surface sensitivities for all patches.
const dictionary & dict() const
Return objective dictionary.
virtual void update_dJdStressMultiplier()
Update dJ/dStress field.
TypeName("objective")
Runtime type information.
scalar weight() const
Return the objective function weight.
bool hasdxdbDirectMult() const
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
void incrementIntegrationTimes(const scalar timeSpan)
Increment integration times.
virtual bool readDict(const dictionary &dict)
bool isWithinIntegrationTime() const
Check whether this is an objective integration time.
virtual scalar J()=0
Return the instantaneous objective function value.
const volScalarField & dJdb()
Contribution to field sensitivities.
virtual void write() const
Write objective function history.
virtual void update()=0
Update objective function derivatives.