Loading...
Searching...
No Matches
objectiveNutSqr.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "objectiveNutSqr.H"
32#include "createZeroField.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
40namespace objectives
41{
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
47(
51);
52
53
54// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
55
56void objectiveNutSqr::populateFieldNames()
57{
58 if (adjointTurbulenceNames_.empty())
59 {
60 const incompressibleAdjointSolver& adjSolver =
64 const wordList& baseNames =
65 adjointRAS().getAdjointTMVariablesBaseNames();
66 forAll(baseNames, nI)
67 {
69 (adjSolver.extendedVariableName(baseNames[nI]));
70 adjointTurbulenceNames_.
71 push_back(adjSolver.extendedVariableName(baseNames[nI]));
72 }
73 }
74}
75
76
77// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78
80(
81 const fvMesh& mesh,
82 const dictionary& dict,
83 const word& adjointSolverName,
84 const word& primalSolverName
85)
86:
87 objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
88 zones_(mesh_.cellZones().indices(dict.get<wordRes>("zones"))),
89 adjointTurbulenceNames_()
90{
91 // Check if cellZones provided include at least one cell
92 checkCellZonesSize(zones_);
93 // Allocate dJdTMvar1Ptr_ and dJdTMvar2Ptr_ if needed
95 // Allocate term to be added to volume-based sensitivity derivatives
96 divDxDbMultPtr_.reset
97 (
99 (
100 IOobject
101 (
102 "divDxDbMult" + objectiveName_,
103 mesh_.time().timeName(),
104 mesh_,
107 ),
108 mesh_,
112 );
113}
114
115
116// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117
118scalar objectiveNutSqr::J()
119{
120 J_ = Zero;
121
123 turbVars = vars_.RASModelVariables();
124 const volScalarField& nut = turbVars->nutRefInst();
125
126 //scalar zoneVol(Zero);
127 for (const label zI : zones_)
128 {
129 const cellZone& zoneI = mesh_.cellZones()[zI];
130 for (const label cellI : zoneI)
131 {
132 J_ += sqr(nut[cellI])*(mesh_.V()[cellI]);
133 //zoneVol += mesh_.V()[cellI];
134 }
135 }
137 //reduce(zoneVol, sumOp<scalar>());
138
139 return J_;
140}
141
142
144{
145 // Add source from possible dependencies of nut on U
147 {
148 const incompressibleAdjointSolver& adjSolver =
152 adjSolver.getAdjointVars().adjointTurbulence();
155 tmp<volScalarField> dnutdUMult = 2*turbVars->nutRef();
156 tmp<volVectorField> dnutdU = adjointRAS->nutJacobianU(dnutdUMult);
157 if (dnutdU)
158 {
159 // If nut depends on U, allocate dJdv and add Ua to the fieldNames.
160 // It should be safe to do this here since objectives are updated
161 // before the first adjoint solution
162 if (!dJdvPtr_)
163 {
164 dJdvPtr_.reset
165 (
167 (
168 mesh_,
169 "dJdv" + type(),
171 )
172 );
173 }
174 if (!fieldNames_.size())
175 {
177 }
178 for (const label zI : zones_)
179 {
180 const cellZone& zoneI = mesh_.cellZones()[zI];
181 for (const label cellI : zoneI)
182 {
183 dJdvPtr_()[cellI] = dnutdU()[cellI];
185 }
186 }
187 }
188}
189
190
192{
193 const autoPtr<incompressible::RASModelVariables>&
194 turbVars = vars_.RASModelVariables();
195 const volScalarField& nut = turbVars->nutRef();
196 volScalarField JacobianMultiplier(2*nut);
197
199 (
202 JacobianMultiplier,
203 zones_
204 );
205}
206
207
209{
211 turbVars = vars_.RASModelVariables();
212 const volScalarField& nut = turbVars->nutRef();
213 volScalarField JacobianMultiplier(2*nut);
214
216 (
219 JacobianMultiplier,
220 zones_
221 );
222}
223
224
226{
228 turbVars = vars_.RASModelVariables();
229 const volScalarField& nut = turbVars->nutRef();
230
231 volScalarField& divDxDbMult = divDxDbMultPtr_();
232
233 for (const label zI : zones_)
234 {
235 const cellZone& zoneI = mesh_.cellZones()[zI];
236 for (const label cellI : zoneI)
237 {
238 divDxDbMult[cellI] = sqr(nut[cellI]);
239 }
240 }
241 divDxDbMult.correctBoundaryConditions();
242}
243
244
246{
247 populateFieldNames();
248 const label fieldI = fieldNames_.find(matrix.psi().name());
249
250 if (fieldI == 0)
251 {
252 matrix += weight()*dJdTMvar1();
253 }
254 if (fieldI == 1)
255 {
256 matrix += weight()*dJdTMvar2();
257 }
258}
259
260
261// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262
263} // End namespace objectives
264} // End namespace Foam
265
266// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void push_back(link *item)
Add at back of list.
Definition DLListBase.C:52
void correctBoundaryConditions()
Correct boundary field.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
void push_back(const T &val)
Append an element at the end of the list.
Definition ListI.H:221
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A subset of mesh cells.
Definition cellZone.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Base class for incompressibleAdjoint solvers.
virtual const incompressibleAdjointVars & getAdjointVars() const
Access to the incompressible adjoint variables set.
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
virtual tmp< volScalarField > nutJacobianTMVar2() const
Jacobian of nut wrt the second turbulence model variable.
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt the first turbulence model variable.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
bool foundObject(const word &name, const bool recursive=false) const
Contains the named Type?
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Abstract base class for objective functions in incompressible flows.
const volScalarField & dJdTMvar1()
Contribution to field adjoint turbulence model variable 1.
autoPtr< volScalarField > dJdTMvar2Ptr_
Second turbulence model variable.
void update_dJdTMvar(autoPtr< volScalarField > &dJdTMvarPtr, tmp< volScalarField >(incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const, const volScalarField &JacobianMultiplier, const labelList &zones)
Compute dJdTMVar{1,2}.
const volScalarField & dJdTMvar2()
Contribution to field adjoint turbulence model variable 2.
void allocatedJdTurbulence()
Allocate fields related to the differentiation of turbulence models, if necessary.
autoPtr< volScalarField > dJdTMvar1Ptr_
First turbulence model variable.
const incompressibleVars & vars_
void checkCellZonesSize(const labelList &zoneIDs) const
Check if cellZones provided include at least one cell.
autoPtr< volVectorField > dJdvPtr_
const fvMesh & mesh_
Definition objective.H:63
const word objectiveName_
Definition objective.H:67
wordList fieldNames_
List of adjoint fields for which this objective will contribute sources to their equations.
Definition objective.H:130
virtual scalar weight() const
Return the objective function weight.
Definition objective.C:353
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition objective.H:194
scalar J_
Objective function value and weight.
Definition objective.H:76
const dictionary & dict() const
Return objective dictionary.
Definition objective.C:91
const word adjointSolverName_
Definition objective.H:65
Objective qualitatively quantifying noise through the integral of the squared turbulent viscosity ove...
virtual void update_dJdTMvar2()
Update field to be added to the second adjoint turbulence model PDE.
objectiveNutSqr(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Construct from components.
virtual void addSource(fvScalarMatrix &matrix)
Add source terms to the adjoint turbulence model equations.
virtual scalar J()
Return the objective function value.
virtual void update_dJdTMvar1()
Update field to be added to the first adjoint turbulence model PDE.
virtual void update_divDxDbMultiplier()
Update field to be added to be added to volume-based sensitivity derivatives, emerging from delta ( d...
virtual void update_dJdv()
Update values to be added to the adjoint outlet velocity.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition polyMesh.H:679
word extendedVariableName(const word &varName) const
Given a variable name, return a name that is possibly appended by the solverName, depending on useSol...
Definition solverI.H:42
A class for managing temporary objects.
Definition tmp.H:75
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
scalar nut
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
autoPtr< GeometricField< Type, fvPatchField, volMesh > > createZeroFieldPtr(const fvMesh &mesh, const word &name, const dimensionSet dims, bool printAllocation=false)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299