Loading...
Searching...
No Matches
adjointEikonalSolver Class Reference

Solver of the adjoint to the eikonal PDE. More...

#include <adjointEikonalSolver.H>

Collaboration diagram for adjointEikonalSolver:

Public Member Functions

 TypeName ("adjointEikonalSolver")
 Runtime type information.
 adjointEikonalSolver (const fvMesh &mesh, const dictionary &dict, adjointSolver &adjointSolver, const labelHashSet &sensitivityPatchIDs)
 Construct from components.
virtual ~adjointEikonalSolver ()=default
virtual bool readDict (const dictionary &dict)
 Read dict if changed.
void accumulateIntegrand (const scalar dt)
 Accumulate source term.
void solve ()
 Calculate the adjoint distance field.
void reset ()
 Reset source term.
boundaryVectorFielddistanceSensitivities ()
 Return the sensitivity term depending on da.
tmp< volTensorFieldgetFISensitivityTerm () const
 Return the volume-based sensitivity term depending on da.
tmp< scalarFieldtopologySensitivities (const word &designVarsName) const
 Return sensitivity contribution to topology optimisation.
const volScalarFieldda ()
 Return the adjoint distance field.
const volScalarFieldd ()
 Return the distance field.
tmp< volVectorFieldgradEikonal ()
 Return the gradient of the eikonal equation.

Protected Member Functions

wordList patchTypes () const
 Return the boundary condition types for da.
tmp< surfaceScalarFieldcomputeYPhi ()
 Compute convecting velocity.
void read ()
 Read options each time a new solution is found.

Protected Attributes

const fvMeshmesh_
dictionary dict_
adjointSolveradjointSolver_
const labelHashSetsensitivityPatchIDs_
label nEikonalIters_
scalar tolerance_
scalar epsilon_
labelHashSet wallPatchIDs_
volScalarField da_
volScalarField source_
autoPtr< boundaryVectorFielddistanceSensPtr_
 Wall face sens w.r.t. (x,y.z).

Detailed Description

Solver of the adjoint to the eikonal PDE.

Reference:

    For the development of the adjoint eikonal PDE and its boundary
    conditions

        Papoutsis-Kiachagias, E. M., & Giannakoglou, K. C. (2014).
        Continuous Adjoint Methods for Turbulent Flows, Applied to Shape
        and Topology Optimization: Industrial Applications.
        Archives of Computational Methods in Engineering, 23(2), 255-299.
        http://doi.org/10.1007/s11831-014-9141-9

To be as consistent as possible, it is recommended to use the advectionDiffusion wallDist method in fvSchemes, instead of the more widely used meshWave

Example of the adjointEikonalSolver specification in optimisationDict:

    optimisation
    {
        sensitivities
        {
            includeDistance true;
            adjointEikonalSolver
            {
                // epsilon should be the same as the one used
                // in fvSchemes/wallDist/advectionDiffusionCoeffs
                epsilon   0.1;
                iters     1000;
                tolerance 1e-6;
            }
        }
    }

Example of the entries in fvSchemes:

    divSchemes
    {
        .
        .
        // avoid bounded schemes since yPhi is not conservative
        div(-yPhi,da) Gauss linearUpwind grad(da);
        .
        .
    }
    laplacianSchemes
    {
        .
        .
        laplacian(yWall,da) Gauss linear corrected;
        .
        .
    }

Also, the solver specification and a relaxation factor for da are required in fvSolution

    da
    {
        solver          PBiCGStab;
        preconditioner  DILU;
        tolerance       1e-9;
        relTol          0.1;
    }

    relaxationFactors
    {
        equations
        {
            .
            .
            da           0.5;
            .
            .
        }
    }
See also
Foam::patchDistMethod::advectionDiffusion Foam::wallDist
Source files

Definition at line 143 of file adjointEikonalSolver.H.

Constructor & Destructor Documentation

◆ adjointEikonalSolver()

adjointEikonalSolver ( const fvMesh & mesh,
const dictionary & dict,
adjointSolver & adjointSolver,
const labelHashSet & sensitivityPatchIDs )

Construct from components.

Definition at line 116 of file adjointEikonalSolver.C.

References adjointSolver_, Foam::createZeroBoundaryPtr(), da_, dict, dict_, Foam::dimLength, distanceSensPtr_, epsilon_, mesh, mesh_, nEikonalIters_, patchTypes(), read(), sensitivityPatchIDs_, source_, timeName, tolerance_, wallPatchIDs_, and Foam::Zero.

Here is the call graph for this function:

◆ ~adjointEikonalSolver()

virtual ~adjointEikonalSolver ( )
virtualdefault

References dict.

Member Function Documentation

◆ patchTypes()

wordList patchTypes ( ) const
protected

Return the boundary condition types for da.

Definition at line 47 of file adjointEikonalSolver.C.

References mesh_, wallPatchIDs_, and fvPatchFieldBase::zeroGradientType().

Referenced by adjointEikonalSolver().

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

◆ computeYPhi()

tmp< surfaceScalarField > computeYPhi ( )
protected

◆ read()

void read ( )
protected

Read options each time a new solution is found.

Definition at line 64 of file adjointEikonalSolver.C.

References dict_, e, epsilon_, mesh_, nEikonalIters_, and tolerance_.

Referenced by adjointEikonalSolver(), readDict(), and solve().

Here is the caller graph for this function:

◆ TypeName()

TypeName ( "adjointEikonalSolver" )

Runtime type information.

References dict, and mesh.

◆ readDict()

bool readDict ( const dictionary & dict)
virtual

Read dict if changed.

Definition at line 174 of file adjointEikonalSolver.C.

References dict, dict_, and read().

Here is the call graph for this function:

◆ accumulateIntegrand()

void accumulateIntegrand ( const scalar dt)

Accumulate source term.

Definition at line 183 of file adjointEikonalSolver.C.

References adjointSolver_, and source_.

◆ solve()

◆ reset()

void reset ( )

Reset source term.

Definition at line 260 of file adjointEikonalSolver.C.

References distanceSensPtr_, source_, and Foam::Zero.

◆ distanceSensitivities()

boundaryVectorField & distanceSensitivities ( )

Return the sensitivity term depending on da.

Definition at line 267 of file adjointEikonalSolver.C.

References adjointSolver_, d(), da_, distanceSensPtr_, Foam::endl(), Foam::Info, mesh_, sensitivityPatchIDs_, and td().

Here is the call graph for this function:

◆ getFISensitivityTerm()

◆ topologySensitivities()

tmp< scalarField > topologySensitivities ( const word & designVarsName) const

Return sensitivity contribution to topology optimisation.

Definition at line 350 of file adjointEikonalSolver.C.

References adjointSolver_, d(), da_, fvOptions, mesh_, options::New(), tmp< T >::New(), sensitivityTopO::postProcessSens(), td(), and Foam::Zero.

Here is the call graph for this function:

◆ da()

const volScalarField & da ( )

Return the adjoint distance field.

Definition at line 371 of file adjointEikonalSolver.C.

References da_.

◆ d()

const volScalarField & d ( )

Return the distance field.

Referenced by computeYPhi(), distanceSensitivities(), getFISensitivityTerm(), gradEikonal(), solve(), and topologySensitivities().

Here is the caller graph for this function:

◆ gradEikonal()

tmp< volVectorField > gradEikonal ( )

Return the gradient of the eikonal equation.

Definition at line 377 of file adjointEikonalSolver.C.

References adjointSolver_, d(), Foam::fvc::grad(), tmp< T >::New(), and td().

Here is the call graph for this function:

Member Data Documentation

◆ mesh_

◆ dict_

dictionary dict_
protected

Definition at line 166 of file adjointEikonalSolver.H.

Referenced by adjointEikonalSolver(), read(), and readDict().

◆ adjointSolver_

◆ sensitivityPatchIDs_

const labelHashSet& sensitivityPatchIDs_
protected

Definition at line 170 of file adjointEikonalSolver.H.

Referenced by adjointEikonalSolver(), and distanceSensitivities().

◆ nEikonalIters_

label nEikonalIters_
protected

Definition at line 172 of file adjointEikonalSolver.H.

Referenced by adjointEikonalSolver(), read(), and solve().

◆ tolerance_

scalar tolerance_
protected

Definition at line 174 of file adjointEikonalSolver.H.

Referenced by adjointEikonalSolver(), read(), and solve().

◆ epsilon_

scalar epsilon_
protected

Definition at line 176 of file adjointEikonalSolver.H.

Referenced by adjointEikonalSolver(), getFISensitivityTerm(), read(), and solve().

◆ wallPatchIDs_

labelHashSet wallPatchIDs_
protected

◆ da_

◆ source_

volScalarField source_
protected

Definition at line 182 of file adjointEikonalSolver.H.

Referenced by accumulateIntegrand(), adjointEikonalSolver(), reset(), and solve().

◆ distanceSensPtr_

autoPtr<boundaryVectorField> distanceSensPtr_
protected

Wall face sens w.r.t. (x,y.z).

Definition at line 187 of file adjointEikonalSolver.H.

Referenced by adjointEikonalSolver(), distanceSensitivities(), and reset().


The documentation for this class was generated from the following files:
  • src/optimisation/adjointOptimisation/adjoint/optimisation/adjointSensitivity/adjointSensitivity/adjointEikonalSolver/adjointEikonalSolver.H
  • src/optimisation/adjointOptimisation/adjoint/optimisation/adjointSensitivity/adjointSensitivity/adjointEikonalSolver/adjointEikonalSolver.C