objective Class Referenceabstract

Abstract base class for objective functions. No point in making this runTime selectable since its children will have different constructors. More...

Inheritance diagram for objective:
[legend]
Collaboration diagram for objective:
[legend]

Public Member Functions

 TypeName ("objective")
 Runtime type information. More...
 
 declareRunTimeNewSelectionTable (autoPtr, objective, objective,(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName),(mesh, dict, adjointSolverName, primalSolverName))
 
 objective (const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
 Construct from components. More...
 
virtual ~objective ()=default
 Destructor. More...
 
virtual bool readDict (const dictionary &dict)
 
virtual scalar J ()=0
 Return the instantaneous objective function value. More...
 
scalar JCycle () const
 
void accumulateJMean (solverControl &solverControl)
 Accumulate contribution for the mean objective value. More...
 
void accumulateJMean ()
 Accumulate contribution for the mean objective value. More...
 
scalar weight () const
 Return the objective function weight. More...
 
bool isWithinIntegrationTime () const
 Check whether this is an objective integration time. More...
 
void incrementIntegrationTimes (const scalar timeSpan)
 Increment integration times. More...
 
const volScalarFielddJdb ()
 Contribution to field sensitivities. More...
 
const fvPatchVectorFieldboundarydJdb (const label)
 Contribution to surface sensitivities for a specific patch. More...
 
const fvPatchVectorFielddSdbMultiplier (const label)
 Multiplier of delta(n dS)/delta b. More...
 
const fvPatchVectorFielddndbMultiplier (const label)
 Multiplier of delta(n dS)/delta b. More...
 
const fvPatchVectorFielddxdbMultiplier (const label)
 Multiplier of delta(x)/delta b. More...
 
const fvPatchVectorFielddxdbDirectMultiplier (const label)
 Multiplier of delta(x)/delta b. More...
 
const vectorFieldboundaryEdgeMultiplier (const label patchI, const label edgeI)
 Multiplier located at patch boundary edges. More...
 
const fvPatchTensorFieldboundarydJdStress (const label)
 Objective partial deriv wrt stress tensor. More...
 
const boundaryVectorFieldboundarydJdb ()
 Contribution to surface sensitivities for all patches. More...
 
const boundaryVectorFielddSdbMultiplier ()
 Multiplier of delta(n dS)/delta b for all patches. More...
 
const boundaryVectorFielddndbMultiplier ()
 Multiplier of delta(n dS)/delta b for all patches. More...
 
const boundaryVectorFielddxdbMultiplier ()
 Multiplier of delta(x)/delta b for all patches. More...
 
const boundaryVectorFielddxdbDirectMultiplier ()
 Multiplier of delta(x)/delta b for all patches. More...
 
const vectorField3boundaryEdgeMultiplier ()
 Multiplier located at patch boundary edges. More...
 
const boundaryTensorFieldboundarydJdStress ()
 Objective partial deriv wrt stress tensor. More...
 
const volScalarFielddivDxDbMultiplier ()
 Multiplier of grad( delta(x)/delta b) for volume-based sensitivities. More...
 
const volTensorFieldgradDxDbMultiplier ()
 Multiplier of grad( delta(x)/delta b) for volume-based sensitivities. More...
 
virtual void update ()=0
 Update objective function derivatives. More...
 
virtual void nullify ()
 Nullify adjoint contributions. More...
 
virtual void updateNormalizationFactor ()
 
virtual void update_boundarydJdb ()
 Update objective function derivative term. More...
 
virtual void update_dSdbMultiplier ()
 Update d (normal dS) / db multiplier. Surface-based sensitivity term. More...
 
virtual void update_dndbMultiplier ()
 Update d (normal) / db multiplier. Surface-based sensitivity term. More...
 
virtual void update_dxdbMultiplier ()
 Update d (x) / db multiplier. Surface-based sensitivity term. More...
 
virtual void update_dxdbDirectMultiplier ()
 
virtual void update_boundaryEdgeContribution ()
 Update boundary edge contributions. More...
 
virtual void update_dJdStressMultiplier ()
 Update dJ/dStress field. More...
 
virtual void update_divDxDbMultiplier ()
 Update div( dx/db multiplier). Volume-based sensitivity term. More...
 
virtual void update_gradDxDbMultiplier ()
 Update grad( dx/db multiplier). Volume-based sensitivity term. More...
 
virtual void write () const
 Write objective function history. More...
 
virtual void writeInstantaneousValue () const
 Write objective function history at each primal solver iteration. More...
 
virtual void writeMeanValue () const
 Write mean objective function history. More...
 
const wordobjectiveName () const
 
bool hasdJdb () const
 
bool hasBoundarydJdb () const
 
bool hasdSdbMult () const
 
bool hasdndbMult () const
 
bool hasdxdbMult () const
 
bool hasdxdbDirectMult () const
 
bool hasBoundaryEdgeContribution () const
 
bool hasBoundarydJdStress () const
 
bool hasDivDxDbMult () const
 
bool hasGradDxDbMult () const
 
bool hasIntegrationStartTime () const
 
bool hasIntegrationEndTime () const
 

Static Public Member Functions

static autoPtr< objectiveNew (const fvMesh &mesh, const dictionary &dict, const word &objectiveType, const word &adjointSolverName, const word &primalSolverName)
 Return a reference to the selected turbulence model. More...
 

Protected Member Functions

const dictionarydict () const
 Return objective dictionary. More...
 
void setObjectiveFilePtr () const
 Set the output file ptr. More...
 
void setInstantValueFilePtr () const
 Set the output file ptr for the instantaneous value. More...
 
void setMeanValueFilePtr () const
 Set the output file ptr for the mean value. More...
 

Protected Attributes

const fvMeshmesh_
 
dictionary dict_
 
const word adjointSolverName_
 
const word primalSolverName_
 
const word objectiveName_
 
bool computeMeanFields_
 
bool nullified_
 
scalar J_
 
scalar JMean_
 
scalar weight_
 
autoPtr< scalar > integrationStartTimePtr_
 
autoPtr< scalar > integrationEndTimePtr_
 
autoPtr< volScalarFielddJdbPtr_
 
autoPtr< boundaryVectorFieldbdJdbPtr_
 Term from material derivative. More...
 
autoPtr< boundaryVectorFieldbdSdbMultPtr_
 Term multiplying delta(n dS)/delta b. More...
 
autoPtr< boundaryVectorFieldbdndbMultPtr_
 Term multiplying delta(n)/delta b. More...
 
autoPtr< boundaryVectorFieldbdxdbMultPtr_
 Term multiplying delta(x)/delta b at the boundary. More...
 
autoPtr< boundaryVectorFieldbdxdbDirectMultPtr_
 
autoPtr< vectorField3bEdgeContribution_
 
autoPtr< boundaryTensorFieldbdJdStressPtr_
 For use with discrete-like sensitivities. More...
 
autoPtr< volScalarFielddivDxDbMultPtr_
 Multiplier of d(Volume)/db. More...
 
autoPtr< volTensorFieldgradDxDbMultPtr_
 Emerging from volume objectives that include spatial derivatives. More...
 
fileName objFunctionFolder_
 Output file variables. More...
 
autoPtr< OFstreamobjFunctionFilePtr_
 File to keep the objective values after the end of the primal solver. More...
 
autoPtr< OFstreaminstantValueFilePtr_
 
autoPtr< OFstreammeanValueFilePtr_
 

Detailed Description

Abstract base class for objective functions. No point in making this runTime selectable since its children will have different constructors.

Source files

Definition at line 58 of file objective.H.

Constructor & Destructor Documentation

◆ objective()

objective ( const fvMesh mesh,
const dictionary dict,
const word adjointSolverName,
const word primalSolverName 
)

Construct from components.

Definition at line 102 of file objective.C.

References dict, dictionary::found(), dictionary::get(), IOobject::NO_WRITE, and IOobject::READ_IF_PRESENT.

Here is the call graph for this function:

◆ ~objective()

virtual ~objective ( )
virtualdefault

Destructor.

Member Function Documentation

◆ dict()

const dictionary & dict ( ) const
protected

Return objective dictionary.

Definition at line 93 of file objective.C.

References objective::dict_.

Referenced by objectivePtLosses::initialize(), and objective::readDict().

Here is the caller graph for this function:

◆ setObjectiveFilePtr()

void setObjectiveFilePtr ( ) const
protected

Set the output file ptr.

Definition at line 58 of file objective.C.

References objective::adjointSolverName_, objective::objectiveName_, objective::objFunctionFilePtr_, and objective::objFunctionFolder_.

Referenced by objectiveForceTarget::write(), objectivePartialVolume::write(), objectivePtLosses::write(), and objective::write().

Here is the caller graph for this function:

◆ setInstantValueFilePtr()

void setInstantValueFilePtr ( ) const
protected

Set the output file ptr for the instantaneous value.

Definition at line 67 of file objective.C.

References objective::adjointSolverName_, objective::instantValueFilePtr_, objective::objectiveName_, and objective::objFunctionFolder_.

Referenced by objective::writeInstantaneousValue().

Here is the caller graph for this function:

◆ setMeanValueFilePtr()

void setMeanValueFilePtr ( ) const
protected

Set the output file ptr for the mean value.

Definition at line 79 of file objective.C.

References objective::adjointSolverName_, objective::meanValueFilePtr_, objective::objectiveName_, and objective::objFunctionFolder_.

Referenced by objective::writeMeanValue().

Here is the caller graph for this function:

◆ TypeName()

TypeName ( "objective"  )

Runtime type information.

◆ declareRunTimeNewSelectionTable()

declareRunTimeNewSelectionTable ( autoPtr  ,
objective  ,
objective  ,
(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)  ,
(mesh, dict, adjointSolverName, primalSolverName)   
)

◆ New()

autoPtr< objective > New ( const fvMesh mesh,
const dictionary dict,
const word objectiveType,
const word adjointSolverName,
const word primalSolverName 
)
static

Return a reference to the selected turbulence model.

Definition at line 182 of file objective.C.

References dict, Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, and mesh.

Referenced by objectiveManager::objectiveManager().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ readDict()

bool readDict ( const dictionary dict)
virtual

Definition at line 218 of file objective.C.

References objective::dict(), and objective::dict_.

Here is the call graph for this function:

◆ J()

virtual scalar J ( )
pure virtual

Return the instantaneous objective function value.

Implemented in objectiveIncompressible, objectiveMoment, objectiveForce, objectivePtLosses, objectiveForceTarget, and objectivePartialVolume.

Referenced by objective::JCycle().

Here is the caller graph for this function:

◆ JCycle()

scalar JCycle ( ) const

Return the mean objective function value, if it exists, otherwise the mean one

Definition at line 225 of file objective.C.

References objective::computeMeanFields_, objective::hasIntegrationEndTime(), objective::hasIntegrationStartTime(), objective::J(), objective::J_, and objective::JMean_.

Here is the call graph for this function:

◆ accumulateJMean() [1/2]

void accumulateJMean ( solverControl solverControl)

Accumulate contribution for the mean objective value.

For steady-state runs

Definition at line 246 of file objective.C.

References solverControl::averageIter(), solverControl::doAverageIter(), objective::J_, objective::JMean_, and Foam::Zero.

Here is the call graph for this function:

◆ accumulateJMean() [2/2]

void accumulateJMean ( )

Accumulate contribution for the mean objective value.

For unsteady runs

Definition at line 263 of file objective.C.

References TimeState::deltaT(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, objective::hasIntegrationEndTime(), objective::hasIntegrationStartTime(), objective::integrationStartTimePtr_, objective::isWithinIntegrationTime(), objective::J_, objective::JMean_, objective::mesh_, fvMesh::time(), and dimensioned< Type >::value().

Here is the call graph for this function:

◆ weight()

scalar weight ( ) const

Return the objective function weight.

Definition at line 285 of file objective.C.

References objective::weight_.

◆ isWithinIntegrationTime()

bool isWithinIntegrationTime ( ) const

Check whether this is an objective integration time.

Definition at line 291 of file objective.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, objective::hasIntegrationEndTime(), objective::hasIntegrationStartTime(), objective::integrationEndTimePtr_, objective::integrationStartTimePtr_, objective::mesh_, fvMesh::time(), and dimensioned< Type >::value().

Referenced by objective::accumulateJMean().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ incrementIntegrationTimes()

void incrementIntegrationTimes ( const scalar  timeSpan)

Increment integration times.

Definition at line 312 of file objective.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, objective::hasIntegrationEndTime(), objective::hasIntegrationStartTime(), objective::integrationEndTimePtr_, and objective::integrationStartTimePtr_.

Here is the call graph for this function:

◆ dJdb()

const volScalarField & dJdb ( )

Contribution to field sensitivities.

Definition at line 328 of file objective.C.

References objective::dJdbPtr_, objective::mesh_, and objective::objectiveName_.

◆ boundarydJdb() [1/2]

const fvPatchVectorField & boundarydJdb ( const label  patchI)

Contribution to surface sensitivities for a specific patch.

Definition at line 348 of file objective.C.

References objective::bdJdbPtr_, and objective::mesh_.

◆ dSdbMultiplier() [1/2]

const fvPatchVectorField & dSdbMultiplier ( const label  patchI)

Multiplier of delta(n dS)/delta b.

Definition at line 358 of file objective.C.

References objective::bdSdbMultPtr_, and objective::mesh_.

◆ dndbMultiplier() [1/2]

const fvPatchVectorField & dndbMultiplier ( const label  patchI)

Multiplier of delta(n dS)/delta b.

Definition at line 368 of file objective.C.

References objective::bdndbMultPtr_, and objective::mesh_.

◆ dxdbMultiplier() [1/2]

const fvPatchVectorField & dxdbMultiplier ( const label  patchI)

Multiplier of delta(x)/delta b.

Definition at line 378 of file objective.C.

References objective::bdxdbMultPtr_, and objective::mesh_.

◆ dxdbDirectMultiplier() [1/2]

const fvPatchVectorField & dxdbDirectMultiplier ( const label  patchI)

Multiplier of delta(x)/delta b.

Definition at line 388 of file objective.C.

References objective::bdxdbDirectMultPtr_, and objective::mesh_.

◆ boundaryEdgeMultiplier() [1/2]

const vectorField & boundaryEdgeMultiplier ( const label  patchI,
const label  edgeI 
)

Multiplier located at patch boundary edges.

Definition at line 399 of file objective.C.

References Foam::exit(), Foam::FatalError, and FatalErrorInFunction.

Here is the call graph for this function:

◆ boundarydJdStress() [1/2]

const fvPatchTensorField & boundarydJdStress ( const label  patchI)

Objective partial deriv wrt stress tensor.

Definition at line 414 of file objective.C.

References objective::bdJdStressPtr_, and objective::mesh_.

◆ boundarydJdb() [2/2]

const boundaryVectorField & boundarydJdb ( )

Contribution to surface sensitivities for all patches.

Definition at line 424 of file objective.C.

References objective::bdJdbPtr_, and objective::mesh_.

◆ dSdbMultiplier() [2/2]

const boundaryVectorField & dSdbMultiplier ( )

Multiplier of delta(n dS)/delta b for all patches.

Definition at line 434 of file objective.C.

References objective::bdSdbMultPtr_, and objective::mesh_.

◆ dndbMultiplier() [2/2]

const boundaryVectorField & dndbMultiplier ( )

Multiplier of delta(n dS)/delta b for all patches.

Definition at line 444 of file objective.C.

References objective::bdndbMultPtr_, and objective::mesh_.

◆ dxdbMultiplier() [2/2]

const boundaryVectorField & dxdbMultiplier ( )

Multiplier of delta(x)/delta b for all patches.

Definition at line 454 of file objective.C.

References objective::bdxdbMultPtr_, and objective::mesh_.

◆ dxdbDirectMultiplier() [2/2]

const boundaryVectorField & dxdbDirectMultiplier ( )

Multiplier of delta(x)/delta b for all patches.

Definition at line 464 of file objective.C.

References objective::bdxdbDirectMultPtr_, and objective::mesh_.

◆ boundaryEdgeMultiplier() [2/2]

const vectorField3 & boundaryEdgeMultiplier ( )

Multiplier located at patch boundary edges.

Definition at line 474 of file objective.C.

References objective::bdxdbDirectMultPtr_, objective::bEdgeContribution_, Foam::endl(), Foam::exit(), Foam::FatalError, and FatalErrorInFunction.

Here is the call graph for this function:

◆ boundarydJdStress() [2/2]

const boundaryTensorField & boundarydJdStress ( )

Objective partial deriv wrt stress tensor.

Definition at line 487 of file objective.C.

References objective::bdJdStressPtr_, and objective::mesh_.

◆ divDxDbMultiplier()

const volScalarField & divDxDbMultiplier ( )

Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.

Definition at line 497 of file objective.C.

References Foam::dimless, objective::divDxDbMultPtr_, objective::mesh_, and objective::objectiveName_.

◆ gradDxDbMultiplier()

const volTensorField & gradDxDbMultiplier ( )

Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.

Definition at line 518 of file objective.C.

References Foam::dimLength, Foam::dimTime, objective::gradDxDbMultPtr_, objective::mesh_, objective::objectiveName_, Foam::pow2(), and Foam::pow3().

Here is the call graph for this function:

◆ update()

virtual void update ( )
pure virtual

Update objective function derivatives.

Implemented in objectiveIncompressible.

◆ nullify()

void nullify ( )
virtual

◆ updateNormalizationFactor()

void updateNormalizationFactor ( )
virtual

Update normalization factors, for objectives in which the factor is not known a priori

Definition at line 240 of file objective.C.

◆ update_boundarydJdb()

virtual void update_boundarydJdb ( )
inlinevirtual

Update objective function derivative term.

Reimplemented in objectiveIncompressible.

Definition at line 315 of file objective.H.

◆ update_dSdbMultiplier()

virtual void update_dSdbMultiplier ( )
inlinevirtual

Update d (normal dS) / db multiplier. Surface-based sensitivity term.

Reimplemented in objectiveIncompressible, objectiveMoment, objectiveForce, and objectivePartialVolume.

Definition at line 319 of file objective.H.

◆ update_dndbMultiplier()

virtual void update_dndbMultiplier ( )
inlinevirtual

Update d (normal) / db multiplier. Surface-based sensitivity term.

Reimplemented in objectiveIncompressible.

Definition at line 323 of file objective.H.

◆ update_dxdbMultiplier()

virtual void update_dxdbMultiplier ( )
inlinevirtual

Update d (x) / db multiplier. Surface-based sensitivity term.

Reimplemented in objectiveIncompressible, objectiveMoment, and objectiveForce.

Definition at line 327 of file objective.H.

◆ update_dxdbDirectMultiplier()

virtual void update_dxdbDirectMultiplier ( )
inlinevirtual

Update d (x) / db multiplier. Surface and volume-based sensitivity term

Reimplemented in objectiveIncompressible, objectiveMoment, and objectivePartialVolume.

Definition at line 332 of file objective.H.

◆ update_boundaryEdgeContribution()

virtual void update_boundaryEdgeContribution ( )
inlinevirtual

Update boundary edge contributions.

Definition at line 336 of file objective.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_dJdStressMultiplier()

virtual void update_dJdStressMultiplier ( )
inlinevirtual

Update dJ/dStress field.

Reimplemented in objectiveForce.

Definition at line 340 of file objective.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_divDxDbMultiplier()

virtual void update_divDxDbMultiplier ( )
inlinevirtual

Update div( dx/db multiplier). Volume-based sensitivity term.

Reimplemented in objectiveIncompressible.

Definition at line 344 of file objective.H.

◆ update_gradDxDbMultiplier()

virtual void update_gradDxDbMultiplier ( )
inlinevirtual

Update grad( dx/db multiplier). Volume-based sensitivity term.

Reimplemented in objectiveIncompressible.

Definition at line 348 of file objective.H.

◆ write()

void write ( ) const
virtual

Write objective function history.

Reimplemented in objectiveIncompressible, objectivePtLosses, objectivePartialVolume, and objectiveForceTarget.

Definition at line 593 of file objective.C.

References Foam::endl(), objective::J_, UPstream::master(), objective::mesh_, objective::objFunctionFilePtr_, objective::setObjectiveFilePtr(), Foam::tab, fvMesh::time(), and dimensioned< Type >::value().

Referenced by objectiveIncompressible::write().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ writeInstantaneousValue()

void writeInstantaneousValue ( ) const
virtual

Write objective function history at each primal solver iteration.

Definition at line 610 of file objective.C.

References Foam::endl(), objective::instantValueFilePtr_, objective::J_, UPstream::master(), objective::mesh_, objective::setInstantValueFilePtr(), Foam::tab, fvMesh::time(), and dimensioned< Type >::value().

Here is the call graph for this function:

◆ writeMeanValue()

void writeMeanValue ( ) const
virtual

Write mean objective function history.

Definition at line 627 of file objective.C.

References dictionary::add(), objective::computeMeanFields_, Foam::endl(), objective::hasIntegrationEndTime(), objective::hasIntegrationStartTime(), objective::JMean_, UPstream::master(), objective::meanValueFilePtr_, objective::mesh_, IOobject::NO_READ, IOobject::NO_WRITE, objective::objectiveName_, objective::setMeanValueFilePtr(), Foam::tab, fvMesh::time(), Time::timeName(), and dimensioned< Type >::value().

Here is the call graph for this function:

◆ objectiveName()

const Foam::word & objectiveName ( ) const
inline

Definition at line 33 of file objectiveI.H.

References objective::objectiveName_.

◆ hasdJdb()

bool hasdJdb ( ) const
inline

Definition at line 39 of file objectiveI.H.

Referenced by objective::nullify().

Here is the caller graph for this function:

◆ hasBoundarydJdb()

bool hasBoundarydJdb ( ) const
inline

Definition at line 45 of file objectiveI.H.

Referenced by objective::nullify().

Here is the caller graph for this function:

◆ hasdSdbMult()

bool hasdSdbMult ( ) const
inline

Definition at line 51 of file objectiveI.H.

Referenced by objective::nullify().

Here is the caller graph for this function:

◆ hasdndbMult()

bool hasdndbMult ( ) const
inline

Definition at line 57 of file objectiveI.H.

Referenced by objective::nullify().

Here is the caller graph for this function:

◆ hasdxdbMult()

bool hasdxdbMult ( ) const
inline

Definition at line 63 of file objectiveI.H.

Referenced by objective::nullify().

Here is the caller graph for this function:

◆ hasdxdbDirectMult()

bool hasdxdbDirectMult ( ) const
inline

Definition at line 69 of file objectiveI.H.

Referenced by objective::nullify().

Here is the caller graph for this function:

◆ hasBoundaryEdgeContribution()

bool hasBoundaryEdgeContribution ( ) const
inline

Definition at line 75 of file objectiveI.H.

Referenced by objective::nullify().

Here is the caller graph for this function:

◆ hasBoundarydJdStress()

bool hasBoundarydJdStress ( ) const
inline

Definition at line 93 of file objectiveI.H.

Referenced by objective::nullify().

Here is the caller graph for this function:

◆ hasDivDxDbMult()

bool hasDivDxDbMult ( ) const
inline

Definition at line 81 of file objectiveI.H.

Referenced by objective::nullify().

Here is the caller graph for this function:

◆ hasGradDxDbMult()

bool hasGradDxDbMult ( ) const
inline

Definition at line 87 of file objectiveI.H.

Referenced by objective::nullify().

Here is the caller graph for this function:

◆ hasIntegrationStartTime()

bool hasIntegrationStartTime ( ) const
inline

Definition at line 99 of file objectiveI.H.

Referenced by objective::accumulateJMean(), objective::incrementIntegrationTimes(), objective::isWithinIntegrationTime(), objective::JCycle(), and objective::writeMeanValue().

Here is the caller graph for this function:

◆ hasIntegrationEndTime()

bool hasIntegrationEndTime ( ) const
inline

Definition at line 105 of file objectiveI.H.

Referenced by objective::accumulateJMean(), objective::incrementIntegrationTimes(), objective::isWithinIntegrationTime(), objective::JCycle(), and objective::writeMeanValue().

Here is the caller graph for this function:

Member Data Documentation

◆ mesh_

const fvMesh& mesh_
protected

Definition at line 64 of file objective.H.

Referenced by objective::accumulateJMean(), objective::boundarydJdb(), objectiveIncompressible::boundarydJdp(), objective::boundarydJdStress(), objectiveIncompressible::boundarydJdT(), objectiveIncompressible::boundarydJdTMvar1(), objectiveIncompressible::boundarydJdTMvar2(), objectiveIncompressible::boundarydJdv(), objectiveIncompressible::boundarydJdvn(), objectiveIncompressible::boundarydJdvt(), objective::divDxDbMultiplier(), objective::dJdb(), objectiveIncompressible::dJdp(), objectiveIncompressible::dJdT(), objectiveIncompressible::dJdTMvar1(), objectiveIncompressible::dJdTMvar2(), objectiveIncompressible::dJdv(), objective::dndbMultiplier(), objective::dSdbMultiplier(), objective::dxdbDirectMultiplier(), objective::dxdbMultiplier(), objective::gradDxDbMultiplier(), objectivePtLosses::initialize(), objective::isWithinIntegrationTime(), objectivePartialVolume::J(), objectivePtLosses::J(), objectiveForce::J(), objectiveMoment::J(), objectivePtLosses::update_boundarydJdp(), objectiveMoment::update_boundarydJdp(), objectivePtLosses::update_boundarydJdv(), objectivePtLosses::update_boundarydJdvn(), objectivePtLosses::update_boundarydJdvt(), objectiveForce::update_dJdStressMultiplier(), objectivePartialVolume::update_dSdbMultiplier(), objectiveMoment::update_dSdbMultiplier(), objectivePartialVolume::update_dxdbDirectMultiplier(), objectiveMoment::update_dxdbDirectMultiplier(), objectiveForce::update_dxdbMultiplier(), objectiveMoment::update_dxdbMultiplier(), objectiveForceTarget::write(), objectivePartialVolume::write(), objectivePtLosses::write(), objective::write(), objective::writeInstantaneousValue(), and objective::writeMeanValue().

◆ dict_

dictionary dict_
protected

Definition at line 65 of file objective.H.

Referenced by objective::dict(), and objective::readDict().

◆ adjointSolverName_

const word adjointSolverName_
protected

◆ primalSolverName_

const word primalSolverName_
protected

Definition at line 67 of file objective.H.

◆ objectiveName_

◆ computeMeanFields_

bool computeMeanFields_
protected

◆ nullified_

bool nullified_
protected

Definition at line 70 of file objective.H.

Referenced by objectiveIncompressible::nullify(), and objective::nullify().

◆ J_

◆ JMean_

scalar JMean_
protected

◆ weight_

scalar weight_
protected

Definition at line 75 of file objective.H.

Referenced by objective::weight().

◆ integrationStartTimePtr_

autoPtr<scalar> integrationStartTimePtr_
protected

◆ integrationEndTimePtr_

autoPtr<scalar> integrationEndTimePtr_
protected

◆ dJdbPtr_

autoPtr<volScalarField> dJdbPtr_
protected

Definition at line 85 of file objective.H.

Referenced by objective::dJdb(), and objective::nullify().

◆ bdJdbPtr_

autoPtr<boundaryVectorField> bdJdbPtr_
protected

Term from material derivative.

Definition at line 91 of file objective.H.

Referenced by objective::boundarydJdb(), and objective::nullify().

◆ bdSdbMultPtr_

◆ bdndbMultPtr_

autoPtr<boundaryVectorField> bdndbMultPtr_
protected

Term multiplying delta(n)/delta b.

Definition at line 97 of file objective.H.

Referenced by objective::dndbMultiplier(), and objective::nullify().

◆ bdxdbMultPtr_

autoPtr<boundaryVectorField> bdxdbMultPtr_
protected

Term multiplying delta(x)/delta b at the boundary.

Definition at line 100 of file objective.H.

Referenced by objective::dxdbMultiplier(), objective::nullify(), objectiveForce::update_dxdbMultiplier(), and objectiveMoment::update_dxdbMultiplier().

◆ bdxdbDirectMultPtr_

autoPtr<boundaryVectorField> bdxdbDirectMultPtr_
protected

Term multiplying delta(x)/delta b at the boundary for objectives that directly depend on x, e.g. moment Needed in both FI and SI computations

Definition at line 105 of file objective.H.

Referenced by objective::boundaryEdgeMultiplier(), objective::dxdbDirectMultiplier(), objective::nullify(), objectivePartialVolume::update_dxdbDirectMultiplier(), and objectiveMoment::update_dxdbDirectMultiplier().

◆ bEdgeContribution_

autoPtr<vectorField3> bEdgeContribution_
protected

Contribution located in specific parts of a patch. Mainly intended for patch boundary edges contributions, e.g. NURBS surfaces G1 continuity

Definition at line 110 of file objective.H.

Referenced by objective::boundaryEdgeMultiplier(), and objective::nullify().

◆ bdJdStressPtr_

autoPtr<boundaryTensorField> bdJdStressPtr_
protected

For use with discrete-like sensitivities.

Definition at line 113 of file objective.H.

Referenced by objective::boundarydJdStress(), objective::nullify(), and objectiveForce::update_dJdStressMultiplier().

◆ divDxDbMultPtr_

autoPtr<volScalarField> divDxDbMultPtr_
protected

Multiplier of d(Volume)/db.

Definition at line 120 of file objective.H.

Referenced by objective::divDxDbMultiplier(), and objective::nullify().

◆ gradDxDbMultPtr_

autoPtr<volTensorField> gradDxDbMultPtr_
protected

Emerging from volume objectives that include spatial derivatives.

Definition at line 123 of file objective.H.

Referenced by objective::gradDxDbMultiplier(), and objective::nullify().

◆ objFunctionFolder_

fileName objFunctionFolder_
protected

Output file variables.

Definition at line 126 of file objective.H.

Referenced by objective::setInstantValueFilePtr(), objective::setMeanValueFilePtr(), and objective::setObjectiveFilePtr().

◆ objFunctionFilePtr_

autoPtr<OFstream> objFunctionFilePtr_
mutableprotected

File to keep the objective values after the end of the primal solver.

Definition at line 129 of file objective.H.

Referenced by objective::setObjectiveFilePtr(), objectiveForceTarget::write(), objectivePartialVolume::write(), objectivePtLosses::write(), and objective::write().

◆ instantValueFilePtr_

autoPtr<OFstream> instantValueFilePtr_
mutableprotected

File to keep the objective values at each iteration of the primal solver

Definition at line 133 of file objective.H.

Referenced by objective::setInstantValueFilePtr(), and objective::writeInstantaneousValue().

◆ meanValueFilePtr_

autoPtr<OFstream> meanValueFilePtr_
mutableprotected

File to keep the average objective values after the end of the primal solver

Definition at line 137 of file objective.H.

Referenced by objective::setMeanValueFilePtr(), and objective::writeMeanValue().


The documentation for this class was generated from the following files: