topODesignVariables.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 
29 #include "localIOdictionary.H"
30 #include "topODesignVariables.H"
31 #include "MeshObject.H"
32 #include "wallFvPatch.H"
33 #include "cutFaceIso.H"
34 #include "cutCellIso.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(topODesignVariables, 0);
43  (
44  designVariables,
45  topODesignVariables,
46  designVariables
47  );
48 }
49 
50 
51 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
52 
54 (
55  const scalarField& correction,
56  const label fluidID
57 )
58 {
59  DebugInfo
60  << "Updating design variables for field " << fluidID << endl;
61  const label n = mesh_.nCells();
62  SubField<scalar> localCorrection(correction, n, fluidID*n);
63  SubField<scalar> field(*this, n, fluidID*n);
64 
65  // Update porosity in adjoint porous cells
67  {
68  forAll(field, cellI)
69  {
70  field[cellI] +=
71  min
72  (
73  max
74  (
75  field[cellI] + localCorrection[cellI],
76  scalar(0)
77  ),
78  1.
79  )
80  - field[cellI];
81  }
82  }
83  else
84  {
85  for (label cellZoneID : zones_.adjointPorousZoneIDs())
86  {
87  const labelList& zoneCells = mesh_.cellZones()[cellZoneID];
88  for (label cellI : zoneCells)
89  {
90  field[cellI] +=
91  min
92  (
93  max
94  (
95  field[cellI] + localCorrection[cellI],
96  scalar(0)
97  ),
98  1.
99  )
100  - field[cellI];
101  }
102  }
103  }
104 }
105 
106 
108 {
109  SubField<scalar> alpha(*this, mesh_.nCells());
110  // Zero porosity in the cells next to IO
111  for (label cellI : zones_.IOCells())
112  {
113  alpha[cellI] = 0.;
114  }
115 
116  // Apply fixed porosity
117  forAll(zones_.fixedPorousZoneIDs(), zI)
118  {
119  const label cellZoneID = zones_.fixedPorousZoneIDs()[zI];
120  const labelList& zoneCells = mesh_.cellZones()[cellZoneID];
121  const scalar alphaValue(zones_.fixedPorousValues()[zI]);
122  for (label cellI : zoneCells)
123  {
124  alpha[cellI] = alphaValue;
125  }
126  }
127 
128  // Apply fixed zero porosity
129  for (label cellZoneID : zones_.fixedZeroPorousZoneIDs())
130  {
131  const labelList& zoneCells = mesh_.cellZones()[cellZoneID];
132  for (label cellI : zoneCells)
133  {
134  alpha[cellI] = 0.;
135  }
136  }
137 }
138 
139 
141 {
142  const scalar maxChange(gMax(mag(correction)));
143  Info<< "maxInitChange/maxChange \t"
144  << maxInitChange_() << "/" << maxChange << endl;
145  const scalar eta(maxInitChange_() / maxChange);
146  Info<< "Setting eta value to " << eta << endl;
147  correction *= eta;
148 
149  return eta;
150 }
151 
152 
154 (
155  const label fluidID,
156  const bool activeIO
157 )
158 {
159  const label offset(fluidID*mesh_.nCells());
160  label varI(activeDesignVariables_.size());
161  activeDesignVariables_.setSize(offset + mesh_.nCells(), -1);
162  // Set active design variables
163  // If specific porosity zones are prescribed, use them directly
164  if (!zones_.adjointPorousZoneIDs().empty())
165  {
166  for (label cellZoneID : zones_.adjointPorousZoneIDs())
167  {
168  for (const label var : mesh_.cellZones()[cellZoneID])
169  {
170  activeDesignVariables_[varI++] = var + offset;
171  }
172  }
173  }
174  // Else, pick up all cells in non-constant porosity zones
175  else
176  {
177  boolList isActiveDV(mesh_.nCells(), true);
178  // Exclude cells with fixed porosity
179  for (label cellZoneID : zones_.fixedPorousZoneIDs())
180  {
181  for (label cellI : mesh_.cellZones()[cellZoneID])
182  {
183  isActiveDV[cellI] = false;
184  }
185  }
186  for (label cellZoneID : zones_.fixedZeroPorousZoneIDs())
187  {
188  for (label cellI : mesh_.cellZones()[cellZoneID])
189  {
190  isActiveDV[cellI] = false;
191  }
192  }
193  if (!activeIO)
194  {
195  for (label cellI : zones_.IOCells())
196  {
197  isActiveDV[cellI] = false;
198  }
199  }
200 
201  // Set active design variables
202  forAll(isActiveDV, cellI)
203  {
204  if (isActiveDV[cellI])
205  {
206  activeDesignVariables_[varI++] = offset + cellI;
207  }
208  }
209  }
210  activeDesignVariables_.setSize(varI);
211 }
212 
213 
215 (
216  const word& name,
217  const label fluidID,
218  const bool setIOValues
219 )
220 {
221  const label offset(fluidID*mesh_.nCells());
223  {
224  SubField<scalar>(*this, mesh_.nCells(), offset) =
225  scalarField(name, *this, mesh_.nCells());
226 
227  /*
228  // Set values next to IO boundaries if needed
229  if (setIOValues)
230  {
231  forAll(mesh_.boundary(), patchI)
232  {
233  const fvPatch& patch = mesh_.boundary()[patchI];
234  if (patch.type() == "patch")
235  {
236  const labelList& faceCells = patch.faceCells();
237  const scalarField& pf = volField.boundaryField()[patchI];
238  forAll(faceCells, fI)
239  {
240  scalarField::operator[](offset + faceCells[fI]) = pf[fI];
241  }
242  }
243  }
244  }
245  */
246  }
247 }
248 
249 
251 {
252  // Set active design variables
253  setActiveDesignVariables();
254 
255  // Read in values from file, if present
256  readField("alpha", 0, true);
257 
258  if (regularisation_.growFromWalls())
259  {
260  scalarField& alpha = *this;
261  for (const fvPatch& patch : mesh_.boundary())
262  {
263  if (isA<wallFvPatch>(patch))
264  {
265  const labelList& faceCells = patch.faceCells();
266  forAll(faceCells, cI)
267  {
268  alpha[faceCells[cI]] = 1.;
269  }
270  }
271  }
272  }
273 
274  // Make sure alpha has fixed values where it should
275  scalarField dummyCorrection(mesh_.nCells(), Zero);
276  update(dummyCorrection);
277 
278  // Read bounds for design variables, if present
280 }
281 
282 
283 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
284 
285 Foam::topODesignVariables::topODesignVariables
286 (
287  fvMesh& mesh,
288  const dictionary& dict
289 )
290 :
291  topODesignVariables(mesh, dict, mesh.nCells())
292 {}
293 
294 
295 Foam::topODesignVariables::topODesignVariables
296 (
297  fvMesh& mesh,
298  const dictionary& dict,
299  const label size
300 )
301 :
303  designVariables(mesh, dict, size),
304  alpha_(SubField<scalar>(*this, mesh.nCells(), 0)),
305  regularisation_
306  (
307  mesh,
308  alpha_,
309  zones_,
310  dict_.subDict("regularisation")
311  ),
312  writeAllFields_
313  (
314  dict.getOrDefaultCompat<bool>
315  (
316  "writeAllFields", {{"writeAllAlphaFields", 2306}}, false
317  )
318  ),
319  addFvOptions_(dict.getOrDefault<bool>("addFvOptions", false))
320 {}
321 
322 
323 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
324 
326 (
327  fvMesh& mesh,
328  const dictionary& dict
329 )
330 {
332 }
333 
334 
335 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
336 
338 {
339  return regularisation_.beta();
340 }
341 
342 
344 (
345  const word& interpolationFieldName
346 ) const
347 {
348  return beta().primitiveField();
349 }
350 
351 
353 (
355  const topOInterpolationFunction& interpolationFunc,
356  const FieldField<Field, scalar>& fluidValues,
357  const scalarField& solidValues,
358  const label fieldi,
359  const word& interpolationFieldName
360 ) const
361 {
362  const scalarField& indicator = interpolationField(interpolationFieldName);
363  scalarField interpolant(indicator.size(), Zero);
364  interpolationFunc.interpolate(indicator, interpolant);
365 
366  // Interpolate field values
367  const scalar diff(solidValues[fieldi] - fluidValues[0][fieldi]);
368  field.primitiveFieldRef() = fluidValues[0][fieldi] + interpolant*diff;
369  field.correctBoundaryConditions();
370 }
371 
372 
374 (
375  scalarField& sens,
376  const topOInterpolationFunction& interpolationFunc,
377  const FieldField<Field, scalar>& fluidValues,
378  const scalarField& solidValues,
379  const label fieldi,
380  const word& designVariablesName,
381  const word& interpolationFieldName
382 ) const
383 {
384  const scalarField& indicator = interpolationField(interpolationFieldName);
385  sens *=
386  (solidValues[fieldi] - fluidValues[0][fieldi])
387  *interpolationFunc.derivative(indicator);
388 }
389 
390 
392 (
394  const topOInterpolationFunction& interpolationFunc
395 ) const
396 {
397  const scalarField& beta = this->beta().primitiveField();
398  scalarField interpolant(beta.size(), Zero);
399  interpolationFunc.interpolate(beta, interpolant);
400  field *= scalar(1) - interpolant;
401 }
402 
403 
405 (
406  scalarField& sens,
407  const topOInterpolationFunction& interpolationFunc,
408  const word& designVariablesName
409 ) const
410 {
411  const scalarField& beta = this->beta().primitiveField();
412  sens *= - interpolationFunc.derivative(beta);
413 }
414 
415 
417 (
418  const word& interpolationFieldName,
419  const topOInterpolationFunction& interpolationFunc
420 ) const
421 {
422  const scalarField& beta = this->beta().primitiveField();
424  interpolationFunc.interpolate(beta, tinterpolant.ref());
425  return tinterpolant;
426 }
427 
428 
430 (
431  const word& interpolationFieldName,
432  const topOInterpolationFunction& interpolationFunc
433 ) const
434 {
435  const scalarField& beta = this->beta().primitiveField();
436  return interpolationFunc.derivative(beta);
437 }
438 
439 
441 {
442  // Update alpha values
443  updateField(correction);
444 
445  // Fix alpha in zones
446  applyFixedValues();
447 
448  // Update the beta field
449  regularisation_.updateBeta();
450 
451  // Though the mesh is kept constant, the distance from wall may change
452  // if the method computing it includes fvOptions that depend on the
453  // indicator field.
454  // Trick wallDist into updating it
455  if (mesh_.foundObject<UpdateableMeshObject<fvMesh>>("wallDist"))
456  {
457  mesh_.lookupObjectRef<UpdateableMeshObject<fvMesh>>("wallDist").
458  movePoints();
459  }
460 
461  // Write the 0.5 beta iso-line to files, as an indication of the
462  // fluid-solid interface
463  labelList changedFaces(mesh_.nFaces(), -1);
464  List<wallPointData<label>> changedFacesInfo(mesh_.nFaces());
465  writeFluidSolidInterface(-beta(), -0.5, changedFaces, changedFacesInfo);
466 }
467 
468 
470 {
471  return true;
472 }
473 
474 
476 (
477  adjointSensitivity& adjointSens
478 )
479 {
480  // Raw sensitivities field.
481  // Does not include the adjoint to the regularisation and projection steps
482  const scalarField& fieldSens = adjointSens.fieldSensPtr()->primitiveField();
483 
484  // Return field (complete sensitivities)
485  auto tobjectiveSens(tmp<scalarField>::New(fieldSens));
486  scalarField& objectiveSens = tobjectiveSens.ref();
487 
488  // Add part due to regularisation and projection
489  regularisation_.postProcessSens(objectiveSens);
490 
491  // Write final sensitivities field
492  if (writeAllFields_ && mesh_.time().writeTime())
493  {
494  volScalarField sens
495  (
496  IOobject
497  (
498  "topOSens" + adjointSens.getAdjointSolver().solverName(),
499  mesh_.time().timeName(),
500  mesh_,
503  ),
504  mesh_,
506  );
507  sens.primitiveFieldRef() = objectiveSens;
508  sens.write();
509  }
510 
511  return tobjectiveSens;
512 }
513 
514 
516 {
517  // Rest of the contrsuctor initialisation
518  initialize();
519 }
520 
521 
523 (
524  const PtrList<primalSolver>& primalSolvers,
526 )
527 {
528  // WIP
529  if (addFvOptions_)
530  {
531  for (const primalSolver& solver : primalSolvers)
532  {
534  }
535  for (const adjointSolverManager& manager : adjointSolverManagers)
536  {
537  const PtrList<adjointSolver>& adjointSolvers =
538  manager.adjointSolvers();
539  for (const adjointSolver& solver : adjointSolvers)
540  {
542  }
543  }
544  }
545 }
546 
547 
549 {
550  if (writeAllFields_ && mesh_.time().writeTime())
551  {
553  (
554  IOobject
555  (
556  "alpha",
557  mesh_.time().timeName(),
558  mesh_,
561  ),
562  mesh_,
564  );
565  alpha.primitiveFieldRef() = alpha_;
566 
567  alpha.write();
568  }
569 }
570 
571 
572 bool Foam::topODesignVariables::writeData(Ostream& os) const
573 {
574  const scalarField& alpha = alpha_;
575  alpha.writeEntry("alpha", os);
576 
577  return true;
578 }
579 
580 
581 // ************************************************************************* //
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
virtual tmp< scalarField > assembleSensitivities(adjointSensitivity &sens)
Assemble sensitivity derivatives, by combining the part related to the primal and adjoint solution wi...
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:373
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual void interpolationSensitivities(scalarField &sens, const topOInterpolationFunction &interpolationFunc, const FieldField< Field, scalar > &fluidValues, const scalarField &solidValues, const label fieldi, const word &designVariablesName, const word &interpolationFieldName="beta") const
Post-processing sensitivities due to interpolations based on the indicator fields.
rDeltaTY field()
Base solver class.
Definition: solver.H:45
virtual tmp< scalarField > penaltySensitivities(const word &interpolationFieldName, const topOInterpolationFunction &interpolationFunc) const
Return the penalty term derivative.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Abstract base class for adjoint-based sensitivities.
Base class for primal solvers.
Definition: primalSolver.H:46
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
virtual void interpolate(volScalarField &field, const topOInterpolationFunction &interpolationFunc, const FieldField< Field, scalar > &fluidValues, const scalarField &solidValues, const label fieldi, const word &interpolationFieldName="beta") const
Interpolate between the given field and solid values.
topOZones zones_
Cell zones useful for defining the constant and changing parts of the domain in topO.
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
A class that holds the data needed to identify things (zones, patches) in a dynamic mesh...
Definition: DynamicID.H:48
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
Design variables for porosity-based topology optimisation (topO) problems.
const adjointSolver & getAdjointSolver() const
Const access to adjoint solver.
virtual const volScalarField & beta() const
Get the indicator field.
Macros for easy insertion into run-time selection tables.
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:104
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual scalar computeEta(scalarField &correction)
Compute eta if not set in the first step.
virtual void addTopOFvOptions() const
Add topO fvOptions.
Definition: solver.C:102
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
void writeEntry(const word &keyword, Ostream &os) const
Write as a dictionary entry with keyword.
virtual const scalarField & interpolationField(const word &interpolationFieldName="beta") const
Return interpolant.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void writeDesignVars()
Write porosity field to file.
virtual void nullifyInSolidSensitivities(scalarField &sens, const topOInterpolationFunction &interpolationFunc, const word &designVariablesName) const
Nullify given field in the solid domain.
static autoPtr< topODesignVariables > New(fvMesh &mesh, const dictionary &dict)
Construct and return the selected design variables.
virtual void setActiveDesignVariables(const label fluidID=0, const bool activeIO=false)
Set active design variables.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
virtual tmp< scalarField > derivative(const scalarField &arg) const =0
Return of function with respect to the argument field.
#define DebugInfo
Report an information message using Foam::Info.
virtual void update(scalarField &correction)
Update design variables based on a given correction.
virtual bool writeData(Ostream &) const
The writeData function required by the regIOobject write operation.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
mesh update()
defineTypeNameAndDebug(combustionModel, 0)
DynamicID< cellZoneMesh > cellZoneID
Foam::cellZoneID.
Definition: ZoneIDs.H:38
void applyFixedValues()
Apply fixed values in certain zones.
virtual void nullifyInSolid(scalarField &field, const topOInterpolationFunction &interpolationFunc) const
Nullify given field in the solid domain.
virtual void setInitialValues()
Set initial values of the design variables.
virtual void updateField(const scalarField &correction, const label fluidID=0)
Update the design variables given their correction.
Abstract base class for defining design variables.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nCells() const noexcept
Number of mesh cells.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
const word & solverName() const
Return the solver name.
Definition: solverI.H:30
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
virtual tmp< scalarField > penalty(const word &interpolationFieldName, const topOInterpolationFunction &interpolationFunc) const
Return the Brinkman penalisation term.
Automatically write from objectRegistry::writeObject()
const labelList & adjointPorousZoneIDs() const
Cell zone IDs in which porosity is allowed to change.
Definition: topOZones.H:193
const std::string patch
OpenFOAM patch number as a std::string.
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
virtual void initialize()
Part of the constructor initialisation.
virtual void interpolate(const scalarField &arg, scalarField &res) const =0
Interpolate argument to result.
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
label n
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
virtual void addFvOptions(const PtrList< primalSolver > &primalSolver, const PtrList< adjointSolverManager > &adjointSolverManagers)
Automatically add fvOptions depending on the design variables to the primal and adjoint solvers...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Base class for all design variables related to topology optimisation (topO). Provides the lookup func...
PtrList< adjointSolverManager > & adjointSolverManagers
Definition: createFields.H:8
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual bool globalSum() const
Whether to use global sum when computing matrix-vector products in update methods.
List< bool > boolList
A List of bools.
Definition: List.H:60
const autoPtr< volScalarField > & fieldSensPtr() const
Get the fieldSensPtr.
Definition: sensitivity.H:152
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
void readField(const word &name, const label fluidID=0, const bool setIOValues=false)
Read field with (path of) the design variables and store input in the design variables list with opti...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127