objectiveIncompressible Class Referenceabstract

Abstract base class for objective functions in incompressible flows. More...

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

Public Member Functions

 TypeName ("incompressible")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, objectiveIncompressible, dictionary,(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName),(mesh, dict, adjointSolverName, primalSolverName))
 
 objectiveIncompressible (const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
 Construct from components. More...
 
virtual ~objectiveIncompressible ()=default
 Destructor. More...
 
virtual scalar J ()=0
 Return the objective function value. More...
 
const volVectorFielddJdv ()
 Contribution to field adjoint momentum eqs. More...
 
const volScalarFielddJdp ()
 Contribution to field adjoint continuity eq. More...
 
const volScalarFielddJdT ()
 Contribution to field adjoint energy eq. More...
 
const volScalarFielddJdTMvar1 ()
 Contribution to field adjoint turbulence model variable 1. More...
 
const volScalarFielddJdTMvar2 ()
 Contribution to field adjoint turbulence model variable 2. More...
 
const fvPatchVectorFieldboundarydJdv (const label)
 Objective partial deriv wrt velocity for a specific patch. More...
 
const fvPatchScalarFieldboundarydJdvn (const label)
 Objective partial deriv wrt normal velocity for a specific patch. More...
 
const fvPatchVectorFieldboundarydJdvt (const label)
 Objective partial deriv wrt tangent velocity for a specific patch. More...
 
const fvPatchVectorFieldboundarydJdp (const label)
 
const fvPatchScalarFieldboundarydJdT (const label)
 Objective partial deriv wrt temperature for a specific patch. More...
 
const fvPatchScalarFieldboundarydJdTMvar1 (const label)
 
const fvPatchScalarFieldboundarydJdTMvar2 (const label)
 
const boundaryVectorFieldboundarydJdv ()
 Objective partial deriv wrt velocity for all patches. More...
 
const boundaryScalarFieldboundarydJdvn ()
 Objective partial deriv wrt normal velocity for all patches. More...
 
const boundaryVectorFieldboundarydJdvt ()
 Objective partial deriv wrt tangent velocity for all patches. More...
 
const boundaryVectorFieldboundarydJdp ()
 Objective partial deriv wrt pressure (times normal) for all patches. More...
 
const boundaryScalarFieldboundarydJdT ()
 Objective partial deriv wrt temperature for all patches. More...
 
const boundaryScalarFieldboundarydJdTMvar1 ()
 Objective partial deriv wrt turbulence model var 1 for all patches. More...
 
const boundaryScalarFieldboundarydJdTMvar2 ()
 Objective partial deriv wrt turbulence model var 2 for all patches. More...
 
virtual void update ()
 Update objective function derivatives. More...
 
virtual void nullify ()
 Update objective function derivatives. More...
 
virtual void update_dJdv ()
 Update vol and boundary fields and derivatives. More...
 
virtual void update_dJdp ()
 
virtual void update_dJdT ()
 
virtual void update_dJdTMvar1 ()
 
virtual void update_dJdTMvar2 ()
 
virtual void update_dJdb ()
 
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 update_boundarydJdv ()
 
virtual void update_boundarydJdvn ()
 
virtual void update_boundarydJdvt ()
 
virtual void update_boundarydJdp ()
 
virtual void update_boundarydJdT ()
 
virtual void update_boundarydJdTMvar1 ()
 
virtual void update_boundarydJdTMvar2 ()
 
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_meanValues ()
 
virtual void write () const
 Write objective function history. More...
 
bool hasdJdv () const
 Inline functions for checking whether pointers are set or not. More...
 
bool hasdJdp () const
 
bool hasdJdT () const
 
bool hasdJdTMVar1 () const
 
bool hasdJdTMVar2 () const
 
bool hasBoundarydJdv () const
 
bool hasBoundarydJdvn () const
 
bool hasBoundarydJdvt () const
 
bool hasBoundarydJdp () const
 
bool hasBoundarydJdT () const
 
bool hasBoundarydJdTMVar1 () const
 
bool hasBoundarydJdTMVar2 () const
 
- Public Member Functions inherited from objective
 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)
 
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 updateNormalizationFactor ()
 
virtual void update_boundaryEdgeContribution ()
 Update boundary edge contributions. More...
 
virtual void update_dJdStressMultiplier ()
 Update dJ/dStress field. 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< objectiveIncompressibleNew (const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
 Return a reference to the selected turbulence model. More...
 
- Static Public Member Functions inherited from objective
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 Attributes

const incompressibleVarsvars_
 
autoPtr< volVectorFielddJdvPtr_
 
autoPtr< volScalarFielddJdpPtr_
 
autoPtr< volScalarFielddJdTPtr_
 
autoPtr< volScalarFielddJdTMvar1Ptr_
 First turbulence model variable. More...
 
autoPtr< volScalarFielddJdTMvar2Ptr_
 Second turbulence model variable. More...
 
autoPtr< boundaryVectorFieldbdJdvPtr_
 
autoPtr< boundaryScalarFieldbdJdvnPtr_
 Adjoint outlet pressure. More...
 
autoPtr< boundaryVectorFieldbdJdvtPtr_
 Adjoint outlet velocity. More...
 
autoPtr< boundaryVectorFieldbdJdpPtr_
 Adjoint (intlet,wall) velocity. More...
 
autoPtr< boundaryScalarFieldbdJdTPtr_
 Adjoint outlet temperature. More...
 
autoPtr< boundaryScalarFieldbdJdTMvar1Ptr_
 Adjoint outlet turbulence model var 1. More...
 
autoPtr< boundaryScalarFieldbdJdTMvar2Ptr_
 Adjoint outlet turbulence model var 2. More...
 
- Protected Attributes inherited from objective
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_
 

Additional Inherited Members

- Protected Member Functions inherited from objective
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...
 

Detailed Description

Abstract base class for objective functions in incompressible flows.

Source files

Definition at line 54 of file objectiveIncompressible.H.

Constructor & Destructor Documentation

◆ objectiveIncompressible()

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

Construct from components.

Definition at line 54 of file objectiveIncompressible.C.

References dict, and dictionary::get().

Here is the call graph for this function:

◆ ~objectiveIncompressible()

virtual ~objectiveIncompressible ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName()

TypeName ( "incompressible"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
objectiveIncompressible  ,
dictionary  ,
(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)  ,
(mesh, dict, adjointSolverName, primalSolverName)   
)

◆ New()

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

Return a reference to the selected turbulence model.

Definition at line 97 of file objectiveIncompressible.C.

References dict, dictionary::dictName(), Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, dictionary::get(), Foam::Info, and mesh.

Here is the call graph for this function:

◆ J()

virtual scalar J ( )
pure virtual

Return the objective function value.

Implements objective.

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

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ dJdv()

const volVectorField & dJdv ( )

Contribution to field adjoint momentum eqs.

Definition at line 131 of file objectiveIncompressible.C.

References objectiveIncompressible::dJdvPtr_, objective::mesh_, and Foam::type().

Here is the call graph for this function:

◆ dJdp()

const volScalarField & dJdp ( )

Contribution to field adjoint continuity eq.

Definition at line 150 of file objectiveIncompressible.C.

References objectiveIncompressible::dJdpPtr_, objective::mesh_, and Foam::type().

Here is the call graph for this function:

◆ dJdT()

const volScalarField & dJdT ( )

Contribution to field adjoint energy eq.

Definition at line 169 of file objectiveIncompressible.C.

References objectiveIncompressible::dJdTPtr_, objective::mesh_, and Foam::type().

Here is the call graph for this function:

◆ dJdTMvar1()

const volScalarField & dJdTMvar1 ( )

Contribution to field adjoint turbulence model variable 1.

Definition at line 188 of file objectiveIncompressible.C.

References objectiveIncompressible::dJdTMvar1Ptr_, objective::mesh_, and Foam::type().

Here is the call graph for this function:

◆ dJdTMvar2()

const volScalarField & dJdTMvar2 ( )

Contribution to field adjoint turbulence model variable 2.

Definition at line 207 of file objectiveIncompressible.C.

References objectiveIncompressible::dJdTMvar2Ptr_, objective::mesh_, and Foam::type().

Here is the call graph for this function:

◆ boundarydJdv() [1/2]

const fvPatchVectorField & boundarydJdv ( const label  patchI)

Objective partial deriv wrt velocity for a specific patch.

Definition at line 227 of file objectiveIncompressible.C.

◆ boundarydJdvn() [1/2]

const fvPatchScalarField & boundarydJdvn ( const label  patchI)

Objective partial deriv wrt normal velocity for a specific patch.

Definition at line 240 of file objectiveIncompressible.C.

◆ boundarydJdvt() [1/2]

const fvPatchVectorField & boundarydJdvt ( const label  patchI)

Objective partial deriv wrt tangent velocity for a specific patch.

Definition at line 253 of file objectiveIncompressible.C.

◆ boundarydJdp() [1/2]

const fvPatchVectorField & boundarydJdp ( const label  patchI)

Objective partial deriv wrt pressure (times normal) for a specific patch

Definition at line 266 of file objectiveIncompressible.C.

◆ boundarydJdT() [1/2]

const fvPatchScalarField & boundarydJdT ( const label  patchI)

Objective partial deriv wrt temperature for a specific patch.

Definition at line 279 of file objectiveIncompressible.C.

◆ boundarydJdTMvar1() [1/2]

const fvPatchScalarField & boundarydJdTMvar1 ( const label  patchI)

Objective partial deriv wrt turbulence model var 1 for a specific patch

Definition at line 292 of file objectiveIncompressible.C.

◆ boundarydJdTMvar2() [1/2]

const fvPatchScalarField & boundarydJdTMvar2 ( const label  patchI)

Objective partial deriv wrt turbulence model var 2 for a specific patch

Definition at line 305 of file objectiveIncompressible.C.

◆ boundarydJdv() [2/2]

const boundaryVectorField & boundarydJdv ( )

Objective partial deriv wrt velocity for all patches.

Definition at line 317 of file objectiveIncompressible.C.

References objectiveIncompressible::bdJdvPtr_, and objective::mesh_.

Referenced by boundaryAdjointContributionIncompressible::velocitySource().

Here is the caller graph for this function:

◆ boundarydJdvn() [2/2]

const boundaryScalarField & boundarydJdvn ( )

Objective partial deriv wrt normal velocity for all patches.

Definition at line 327 of file objectiveIncompressible.C.

References objectiveIncompressible::bdJdvnPtr_, and objective::mesh_.

Referenced by boundaryAdjointContributionIncompressible::pressureSource().

Here is the caller graph for this function:

◆ boundarydJdvt() [2/2]

const boundaryVectorField & boundarydJdvt ( )

Objective partial deriv wrt tangent velocity for all patches.

Definition at line 337 of file objectiveIncompressible.C.

References objectiveIncompressible::bdJdvtPtr_, and objective::mesh_.

Referenced by boundaryAdjointContributionIncompressible::tangentVelocitySource().

Here is the caller graph for this function:

◆ boundarydJdp() [2/2]

const boundaryVectorField & boundarydJdp ( )

Objective partial deriv wrt pressure (times normal) for all patches.

Definition at line 347 of file objectiveIncompressible.C.

References objectiveIncompressible::bdJdpPtr_, and objective::mesh_.

Referenced by boundaryAdjointContributionIncompressible::normalVelocitySource().

Here is the caller graph for this function:

◆ boundarydJdT() [2/2]

const boundaryScalarField & boundarydJdT ( )

Objective partial deriv wrt temperature for all patches.

Definition at line 357 of file objectiveIncompressible.C.

References objectiveIncompressible::bdJdTPtr_, and objective::mesh_.

Referenced by boundaryAdjointContributionIncompressible::energySource().

Here is the caller graph for this function:

◆ boundarydJdTMvar1() [2/2]

const boundaryScalarField & boundarydJdTMvar1 ( )

Objective partial deriv wrt turbulence model var 1 for all patches.

Definition at line 367 of file objectiveIncompressible.C.

References objectiveIncompressible::bdJdTMvar1Ptr_, and objective::mesh_.

Referenced by boundaryAdjointContributionIncompressible::adjointTMVariable1Source().

Here is the caller graph for this function:

◆ boundarydJdTMvar2() [2/2]

const boundaryScalarField & boundarydJdTMvar2 ( )

Objective partial deriv wrt turbulence model var 2 for all patches.

Definition at line 377 of file objectiveIncompressible.C.

References objectiveIncompressible::bdJdTMvar2Ptr_, and objective::mesh_.

Referenced by boundaryAdjointContributionIncompressible::adjointTMVariable2Source().

Here is the caller graph for this function:

◆ update()

◆ nullify()

◆ update_dJdv()

virtual void update_dJdv ( )
inlinevirtual

Update vol and boundary fields and derivatives.

Do nothing in the base. The relevant ones should be overwritten in the child objective functions

Definition at line 236 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_dJdp()

virtual void update_dJdp ( )
inlinevirtual

Definition at line 239 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_dJdT()

virtual void update_dJdT ( )
inlinevirtual

Definition at line 242 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_dJdTMvar1()

virtual void update_dJdTMvar1 ( )
inlinevirtual

Definition at line 245 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_dJdTMvar2()

virtual void update_dJdTMvar2 ( )
inlinevirtual

Definition at line 248 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_dJdb()

virtual void update_dJdb ( )
inlinevirtual

Definition at line 251 of file objectiveIncompressible.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 from objective.

Definition at line 254 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_gradDxDbMultiplier()

virtual void update_gradDxDbMultiplier ( )
inlinevirtual

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

Reimplemented from objective.

Definition at line 257 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_boundarydJdv()

virtual void update_boundarydJdv ( )
inlinevirtual

Reimplemented in objectivePtLosses.

Definition at line 260 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_boundarydJdvn()

virtual void update_boundarydJdvn ( )
inlinevirtual

Reimplemented in objectivePtLosses.

Definition at line 263 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_boundarydJdvt()

virtual void update_boundarydJdvt ( )
inlinevirtual

Reimplemented in objectivePtLosses.

Definition at line 266 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_boundarydJdp()

virtual void update_boundarydJdp ( )
inlinevirtual

Reimplemented in objectiveMoment, objectiveForce, and objectivePtLosses.

Definition at line 269 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_boundarydJdT()

virtual void update_boundarydJdT ( )
inlinevirtual

Definition at line 272 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_boundarydJdTMvar1()

virtual void update_boundarydJdTMvar1 ( )
inlinevirtual

Definition at line 275 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_boundarydJdTMvar2()

virtual void update_boundarydJdTMvar2 ( )
inlinevirtual

Definition at line 278 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_boundarydJdb()

virtual void update_boundarydJdb ( )
inlinevirtual

Update objective function derivative term.

Reimplemented from objective.

Definition at line 281 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_dSdbMultiplier()

virtual void update_dSdbMultiplier ( )
inlinevirtual

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

Reimplemented from objective.

Reimplemented in objectiveMoment, objectiveForce, and objectivePartialVolume.

Definition at line 284 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_dndbMultiplier()

virtual void update_dndbMultiplier ( )
inlinevirtual

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

Reimplemented from objective.

Definition at line 287 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_dxdbMultiplier()

virtual void update_dxdbMultiplier ( )
inlinevirtual

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

Reimplemented from objective.

Reimplemented in objectiveMoment, and objectiveForce.

Definition at line 290 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_dxdbDirectMultiplier()

virtual void update_dxdbDirectMultiplier ( )
inlinevirtual

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

Reimplemented from objective.

Reimplemented in objectiveMoment, and objectivePartialVolume.

Definition at line 293 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ update_meanValues()

virtual void update_meanValues ( )
inlinevirtual

Some objectives need to store some auxiliary values. If averaging is enabled, update these mean values here.

By convention, the mean values (eg mean drag) refer to these flow values computed using the mean fields, rather than averaging the values themselves

Reimplemented in objectiveMoment.

Definition at line 301 of file objectiveIncompressible.H.

Referenced by objectiveIncompressible::update().

Here is the caller graph for this function:

◆ write()

void write ( ) const
virtual

Write objective function history.

Reimplemented from objective.

Reimplemented in objectivePtLosses, objectivePartialVolume, and objectiveForceTarget.

Definition at line 485 of file objectiveIncompressible.C.

References objective::write().

Here is the call graph for this function:

◆ hasdJdv()

bool hasdJdv ( ) const
inline

Inline functions for checking whether pointers are set or not.

Definition at line 33 of file objectiveIncompressibleI.H.

References objectiveIncompressible::dJdvPtr_.

Referenced by objectiveIncompressible::nullify().

Here is the caller graph for this function:

◆ hasdJdp()

bool hasdJdp ( ) const
inline

Definition at line 39 of file objectiveIncompressibleI.H.

Referenced by objectiveIncompressible::nullify().

Here is the caller graph for this function:

◆ hasdJdT()

bool hasdJdT ( ) const
inline

Definition at line 45 of file objectiveIncompressibleI.H.

Referenced by objectiveIncompressible::nullify().

Here is the caller graph for this function:

◆ hasdJdTMVar1()

bool hasdJdTMVar1 ( ) const
inline

Definition at line 51 of file objectiveIncompressibleI.H.

Referenced by objectiveIncompressible::nullify().

Here is the caller graph for this function:

◆ hasdJdTMVar2()

bool hasdJdTMVar2 ( ) const
inline

Definition at line 57 of file objectiveIncompressibleI.H.

Referenced by objectiveIncompressible::nullify().

Here is the caller graph for this function:

◆ hasBoundarydJdv()

bool hasBoundarydJdv ( ) const
inline

Definition at line 63 of file objectiveIncompressibleI.H.

Referenced by objectiveIncompressible::nullify().

Here is the caller graph for this function:

◆ hasBoundarydJdvn()

bool hasBoundarydJdvn ( ) const
inline

Definition at line 69 of file objectiveIncompressibleI.H.

Referenced by objectiveIncompressible::nullify().

Here is the caller graph for this function:

◆ hasBoundarydJdvt()

bool hasBoundarydJdvt ( ) const
inline

Definition at line 75 of file objectiveIncompressibleI.H.

Referenced by objectiveIncompressible::nullify().

Here is the caller graph for this function:

◆ hasBoundarydJdp()

bool hasBoundarydJdp ( ) const
inline

Definition at line 81 of file objectiveIncompressibleI.H.

Referenced by objectiveIncompressible::nullify().

Here is the caller graph for this function:

◆ hasBoundarydJdT()

bool hasBoundarydJdT ( ) const
inline

Definition at line 87 of file objectiveIncompressibleI.H.

Referenced by objectiveIncompressible::nullify().

Here is the caller graph for this function:

◆ hasBoundarydJdTMVar1()

bool hasBoundarydJdTMVar1 ( ) const
inline

Definition at line 93 of file objectiveIncompressibleI.H.

Referenced by objectiveIncompressible::nullify().

Here is the caller graph for this function:

◆ hasBoundarydJdTMVar2()

bool hasBoundarydJdTMVar2 ( ) const
inline

Definition at line 99 of file objectiveIncompressibleI.H.

Referenced by objectiveIncompressible::nullify().

Here is the caller graph for this function:

Member Data Documentation

◆ vars_

◆ dJdvPtr_

◆ dJdpPtr_

◆ dJdTPtr_

◆ dJdTMvar1Ptr_

autoPtr<volScalarField> dJdTMvar1Ptr_
protected

First turbulence model variable.

Definition at line 72 of file objectiveIncompressible.H.

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

◆ dJdTMvar2Ptr_

autoPtr<volScalarField> dJdTMvar2Ptr_
protected

Second turbulence model variable.

Definition at line 75 of file objectiveIncompressible.H.

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

◆ bdJdvPtr_

◆ bdJdvnPtr_

◆ bdJdvtPtr_

◆ bdJdpPtr_

◆ bdJdTPtr_

autoPtr<boundaryScalarField> bdJdTPtr_
protected

Adjoint outlet temperature.

Definition at line 91 of file objectiveIncompressible.H.

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

◆ bdJdTMvar1Ptr_

autoPtr<boundaryScalarField> bdJdTMvar1Ptr_
protected

Adjoint outlet turbulence model var 1.

Definition at line 94 of file objectiveIncompressible.H.

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

◆ bdJdTMvar2Ptr_

autoPtr<boundaryScalarField> bdJdTMvar2Ptr_
protected

Adjoint outlet turbulence model var 2.

Definition at line 97 of file objectiveIncompressible.H.

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


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