ENH: Added the nutSqr surrogate noise objective

which qualitatively quantifies noise through a volume integral of the squared
turbulent viscosity.
This commit is contained in:
Vaggelis Papoutsis
2020-04-01 13:40:26 +03:00
committed by Andrew Heather
parent ad76df43c3
commit dcc039ce80
3 changed files with 297 additions and 0 deletions

View File

@ -55,6 +55,7 @@ objectives/incompressible/objectiveMoment/objectiveMoment.C
objectives/incompressible/objectivePtLosses/objectivePtLosses.C
objectives/incompressible/objectiveForceTarget/objectiveForceTarget.C
objectives/incompressible/objectivePartialVolume/objectivePartialVolume.C
objectives/incompressible/objectiveNutSqr/objectiveNutSqr.C
/* OBJECTIVE MANAGER*/
objectiveManager/objectiveManager/objectiveManager.C

View File

@ -0,0 +1,169 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2007-2020 PCOpt/NTUA
Copyright (C) 2013-2020 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/>.
\*---------------------------------------------------------------------------*/
#include "objectiveNutSqr.H"
#include "createZeroField.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace objectives
{
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
defineTypeNameAndDebug(objectiveNutSqr, 1);
addToRunTimeSelectionTable
(
objectiveIncompressible,
objectiveNutSqr,
dictionary
);
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
objectiveNutSqr::objectiveNutSqr
(
const fvMesh& mesh,
const dictionary& dict,
const word& adjointSolverName,
const word& primalSolverName
)
:
objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
zones_
(
mesh_.cellZones().indices(this->dict().get<wordRes>("zones"))
)
{
//- Allocate source term for the adjoint turbulence model
dJdTMvar1Ptr_.reset
(
createZeroFieldPtr<scalar>
(
mesh_,
"ATMSource",
(dimless/dimTime/dimTime)
)
);
//- Allocate term to be added to volume-based sensitivity derivatives
divDxDbMultPtr_.reset
(
createZeroFieldPtr<scalar>
(
mesh_,
("divDxdbMult"+type()) ,
//variable dimensions!!
//Dummy zeroes. Only the internalField will be used
dimensionSet(0,0,0,0,0,0,0)
)
);
//- set file pointer
//objFunctionFilePtr_ = new OFstream(objFunctionFolder_/type());
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
scalar objectiveNutSqr::J()
{
J_ = Zero;
const autoPtr<incompressible::RASModelVariables>&
turbVars = vars_.RASModelVariables();
const volScalarField& nut = turbVars->nutRefInst();
//scalar zoneVol(Zero);
for (const label zI : zones_)
{
const cellZone& zoneI = mesh_.cellZones()[zI];
for (const label cellI : zoneI)
{
J_ += sqr(nut[cellI])*(mesh_.V()[cellI]);
//zoneVol += mesh_.V()[cellI];
}
}
reduce(J_, sumOp<scalar>());
//reduce(zoneVol, sumOp<scalar>());
return J_;
}
void objectiveNutSqr::update_dJdTMvar1()
{
const autoPtr<incompressible::RASModelVariables>&
turbVars = vars_.RASModelVariables();
const singlePhaseTransportModel& lamTransp = vars_.laminarTransport();
const volScalarField& nut = turbVars->nutRef();
tmp<volScalarField> tnutJacobian = turbVars->nutJacobianVar1(lamTransp);
const volScalarField& nutJacobian = tnutJacobian();
volScalarField& dJdTMvar1 = dJdTMvar1Ptr_();
for (const label zI : zones_)
{
const cellZone& zoneI = mesh_.cellZones()[zI];
for (const label cellI : zoneI)
{
dJdTMvar1[cellI] = 2.*nut[cellI]*nutJacobian[cellI];
}
}
}
void objectiveNutSqr::update_divDxDbMultiplier()
{
const autoPtr<incompressible::RASModelVariables>&
turbVars = vars_.RASModelVariables();
const volScalarField& nut = turbVars->nutRef();
volScalarField& divDxDbMult = divDxDbMultPtr_();
for (const label zI : zones_)
{
const cellZone& zoneI = mesh_.cellZones()[zI];
for (const label cellI : zoneI)
{
divDxDbMult[cellI] = sqr(nut[cellI]);
}
}
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace objectives
} // End namespace Foam
// ************************************************************************* //

View File

@ -0,0 +1,127 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2007-2020 PCOpt/NTUA
Copyright (C) 2013-2020 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/>.
Class
Foam::objectives::objectiveNutSqr
Description
Objective qualitatively quantifying noise through the integral of the
squared turbulent viscosity over specified cellZones. Requires the adjoint
to the turbulence model to be incorporated into the optimisation loop.
Reference:
\verbatim
Objective function initially presented in
Papoutsis-Kiachagias, E. M., Magoulas, N., Mueller, J., Othmer, C.,
& Giannakoglou, K. C. (2015).
Noise reduction in car aerodynamics using a surrogate objective
function and the continuous adjoint method with wall functions.
Computers & Fluids, 122(20), 223-232.
https://doi.org/10.1016/j.compfluid.2015.09.002
\endverbatim
SourceFiles
objectiveNutSqrNoise.C
\*---------------------------------------------------------------------------*/
#ifndef objectiveNutSqr_H
#define objectiveNutSqr_H
#include "objectiveIncompressible.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace objectives
{
/*---------------------------------------------------------------------------*\
Class objectiveNutSqr Declaration
\*---------------------------------------------------------------------------*/
class objectiveNutSqr
:
public objectiveIncompressible
{
// Private data
//- Where to define the objective
labelList zones_;
public:
//- Runtime type information
TypeName("nutSqr");
// Constructors
//- from components
objectiveNutSqr
(
const fvMesh& mesh,
const dictionary& dict,
const word& adjointSolverName,
const word& primalSolverName
);
//- Destructor
virtual ~objectiveNutSqr() = default;
// Member Functions
//- Return the objective function value
scalar J();
//- Update field to be added to the adjoint turbulence model PDE
// The equivalent for the second turbulence model variable should
// be implemented for turbulence models other than SA
void update_dJdTMvar1();
//- Update field to be added to be added to volume-based
// sensitivity derivatives, emerging from delta ( dV ) / delta b
void update_divDxDbMultiplier();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace objectiveNutSqrNoise
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //