sensitivitySurface Class Reference

Calculation of adjoint based sensitivities at wall faces. More...

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

Public Member Functions

 TypeName ("surface")
 Runtime type information. More...
 
 sensitivitySurface (const fvMesh &mesh, const dictionary &dict, incompressibleAdjointSolver &adjointSolver)
 Construct from components. More...
 
virtual ~sensitivitySurface ()=default
 Destructor. More...
 
void read ()
 Read controls and update solver pointers if necessary. More...
 
virtual bool readDict (const dictionary &dict)
 Read dict if changed. More...
 
void computeDerivativesSize ()
 Compute the number of faces on sensitivityPatchIDs_. More...
 
virtual void accumulateIntegrand (const scalar dt)
 Accumulate sensitivity integrands. More...
 
virtual void assembleSensitivities ()
 Assemble sensitivities. More...
 
virtual void clearSensitivities ()
 Zero sensitivity fields and their constituents. More...
 
autoPtr< adjointEikonalSolver > & getAdjointEikonalSolver ()
 Get adjoint eikonal solver. More...
 
virtual void write (const word &baseName=word::null)
 Write sensitivity maps. More...
 
bool getIncludeObjective () const
 Get access to the includeObjective bool. More...
 
bool getIncludeSurfaceArea () const
 Get access to the includeSurfaceArea bool. More...
 
void setIncludeObjective (const bool includeObjective)
 Set includeObjective bool. More...
 
void setIncludeSurfaceArea (const bool includeSurfaceArea)
 Set includeSurfaceArea bool. More...
 
- Public Member Functions inherited from adjointSensitivity
 TypeName ("adjointSensitivity")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, adjointSensitivity, dictionary,(const fvMesh &mesh, const dictionary &dict, incompressibleAdjointSolver &adjointSolver),(mesh, dict, adjointSolver))
 
 adjointSensitivity (const fvMesh &mesh, const dictionary &dict, incompressibleAdjointSolver &adjointSolver)
 Construct from components. More...
 
virtual ~adjointSensitivity ()=default
 Destructor. More...
 
const incompressibleVarsprimalVars () const
 Get primal variables. More...
 
const incompressibleAdjointVarsadjointVars () const
 Get adjoint variables. More...
 
const incompressibleAdjointSolveradjointSolver () const
 Get adjoint solver. More...
 
virtual void accumulateIntegrand (const scalar dt)=0
 Accumulate sensitivity integrands. More...
 
virtual void assembleSensitivities ()=0
 Assemble sensitivities. More...
 
virtual const scalarFieldcalculateSensitivities ()
 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...
 
tmp< volTensorFieldcomputeGradDxDbMultiplier ()
 
tmp< volVectorFieldadjointMeshMovementSource ()
 Compute source term for adjoint mesh movement equation. More...
 
- 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 dictionarydict () const
 Return the construction dictionary. More...
 
virtual bool readDict (const dictionary &dict)
 Read dictionary if changed. More...
 
virtual void computeDerivativesSize ()
 Compute design variables number. Does nothing in the base. More...
 
virtual const scalarFieldcalculateSensitivities ()=0
 Calculates and returns sensitivity fields. More...
 
virtual void write (const word &baseName=word::null)
 Write sensitivity fields. More...
 

Protected Member Functions

void addGeometricSens ()
 
void setSuffixName ()
 Set suffix name for sensitivity fields. More...
 
void smoothSensitivities ()
 
scalar computeRadius (const faMesh &aMesh)
 

Protected Attributes

bool includeSurfaceArea_
 Include surface area in sens computation. More...
 
bool includePressureTerm_
 Include the adjoint pressure term in sens computation. More...
 
bool includeGradStressTerm_
 Include the term containing the grad of the stress at the boundary. More...
 
bool includeTransposeStresses_
 Include the transpose part of the adjoint stresses. More...
 
bool useSnGradInTranposeStresses_
 Use snGrad in the transpose part of the adjoint stresses. More...
 
bool includeDivTerm_
 Include the term from the deviatoric part of the stresses. More...
 
bool includeDistance_
 Include distance variation in sens computation. More...
 
bool includeMeshMovement_
 Include mesh movement variation in sens computation. More...
 
bool includeObjective_
 Include terms directly emerging from the objective function. More...
 
bool writeGeometricInfo_
 Write geometric info for use by external programs. More...
 
bool smoothSensitivities_
 
autoPtr< adjointEikonalSolvereikonalSolver_
 
autoPtr< adjointMeshMovementSolvermeshMovementSolver_
 
autoPtr< volVectorFieldnfOnPatchPtr_
 
autoPtr< volVectorFieldSfOnPatchPtr_
 
autoPtr< volVectorFieldCfOnPatchPtr_
 
- Protected Attributes inherited from adjointSensitivity
scalarField derivatives_
 
incompressibleAdjointSolveradjointSolver_
 
const incompressibleVarsprimalVars_
 
incompressibleAdjointVarsadjointVars_
 
objectiveManagerobjectiveManager_
 
- Protected Attributes inherited from sensitivity
const fvMeshmesh_
 
dictionary dict_
 
autoPtr< volScalarFieldfieldSensPtr_
 

Additional Inherited Members

- Static Public Member Functions inherited from adjointSensitivity
static autoPtr< adjointSensitivityNew (const fvMesh &mesh, const dictionary &dict, incompressibleAdjointSolver &adjointSolver)
 Return a reference to the selected turbulence model. More...
 

Detailed Description

Calculation of adjoint based sensitivities at wall faces.

Reference:

    The computation of the sensitivity derivatives  follows the (E-)SI
    formulation of

        Kavvadias, I. S., Papoutsis-Kiachagias, E. M., 
        Giannakoglou, K. C. (2015).
        On the proper treatment of grid sensitivities in continuous adjoint
        methods for shape optimization.
        Journal of computational physics, 301, 1-18.
        https://doi.org/10.1016/j.jcp.2015.08.012

    whereas their smoothing based on the computation of the 'Sobolev
    gradient' is derived from

        Vassberg J. C., Jameson A. (2006).
        Aerodynamic Shape Optimization Part I: Theoretical Background.
        VKI Lecture Series, Introduction to Optimization and
        Multidisciplinary Design, Brussels, Belgium, 8 March, 2006.
Source files

Definition at line 81 of file sensitivitySurfaceIncompressible.H.

Constructor & Destructor Documentation

◆ sensitivitySurface()

sensitivitySurface ( const fvMesh mesh,
const dictionary dict,
incompressibleAdjointSolver adjointSolver 
)

◆ ~sensitivitySurface()

virtual ~sensitivitySurface ( )
virtualdefault

Destructor.

Member Function Documentation

◆ addGeometricSens()

void addGeometricSens ( )
protected

Add sensitivities from dSd/db and dnf/db computed at points and mapped to faces

Definition at line 62 of file sensitivitySurfaceIncompressible.C.

References fvMesh::boundary(), polyMesh::boundaryMesh(), polyMesh::faces(), forAll, objectiveManager::getObjectiveFunctions(), Tensor< scalar >::I, sensitivitySurface::includeObjective_, deltaBoundary::makeFaceCentresAndAreas_d(), sensitivity::mesh_, primitiveMesh::nPoints(), adjointSensitivity::objectiveManager_, p, polyMesh::points(), face::points(), PrimitivePatchInterpolation< Patch >::pointToFaceInterpolate(), UList< T >::size(), syncTools::syncPointList(), and Foam::Zero.

Referenced by sensitivitySurface::accumulateIntegrand().

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

◆ setSuffixName()

void setSuffixName ( )
protected

Set suffix name for sensitivity fields.

Definition at line 200 of file sensitivitySurfaceIncompressible.C.

References adjointSensitivity::adjointVars_, sensitivity::dict(), sensitivitySurface::includeMeshMovement_, shapeSensitivitiesBase::setSuffix(), and variablesSet::solverName().

Referenced by sensitivitySurface::sensitivitySurface(), and sensitivitySurface::write().

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

◆ smoothSensitivities()

◆ computeRadius()

scalar computeRadius ( const faMesh aMesh)
protected

Compute the physical smoothing radius based on the average boundary face 'length'

Definition at line 385 of file sensitivitySurfaceIncompressible.C.

References aMesh(), polyMesh::bounds(), sensitivity::dict(), DimensionedField< Type, GeoMesh >::field(), forAll, Foam::gAverage(), polyMesh::geometricD(), dictionary::getOrDefault(), sensitivity::mesh_, polyMesh::nGeometricD(), Foam::pow(), faMesh::S(), and boundBox::span().

Referenced by sensitivitySurface::smoothSensitivities().

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

◆ TypeName()

TypeName ( "surface"  )

Runtime type information.

◆ read()

◆ readDict()

bool readDict ( const dictionary dict)
virtual

Read dict if changed.

Reimplemented from sensitivity.

Definition at line 562 of file sensitivitySurfaceIncompressible.C.

References sensitivity::dict(), sensitivitySurface::eikonalSolver_, and sensitivitySurface::meshMovementSolver_.

Referenced by SIBase::readDict().

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

◆ computeDerivativesSize()

void computeDerivativesSize ( )
virtual

Compute the number of faces on sensitivityPatchIDs_.

Reimplemented from sensitivity.

Definition at line 583 of file sensitivitySurfaceIncompressible.C.

References fvMesh::boundary(), adjointSensitivity::derivatives_, sensitivity::mesh_, nFaces(), List< T >::setSize(), and UPtrList< T >::size().

Referenced by sensitivitySurface::sensitivitySurface().

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

◆ accumulateIntegrand()

void accumulateIntegrand ( const scalar  dt)
virtual

Accumulate sensitivity integrands.

Implements adjointSensitivity.

Definition at line 594 of file sensitivitySurfaceIncompressible.C.

References sensitivitySurface::addGeometricSens(), incompressibleAdjointSolver::additionalSensitivityMapTerms(), adjointSensitivity::adjointSolver_, incompressibleAdjointVars::adjointTurbulence(), adjointSensitivity::adjointVars_, fvMesh::boundary(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), tmp< T >::clear(), Field< Type >::component(), tmp< T >::cref(), DebugInfo, Foam::dimVelocity, sensitivitySurface::eikonalSolver_, Foam::endl(), forAll, objectiveManager::getObjectiveFunctions(), Foam::fvc::grad(), sensitivitySurface::includeDistance_, sensitivitySurface::includeDivTerm_, sensitivitySurface::includeGradStressTerm_, sensitivitySurface::includeMeshMovement_, sensitivitySurface::includeObjective_, sensitivitySurface::includePressureTerm_, sensitivitySurface::includeTransposeStresses_, sensitivity::mesh_, sensitivitySurface::meshMovementSolver_, Time::New(), adjointSensitivity::objectiveManager_, p, incompressibleVars::p(), incompressibleAdjointMeanFlowVars::pa(), adjointSensitivity::primalVars_, autoPtr< T >::ref(), tmp< T >::ref(), UList< T >::size(), Foam::sqr(), Foam::twoSymm(), U, incompressibleVars::U(), incompressibleAdjointMeanFlowVars::Ua(), Foam::unzipRow(), sensitivitySurface::useSnGradInTranposeStresses_, and Foam::Zero.

Referenced by SIBase::accumulateIntegrand().

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

◆ assembleSensitivities()

void assembleSensitivities ( )
virtual

◆ clearSensitivities()

void clearSensitivities ( )
virtual

Zero sensitivity fields and their constituents.

Reimplemented from adjointSensitivity.

Definition at line 903 of file sensitivitySurfaceIncompressible.C.

References adjointSensitivity::clearSensitivities(), shapeSensitivitiesBase::clearSensitivities(), sensitivitySurface::eikonalSolver_, sensitivitySurface::includeDistance_, sensitivitySurface::includeMeshMovement_, and sensitivitySurface::meshMovementSolver_.

Referenced by SIBase::clearSensitivities().

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

◆ getAdjointEikonalSolver()

autoPtr< adjointEikonalSolver > & getAdjointEikonalSolver ( )

Get adjoint eikonal solver.

Definition at line 920 of file sensitivitySurfaceIncompressible.C.

References sensitivitySurface::eikonalSolver_.

◆ write()

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

Write sensitivity maps.

Reimplemented from adjointSensitivity.

Definition at line 926 of file sensitivitySurfaceIncompressible.C.

References sensitivitySurface::CfOnPatchPtr_, sensitivitySurface::nfOnPatchPtr_, sensitivitySurface::setSuffixName(), sensitivitySurface::SfOnPatchPtr_, ObukhovLength::write(), shapeSensitivitiesBase::write(), and sensitivitySurface::writeGeometricInfo_.

Referenced by SIBase::write().

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

◆ getIncludeObjective()

bool getIncludeObjective ( ) const
inline

Get access to the includeObjective bool.

Definition at line 38 of file sensitivitySurfaceIncompressibleI.H.

References sensitivitySurface::includeObjective_.

Referenced by SIBase::read().

Here is the caller graph for this function:

◆ getIncludeSurfaceArea()

bool getIncludeSurfaceArea ( ) const
inline

Get access to the includeSurfaceArea bool.

Definition at line 44 of file sensitivitySurfaceIncompressibleI.H.

References sensitivitySurface::includeSurfaceArea_.

◆ setIncludeObjective()

void setIncludeObjective ( const bool  includeObjective)
inline

Set includeObjective bool.

Definition at line 50 of file sensitivitySurfaceIncompressibleI.H.

References sensitivitySurface::includeObjective_.

Referenced by SIBase::read().

Here is the caller graph for this function:

◆ setIncludeSurfaceArea()

void setIncludeSurfaceArea ( const bool  includeSurfaceArea)
inline

Set includeSurfaceArea bool.

Definition at line 59 of file sensitivitySurfaceIncompressibleI.H.

References sensitivitySurface::includeSurfaceArea_.

Referenced by SIBase::read().

Here is the caller graph for this function:

Member Data Documentation

◆ includeSurfaceArea_

bool includeSurfaceArea_
protected

◆ includePressureTerm_

bool includePressureTerm_
protected

Include the adjoint pressure term in sens computation.

Definition at line 95 of file sensitivitySurfaceIncompressible.H.

Referenced by sensitivitySurface::accumulateIntegrand(), and sensitivitySurface::read().

◆ includeGradStressTerm_

bool includeGradStressTerm_
protected

Include the term containing the grad of the stress at the boundary.

Definition at line 98 of file sensitivitySurfaceIncompressible.H.

Referenced by sensitivitySurface::accumulateIntegrand(), and sensitivitySurface::read().

◆ includeTransposeStresses_

bool includeTransposeStresses_
protected

Include the transpose part of the adjoint stresses.

Definition at line 101 of file sensitivitySurfaceIncompressible.H.

Referenced by sensitivitySurface::accumulateIntegrand(), and sensitivitySurface::read().

◆ useSnGradInTranposeStresses_

bool useSnGradInTranposeStresses_
protected

Use snGrad in the transpose part of the adjoint stresses.

Definition at line 104 of file sensitivitySurfaceIncompressible.H.

Referenced by sensitivitySurface::accumulateIntegrand(), and sensitivitySurface::read().

◆ includeDivTerm_

bool includeDivTerm_
protected

Include the term from the deviatoric part of the stresses.

Definition at line 107 of file sensitivitySurfaceIncompressible.H.

Referenced by sensitivitySurface::accumulateIntegrand(), and sensitivitySurface::read().

◆ includeDistance_

bool includeDistance_
protected

◆ includeMeshMovement_

◆ includeObjective_

◆ writeGeometricInfo_

bool writeGeometricInfo_
protected

◆ smoothSensitivities_

bool smoothSensitivities_
protected

Smooth sensitivity derivatives based on the computation of the 'Sobolev gradient'

Definition at line 123 of file sensitivitySurfaceIncompressible.H.

Referenced by sensitivitySurface::assembleSensitivities(), and sensitivitySurface::read().

◆ eikonalSolver_

◆ meshMovementSolver_

◆ nfOnPatchPtr_

◆ SfOnPatchPtr_

◆ CfOnPatchPtr_


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