201 void operator=(
const objective&) =
delete;
223 const word& adjointSolverName,
224 const word& primalSolverName
237 const word& adjointSolverName,
238 const word& primalSolverName
249 const word& objectiveType,
250 const word& adjointSolverName,
251 const word& primalSolverName
264 virtual scalar
J() = 0;
281 virtual scalar
weight()
const;
351 virtual void update() = 0;
410 virtual bool write(
const bool valid =
true)
const;
477 #include "objectiveI.H" void setInstantValueFilePtr() const
Set the output file ptr for the instantaneous value.
virtual void update_dSdbMultiplier()
Update d (normal dS) / db multiplier. Surface-based sensitivity term.
virtual void updateNormalizationFactor()
Update normalization factors, for objectives in which the factor is not known a priori.
const scalarField & dJdbField() const
Contribution to sensitivities with a random number of designVars.
autoPtr< boundaryVectorField > bdxdbMultPtr_
Term multiplying delta(x)/delta b at the boundary.
fvMatrix< scalar > fvScalarMatrix
autoPtr< scalar > target_
Target value, in case the objective is used as a constraint.
A class for handling file names.
virtual void incrementIntegrationTimes(const scalar timeSpan)
Increment integration times.
virtual void update_dJdb()
const word adjointSolverName_
virtual void addSource(fvVectorMatrix &matrix)
Manipulate fvVectorMatrix through the objective.
const dictionary & dict() const
Return objective dictionary.
virtual void writeMeanValue() const
Write mean objective function history.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
const vectorField3 & boundaryEdgeMultiplier() const
Multiplier located at patch boundary edges.
bool hasdxdbDirectMult() const noexcept
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
bool hasBoundarydJdb() const noexcept
bool hasdJdbField() const noexcept
autoPtr< scalar > normFactor_
Normalization factor.
virtual void update_dxdbMultiplier()
Update d (x) / db multiplier. Surface-based sensitivity term.
autoPtr< scalar > targetLeft_
Target on the left hand-side of a double inequality, for double sided constraints.
void setWrite(const bool shouldWrite)
Should the objective be written to file upon calling write()?
bool hasdndbMult() const noexcept
virtual void doNormalization()
Normalize all fields allocated by the objective.
virtual scalar J()=0
Return the instantaneous objective function value.
virtual void accumulateJMean()
Accumulate contribution for the mean objective value.
autoPtr< vectorField3 > bEdgeContribution_
Contribution located in specific parts of a patch. Mainly intended for patch boundary edges contribut...
autoPtr< boundaryVectorField > bdJdbPtr_
Term from material derivative.
const boundaryVectorField & dndbMultiplier() const
Multiplier of delta(n dS)/delta b for all patches.
bool hasBoundaryEdgeContribution() const noexcept
const word objectiveName_
virtual void update_gradDxDbMultiplier()
Update grad( dx/db multiplier). Volume-based sensitivity term.
wordList fieldNames_
List of adjoint fields for which this objective will contribute sources to their equations.
localIOdictionary is derived from IOdictionary but excludes parallel master reading.
dictionary()
Default construct, a top-level empty dictionary.
bool hasGradDxDbMult() const noexcept
virtual void update_divDxDbMultiplier()
Update div( dx/db multiplier). Volume-based sensitivity term.
virtual bool write(const bool valid=true) const
Write objective function history.
bool hasdSdbMult() const noexcept
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.
autoPtr< scalar > integrationEndTimePtr_
const boundaryVectorField & dxdbMultiplier() const
Multiplier of delta(x)/delta b for all patches.
const volScalarField & dJdb() const
Contribution to field sensitivities.
void setObjectiveFilePtr() const
Set the output file ptr.
A class for handling words, derived from Foam::string.
void setMeanValueFilePtr() const
Set the output file ptr for the mean value.
virtual void addHeaderColumns() const
Write headers for additional columns.
virtual void update_dndbMultiplier()
Update d (normal) / db multiplier. Surface-based sensitivity term.
virtual ~objective()=default
Destructor.
Base class for solver control classes.
bool computed_
Whether the objective is computed or not.
bool hasdJdb() const noexcept
virtual void update_boundarydJdb()
Update objective function derivative term.
const boundaryVectorField & boundarydJdb() const
Contribution to surface sensitivities for all patches.
void negate(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1)
scalar JMean_
Average objective value.
scalar J_
Objective function value and weight.
TypeName("objective")
Runtime type information.
const volScalarField & divDxDbMultiplier() const
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
virtual scalar weight() const
Return the objective function weight.
virtual void update_dJdbField()
virtual scalar JCycle(bool negate=false) const
Return the objective function of the optimisation cycle.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
autoPtr< boundaryVectorField > bdndbMultPtr_
Term multiplying delta(n)/delta b.
virtual bool writeData(Ostream &os) const
Write averaged objective for continuation.
virtual void writeInstantaneousValue() const
Write objective function history at each primal solver iteration.
OBJstream os(runTime.globalPath()/outputName)
const word primalSolverName_
const boundaryVectorField & dSdbMultiplier() const
Multiplier of delta(n dS)/delta b for all patches.
const boundaryVectorField & dxdbDirectMultiplier() const
Multiplier of delta(x)/delta b for all patches.
Useful typenames for fields defined only at the boundaries.
virtual bool isWithinIntegrationTime() const
Check whether this is an objective integration time.
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
virtual void addColumnValues() const
Write information to additional columns.
autoPtr< OFstream > meanValueFilePtr_
File to keep the average objective values after the end of the primal solver.
unsigned int width_
Default width of entries when writing in the objective files.
virtual void update_boundaryEdgeContribution()
Update boundary edge contributions.
bool hasIntegrationStartTime() const noexcept
virtual void nullify()
Nullify adjoint contributions.
autoPtr< scalarField > dJdbFieldPtr_
Contribution to field sensitivity derivatives.
virtual void writeInstantaneousSeparator() const
Append a blank line after the end of optimisation cycle to the file holding the instantaneous objecti...
bool shouldWrite() const
Should the objective be written to file upon calling write()?
virtual bool normalize() const
Is the objective normalized.
scalar weight_
Objective weight.
Mesh data needed to do the Finite Volume discretisation.
fvMatrix< vector > fvVectorMatrix
virtual void addHeaderInfo() const
Write any information that needs to go the header of the file.
virtual bool readDict(const dictionary &dict)
bool hasDivDxDbMult() const noexcept
void setComputed(const bool isComputed) noexcept
Set the computed status of the objective.
const volTensorField & gradDxDbMultiplier() const
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
bool hasdxdbMult() const noexcept
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
autoPtr< scalar > integrationStartTimePtr_
Objective integration start and end times (for unsteady flows)
declareRunTimeNewSelectionTable(autoPtr, objective, objective,(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName),(mesh, dict, adjointSolverName, primalSolverName))
virtual void update()=0
Update objective function derivatives.
autoPtr< OFstream > instantValueFilePtr_
File to keep the objective values at each iteration of the primal solver.
Macros to ease declaration of run-time selection tables.
autoPtr< boundaryVectorField > bdxdbDirectMultPtr_
Term multiplying delta(x)/delta b at the boundary for objectives that directly depend on x...
virtual void update_dxdbDirectMultiplier()
Update d (x) / db multiplier. Surface and volume-based sensitivity term.
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
fileName objFunctionFolder_
Output file variables.
bool hasIntegrationEndTime() const noexcept
const word & objectiveName() const
Return the objective name.