shapeOptimisationIncompressible.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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 namespace incompressible
38 {
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 defineTypeNameAndDebug(shapeOptimisation, 0);
44 (
45  optimisationType,
46  shapeOptimisation,
47  dictionary
48 );
49 
50 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
51 
53 (
55 )
56 {
57  if (!updateMethod_->initialEtaSet())
58  {
59  // In the unlikely event that eta is not set and the line search step
60  // is not 1, multiply with it
61  // if (lineSearch_.valid()) correction *= lineSearch_->step();
62 
63  // Compute eta based on desirable mesh movement size
64  scalar eta = optMeshMovement_->computeEta(correction);
65  correction *= eta;
66 
67  // Update eta known by the optimisation method and inform it that is
68  // has been set
69  updateMethod_->setStep(eta);
70  updateMethod_->initialEtaSet() = true;
71 
72  // If a backtracking should be made at the first optimisation cycle,
73  // the direction of the subsequent line searches of the same cycle
74  // should also be scaled with the newly computed eta. We do this by
75  // changing the line search step. This will happen only at the first
76  // optimisation cycle since the updated value of eta will be included
77  // in the line search direction in all subsequent optimisation cycles
78  //correction *= eta;
79  }
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
85 shapeOptimisation::shapeOptimisation
86 (
87  fvMesh& mesh,
88  const dictionary& dict,
90 )
91 :
93  optMeshMovement_(nullptr),
94  writeEachMesh_
95  (
96  dict.subDict("optimisationType").
97  lookupOrDefault<bool>("writeEachMesh", false)
98  ),
99  updateGeometry_
100  (
101  dict.subDict("optimisationType").
102  lookupOrDefault<bool>("updateGeometry", true)
103  )
104 {
105  // Note: to be updated
107  (
108  mesh_.boundaryMesh().patchSet
109  (
110  dict_.subDict("sensitivities").get<wordRes>("patches")
111  )
112  );
113  if (patches.empty())
114  {
116  << "There is no patch on which to compute sensitivities. "
117  << "Check optimisationDict \n"
118  << endl;
119  }
120  labelList sensitivityPatchIDs = patches.toc();
121  optMeshMovement_.reset
122  (
124  (
125  mesh_,
126  dict_.subDict("meshMovement"),
127  sensitivityPatchIDs
128  ).ptr()
129  );
130 
131  // Sanity checks: at least one of eta or maxAllowedDisplacement must be set
132  if
133  (
134  !updateMethod_->initialEtaSet()
135  && !optMeshMovement_().maxAllowedDisplacementSet()
136  )
137  {
139  << "Neither eta (updateMethod) "
140  << "nor maxAllowedDisplacement (meshMovement) have been set"
141  << nl
142  << exit(FatalError);
143  }
144 }
145 
146 
147 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
148 
150 {
151  Info<< nl << "Moving Mesh..." << nl << endl;
152 
153  // Sum contributions
154  scalarField objectiveSens(0);
155  PtrList<scalarField> constraintSens(0);
156  scalar objectiveValue(Zero);
157  scalarField constraintValues(0);
159  {
160  adjointSolverManager& adjSolvManager(adjointSolvManagers_[amI]);
161  const scalar opWeight(adjSolvManager.operatingPointWeight());
162 
163  // Allocate objective sens size if necessary
164  tmp<scalarField> tadjointSolverManagerSens =
165  adjSolvManager.aggregateSensitivities();
166 
167  if (objectiveSens.empty())
168  {
169  objectiveSens.setSize(tadjointSolverManagerSens().size(), Zero);
170  }
171  objectiveSens += opWeight*tadjointSolverManagerSens();
172  objectiveValue += opWeight*adjSolvManager.objectiveValue();
173 
174  // Allocate constraint sens size if necessary
175  PtrList<scalarField> adjointSolverManagerConstSens
176  = adjSolvManager.constraintSensitivities();
177  tmp<scalarField> cValues = adjSolvManager.constraintValues();
178  if (constraintSens.empty())
179  {
180  constraintSens.setSize(adjointSolverManagerConstSens.size());
181  forAll(constraintSens, cI)
182  {
183  constraintSens.set
184  (
185  cI,
186  new scalarField
187  (
188  adjointSolverManagerConstSens[cI].size(),
189  Zero
190  )
191  );
192  constraintValues.setSize(cValues().size());
193  constraintValues = Zero;
194  }
195  }
196 
197  forAll(constraintSens, cI)
198  {
199  constraintSens[cI] += opWeight*adjointSolverManagerConstSens[cI];
200  }
201  constraintValues += opWeight*cValues();
202  }
203 
204  // Based on the sensitivities, return design variables correction
205  updateMethod_->setObjectiveDeriv(objectiveSens);
206  updateMethod_->setConstraintDeriv(constraintSens);
207  updateMethod_->setObjectiveValue(objectiveValue);
208  updateMethod_->setConstraintValues(constraintValues);
209  scalarField& correction = updateMethod_->returnCorrection();
210 
211  // Computed eta if needed
213  updateMethod_->writeCorrection();
214 
215  // Communicate the movement to optMeshMovement
216  optMeshMovement_->setCorrection(correction);
217  if (updateGeometry_)
218  {
219  optMeshMovement_->moveMesh();
220 
221  if (writeEachMesh_)
222  {
223  Info<< " Writing new mesh points " << endl;
225  (
226  IOobject
227  (
228  "points",
231  mesh_,
234  false
235  ),
236  mesh_.points()
237  );
238  points.write();
239  }
240  }
241 }
242 
243 
245 {
246  // Computed eta if needed
248 
249  // Multiply with line search step, if necessary
251  if (lineSearch_.valid())
252  {
253  correction *= lineSearch_->step();
254  }
255 
256  // Communicate the movement to optMeshMovement
257  optMeshMovement_->setCorrection(correction);
258 
259  if (updateGeometry_)
260  {
261  // Update the mesh
262  optMeshMovement_->moveMesh();
263 
264  if (writeEachMesh_)
265  {
266  Info<< " Writing new mesh points " << endl;
268  (
269  IOobject
270  (
271  "points",
274  mesh_,
277  false
278  ),
279  mesh_.points()
280  );
281  points.write();
282  }
283  }
284 }
285 
286 
288 {
289  optMeshMovement_->storeDesignVariables();
290 }
291 
292 
294 {
295  optMeshMovement_->resetDesignVariables();
296 }
297 
298 
300 {
302  updateMethod_->writeCorrection();
303 }
304 
305 
306 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
307 
308 } // End namespace incompressible
309 } // End namespace Foam
310 
311 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
Foam::incompressible::shapeOptimisation::storeDesignVariables
virtual void storeDesignVariables()
Store design variables, as the starting point for line search.
Definition: shapeOptimisationIncompressible.C:287
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::adjointSolverManager::operatingPointWeight
scalar operatingPointWeight() const
Const access to adjoint solvers.
Definition: adjointSolverManager.C:172
Foam::adjointSolverManager::constraintSensitivities
virtual PtrList< scalarField > constraintSensitivities()
Get constraint sensitivities. One scalarField per constraint.
Definition: adjointSolverManager.C:232
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::incompressible::defineTypeNameAndDebug
defineTypeNameAndDebug(adjointEikonalSolver, 0)
Foam::incompressible::optimisationType::write
virtual void write()
Write useful quantities to files.
Definition: optimisationTypeIncompressible.C:271
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::optMeshMovement::New
static autoPtr< optMeshMovement > New(fvMesh &mesh, const dictionary &dict, const labelList &patchIDs)
Definition: optMeshMovement.C:93
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:315
Foam::incompressible::optimisationType
Abstract base class for optimisation methods.
Definition: optimisationTypeIncompressible.H:58
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::adjointSolverManager::aggregateSensitivities
virtual tmp< scalarField > aggregateSensitivities()
Aggregate sensitivities from various adjoint solvers.
Definition: adjointSolverManager.C:208
Foam::correction
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
Foam::HashSet< label, Hash< label > >
Foam::incompressible::shapeOptimisation::optMeshMovement_
autoPtr< optMeshMovement > optMeshMovement_
Definition: shapeOptimisationIncompressible.H:70
Foam::incompressible::shapeOptimisation::updateGeometry_
bool updateGeometry_
Definition: shapeOptimisationIncompressible.H:73
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:815
Foam::incompressible::shapeOptimisation::write
virtual void write()
Write useful quantities to files.
Definition: shapeOptimisationIncompressible.C:299
Foam::incompressible::addToRunTimeSelectionTable
addToRunTimeSelectionTable(adjointSensitivity, sensitivityBezier, dictionary)
Foam::adjointSolverManager::constraintValues
virtual tmp< scalarField > constraintValues()
Get constraint values.
Definition: adjointSolverManager.C:281
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::PtrList::setSize
void setSize(const label newLen)
Same as resize()
Definition: PtrListI.H:108
adjointSolverManagers
PtrList< adjointSolverManager > & adjointSolverManagers
Definition: createFields.H:8
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
Foam::incompressible::optimisationType::mesh_
fvMesh & mesh_
Definition: optimisationTypeIncompressible.H:64
Foam::incompressible::shapeOptimisation::resetDesignVariables
virtual void resetDesignVariables()
Store design variables, as the starting point for line search.
Definition: shapeOptimisationIncompressible.C:293
Foam::adjointSolverManager
Class for managing adjoint solvers, which may be more than one per operating point.
Definition: adjointSolverManager.H:54
Foam::incompressible::optimisationType::adjointSolvManagers_
PtrList< adjointSolverManager > & adjointSolvManagers_
Definition: optimisationTypeIncompressible.H:66
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::PtrList::set
const T * set(const label i) const
Return const pointer to element (if set) or nullptr.
Definition: PtrListI.H:143
Foam::incompressible::shapeOptimisation::writeEachMesh_
bool writeEachMesh_
Definition: shapeOptimisationIncompressible.H:72
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
shapeOptimisationIncompressible.H
Shape optimisation support library.
Foam::adjointSolverManager::objectiveValue
scalar objectiveValue()
Get objective value.
Definition: adjointSolverManager.C:267
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::incompressible::optimisationType::updateMethod_
autoPtr< updateMethod > updateMethod_
Definition: optimisationTypeIncompressible.H:67
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::incompressible::shapeOptimisation::update
void update()
Master function. Calls all the others.
Definition: shapeOptimisationIncompressible.C:149
Foam::List< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::direction
uint8_t direction
Definition: direction.H:47
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::incompressible::shapeOptimisation::computeEta
virtual void computeEta(scalarField &correction)
Definition: shapeOptimisationIncompressible.C:53
Foam::IOobject::NO_READ
Definition: IOobject.H:123
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::incompressible::optimisationType::lineSearch_
autoPtr< lineSearch > lineSearch_
Definition: optimisationTypeIncompressible.H:69