levelSetDesignVariables.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) 2022-2023 PCOpt/NTUA
9  Copyright (C) 2022-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 
30 #include "wallDist.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(levelSetDesignVariables, 1);
43 (
47 );
48 
49 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
50 
52 {
53  scalarField& vars = *this;
54  if (localIOdictionary::found("alpha"))
55  {
56  vars = (scalarField("alpha", *this, vars.size()));
57  }
58  else
59  {
60  // Initialise as the distance from the wall patches of the initial
61  // domain
62  const labelHashSet wallPatchIDs =
65  (
66  IOobject
67  (
68  "yLevelSet",
69  mesh_.time().timeName(),
70  mesh_,
73  ),
74  mesh_,
76  patchDistMethod::patchTypes<scalar>(mesh_, wallPatchIDs)
77  );
79  (
80  dict_.subDict("initialisation"),
81  mesh_,
82  wallPatchIDs
83  )->correct(y);
84  vars = y.primitiveField();
85 
86  if (debug)
87  {
89  }
90  }
91 }
92 
93 
95 {
97 
98  // Safety - set beta equal to zero next to the IO cells
99  for (const label IOcell : zones_.IOCells())
100  {
101  betaIf[IOcell] = Zero;
102  }
103 
104  const labelList& fixedZeroPorousZones = zones_.fixedZeroPorousZoneIDs();
105  const labelList& fixedPorousZones = zones_.fixedPorousZoneIDs();
106  const scalarList& fixedPorousValues = zones_.fixedPorousValues();
107 
108  // Apply fixed porosity
109  for (const label zoneID : fixedZeroPorousZones)
110  {
111  const labelList& zone = mesh_.cellZones()[zoneID];
112  for (const label cellI : zone)
113  {
114  betaIf[cellI] = Zero;
115  }
116  }
117 
118  // Apply fixed porosity
119  forAll(fixedPorousZones, zI)
120  {
121  const label zoneID = fixedPorousZones[zI];
122  const scalar value = fixedPorousValues[zI];
123  const labelList& zone = mesh_.cellZones()[zoneID];
124  for (const label cellI : zone)
125  {
126  betaIf[cellI] = value >= 0 ? 0 : 1;
127  }
128  }
129 
131 }
132 
133 
135 {
138  {
139  boolList isActiveVar(mesh_.nCells(), true);
140 
141  const labelList& fixedZeroPorousZones =
143  for (const label zoneID : fixedZeroPorousZones)
144  {
145  const labelList& zone = mesh_.cellZones()[zoneID];
146  for (const label cellI : zone)
147  {
148  isActiveVar[cellI] = false;
149  }
150  }
151 
152  const labelList& fixedPorousZones = zones_.fixedPorousZoneIDs();
153  for (const label zoneID : fixedPorousZones)
154  {
155  const labelList& zone = mesh_.cellZones()[zoneID];
156  for (const label cellI : zone)
157  {
158  isActiveVar[cellI] = false;
159  }
160  }
161 
162  if (!activeIO)
163  {
164  for (label cellI : zones_.IOCells())
165  {
166  isActiveVar[cellI] = false;
167  }
168  }
169 
170  label iVar(0);
171  forAll(isActiveVar, cI)
172  {
173  if (isActiveVar[cI])
174  {
175  activeDesignVariables_[iVar++] = cI;
176  }
177  }
179  }
180  else
181  {
182  const labelList& adjointPorousZoneIDs = zones_.adjointPorousZoneIDs();
183 
184  label iVar(0);
185  for (const label cellZoneID : adjointPorousZoneIDs)
186  {
187  const labelList& zone = mesh_.cellZones()[cellZoneID];
188 
189  for (const label cellI : zone)
190  {
191  activeDesignVariables_[iVar] = cellI;
192  iVar++;
193  }
194  }
196  }
197 }
198 
199 
201 {
202  // Compute the beta field by passing the distance field through
203  // a Heaviside function
206  beta = 1 - beta;
207  // Apply fixed values if necessary
209 
211 }
212 
213 
215 {
216  Info<< "Re-initilising the level-set distance field" << nl << endl;
217 
219  (
220  IOobject
221  (
222  "yLevelSet",
223  mesh_.time().timeName(),
224  mesh_,
227  ),
228  mesh_,
230  wordList
231  (
232  mesh_.boundary().size(),
234  )
235  );
236  y.primitiveFieldRef() = aTilda_.primitiveFieldRef();
237  y.correctBoundaryConditions();
238 
239  changedFaces_.clear();
240  changedFaces_.setSize(mesh_.nFaces(), -1);
241 
242  changedFacesInfo_.clear();
243  changedFacesInfo_.setSize(mesh_.nFaces());
244 
245  writeFluidSolidInterface(aTilda_, 0, changedFaces_, changedFacesInfo_);
246 
247  List<wallPointData<label>> allFaceInfo(mesh_.nFaces());
248  allCellInfo_.clear();
249  allCellInfo_.setSize(mesh_.nCells());
250 
252  (
253  mesh_,
254  changedFaces_,
255  changedFacesInfo_,
256  allFaceInfo,
257  allCellInfo_,
258  mesh_.globalData().nTotalCells() + 1
259  );
260 
261  // Transfer the distance from cellInfo to the alphaTilda field
262  forAll(allCellInfo_, celli)
263  {
264  if (allCellInfo_[celli].valid(wave.data()))
265  {
266  signedDistances_[celli] =
267  sign(aTilda_[celli])*Foam::sqrt(allCellInfo_[celli].distSqr());
268  }
269  }
270  signedDistances_.correctBoundaryConditions();
271 }
272 
273 
274 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
275 
276 levelSetDesignVariables::levelSetDesignVariables
277 (
278  fvMesh& mesh,
279  const dictionary& dict
280 )
281 :
283  designVariables(mesh, dict, mesh.nCells()),
284  radius_
285  (regularisationRadius::New(mesh, dict.subDict("regularisation"), false)),
286  regularisation_
287  (regularisationPDE::New(mesh, dict.subDict("regularisation"), zones_)),
288  aTilda_
289  (
290  IOobject
291  (
292  "aTilda",
293  mesh_.time().timeName(),
294  mesh_,
295  IOobject::READ_IF_PRESENT,
296  IOobject::AUTO_WRITE
297  ),
298  mesh_,
300  zeroGradientFvPatchField<scalar>::typeName
301  ),
302  signedDistances_
303  (
304  IOobject
305  (
306  "signedDistances",
307  mesh_.time().timeName(),
308  mesh_,
309  IOobject::READ_IF_PRESENT,
310  IOobject::AUTO_WRITE
311  ),
312  mesh_,
314  zeroGradientFvPatchField<scalar>::typeName
315  ),
316  interpolation_
317  (
318  topOInterpolationFunction::New(mesh_, dict_.subDict("interpolation"))
319  ),
320  beta_
321  (
322  IOobject
323  (
324  "beta",
325  mesh_.time().timeName(),
326  mesh_,
327  IOobject::READ_IF_PRESENT,
328  IOobject::AUTO_WRITE
329  ),
330  mesh_,
332  ),
333  fixATildaValues_(dict.getOrDefault<bool>("fixATildaValues", true)),
334  writeAllDistanceFields_
335  (dict.getOrDefault<bool>("writeAllDistanceFields", false)),
336  changedFaces_(),
337  changedFacesInfo_(),
338  allCellInfo_()
339 {
340  // Read the alpha field if present, or set it based on the distance field
341  readField();
342 
343  // Read bounds of design variables if present.
344  // If not, use the maximum distance field in the fluid to set them.
345  scalar maxDist = gMax(*this);
346  scalar lowerBound =
347  localIOdictionary::getOrDefault<scalar>("lowerBound", -maxDist - SMALL);
348  scalar upperBound =
349  localIOdictionary::getOrDefault<scalar>("upperBound", maxDist + SMALL);
350  readBounds
351  (
352  autoPtr<scalar>(new scalar(lowerBound)),
353  autoPtr<scalar>(new scalar(upperBound))
354  );
355  DebugInfo
356  << "Using lower/upper bounds "
357  << lowerBounds_()[0] << "/" << upperBounds_()[0]
358  << endl;
359 
360  // Update the beta field based on the initial design vars
361  scalarField zeroUpdate(scalarField::size(), Zero);
362  update(zeroUpdate);
363 
364  // Determine which design variables are active
366 }
367 
368 
369 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
370 
372 (
373  fvMesh& mesh,
374  const dictionary& dict
375 )
376 {
378  (
380  );
381 }
382 
383 
384 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
387 {
388  return beta_;
389 }
390 
391 
393 {
395 
396  // Compute the regularised design variables
397  regularisation_->regularise
398  (
400  true, radius_(), upperBounds_()[0], fixATildaValues_
401  );
403 
405  {
406  writeDesignVars();
407  aTilda_.write();
408  }
409 
410  // Compute signed distances based on aTilda
412 
413  // Set beta based on the signed distances
414  updateBeta();
415 
417  {
419  beta_.write();
420  }
421 
422  // Though the mesh is kept constant, the distance from wall may change
423  // due to fvOptions depending on beta. Trick wallDist into updating it
424  if (mesh_.foundObject<UpdateableMeshObject<fvMesh>>("wallDist"))
425  {
427  movePoints();
428  }
429 }
430 
431 
433 {
434  const scalar maxChange(gMax(mag(correction)));
435  Info<< "maxInitChange/maxChange \t"
436  << maxInitChange_() << "/" << maxChange << endl;
437  const scalar eta(maxInitChange_() / maxChange);
438  Info<< "Setting eta value to " << eta << endl;
439  correction *= eta;
440 
441  return eta;
442 }
443 
444 
446 {
447  return true;
448 }
449 
450 
452 (
453  adjointSensitivity& adjointSens
454 )
455 {
456  // Raw sensitivities field
457  const scalarField& fieldSens = adjointSens.fieldSensPtr()->primitiveField();
458 
459  // Return field
460  auto tobjectiveSens(tmp<scalarField>::New(fieldSens));
461  scalarField& objectiveSens = tobjectiveSens.ref();
462 
463  // Multiply with dBetadAtilda
464  objectiveSens *=
466 
467  // Solve the adjoint to the regularisation equation
469  regularise(aTilda_, objectiveSens, objectiveSens, false, radius_());
470 
471  objectiveSens *= mesh_.V();
472 
474  {
475  volScalarField sens
476  (
477  IOobject
478  (
479  "sens" + adjointSens.getAdjointSolver().solverName(),
480  mesh_.time().timeName(),
481  mesh_,
484  ),
485  mesh_,
487  );
488  sens.primitiveFieldRef() = objectiveSens;
489  sens.write();
490  }
491 
492  return tobjectiveSens;
493 }
494 
495 
497 {
499  {
501  (
502  IOobject
503  (
504  "alpha",
505  mesh_.time().timeName(),
506  mesh_,
509  ),
510  mesh_,
512  );
513  alpha.primitiveFieldRef() = *this;
514  alpha.correctBoundaryConditions();
515 
516  alpha.write();
517  }
518 }
519 
520 
521 bool levelSetDesignVariables::writeData(Ostream& os) const
522 {
523  // Lower and upper bound values might be computed based on the initial
524  // distance field. Write to file to enable continuation
525  os.writeEntry("lowerBound", lowerBounds_()[0]);
526  os.writeEntry("upperBound", upperBounds_()[0]);
527 
528  scalarField::writeEntry("alpha", os);
529 
530  return true;
531 }
532 
533 
534 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
535 
536 } // End namespace Foam
537 
538 // ************************************************************************* //
Base class for selecting the regulatisation radius.
labelList activeDesignVariables_
Which of the design variables will be updated.
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
dimensionedScalar sign(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Lookup and return non-const reference to the object of the given Type. Fatal if not found or the wron...
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.
autoPtr< regularisationPDE > regularisation_
Regularisation mechanism.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
virtual bool globalSum() const
Whether to use global sum when computing matrix-vector products in update methods.
volScalarField aTilda_
The regularised field.
void updateSignedDistances()
Make aTilda a signed distance field based on the zero iso-surface of the current aTilda field...
topOZones zones_
Cell zones useful for defining the constant and changing parts of the domain in topO.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
void applyFixedPorosityValues()
Apply fixed values in certain zones.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const scalarList & fixedPorousValues() const
Values of alpha for each fixed porous zone.
Definition: topOZones.H:177
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:206
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:720
Base class for selecting the regulatisation PDE.
autoPtr< regularisationRadius > radius_
The regularisation radius.
Ignore writing from objectRegistry::writeObject()
bool writeTime() const noexcept
True if this is a write interval.
Definition: TimeStateI.H:73
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual void update(scalarField &correction)
Update design variables based on a given correction.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const adjointSolver & getAdjointSolver() const
Const access to adjoint solver.
volScalarField signedDistances_
The signed distances field.
const cellZone & IOCells() const
Cells next to IO boundaries.
Definition: topOZones.H:201
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
Macros for easy insertion into run-time selection tables.
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
word timeName
Definition: getTimeIndex.H:3
scalar y
static autoPtr< levelSetDesignVariables > New(fvMesh &mesh, const dictionary &dict)
Return a reference to the selected turbulence model.
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:46
const labelList & fixedZeroPorousZoneIDs() const
Cell zone IDs with fixed zero porosity values.
Definition: topOZones.H:185
autoPtr< topOInterpolationFunction > interpolation_
Function to transorm signed distances to the indicator field beta_.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
virtual const volScalarField & beta() const
Const reference to the beta field.
Signed distance field design variables for level-set based topology optimization (topO).
autoPtr< scalar > maxInitChange_
Maximum design variables&#39; change in the first optimisation cycle.
void readBounds(autoPtr< scalar > lowerBoundPtr=nullptr, autoPtr< scalar > upperBoundPtr=nullptr)
Read bounds for design variables, if present.
void readField()
Read the design variables field.
bool writeAllDistanceFields_
Write all fields related to the distance calculation (debugging)
#define DebugInfo
Report an information message using Foam::Info.
volScalarField beta_
The indicator field.
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
label size() const noexcept
The number of elements in the container.
Definition: UList.H:671
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
DynamicID< cellZoneMesh > cellZoneID
Foam::cellZoneID.
Definition: ZoneIDs.H:38
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
List< word > wordList
List of word.
Definition: fileName.H:59
autoPtr< scalarField > lowerBounds_
Lower bounds of the design variables.
virtual tmp< scalarField > assembleSensitivities(adjointSensitivity &adjointSens)
Add contributions from the adjoint to the regularization PDE, the derivative of the interpolation fun...
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Abstract base class for defining design variables.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nCells() const noexcept
Number of mesh cells.
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.
Automatically write from objectRegistry::writeObject()
autoPtr< scalarField > upperBounds_
Upper bounds of the design variables.
const labelList & adjointPorousZoneIDs() const
Cell zone IDs in which porosity is allowed to change.
Definition: topOZones.H:193
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:678
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
void updateBeta()
Update the beta field.
void operator+=(const UList< scalar > &)
Definition: Field.C:799
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
virtual scalar computeEta(scalarField &correction)
Compute eta if not set in the first step.
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void setActiveDesignVariables(bool activeIO=false)
Determine which design variables are active.
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
static autoPtr< patchDistMethod > New(const dictionary &dict, const fvMesh &mesh, const labelHashSet &patchIDs, const word &defaultPatchDistMethod=word::null)
List< bool > boolList
A List of bools.
Definition: List.H:60
const labelList & fixedPorousZoneIDs() const
Cell zone IDs with fixed porosity values.
Definition: topOZones.H:169
labelHashSet findPatchIDs() const
Find patch indices for a given polyPatch type.
bool fixATildaValues_
Fix aTilda values in fixed{Zero}PorousZones and IOcells.
const autoPtr< volScalarField > & fieldSensPtr() const
Get the fieldSensPtr.
Definition: sensitivity.H:152
virtual void writeDesignVars()
Write useful quantities to files.
virtual bool writeData(Ostream &) const
The writeData function required by the regIOobject write operation.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127