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 -------------------------------------------------------------------------------
11 License
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 
38 namespace Foam
39 {
40 
41 namespace objectives
42 {
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 defineTypeNameAndDebug(objectivePowerDissipation, 0);
48 (
49  objectiveIncompressible,
50  objectivePowerDissipation,
51  dictionary
52 );
53 
54 
55 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
56 
57 void objectivePowerDissipation::populateFieldNames()
58 {
59  if (fieldNames_.size() == 1)
60  {
61  const incompressibleAdjointSolver& adjSolver =
62  mesh_.lookupObject<incompressibleAdjointSolver>(adjointSolverName_);
63  const autoPtr<incompressibleAdjoint::adjointRASModel>& adjointRAS =
64  adjSolver.getAdjointVars().adjointTurbulence();
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  (
106  createZeroFieldPtr<vector>
107  (
108  mesh_,
109  "dJdv" + type(),
111  )
112  );
113  // Allocate terms to be added to volume-based sensitivity derivatives
114  divDxDbMultPtr_.reset
115  (
116  new volScalarField
117  (
118  IOobject
119  (
120  "divDxDbMult" + objectiveName_,
121  mesh_.time().timeName(),
122  mesh_,
125  ),
126  mesh_,
129  )
130  );
131  gradDxDbMultPtr_.reset
132  (
133  createZeroFieldPtr<tensor>
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 =
169  mesh_.lookupObject<topOVariablesBase>("topOVars");
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>());
178  }
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  }
332  }
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 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
objectivePowerDissipation(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Construct from components.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
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...
const word adjointSolverName_
Definition: objective.H:65
Base solver class.
Definition: solver.H:45
scalar getBetaMax() const
Get betaMax value.
autoPtr< volScalarField > dJdTMvar2Ptr_
Second turbulence model variable.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:86
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
void update_dJdTMvar(autoPtr< volScalarField > &dJdTMvarPtr, tmp< volScalarField >(incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const, const volScalarField &JacobianMultiplier, const labelList &zones)
Compute dJdTMVar{1,2}.
compressible::turbulenceModel & turb
void checkCellZonesSize(const labelList &zoneIDs) const
Check if cellZones provided include at least one cell.
virtual void update_dJdv()
Update values to be added to the adjoint outlet velocity.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
void push_back(const T &val)
Append an element at the end of the list.
Definition: ListI.H:227
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Macros for easy insertion into run-time selection tables.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const word objectiveName_
Definition: objective.H:67
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
wordList fieldNames_
List of adjoint fields for which this objective will contribute sources to their equations.
Definition: objective.H:130
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual scalar J()
Return the objective function value.
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
const volScalarField & dJdb() const
Contribution to field sensitivities.
Definition: objectiveI.H:44
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
autoPtr< volVectorField > dJdvPtr_
defineTypeNameAndDebug(objectivePartialVolume, 1)
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
const fvMesh & mesh_
Definition: objective.H:63
virtual void update_dJdTMvar1()
Update field to be added to the first adjoint turbulence model PDE.
virtual const volScalarField & beta() const =0
Get field used for physical interpolations.
virtual void addSource(fvScalarMatrix &matrix)
Add source terms to the adjoint turbulence model equations.
scalar J_
Objective function value and weight.
Definition: objective.H:76
label find(const T &val) const
Find index of the first occurrence of the value.
Definition: UList.C:173
const volVectorField & U() const
Return const reference to velocity.
virtual scalar weight() const
Return the objective function weight.
Definition: objective.C:347
Istream and Ostream manipulators taking arguments.
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
Base class for selecting the betaMax value, i.e. the value multiplying the Brinkman penalisation term...
Definition: betaMax.H:49
virtual void update_divDxDbMultiplier()
Update div(dx/db multiplier). Volume-based sensitivity term.
virtual void update_dJdTMvar2()
Update field to be added to the second adjoint turbulence model PDE.
A subset of mesh cells.
Definition: cellZone.H:58
addToRunTimeSelectionTable(objectiveGeometric, objectivePartialVolume, dictionary)
void allocatedJdTurbulence()
Allocate fields related to the differentiation of turbulence models, if necessary.
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition: objective.H:194
virtual void update_gradDxDbMultiplier()
Update grad(dx/db multiplier). Volume-based sensitivity term.
const incompressibleVars & vars_
List< word > wordList
List of word.
Definition: fileName.H:59
dimensionedScalar pow3(const dimensionedScalar &ds)
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
U
Definition: pEqn.H:72
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual void update_dJdb()
Contribution to field sensitivities.
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:678
void correctBoundaryConditions()
Correct boundary field.
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
Definition: objective.H:199
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
const volVectorField & UInst() const
Return const reference to velocity.
autoPtr< volScalarField > dJdTMvar1Ptr_
First turbulence model variable.
Base class for all design variables related to topology optimisation (topO). Provides the lookup func...
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual tmp< volScalarField > nutJacobianTMVar2() const
Jacobian of nut wrt the second turbulence model variable.
Abstract base class for objective functions in incompressible flows.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt the first turbulence model variable.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127