openfoam/applications/utilities/postProcessing/optimisation/computeSensitivities/computeSensitivities.C
Vaggelis Papoutsis ecc1fb5efb CONTRIB: New adjoint optimisation and tools
A set of libraries and executables creating a workflow for performing
gradient-based optimisation loops. The main executable (adjointOptimisationFoam)
solves the flow (primal) equations, followed by the adjoint equations and,
eventually, the computation of sensitivity derivatives.

Current functionality supports the solution of the adjoint equations for
incompressible turbulent flows, including the adjoint to the Spalart-Allmaras
turbulence model and the adjoint to the nutUSpaldingWallFunction, [1], [2].

Sensitivity derivatives are computed with respect to the normal displacement of
boundary wall nodes/faces (the so-called sensitivity maps) following the
Enhanced Surface Integrals (E-SI) formulation, [3].

The software was developed by PCOpt/NTUA and FOSS GP, with contributions from

Dr. Evangelos Papoutsis-Kiachagias,
Konstantinos Gkaragounis,
Professor Kyriakos Giannakoglou,
Andy Heather

and contributions in earlier version from

Dr. Ioannis Kavvadias,
Dr. Alexandros Zymaris,
Dr. Dimitrios Papadimitriou

[1] A.S. Zymaris, D.I. Papadimitriou, K.C. Giannakoglou, and C. Othmer.
Continuous adjoint approach to the Spalart-Allmaras turbulence model for
incompressible flows. Computers & Fluids, 38(8):1528–1538, 2009.

[2] E.M. Papoutsis-Kiachagias and K.C. Giannakoglou. Continuous adjoint methods
for turbulent flows, applied to shape and topology optimization: Industrial
applications. 23(2):255–299, 2016.

[3] I.S. Kavvadias, E.M. Papoutsis-Kiachagias, and K.C. Giannakoglou. On the
proper treatment of grid sensitivities in continuous adjoint methods for shape
optimization. Journal of Computational Physics, 301:1–18, 2015.

Integration into the official OpenFOAM release by OpenCFD
2019-06-17 12:59:11 +01:00

77 lines
2.6 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2019 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
| Copyright (C) 2007-2019 PCOpt/NTUA
| Copyright (C) 2013-2019 FOSS GP
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Application
computeSensitivities
Description
Computes the sensitivities wrt what is defined in the optimisationDict
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "optimisationManager.H"
#include "primalSolver.H"
#include "adjointSolver.H"
#include "incompressibleVars.H"
#include "incompressibleAdjointVars.H"
#include "adjointBoundaryCondition.H"
#include "adjointSolverManager.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
#include "setRootCase.H"
#include "createTime.H"
#include "createMesh.H"
#include "createFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
forAll(adjointSolverManagers, amI)
{
PtrList<adjointSolver>& adjointSolvers =
adjointSolverManagers[amI].adjointSolvers();
forAll(adjointSolvers, asI)
{
adjointSolvers[asI].getObjectiveManager().updateAndWrite();
adjointSolvers[asI].computeObjectiveSensitivities();
}
}
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
Info<< "End" << endl;
return(0);
}
// ************************************************************************* //