Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
29#include "localIOdictionary.H"
30#include "topODesignVariables.H"
31#include "wallDist.H"
32#include "wallFvPatch.H"
33#include "cutFaceIso.H"
34#include "cutCellIso.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
43 (
46 designVariables
47 );
48}
49
50
51// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
52
54(
56 const label fluidID
57)
58{
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)
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;
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 else
228 {
229 IOobject header
230 (
231 name,
232 mesh_.time().timeName(),
233 mesh_,
236 );
237 if (header.typeHeaderOk<volScalarField>())
238 {
239 Info<< "Setting design variables based on the alpha field "
240 << nl << endl;
241 volScalarField volField
242 (
243 header,
244 mesh_
245 );
246 const scalarField& field = volField.primitiveField();
247 forAll(field, cI)
248 {
249 scalarField::operator[](offset + cI) = field[cI];
250 }
251 }
252 }
253}
254
255
257{
258 // Set active design variables
259 setActiveDesignVariables();
260
261 // Read in values from file, if present
262 readField("alpha", 0, true);
263
264 if (regularisation_.growFromWalls())
265 {
266 scalarField& alpha = *this;
267 for (const fvPatch& patch : mesh_.boundary())
268 {
270 {
271 UIndirectList<scalar>(alpha, patch.faceCells()) = 1;
272 }
273 }
274 }
275
276 // Make sure alpha has fixed values where it should
277 scalarField dummyCorrection(mesh_.nCells(), Zero);
278 update(dummyCorrection);
279
280 // Read bounds for design variables, if present
282}
283
284
285// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
286
287Foam::topODesignVariables::topODesignVariables
288(
289 fvMesh& mesh,
291)
292:
293 topODesignVariables(mesh, dict, mesh.nCells())
294{}
295
296
297Foam::topODesignVariables::topODesignVariables
298(
299 fvMesh& mesh,
300 const dictionary& dict,
301 const label size
302)
303:
305 designVariables(mesh, dict, size),
306 alpha_(SubField<scalar>(*this, mesh.nCells(), 0)),
307 regularisation_
308 (
309 mesh,
310 alpha_,
311 zones_,
312 dict_.subDict("regularisation")
313 ),
314 writeAllFields_
315 (
316 dict.getOrDefaultCompat<bool>
317 (
318 "writeAllFields", {{"writeAllAlphaFields", 2306}}, false
319 )
320 ),
321 addFvOptions_(dict.getOrDefault<bool>("addFvOptions", false))
322{}
323
324
325// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
326
328(
329 fvMesh& mesh,
330 const dictionary& dict
331)
336
337// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
340{
341 return regularisation_.beta();
342}
343
344
346(
347 const word& interpolationFieldName
348) const
349{
350 return beta().primitiveField();
351}
352
353
355(
357 const topOInterpolationFunction& interpolationFunc,
358 const FieldField<Field, scalar>& fluidValues,
359 const scalarField& solidValues,
360 const label fieldi,
361 const word& interpolationFieldName
362) const
363{
364 const scalarField& indicator = interpolationField(interpolationFieldName);
365 scalarField interpolant(indicator.size(), Zero);
366 interpolationFunc.interpolate(indicator, interpolant);
367
368 // Interpolate field values
369 const scalar diff(solidValues[fieldi] - fluidValues[0][fieldi]);
370 field.primitiveFieldRef() = fluidValues[0][fieldi] + interpolant*diff;
371 field.correctBoundaryConditions();
372}
373
374
376(
377 scalarField& sens,
378 const topOInterpolationFunction& interpolationFunc,
379 const FieldField<Field, scalar>& fluidValues,
380 const scalarField& solidValues,
381 const label fieldi,
382 const word& designVariablesName,
383 const word& interpolationFieldName
384) const
385{
386 const scalarField& indicator = interpolationField(interpolationFieldName);
387 sens *=
388 (solidValues[fieldi] - fluidValues[0][fieldi])
389 *interpolationFunc.derivative(indicator);
390}
391
392
394(
396 const topOInterpolationFunction& interpolationFunc
397) const
398{
399 const scalarField& beta = this->beta().primitiveField();
400 scalarField interpolant(beta.size(), Zero);
401 interpolationFunc.interpolate(beta, interpolant);
402 field *= scalar(1) - interpolant;
403}
404
405
407(
408 scalarField& sens,
409 const topOInterpolationFunction& interpolationFunc,
410 const word& designVariablesName
411) const
412{
413 const scalarField& beta = this->beta().primitiveField();
414 sens *= - interpolationFunc.derivative(beta);
415}
416
417
419(
420 const word& interpolationFieldName,
421 const topOInterpolationFunction& interpolationFunc
422) const
423{
424 const scalarField& beta = this->beta().primitiveField();
425 tmp<scalarField> tinterpolant(tmp<scalarField>::New(beta.size(), Zero));
426 interpolationFunc.interpolate(beta, tinterpolant.ref());
427 return tinterpolant;
428}
429
430
432(
433 const word& interpolationFieldName,
434 const topOInterpolationFunction& interpolationFunc
435) const
436{
437 const scalarField& beta = this->beta().primitiveField();
438 return interpolationFunc.derivative(beta);
439}
440
441
443{
444 // Update alpha values
445 updateField(correction);
446
447 // Fix alpha in zones
448 applyFixedValues();
449
450 // Update the beta field
451 regularisation_.updateBeta();
452
453 // Though the mesh is kept constant, the distance from wall may change
454 // if the method computing it includes fvOptions that depend on the
455 // indicator field.
456 // Trick wallDist into updating it
457
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
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 (
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 (
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
572bool 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// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
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.
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ 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
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info. A void type suppresses trait and t...
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
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition SubField.H:58
A List with indirect addressing. Like IndirectList but does not store addressing.
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Definition UList.H:118
scalar & operator[](const label i)
Abstract base class for adjoint-based sensitivities.
const adjointSolver & getAdjointSolver() const
Const access to adjoint solver.
Class for managing adjoint solvers, which may be more than one per operating point.
Base class for adjoint solvers.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
Abstract base class for defining 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.
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 getOrDefaultCompat(const word &keyword, std::initializer_list< std::pair< const char *, int > > compat, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value using any compatibility names if needed.
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
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Base class for primal solvers.
label nCells() const noexcept
Number of mesh cells.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
const autoPtr< volScalarField > & fieldSensPtr() const
Get the fieldSensPtr.
Base solver class.
Definition solver.H:48
virtual void addTopOFvOptions() const
Add topO fvOptions.
Definition solver.C:102
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
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Design variables for porosity-based topology optimisation (topO) problems.
virtual void updateField(const scalarField &correction, const label fluidID=0)
Update the design variables given their correction.
virtual void update(scalarField &correction)
Update design variables based on a given correction.
virtual const scalarField & interpolationField(const word &interpolationFieldName="beta") const
Return interpolant.
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...
fieldRegularisation regularisation_
Mechanism to regularise the field of design variables.
virtual scalar computeEta(scalarField &correction)
Compute eta if not set in the first step.
virtual bool writeData(Ostream &) const
The writeData function required by the regIOobject write operation.
virtual void initialize()
Part of the constructor initialisation.
virtual tmp< scalarField > assembleSensitivities(adjointSensitivity &sens)
Assemble sensitivity derivatives, by combining the part related to the primal and adjoint solution wi...
virtual tmp< scalarField > penaltySensitivities(const word &interpolationFieldName, const topOInterpolationFunction &interpolationFunc) const
Return the penalty term derivative.
bool addFvOptions_
Add the fvOptions necessary for topO automatically.
virtual void nullifyInSolidSensitivities(scalarField &sens, const topOInterpolationFunction &interpolationFunc, const word &designVariablesName) const
Nullify given field in the solid domain.
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.
SubField< scalar > alpha_
A subfield of the design variables corresponding to the porosity field.
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.
bool writeAllFields_
Write all fields related to imposition of the Brinkman penalisation (i.e. design variables,...
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.
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.
virtual const volScalarField & beta() const
Get the indicator field.
virtual void setInitialValues()
Set initial values of the design variables.
virtual bool globalSum() const
Whether to use global sum when computing matrix-vector products in update methods.
void applyFixedValues()
Apply fixed values in certain zones.
virtual tmp< scalarField > penalty(const word &interpolationFieldName, const topOInterpolationFunction &interpolationFunc) const
Return the Brinkman penalisation term.
virtual void nullifyInSolid(scalarField &field, const topOInterpolationFunction &interpolationFunc) const
Nullify given field in the solid domain.
virtual void writeDesignVars()
Write porosity field to file.
virtual tmp< scalarField > derivative(const scalarField &arg) const =0
Return of function with respect to the argument field.
virtual void interpolate(const scalarField &arg, scalarField &res) const =0
Interpolate argument to result.
Base class for all design variables related to topology optimisation (topO). Provides the lookup func...
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.
const labelList & adjointPorousZoneIDs() const
Cell zone IDs in which porosity is allowed to change.
Definition topOZones.H:193
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
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
rDeltaTY field()
mesh update()
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define DebugInfo
Report an information message using Foam::Info.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
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
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition triad.C:373
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
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...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dictionary dict
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
PtrList< adjointSolverManager > & adjointSolverManagers
Definition createFields.H:8