adjointSensitivity Class Referenceabstract

Abstract base class for adjoint-based sensitivities. More...

Inheritance diagram for adjointSensitivity:
Collaboration diagram for adjointSensitivity:

Public Member Functions

 TypeName ("adjointSensitivity")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, adjointSensitivity, dictionary,(const fvMesh &mesh, const dictionary &dict, adjointSolver &adjointSolver),(mesh, dict, adjointSolver))
 
 adjointSensitivity (const fvMesh &mesh, const dictionary &dict, adjointSolver &adjointSolver)
 Construct from components. More...
 
virtual ~adjointSensitivity ()=default
 Destructor. More...
 
virtual bool readDict (const dictionary &dict)
 Read dictionary if changed. More...
 
const adjointSolvergetAdjointSolver () const
 Const access to adjoint solver. More...
 
adjointSolvergetAdjointSolver ()
 Non-const access to adjoint solver. More...
 
bool includeDistance () const
 Should the adjoint eikonal PDE should be solved. More...
 
autoPtr< adjointEikonalSolver > & getAdjointEikonalSolver ()
 Return the adjoint eikonal solver. More...
 
autoPtr< adjointMeshMovementSolver > & getAdjointMeshMovementSolver ()
 Return the adjoint eikonal solver. More...
 
void setSuffix (const word &suffix)
 Set suffix. More...
 
const wordgetSuffix () const
 Get suffix. More...
 
virtual bool computeDxDbInternalField () const
 Should the parameterization compute the internalField of dxdb. More...
 
virtual void accumulateIntegrand (const scalar dt)=0
 Accumulate sensitivity integrands. More...
 
virtual void assembleSensitivities (autoPtr< designVariables > &designVars)
 Assemble sensitivities. More...
 
virtual const scalarFieldcalculateSensitivities (autoPtr< designVariables > &designVars)
 Calculates and returns sensitivity fields. More...
 
const scalarFieldgetSensitivities () const
 Returns the sensitivity fields. More...
 
virtual void clearSensitivities ()
 Zero sensitivity fields and their constituents. More...
 
virtual void write (const word &baseName=word::null)
 Write sensitivity fields. More...
 
const autoPtr< volTensorField > & gradDxDbMult () const
 
autoPtr< volTensorField > & gradDxDbMult ()
 
const autoPtr< scalarField > & divDxDbMult () const
 
const autoPtr< boundaryVectorField > & dxdbMult () const
 
const autoPtr< boundaryVectorField > & dSfdbMult () const
 
const autoPtr< boundaryVectorField > & dnfdbMult () const
 
const autoPtr< boundaryVectorField > & dxdbDirectMult () const
 
const autoPtr< pointBoundaryVectorField > & pointDxDbDirectMult () const
 
const autoPtr< boundaryVectorField > & bcDxDbMult () const
 
const autoPtr< vectorField > & optionsDxDbMult () const
 
- Public Member Functions inherited from sensitivity
 TypeName ("sensitivity")
 Runtime type information. More...
 
 sensitivity (const fvMesh &mesh, const dictionary &dict)
 Construct from components. More...
 
virtual ~sensitivity ()=default
 Destructor. More...
 
const fvMeshmesh () const
 Return reference to mesh. More...
 
const dictionarydict () const
 Return the construction dictionary. More...
 
const autoPtr< volScalarField > & fieldSensPtr () const
 Get the fieldSensPtr. More...
 

Static Public Member Functions

static autoPtr< adjointSensitivityNew (const fvMesh &mesh, const dictionary &dict, adjointSolver &adjointSolver)
 Return a reference to the selected turbulence model. More...
 

Protected Attributes

adjointSolveradjointSolver_
 Reference to the underlaying adjoint solver. More...
 
scalarField derivatives_
 The sensitivity derivative values. More...
 
word suffix_
 Append this word to files related to the sensitivities. More...
 
bool includeDistance_
 Include distance variation in sensitivity computations. More...
 
autoPtr< adjointEikonalSolvereikonalSolver_
 Adjoint eikonal equation solver. More...
 
autoPtr< adjointMeshMovementSolveradjointMeshMovementSolver_
 Adjoint grid displacement solver. More...
 
autoPtr< volTensorFieldgradDxDbMult_
 Multiplier of grad(dx/b) More...
 
autoPtr< scalarFielddivDxDbMult_
 Multiplier of div(dx/db) More...
 
autoPtr< boundaryVectorFielddxdbMult_
 Multiplier of face dx/db. More...
 
autoPtr< boundaryVectorFielddSfdbMult_
 Multiplier of dSf/db. More...
 
autoPtr< boundaryVectorFielddnfdbMult_
 Multiplier of dnf/db. More...
 
autoPtr< boundaryVectorFielddxdbDirectMult_
 Multiplier of dCf/db, found in the objective function. More...
 
autoPtr< pointBoundaryVectorFieldpointDxDbDirectMult_
 Multiplier of dx/db computed at points, found in the objective function. More...
 
autoPtr< boundaryVectorFieldbcDxDbMult_
 Multiplier of dx/db, coming from boundary conditions that depend on the geometry, like rotatingWallVelocity. More...
 
autoPtr< vectorFieldoptionsDxDbMult_
 dx/db multiplier coming from fvOptions More...
 
- Protected Attributes inherited from sensitivity
const fvMeshmesh_
 
dictionary dict_
 
bool writeFieldSens_
 
autoPtr< volScalarFieldfieldSensPtr_
 

Detailed Description

Abstract base class for adjoint-based sensitivities.

Source files

Definition at line 56 of file adjointSensitivity.H.

Constructor & Destructor Documentation

◆ adjointSensitivity()

adjointSensitivity ( const fvMesh mesh,
const dictionary dict,
adjointSolver adjointSolver 
)

Construct from components.

Definition at line 46 of file adjointSensitivity.C.

◆ ~adjointSensitivity()

virtual ~adjointSensitivity ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName()

TypeName ( "adjointSensitivity"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
adjointSensitivity  ,
dictionary  ,
(const fvMesh &mesh, const dictionary &dict, adjointSolver &adjointSolver ,
(mesh, dict, adjointSolver  
)

◆ New()

autoPtr< adjointSensitivity > New ( const fvMesh mesh,
const dictionary dict,
adjointSolver adjointSolver 
)
static

Return a reference to the selected turbulence model.

Definition at line 80 of file adjointSensitivity.C.

References dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, Foam::Info, mesh, fvMesh::name(), and dictionary::optionalSubDict().

Referenced by adjointSolver::allocateSensitivities(), and sensitivityMultiple::sensitivityMultiple().

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

◆ readDict()

bool readDict ( const dictionary dict)
virtual

◆ getAdjointSolver() [1/2]

◆ getAdjointSolver() [2/2]

adjointSolver& getAdjointSolver ( )
inline

Non-const access to adjoint solver.

Definition at line 242 of file adjointSensitivity.H.

References adjointSensitivity::adjointSolver_.

◆ includeDistance()

bool includeDistance ( ) const
inline

Should the adjoint eikonal PDE should be solved.

Definition at line 250 of file adjointSensitivity.H.

References adjointSensitivity::includeDistance_.

Referenced by shapeDesignVariables::assembleSensitivities().

Here is the caller graph for this function:

◆ getAdjointEikonalSolver()

autoPtr<adjointEikonalSolver>& getAdjointEikonalSolver ( )
inline

Return the adjoint eikonal solver.

Definition at line 258 of file adjointSensitivity.H.

References adjointSensitivity::eikonalSolver_.

Referenced by shapeDesignVariables::assembleSensitivities(), and adjointMeshMovementSolver::setSource().

Here is the caller graph for this function:

◆ getAdjointMeshMovementSolver()

autoPtr<adjointMeshMovementSolver>& getAdjointMeshMovementSolver ( )
inline

Return the adjoint eikonal solver.

Definition at line 267 of file adjointSensitivity.H.

References adjointSensitivity::adjointMeshMovementSolver_.

◆ setSuffix()

void setSuffix ( const word suffix)
inline

Set suffix.

Definition at line 275 of file adjointSensitivity.H.

References adjointSensitivity::suffix_.

Referenced by sensitivitySurfacePoints::setSuffixName().

Here is the caller graph for this function:

◆ getSuffix()

const word& getSuffix ( ) const
inline

Get suffix.

Definition at line 283 of file adjointSensitivity.H.

References adjointSensitivity::suffix_.

Referenced by shapeDesignVariables::writeSensitivities().

Here is the caller graph for this function:

◆ computeDxDbInternalField()

bool computeDxDbInternalField ( ) const
virtual

Should the parameterization compute the internalField of dxdb.

Reimplemented in sensitivityShapeFI.

Definition at line 133 of file adjointSensitivity.C.

Referenced by shapeDesignVariables::assembleSensitivities().

Here is the caller graph for this function:

◆ accumulateIntegrand()

virtual void accumulateIntegrand ( const scalar  dt)
pure virtual

Accumulate sensitivity integrands.

Corresponds to the flow and adjoint part of the sensitivities

Implemented in ShapeSensitivitiesBase, sensitivityTopO, and sensitivityMultiple.

◆ assembleSensitivities()

void assembleSensitivities ( autoPtr< designVariables > &  designVars)
virtual

Assemble sensitivities.

Adds the geometric part of the sensitivities

Reimplemented in sensitivitySurfacePoints, sensitivitySurface, sensitivityTopO, sensitivityShapeESI, sensitivityMultiple, and sensitivityShapeFI.

Definition at line 140 of file adjointSensitivity.C.

References designVariables::assembleSensitivities().

Referenced by sensitivityShapeFI::assembleSensitivities(), sensitivityShapeESI::assembleSensitivities(), sensitivityTopO::assembleSensitivities(), and sensitivitySurfacePoints::assembleSensitivities().

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

◆ calculateSensitivities()

const scalarField & calculateSensitivities ( autoPtr< designVariables > &  designVars)
virtual

Calculates and returns sensitivity fields.

Used with optimisation libraries

Implements sensitivity.

Reimplemented in sensitivityMultiple.

Definition at line 149 of file adjointSensitivity.C.

References Foam::type(), and Foam::vtk::write().

Here is the call graph for this function:

◆ getSensitivities()

const scalarField & getSensitivities ( ) const

Returns the sensitivity fields.

Assumes it has already been updated/computed

Definition at line 159 of file adjointSensitivity.C.

References adjointSensitivity::derivatives_.

◆ clearSensitivities()

void clearSensitivities ( )
virtual

Zero sensitivity fields and their constituents.

Reimplemented in ShapeSensitivitiesBase, and sensitivityMultiple.

Definition at line 165 of file adjointSensitivity.C.

References adjointSensitivity::adjointMeshMovementSolver_, adjointSensitivity::derivatives_, adjointSensitivity::eikonalSolver_, sensitivity::fieldSensPtr_, and Foam::Zero.

Referenced by ShapeSensitivitiesBase::clearSensitivities().

Here is the caller graph for this function:

◆ write()

void write ( const word baseName = word::null)
virtual

Write sensitivity fields.

If valid, copies boundaryFields to volFields and writes them. Virtual to be reimplemented by control points-based methods (Bezier, RBF) which do not need to write fields

Reimplemented from sensitivity.

Reimplemented in ShapeSensitivitiesBase, sensitivitySurfacePoints, and sensitivityMultiple.

Definition at line 183 of file adjointSensitivity.C.

References sensitivity::write().

Referenced by sensitivitySurfacePoints::write(), and ShapeSensitivitiesBase::write().

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

◆ gradDxDbMult() [1/2]

const Foam::autoPtr< Foam::volTensorField > & gradDxDbMult ( ) const
inline

Definition at line 25 of file adjointSensitivityI.H.

References adjointSensitivity::gradDxDbMult_.

Referenced by shapeDesignVariables::assembleSensitivities(), and adjointMeshMovementSolver::setSource().

Here is the caller graph for this function:

◆ gradDxDbMult() [2/2]

Foam::autoPtr< Foam::volTensorField > & gradDxDbMult ( )
inline

Definition at line 32 of file adjointSensitivityI.H.

◆ divDxDbMult()

const Foam::autoPtr< Foam::scalarField > & divDxDbMult ( ) const
inline

Definition at line 38 of file adjointSensitivityI.H.

Referenced by shapeDesignVariables::assembleSensitivities().

Here is the caller graph for this function:

◆ dxdbMult()

const Foam::autoPtr< Foam::boundaryVectorField > & dxdbMult ( ) const
inline

Definition at line 45 of file adjointSensitivityI.H.

Referenced by shapeDesignVariables::assembleSensitivities().

Here is the caller graph for this function:

◆ dSfdbMult()

const Foam::autoPtr< Foam::boundaryVectorField > & dSfdbMult ( ) const
inline

Definition at line 52 of file adjointSensitivityI.H.

Referenced by shapeDesignVariables::assembleSensitivities().

Here is the caller graph for this function:

◆ dnfdbMult()

const Foam::autoPtr< Foam::boundaryVectorField > & dnfdbMult ( ) const
inline

Definition at line 58 of file adjointSensitivityI.H.

Referenced by shapeDesignVariables::assembleSensitivities().

Here is the caller graph for this function:

◆ dxdbDirectMult()

const Foam::autoPtr< Foam::boundaryVectorField > & dxdbDirectMult ( ) const
inline

Definition at line 64 of file adjointSensitivityI.H.

Referenced by shapeDesignVariables::assembleSensitivities().

Here is the caller graph for this function:

◆ pointDxDbDirectMult()

const Foam::autoPtr< Foam::pointBoundaryVectorField > & pointDxDbDirectMult ( ) const
inline

Definition at line 70 of file adjointSensitivityI.H.

◆ bcDxDbMult()

const Foam::autoPtr< Foam::boundaryVectorField > & bcDxDbMult ( ) const
inline

Definition at line 76 of file adjointSensitivityI.H.

Referenced by shapeDesignVariables::assembleSensitivities().

Here is the caller graph for this function:

◆ optionsDxDbMult()

const Foam::autoPtr< Foam::vectorField > & optionsDxDbMult ( ) const
inline

Definition at line 82 of file adjointSensitivityI.H.

Referenced by shapeDesignVariables::assembleSensitivities(), and adjointMeshMovementSolver::setSource().

Here is the caller graph for this function:

Member Data Documentation

◆ adjointSolver_

◆ derivatives_

scalarField derivatives_
protected

◆ suffix_

word suffix_
protected

Append this word to files related to the sensitivities.

Definition at line 78 of file adjointSensitivity.H.

Referenced by adjointSensitivity::getSuffix(), adjointSensitivity::setSuffix(), and sensitivitySurface::smoothSensitivities().

◆ includeDistance_

bool includeDistance_
protected

Include distance variation in sensitivity computations.

Definition at line 83 of file adjointSensitivity.H.

Referenced by ShapeSensitivitiesBase::allocateEikonalSolver(), adjointSensitivity::includeDistance(), sensitivityTopO::readDict(), and adjointSensitivity::readDict().

◆ eikonalSolver_

◆ adjointMeshMovementSolver_

◆ gradDxDbMult_

autoPtr<volTensorField> gradDxDbMult_
protected

Multiplier of grad(dx/b)

Definition at line 102 of file adjointSensitivity.H.

Referenced by sensitivityShapeESI::computeDxDbMult(), and adjointSensitivity::gradDxDbMult().

◆ divDxDbMult_

autoPtr<scalarField> divDxDbMult_
protected

Multiplier of div(dx/db)

Definition at line 107 of file adjointSensitivity.H.

◆ dxdbMult_

autoPtr<boundaryVectorField> dxdbMult_
protected

Multiplier of face dx/db.

The term that multiplies the adjoint-related part of the sensitivities in the (E)SI approach

Definition at line 115 of file adjointSensitivity.H.

Referenced by sensitivityShapeESI::computeDxDbMult(), and sensitivitySurfacePoints::finalisePointSensitivities().

◆ dSfdbMult_

autoPtr<boundaryVectorField> dSfdbMult_
protected

Multiplier of dSf/db.

Definition at line 120 of file adjointSensitivity.H.

Referenced by sensitivitySurfacePoints::finalisePointSensitivities().

◆ dnfdbMult_

autoPtr<boundaryVectorField> dnfdbMult_
protected

Multiplier of dnf/db.

Definition at line 125 of file adjointSensitivity.H.

Referenced by sensitivitySurfacePoints::finalisePointSensitivities().

◆ dxdbDirectMult_

autoPtr<boundaryVectorField> dxdbDirectMult_
protected

Multiplier of dCf/db, found in the objective function.

Definition at line 130 of file adjointSensitivity.H.

Referenced by sensitivityShapeESI::computeDxDbMult(), and sensitivitySurfacePoints::finalisePointSensitivities().

◆ pointDxDbDirectMult_

autoPtr<pointBoundaryVectorField> pointDxDbDirectMult_
protected

Multiplier of dx/db computed at points, found in the objective function.

Definition at line 136 of file adjointSensitivity.H.

◆ bcDxDbMult_

autoPtr<boundaryVectorField> bcDxDbMult_
protected

Multiplier of dx/db, coming from boundary conditions that depend on the geometry, like rotatingWallVelocity.

Definition at line 142 of file adjointSensitivity.H.

Referenced by sensitivitySurfacePoints::finalisePointSensitivities().

◆ optionsDxDbMult_

autoPtr<vectorField> optionsDxDbMult_
protected

dx/db multiplier coming from fvOptions

Definition at line 147 of file adjointSensitivity.H.


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