Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
30#include "wallDist.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36namespace Foam
37{
38
39// * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * * //
40
43(
44 designVariables,
45 levelSetDesignVariables,
46 designVariables
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 (
67 (
68 "yLevelSet",
70 mesh_,
73 ),
74 mesh_,
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{
96 scalarField& betaIf = beta_.primitiveFieldRef();
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;
128 }
129
131}
132
133
135{
136 activeDesignVariables_.setSize(mesh_.nCells(), -1);
137 if (zones_.adjointPorousZoneIDs().empty())
138 {
139 boolList isActiveVar(mesh_.nCells(), true);
140
141 const labelList& fixedZeroPorousZones =
142 zones_.fixedZeroPorousZoneIDs();
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 }
178 activeDesignVariables_.setSize(iVar);
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++;
194 }
196 }
197}
198
199
201{
202 // Compute the beta field by passing the distance field through
203 // a Heaviside function
204 scalarField& beta = beta_.primitiveFieldRef();
205 interpolation_->interpolate(signedDistances_.primitiveField(), beta);
206 beta = 1 - beta;
207 // Apply fixed values if necessary
209
210 beta_.correctBoundaryConditions();
211}
212
213
215{
216 Info<< "Re-initilising the level-set distance field" << nl << endl;
217
219 (
221 (
222 "yLevelSet",
223 mesh_.time().timeName(),
224 mesh_,
227 ),
228 mesh_,
231 );
232 y.primitiveFieldRef() = aTilda_.primitiveFieldRef();
233 y.correctBoundaryConditions();
234
235 changedFaces_.clear();
236 changedFaces_.setSize(mesh_.nFaces(), -1);
237
238 changedFacesInfo_.clear();
239 changedFacesInfo_.setSize(mesh_.nFaces());
240
241 writeFluidSolidInterface(aTilda_, 0, changedFaces_, changedFacesInfo_);
242
243 List<wallPointData<label>> allFaceInfo(mesh_.nFaces());
244 allCellInfo_.clear();
245 allCellInfo_.setSize(mesh_.nCells());
246
248 (
249 mesh_,
250 changedFaces_,
251 changedFacesInfo_,
253 allCellInfo_,
254 mesh_.globalData().nTotalCells() + 1
255 );
256
257 // Transfer the distance from cellInfo to the alphaTilda field
258 forAll(allCellInfo_, celli)
259 {
260 if (allCellInfo_[celli].valid(wave.data()))
261 {
262 signedDistances_[celli] =
263 sign(aTilda_[celli])*Foam::sqrt(allCellInfo_[celli].distSqr());
264 }
266 signedDistances_.correctBoundaryConditions();
267}
268
269
270// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
271
272levelSetDesignVariables::levelSetDesignVariables
273(
274 fvMesh& mesh,
275 const dictionary& dict
276)
277:
278 topOVariablesBase(mesh, dict),
279 designVariables(mesh, dict, mesh.nCells()),
280 radius_
281 (regularisationRadius::New(mesh, dict.subDict("regularisation"), false)),
282 regularisation_
283 (regularisationPDE::New(mesh, dict.subDict("regularisation"), zones_)),
284 aTilda_
285 (
286 IOobject
287 (
288 "aTilda",
289 mesh_.time().timeName(),
290 mesh_,
291 IOobject::READ_IF_PRESENT,
292 IOobject::AUTO_WRITE
293 ),
294 mesh_,
296 fvPatchFieldBase::zeroGradientType()
297 ),
298 signedDistances_
299 (
300 IOobject
301 (
302 "signedDistances",
303 mesh_.time().timeName(),
304 mesh_,
305 IOobject::READ_IF_PRESENT,
306 IOobject::AUTO_WRITE
307 ),
308 mesh_,
310 fvPatchFieldBase::zeroGradientType()
311 ),
312 interpolation_
313 (
314 topOInterpolationFunction::New(mesh_, dict_.subDict("interpolation"))
315 ),
316 beta_
317 (
318 IOobject
319 (
320 "beta",
321 mesh_.time().timeName(),
322 mesh_,
323 IOobject::READ_IF_PRESENT,
324 IOobject::AUTO_WRITE
325 ),
326 mesh_,
328 ),
329 fixATildaValues_(dict.getOrDefault<bool>("fixATildaValues", true)),
330 writeAllDistanceFields_
331 (dict.getOrDefault<bool>("writeAllDistanceFields", false)),
332 changedFaces_(),
333 changedFacesInfo_(),
334 allCellInfo_()
335{
336 // Read the alpha field if present, or set it based on the distance field
337 readField();
338
339 // Read bounds of design variables if present.
340 // If not, use the maximum distance field in the fluid to set them.
341 scalar maxDist = gMax(*this);
342 scalar lowerBound =
343 localIOdictionary::getOrDefault<scalar>("lowerBound", -maxDist - SMALL);
344 scalar upperBound =
345 localIOdictionary::getOrDefault<scalar>("upperBound", maxDist + SMALL);
347 (
348 autoPtr<scalar>(new scalar(lowerBound)),
349 autoPtr<scalar>(new scalar(upperBound))
350 );
352 << "Using lower/upper bounds "
353 << lowerBounds_()[0] << "/" << upperBounds_()[0]
354 << endl;
355
356 // Update the beta field based on the initial design vars
357 scalarField zeroUpdate(scalarField::size(), Zero);
358 update(zeroUpdate);
359
360 // Determine which design variables are active
362}
363
364
365// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
366
368(
369 fvMesh& mesh,
370 const dictionary& dict
371)
372{
374 (
375 new levelSetDesignVariables(mesh, dict)
376 );
377}
378
379
380// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
383{
384 return beta_;
385}
386
387
389{
391
392 // Compute the regularised design variables
393 regularisation_->regularise
394 (
397 );
399
401 {
403 aTilda_.write();
404 }
405
406 // Compute signed distances based on aTilda
408
409 // Set beta based on the signed distances
410 updateBeta();
411
413 {
415 beta_.write();
416 }
417
418 // Though the mesh is kept constant, the distance from wall may change
419 // due to fvOptions depending on beta. Trick wallDist into updating it
420
422}
423
424
426{
427 const scalar maxChange(gMax(mag(correction)));
428 Info<< "maxInitChange/maxChange \t"
429 << maxInitChange_() << "/" << maxChange << endl;
430 const scalar eta(maxInitChange_() / maxChange);
431 Info<< "Setting eta value to " << eta << endl;
432 correction *= eta;
433
434 return eta;
435}
436
439{
440 return true;
441}
442
443
445(
446 adjointSensitivity& adjointSens
447)
448{
449 // Raw sensitivities field
450 const scalarField& fieldSens = adjointSens.fieldSensPtr()->primitiveField();
451
452 // Return field
453 auto tobjectiveSens(tmp<scalarField>::New(fieldSens));
454 scalarField& objectiveSens = tobjectiveSens.ref();
455
456 // Multiply with dBetadAtilda
457 objectiveSens *=
459
460 // Solve the adjoint to the regularisation equation
462 regularise(aTilda_, objectiveSens, objectiveSens, false, radius_());
463
464 objectiveSens *= mesh_.V();
465
467 {
468 volScalarField sens
469 (
471 (
472 "sens" + adjointSens.getAdjointSolver().solverName(),
473 mesh_.time().timeName(),
474 mesh_,
477 ),
478 mesh_,
480 );
481 sens.primitiveFieldRef() = objectiveSens;
482 sens.write();
483 }
484
485 return tobjectiveSens;
486}
487
488
490{
492 {
494 (
495 IOobject
496 (
497 "alpha",
498 mesh_.time().timeName(),
499 mesh_,
502 ),
503 mesh_,
505 );
506 alpha.primitiveFieldRef() = *this;
507 alpha.correctBoundaryConditions();
508
509 alpha.write();
510 }
511}
512
513
514bool levelSetDesignVariables::writeData(Ostream& os) const
515{
516 // Lower and upper bound values might be computed based on the initial
517 // distance field. Write to file to enable continuation
518 os.writeEntry("lowerBound", lowerBounds_()[0]);
519 os.writeEntry("upperBound", upperBounds_()[0]);
520
521 scalarField::writeEntry("alpha", os);
522
523 return true;
524}
525
526
527// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
528
529} // End namespace Foam
530
531// ************************************************************************* //
scalar y
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Wave propagation of information through grid. Every iteration information goes through one layer of c...
const TrackingData & data() const noexcept
Additional data to be passed into container.
void operator+=(const UList< scalar > &)
void writeEntry(const word &keyword, Ostream &os) const
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.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
bool writeTime() const noexcept
True if this is a write interval.
Definition TimeStateI.H:73
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
Definition UList.H:706
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Abstract base class for adjoint-based sensitivities.
const adjointSolver & getAdjointSolver() const
Const access to adjoint solver.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Abstract base class for defining design variables.
autoPtr< scalarField > upperBounds_
Upper bounds of the design variables.
autoPtr< scalar > maxInitChange_
Maximum design variables' change in the first optimisation cycle.
void readBounds(autoPtr< scalar > lowerBoundPtr=nullptr, autoPtr< scalar > upperBoundPtr=nullptr)
Read bounds for design variables, if present.
autoPtr< scalarField > lowerBounds_
Lower bounds of the design variables.
labelList activeDesignVariables_
Which of the design variables will be updated.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
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...
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
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
Template invariant parts for fvPatchField.
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Signed distance field design variables for level-set based topology optimization (topO).
virtual void update(scalarField &correction)
Update design variables based on a given correction.
void applyFixedPorosityValues()
Apply fixed values in certain zones.
autoPtr< regularisationRadius > radius_
The regularisation radius.
virtual bool writeData(Ostream &) const
The writeData function required by the regIOobject write operation.
void updateBeta()
Update the beta field.
virtual tmp< scalarField > assembleSensitivities(adjointSensitivity &adjointSens)
Add contributions from the adjoint to the regularization PDE, the derivative of the interpolation fun...
List< wallPointData< label > > changedFacesInfo_
Seed distances to MeshWave and cell distances.
bool fixATildaValues_
Fix aTilda values in fixed{Zero}PorousZones and IOcells.
void setActiveDesignVariables(bool activeIO=false)
Determine which design variables are active.
volScalarField beta_
The indicator field.
List< wallPointData< label > > allCellInfo_
autoPtr< topOInterpolationFunction > interpolation_
Function to transorm signed distances to the indicator field beta_.
void readField()
Read the design variables field.
virtual bool globalSum() const
Whether to use global sum when computing matrix-vector products in update methods.
void updateSignedDistances()
Make aTilda a signed distance field based on the zero iso-surface of the current aTilda field.
autoPtr< regularisationPDE > regularisation_
Regularisation mechanism.
virtual const volScalarField & beta() const
Const reference to the beta field.
labelList changedFaces_
Mesh faces acting as the source of MeshWave.
static autoPtr< levelSetDesignVariables > New(fvMesh &mesh, const dictionary &dict)
Return a reference to the selected turbulence model.
bool writeAllDistanceFields_
Write all fields related to the distance calculation (debugging).
virtual scalar computeEta(scalarField &correction)
Compute eta if not set in the first step.
volScalarField aTilda_
The regularised field.
volScalarField signedDistances_
The signed distances field.
virtual void writeDesignVars()
Write useful quantities to files.
static autoPtr< patchDistMethod > New(const dictionary &dict, const fvMesh &mesh, const labelHashSet &patchIDs, const word &defaultPatchDistMethod=word::null)
static wordList patchTypes(const fvMesh &mesh, const labelHashSet &patchIDs)
Return the patch types for y and n.
labelHashSet findPatchIDs() const
Find patch indices for a given polyPatch type (uses isA test).
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Base class for selecting the regulatisation PDE.
Base class for selecting the regulatisation radius.
const autoPtr< volScalarField > & fieldSensPtr() const
Get the fieldSensPtr.
const word & solverName() const
Return the solver name.
Definition solverI.H:30
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
void writeFluidSolidInterface(const volScalarField &indicator, const scalar isoValue, labelList &changedFaces, List< wallPointData< label > > &changedFacesInfo)
Write the fluid-solid interface to files.
topOZones zones_
Cell zones useful for defining the constant and changing parts of the domain in topO.
static bool try_movePoints(const fvMesh &mesh)
Trigger update of y-field for the "wallDist" MeshObject on the given mesh. A no-op if the wallDist is...
Definition wallDist.C:152
Foam::wallPolyPatch.
Base class for mesh zones.
Definition zone.H:63
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
Abstract base class for interpolation functions used in topology optimisation.
mesh update()
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
List< wallPoints > allFaceInfo(mesh_.nFaces())
word timeName
Definition getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
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.
const dimensionSet dimless
Dimensionless.
dimensionedScalar sign(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
DynamicID< cellZoneMesh > cellZoneID
Foam::cellZoneID.
Definition ZoneIDs.H:38
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299