Loading...
Searching...
No Matches
objectivePowerDissipation.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
31#include "createZeroField.H"
32#include "topOVariablesBase.H"
33#include "IOmanip.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
41namespace objectives
42{
43
44// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45
48(
52);
53
54
55// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
56
57void objectivePowerDissipation::populateFieldNames()
58{
59 if (fieldNames_.size() == 1)
60 {
61 const incompressibleAdjointSolver& adjSolver =
65 const wordList& baseNames =
66 adjointRAS().getAdjointTMVariablesBaseNames();
67 forAll(baseNames, nI)
68 {
70 (adjSolver.extendedVariableName(baseNames[nI]));
71 }
72 }
73}
74
75
76// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
77
79(
80 const fvMesh& mesh,
81 const dictionary& dict,
82 const word& adjointSolverName,
83 const word& primalSolverName
84)
85:
86 objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
87 zones_(mesh_.cellZones().indices(dict.get<wordRes>("zones")))
88{
89 // Append Ua name to fieldNames
91 (
92 1,
94 extendedVariableName("Ua")
95 );
96
97 // Check if cellZones provided include at least one cell
98 checkCellZonesSize(zones_);
99
100 // Allocate dJdTMvar1Ptr_ and dJdTMvar2Ptr_ if needed
102
103 // Allocate source term to the adjoint momentum equations
104 dJdvPtr_.reset
105 (
107 (
108 mesh_,
109 "dJdv" + type(),
111 )
112 );
113 // Allocate terms to be added to volume-based sensitivity derivatives
114 divDxDbMultPtr_.reset
115 (
117 (
118 IOobject
119 (
120 "divDxDbMult" + objectiveName_,
121 mesh_.time().timeName(),
122 mesh_,
125 ),
126 mesh_,
129 )
130 );
131 gradDxDbMultPtr_.reset
132 (
134 (
135 mesh_,
136 ("gradDxdbMult" + type()),
138 )
139 );
140
141 // Allocate direct sensitivity contributions for topology optimisation
142 dJdbPtr_.reset(createZeroFieldPtr<scalar>(mesh_, "dJdb", dimless));
143}
144
145
146// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
147
149{
150 J_ = Zero;
151
152 // References
153 const volVectorField& U = vars_.UInst();
155 const scalarField& V = mesh_.V().field();
156
157 volScalarField integrand(turb->nuEff()*magSqr(twoSymm(fvc::grad(U))));
158
159 for (const label zI : zones_)
160 {
161 const cellZone& zoneI = mesh_.cellZones()[zI];
162 scalarField VZone(V, zoneI);
163 scalarField integrandZone(integrand.primitiveField(), zoneI);
164
165 J_ += 0.5*gSum(integrandZone*VZone);
166 if (mesh_.foundObject<topOVariablesBase>("topOVars"))
167 {
168 const topOVariablesBase& vars =
170 const volScalarField& beta = vars.beta();
171 scalar porosityContr = Zero;
172 for (const label cellI : zoneI)
173 {
174 porosityContr += beta[cellI]*magSqr(U[cellI])*V[cellI];
175 }
176 porosityContr *= vars.getBetaMax();
177 J_ += returnReduce(porosityContr, sumOp<scalar>());
179 }
180
181 return J_;
182}
183
184
186{
187 dJdvPtr_().primitiveFieldRef() = Zero;
188
189 const volVectorField& U = vars_.U();
190 const autoPtr<incompressible::turbulenceModel>& turb = vars_.turbulence();
191 tmp<volSymmTensorField> S = twoSymm(fvc::grad(U));
192
193 // Add source from possible dependencies of nut on U
194 if (mesh_.foundObject<incompressibleAdjointSolver>(adjointSolverName_))
195 {
196 const incompressibleAdjointSolver& adjSolver =
197 mesh_.lookupObject<incompressibleAdjointSolver>
199 const autoPtr<incompressibleAdjoint::adjointRASModel>& adjointRAS =
200 adjSolver.getAdjointVars().adjointTurbulence();
201 tmp<volScalarField> dnutdUMult = 0.5*magSqr(S());
202 tmp<volVectorField> dnutdU = adjointRAS->nutJacobianU(dnutdUMult);
203 if (dnutdU)
204 {
205 for (const label zI : zones_)
206 {
207 const cellZone& zoneI = mesh_.cellZones()[zI];
208 for (const label cellI : zoneI)
209 {
210 dJdvPtr_()[cellI] = dnutdU()[cellI];
211 }
212 }
213 }
214 }
215
216 // Add source from the strain rate magnitude
217 volVectorField integrand(-2.0*fvc::div(turb->nuEff()*S));
218 for (const label zI : zones_)
219 {
220 const cellZone& zoneI = mesh_.cellZones()[zI];
221 for (const label cellI : zoneI)
222 {
223 dJdvPtr_()[cellI] += integrand[cellI];
224 }
225 }
226
227 // Add source from porosity dependencies
228 if (mesh_.foundObject<topOVariablesBase>("topOVars"))
229 {
230 const topOVariablesBase& vars =
231 mesh_.lookupObject<topOVariablesBase>("topOVars");
232 const volScalarField& beta = vars.beta();
233 const scalar betaMax = vars.getBetaMax();
234 for (const label zI : zones_)
235 {
236 const cellZone& zoneI = mesh_.cellZones()[zI];
237 for (const label cellI : zoneI)
238 {
239 dJdvPtr_()[cellI] += 2*betaMax*beta[cellI]*U[cellI];
240 }
241 }
242 }
243}
244
245
247{
248 const volVectorField& U = vars_.U();
249 volScalarField JacobianMultiplier(0.5*magSqr(twoSymm(fvc::grad(U))));
251 (
254 JacobianMultiplier,
255 zones_
256 );
257}
258
259
261{
262 const volVectorField& U = vars_.U();
263 volScalarField JacobianMultiplier(0.5*magSqr(twoSymm(fvc::grad(U))));
265 (
268 JacobianMultiplier,
269 zones_
270 );
271}
272
273
275{
276 // References
277 volScalarField& divDxDbMult = divDxDbMultPtr_();
278 const volVectorField& U = vars_.U();
280 volScalarField integrand(0.5*turb->nuEff()*magSqr(twoSymm(fvc::grad(U))));
281
282 for (const label zI : zones_)
283 {
284 const cellZone& zoneI = mesh_.cellZones()[zI];
285 for (const label cellI : zoneI)
286 {
287 divDxDbMult[cellI] = integrand[cellI];
288 }
289 }
290 divDxDbMult.correctBoundaryConditions();
291}
292
293
295{
296 // References
297 volTensorField& gradDxDbMult = gradDxDbMultPtr_();
298 const volVectorField& U = vars_.U();
299 const autoPtr<incompressible::turbulenceModel>& turb = vars_.turbulence();
300 tmp<volTensorField> gradU = fvc::grad(U);
301 volTensorField integrand(-2.0*turb->nuEff()*(gradU() & twoSymm(gradU())));
302
303 for (const label zI : zones_)
304 {
305 const cellZone& zoneI = mesh_.cellZones()[zI];
306 for (const label cellI : zoneI)
307 {
308 gradDxDbMult[cellI] = integrand[cellI];
309 }
310 }
311 gradDxDbMult.correctBoundaryConditions();
312}
313
314
316{
317 if (mesh_.foundObject<topOVariablesBase>("topOVars"))
318 {
319 scalarField& dJdb = dJdbPtr_().primitiveFieldRef();
320 dJdb = Zero;
321 const volVectorField& U = vars_.UInst();
322 const topOVariablesBase& vars =
323 mesh_.lookupObject<topOVariablesBase>("topOVars");
324 const scalar betaMax = vars.getBetaMax();
325 for (const label zI : zones_)
326 {
327 const cellZone& zoneI = mesh_.cellZones()[zI];
328 for (const label cellI : zoneI)
329 {
330 dJdb[cellI] += betaMax*magSqr(U[cellI]);
331 }
333 }
334}
335
336
337
338
340{
341 populateFieldNames();
342 const label fieldI = fieldNames_.find(matrix.psi().name());
343 if (fieldI == 1)
344 {
345 matrix += weight()*dJdTMvar1Ptr_();
346 }
347 if (fieldI == 2)
348 {
349 matrix += weight()*dJdTMvar2Ptr_();
350 }
351}
352
353
354// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
355
356} // End namespace objectives
357} // End namespace Foam
358
359// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
compressible::turbulenceModel & turb
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
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
void setSize(label n)
Alias for resize().
Definition List.H:536
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
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
Base class for selecting the betaMax value, i.e. the value multiplying the Brinkman penalisation term...
Definition betaMax.H:50
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 volVectorField & UInst() const
Return const reference to velocity.
const volVectorField & U() const
Return const reference to velocity.
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
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.
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}.
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
const volScalarField & dJdb() const
Contribution to field sensitivities.
Definition objectiveI.H:44
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
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
Definition objective.H:199
Computes and minimizes the power dissipation within given cellZones. In the absence of significant vi...
virtual void update_dJdb()
Contribution to field sensitivities.
virtual void update_gradDxDbMultiplier()
Update grad(dx/db multiplier). Volume-based sensitivity term.
objectivePowerDissipation(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Construct from components.
virtual void update_dJdTMvar2()
Update field to be added to the second adjoint turbulence model PDE.
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 div(dx/db multiplier). Volume-based sensitivity term.
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
Base solver class.
Definition solver.H:48
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
Base class for all design variables related to topology optimisation (topO). Provides the lookup func...
scalar getBetaMax() const
Get betaMax value.
virtual const volScalarField & beta() const =0
Get field used for physical interpolations.
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
U
Definition pEqn.H:72
dynamicFvMesh & mesh
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
Type gSum(const FieldField< Field, Type > &f)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
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
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
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.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dictionary dict
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299