Solver of the adjoint to the eikonal PDE. More...
#include <adjointEikonalSolver.H>

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. | |
| boundaryVectorField & | distanceSensitivities () |
| Return the sensitivity term depending on da. | |
| tmp< volTensorField > | getFISensitivityTerm () const |
| Return the volume-based sensitivity term depending on da. | |
| tmp< scalarField > | topologySensitivities (const word &designVarsName) const |
| Return sensitivity contribution to topology optimisation. | |
| const volScalarField & | da () |
| Return the adjoint distance field. | |
| const volScalarField & | d () |
| Return the distance field. | |
| tmp< volVectorField > | gradEikonal () |
| Return the gradient of the eikonal equation. | |
Protected Member Functions | |
| wordList | patchTypes () const |
| Return the boundary condition types for da. | |
| tmp< surfaceScalarField > | computeYPhi () |
| Compute convecting velocity. | |
| void | read () |
| Read options each time a new solution is found. | |
Protected Attributes | |
| const fvMesh & | mesh_ |
| dictionary | dict_ |
| adjointSolver & | adjointSolver_ |
| const labelHashSet & | sensitivityPatchIDs_ |
| label | nEikonalIters_ |
| scalar | tolerance_ |
| scalar | epsilon_ |
| labelHashSet | wallPatchIDs_ |
| volScalarField | da_ |
| volScalarField | source_ |
| autoPtr< boundaryVectorField > | distanceSensPtr_ |
| Wall face sens w.r.t. (x,y.z). | |
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;
.
.
}
}
Definition at line 143 of file adjointEikonalSolver.H.
| 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.

|
virtualdefault |
References dict.
|
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().


|
protected |
Compute convecting velocity.
Definition at line 76 of file adjointEikonalSolver.C.
References adjointSolver_, GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), d(), Foam::dimless, Foam::fvc::grad(), Foam::fvc::interpolate(), mesh_, tmp< T >::New(), IOobjectOption::NO_READ, IOobjectOption::NO_REGISTER, IOobjectOption::NO_WRITE, patches, patchDistMethod::patchTypes(), td(), wallPatchIDs_, and Foam::Zero.
Referenced by solve().


|
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().

|
virtual |
Read dict if changed.
Definition at line 174 of file adjointEikonalSolver.C.
References dict, dict_, and read().

| void accumulateIntegrand | ( | const scalar | dt | ) |
Accumulate source term.
Definition at line 183 of file adjointEikonalSolver.C.
References adjointSolver_, and source_.
| void solve | ( | ) |
Calculate the adjoint distance field.
Definition at line 190 of file adjointEikonalSolver.C.
References adjointSolver_, computeYPhi(), d(), da_, Foam::dimLength, Foam::dimTime, Foam::fvm::div(), Foam::endl(), epsilon_, fvOptions, Foam::gMax(), Foam::Info, SolverPerformance< Type >::initialResidual(), Foam::fvc::laplacian(), Foam::fvm::laplacian(), Foam::mag(), mesh_, nEikonalIters_, options::New(), IOobjectOption::NO_READ, IOobjectOption::NO_REGISTER, IOobjectOption::NO_WRITE, read(), fvMatrix< Type >::relax(), fvMatrix< Type >::solve(), source_, Foam::fvm::SuSp(), td(), and tolerance_.

| void reset | ( | ) |
Reset source term.
Definition at line 260 of file adjointEikonalSolver.C.
References distanceSensPtr_, source_, and Foam::Zero.
| 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().

| tmp< volTensorField > getFISensitivityTerm | ( | ) | const |
Return the volume-based sensitivity term depending on da.
Definition at line 290 of file adjointEikonalSolver.C.
References adjointSolver_, IOobjectOption::AUTO_WRITE, GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), d(), da_, Foam::dimLength, Foam::fvc::div(), Foam::endl(), epsilon_, Foam::fvc::grad(), Foam::Info, Foam::fvc::interpolate(), mesh_, tmp< T >::New(), IOobjectOption::NO_READ, patchDistMethod::patchTypes(), td(), wallPatchIDs_, Foam::Zero, and fvPatchFieldBase::zeroGradientType().

| 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.

| const volScalarField & da | ( | ) |
Return the adjoint distance field.
Definition at line 371 of file adjointEikonalSolver.C.
References da_.
| const volScalarField & d | ( | ) |
Return the distance field.
Referenced by computeYPhi(), distanceSensitivities(), getFISensitivityTerm(), gradEikonal(), solve(), and topologySensitivities().

| 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().

|
protected |
Definition at line 164 of file adjointEikonalSolver.H.
Referenced by adjointEikonalSolver(), computeYPhi(), distanceSensitivities(), getFISensitivityTerm(), patchTypes(), read(), solve(), and topologySensitivities().
|
protected |
Definition at line 166 of file adjointEikonalSolver.H.
Referenced by adjointEikonalSolver(), read(), and readDict().
|
protected |
Definition at line 168 of file adjointEikonalSolver.H.
Referenced by accumulateIntegrand(), adjointEikonalSolver(), computeYPhi(), distanceSensitivities(), getFISensitivityTerm(), gradEikonal(), solve(), and topologySensitivities().
|
protected |
Definition at line 170 of file adjointEikonalSolver.H.
Referenced by adjointEikonalSolver(), and distanceSensitivities().
|
protected |
Definition at line 172 of file adjointEikonalSolver.H.
Referenced by adjointEikonalSolver(), read(), and solve().
|
protected |
Definition at line 174 of file adjointEikonalSolver.H.
Referenced by adjointEikonalSolver(), read(), and solve().
|
protected |
Definition at line 176 of file adjointEikonalSolver.H.
Referenced by adjointEikonalSolver(), getFISensitivityTerm(), read(), and solve().
|
protected |
Definition at line 178 of file adjointEikonalSolver.H.
Referenced by adjointEikonalSolver(), computeYPhi(), getFISensitivityTerm(), and patchTypes().
|
protected |
Definition at line 180 of file adjointEikonalSolver.H.
Referenced by adjointEikonalSolver(), da(), distanceSensitivities(), getFISensitivityTerm(), solve(), and topologySensitivities().
|
protected |
Definition at line 182 of file adjointEikonalSolver.H.
Referenced by accumulateIntegrand(), adjointEikonalSolver(), reset(), and solve().
|
protected |
Wall face sens w.r.t. (x,y.z).
Definition at line 187 of file adjointEikonalSolver.H.
Referenced by adjointEikonalSolver(), distanceSensitivities(), and reset().