sensitivityTopO.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 "adjointSensitivity.H"
30 #include "sensitivityTopO.H"
31 #include "adjointSolver.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
43 
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 {
49  const labelList& adjointPorousIDs = zones_.adjointPorousZoneIDs();
50  if (adjointPorousIDs.empty())
51  {
52  for (label cellZoneID : zones_.fixedPorousZoneIDs())
53  {
54  const labelList& zoneCells = mesh_.cellZones()[cellZoneID];
55  for (label cellI : zoneCells)
56  {
57  sens[cellI] = 0.;
58  }
59  }
61  {
62  const labelList& zoneCells = mesh_.cellZones()[cellZoneID];
63  for (label cellI : zoneCells)
64  {
65  sens[cellI] = 0.;
66  }
67  }
68  for (label cellI : zones_.IOCells())
69  {
70  sens[cellI] = 0.;
71  }
72  }
73  else
74  {
75  // if adjointPorousZones are not empty, zero sensitivities in all cells
76  // not within them
77  scalarField mask(sens.size(), Zero);
78  for (label cellZoneID : adjointPorousIDs)
79  {
80  const labelList& zoneCells = mesh_.cellZones()[cellZoneID];
81  for (label cellI : zoneCells)
82  {
83  mask[cellI] = 1.;
84  }
85  }
86  sens *= mask;
87  }
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
93 sensitivityTopO::sensitivityTopO
94 (
95  const fvMesh& mesh,
96  const dictionary& dict,
97  adjointSolver& adjointSolver
98 )
99 :
100  adjointSensitivity(mesh, dict, adjointSolver),
101  zones_(mesh, dict.parent()),
102  designVariablesName_("beta")
103 {
104  if (includeDistance_)
105  {
106  eikonalSolver_.reset
107  (
108  new adjointEikonalSolver
109  (
110  mesh_,
111  dict_,
112  adjointSolver,
113  labelHashSet(0)
114  )
115  );
116  }
117  // Allocate sensitivities field
118  fieldSensPtr_.reset
119  (
120  createZeroFieldPtr<scalar>
121  (
122  mesh_,
123  "topologySens" + adjointSolver.solverName(),
125  )
126  );
127 
128  // Set return field size
130 }
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
136 {
138  {
139  if (includeDistance_)
140  {
141  if (eikonalSolver_)
142  {
143  eikonalSolver_->readDict(dict);
144  }
145  else
146  {
147  eikonalSolver_.reset
148  (
150  (
151  mesh_,
152  dict_,
154  labelHashSet(0)
155  )
156  );
157  }
158  }
159 
160  return true;
161  }
162 
163  return false;
164 }
165 
166 
167 void sensitivityTopO::accumulateIntegrand(const scalar dt)
168 {
169  // Accumulate source for additional post-processing PDEs, if necessary
170  if (eikonalSolver_)
171  {
172  eikonalSolver_->accumulateIntegrand(dt);
173  }
177 }
178 
179 
181 (
182  autoPtr<designVariables>& designVars
183 )
184 {
185  scalarField& sens = fieldSensPtr_().primitiveFieldRef();
186  if (eikonalSolver_)
187  {
188  eikonalSolver_->solve();
189  sens += eikonalSolver_->topologySensitivities(designVariablesName_);
190  }
192 
194 }
195 
196 
198 (
199  scalarField& sens,
200  scalarField& auxSens,
202  const word& fieldName,
203  const word& designVariablesName
204 )
205 {
206  if (fvOptions.appliesToField(fieldName))
207  {
208  DebugInfo
209  << "Computing SD contributions from the interpolation of "
210  << fieldName << endl;
211  fvOptions.postProcessSens(auxSens, fieldName, designVariablesName);
212  sens += auxSens;
213  }
214 }
215 
216 
218 (
219  scalarField& sens,
220  scalarField& auxSens,
221  const word& fieldName
222 )
223 {
224  fv::options& fvOptions(fv::options::New(this->mesh_));
225  postProcessSens(sens, auxSens, fvOptions, fieldName, designVariablesName_);
226 }
227 
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 } // End namespace Foam
232 
233 // ************************************************************************* //
const fvMesh & mesh_
Definition: sensitivity.H:68
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
adjointSolver & adjointSolver_
Reference to the underlaying adjoint solver.
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.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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
Finite-volume options.
Definition: fvOptions.H:51
const cellZone & IOCells() const
Cells next to IO boundaries.
Definition: topOZones.H:201
word designVariablesName_
Name used as the argument for the post-processing of the sensitivities through fvOptions.
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
dimensionedScalar pow5(const dimensionedScalar &ds)
conserve primitiveFieldRef()+
Macros for easy insertion into run-time selection tables.
fv::options & fvOptions
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
virtual void accumulateIntegrand(const scalar dt)
Accumulate sensitivity integrands.
topOZones zones_
Zones related to topology optimisation.
dynamicFvMesh & mesh
virtual void topOSensMultiplier(scalarField &betaMult, const word &designVariablesName, const scalar dt)
Compute the multiplier of beta.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const labelList & fixedZeroPorousZoneIDs() const
Cell zone IDs with fixed zero porosity values.
Definition: topOZones.H:185
Calculation of adjoint based sensitivities for topology optimisation. This returns just the field par...
autoPtr< volScalarField > fieldSensPtr_
Definition: sensitivity.H:74
#define DebugInfo
Report an information message using Foam::Info.
virtual void assembleSensitivities(autoPtr< designVariables > &designVars)
Assemble sensitivities.
defineTypeNameAndDebug(combustionModel, 0)
DynamicID< cellZoneMesh > cellZoneID
Foam::cellZoneID.
Definition: ZoneIDs.H:38
static void postProcessSens(scalarField &sens, scalarField &auxSens, fv::options &fvOptions, const word &fieldName, const word &designVariablesName)
Add part of the sensitivities coming from fvOptions.
const dictionary & dict() const
Return the construction dictionary.
Definition: sensitivity.H:129
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
bool includeDistance_
Include distance variation in sensitivity computations.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
label nCells() const noexcept
Number of mesh cells.
scalarField derivatives_
The sensitivity derivative values.
const labelList & adjointPorousZoneIDs() const
Cell zone IDs in which porosity is allowed to change.
Definition: topOZones.H:193
void zeroSensInFixedPorousZones(scalarField &sens)
Zero sensitivities in fixed porous zones.
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
autoPtr< adjointEikonalSolver > eikonalSolver_
Adjoint eikonal equation solver.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
dictionary dict_
Definition: sensitivity.H:69
List< label > labelList
A List of labels.
Definition: List.H:62
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
Solver of the adjoint to the eikonal PDE.
virtual void assembleSensitivities(autoPtr< designVariables > &designVars)
Assemble sensitivities.
const labelList & fixedPorousZoneIDs() const
Cell zone IDs with fixed porosity values.
Definition: topOZones.H:169
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127