objective.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  Copyright (C) 2019-2021 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 
30 #include "objective.H"
31 #include "createZeroField.H"
32 #include "IOmanip.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(objective, 0);
42 defineRunTimeSelectionTable(objective, objective);
43 
44 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
45 
46 void objective::makeFolder()
47 {
48  if (Pstream::master())
49  {
50  const Time& time = mesh_.time();
52  time.globalPath()/"optimisation"/type()/time.timeName();
54  {
56  }
57 
59  }
60 }
61 
62 
64 {
66  (
68  );
69 }
70 
71 
73 {
75  (
76  new OFstream
77  (
79  )
80  );
81 }
82 
83 
85 {
86  meanValueFilePtr_.reset
87  (
88  new OFstream
89  (
91  )
92  );
93 }
94 
95 
96 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
97 
98 const dictionary& objective::dict() const
99 {
100  return dict_;
101 }
102 
103 
104 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
105 
106 objective::objective
107 (
108  const fvMesh& mesh,
109  const dictionary& dict,
110  const word& adjointSolverName,
111  const word& primalSolverName
112 )
113 :
115  (
116  IOobject
117  (
118  adjointSolverName + "_" + dict.dictName(),
119  mesh.time().timeName(),
120  fileName("uniform")/fileName("objectives")/adjointSolverName,
121  mesh,
122  IOobject::READ_IF_PRESENT,
123  IOobject::AUTO_WRITE
124  ),
125  // avoid type checking since dictionary is read using the
126  // derived type name and type() will result in "objective"
127  // here
128  word::null
129  ),
130  mesh_(mesh),
131  dict_(dict),
132  adjointSolverName_(adjointSolverName),
133  primalSolverName_(primalSolverName),
134  objectiveName_(dict.dictName()),
135  computeMeanFields_(false), // is reset in derived classes
136  nullified_(false),
137  normalize_(dict.getOrDefault<bool>("normalize", false)),
138  shouldWrite_(true),
139 
140  J_(Zero),
141  JMean_(this->getOrDefault<scalar>("JMean", Zero)),
142  weight_(dict.get<scalar>("weight")),
143  computed_(false),
144  normFactor_(nullptr),
145  target_
146  (
147  dict.found("target") ?
148  autoPtr<scalar>::New(dict.get<scalar>("target")) :
149  nullptr
150  ),
151  targetLeft_
152  (
153  dict.found("targetLeft") ?
154  autoPtr<scalar>::New(dict.get<scalar>("targetLeft")) :
155  nullptr
156  ),
157  integrationStartTimePtr_(nullptr),
158  integrationEndTimePtr_(nullptr),
159  fieldNames_(),
160 
161  // Initialize pointers to nullptr.
162  // Not all of them are required for each objective function.
163  // Each child should allocate whatever is needed.
164 
165  dJdbPtr_(nullptr),
166  dJdbFieldPtr_(nullptr),
167  bdJdbPtr_(nullptr),
168  bdSdbMultPtr_(nullptr),
169  bdndbMultPtr_(nullptr),
170  bdxdbMultPtr_(nullptr),
171  bdxdbDirectMultPtr_(nullptr),
172  bEdgeContribution_(nullptr),
173  divDxDbMultPtr_(nullptr),
174  gradDxDbMultPtr_(nullptr),
175 
176  objFunctionFolder_("word"),
177  objFunctionFilePtr_(nullptr),
178  instantValueFilePtr_(nullptr),
179  meanValueFilePtr_(nullptr),
180  width_(IOstream::defaultPrecision() + 5)
181 {
182  makeFolder();
183  // Read integration start and end times, if present.
184  // For unsteady runs only
185  if (dict.found("integrationStartTime"))
186  {
188  (
189  new scalar(dict.get<scalar>("integrationStartTime"))
190  );
191  }
192  if (dict.found("integrationEndTime"))
193  {
195  (
196  new scalar(dict.get<scalar>("integrationEndTime"))
197  );
198  }
199 
200  // Set normalization factor, if present
201  if (normalize_)
202  {
203  scalar normFactor(Zero);
204  if (dict.readIfPresent("normFactor", normFactor))
205  {
206  normFactor_.reset(new scalar(normFactor));
207  }
208  else if (this->readIfPresent("normFactor", normFactor))
209  {
210  normFactor_.reset(new scalar(normFactor));
211  }
212  }
213 
214  // Set the weight factor in case of continuation
215  this->readIfPresent("weight", weight_);
216 }
217 
218 
219 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
220 
222 (
223  const fvMesh& mesh,
224  const dictionary& dict,
225  const word& objectiveType,
226  const word& adjointSolverName,
227  const word& primalSolverName
228 )
229 {
230  auto* ctorPtr = objectiveConstructorTable(objectiveType);
231 
232  if (!ctorPtr)
233  {
235  (
236  dict,
237  "objective",
238  objectiveType,
239  *objectiveConstructorTablePtr_
240  ) << exit(FatalIOError);
241  }
242 
243  return autoPtr<objective>
244  (
245  ctorPtr
246  (
247  mesh,
248  dict,
249  adjointSolverName,
250  primalSolverName
251  )
252  );
253 }
254 
255 
256 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
257 
259 {
260  dict_ = dict;
261  return true;
262 }
263 
264 
265 scalar objective::JCycle(bool negate) const
266 {
267  scalar J(J_);
268  if
269  (
272  )
273  {
274  J = JMean_;
275  }
276 
277  // Subtract target, in case the objective is used as a constraint
278  if (target_.valid())
279  {
280  if (negate)
281  {
282  J = - J + targetLeft_();
283  }
284  else
285  {
286  J -= target_();
287  }
288  }
289 
290  // Normalize here, in order to get the correct value for line search
291  if (normalize_ && normFactor_)
292  {
293  J /= normFactor_();
294  }
295  J *= weight_;
296 
297  return J;
298 }
299 
300 
302 {
303  if (normalize_ && !normFactor_)
304  {
305  scalar J(JCycle()/weight_);
306  normFactor_.reset(new scalar(J));
307  DebugInfo
308  << "objective " << name() << ":: updating norm factor "
309  << "to " << normFactor_()
310  << " for time = " << mesh_.time().timeName() << endl;
311  }
312 }
313 
314 
315 void objective::accumulateJMean(solverControl& solverControl)
316 {
317  if (solverControl.doAverageIter())
318  {
319  const label iAverageIter = solverControl.averageIter();
320  if (iAverageIter == 0)
321  {
322  JMean_ = Zero;
323  }
324  scalar avIter(iAverageIter);
325  scalar oneOverItP1 = 1./(avIter + 1);
326  scalar mult = avIter*oneOverItP1;
327  JMean_ = JMean_*mult + J_*oneOverItP1;
328  }
329 }
330 
331 
333 {
335  {
336  const scalar time = mesh_.time().value();
338  {
339  const scalar dt = mesh_.time().deltaTValue();
340  const scalar elapsedTime = time - integrationStartTimePtr_();
341  const scalar denom = elapsedTime + dt;
342  JMean_ = (JMean_*elapsedTime + J_*dt)/denom;
343  }
344  }
345  else
346  {
348  << "Unallocated integration start or end time"
349  << exit(FatalError);
350  }
351 }
352 
354 scalar objective::weight() const
355 {
356  return weight_;
357 }
358 
360 bool objective::normalize() const
361 {
362  return normalize_;
363 }
364 
365 
367 {
368  if (normalize_ && normFactor_)
369  {
370  const scalar oneOverNorm(1./normFactor_());
371 
372  if (hasdJdb())
373  {
374  dJdbPtr_().primitiveFieldRef() *= oneOverNorm;
375  }
376  if (hasdJdbField())
377  {
378  dJdbFieldPtr_() *= oneOverNorm;
379  }
380  if (hasBoundarydJdb())
381  {
382  bdJdbPtr_() *= oneOverNorm;
383  }
384  if (hasdSdbMult())
385  {
386  bdSdbMultPtr_() *= oneOverNorm;
387  }
388  if (hasdndbMult())
389  {
390  bdndbMultPtr_() *= oneOverNorm;
391  }
392  if (hasdxdbMult())
393  {
394  bdxdbMultPtr_() *= oneOverNorm;
395  }
396  if (hasdxdbDirectMult())
397  {
398  bdxdbDirectMultPtr_() *= oneOverNorm;
399  }
401  {
402  bEdgeContribution_() *= oneOverNorm;
403  }
404  if (hasDivDxDbMult())
405  {
406  divDxDbMultPtr_() *= oneOverNorm;
407  }
408  if (hasGradDxDbMult())
409  {
410  gradDxDbMultPtr_() *= oneOverNorm;
411  }
412  }
413 }
414 
415 
417 {
419  {
420  const scalar time = mesh_.time().value();
421  return
422  (
425  );
426  }
427  else
428  {
430  << "Unallocated integration start or end time for objective '"
431  << objectiveName_ << "'"
432  << exit(FatalError);
433  }
434  return false;
435 }
436 
437 
438 void objective::incrementIntegrationTimes(const scalar timeSpan)
439 {
441  {
442  integrationStartTimePtr_() += timeSpan;
443  integrationEndTimePtr_() += timeSpan;
444  }
445  else
446  {
448  << "Unallocated integration start or end time"
449  << exit(FatalError);
450  }
451 }
452 
453 
454 void objective::update()
455 {
456  // Objective function value
457  J();
458 
459  // volFields
460  update_dJdb();
464 
465  // boundaryFields
472 }
473 
474 
475 void objective::nullify()
476 {
477  if (!nullified_)
478  {
479  if (hasdJdb())
480  {
481  dJdbPtr_() == dimensionedScalar(dJdbPtr_().dimensions(), Zero);
482  }
483  if (hasdJdbField())
484  {
485  dJdbFieldPtr_() = Zero;
486  }
487  if (hasBoundarydJdb())
488  {
489  bdJdbPtr_() == vector::zero;
490  }
491  if (hasdSdbMult())
492  {
493  bdSdbMultPtr_() == vector::zero;
494  }
495  if (hasdndbMult())
496  {
497  bdndbMultPtr_() == vector::zero;
498  }
499  if (hasdxdbMult())
500  {
501  bdxdbMultPtr_() == vector::zero;
502  }
503  if (hasdxdbDirectMult())
504  {
505  bdxdbDirectMultPtr_() == vector::zero;
506  }
508  {
509  for (Field<vectorField>& field : bEdgeContribution_())
510  {
511  field = vector::zero;
512  }
513  }
514  if (hasDivDxDbMult())
515  {
516  divDxDbMultPtr_() ==
517  dimensionedScalar(divDxDbMultPtr_().dimensions(), Zero);
518  }
519  if (hasGradDxDbMult())
520  {
521  gradDxDbMultPtr_() ==
522  dimensionedTensor(gradDxDbMultPtr_().dimensions(), Zero);
523  }
524 
525  nullified_ = true;
526  }
527 }
528 
529 
530 bool objective::write(const bool valid) const
531 {
532  if (Pstream::master())
533  {
534  // File is opened only upon invocation of the write function
535  // in order to avoid various instantiations of the same objective
536  // opening the same file
537  if (!objFunctionFilePtr_)
538  {
540  OFstream& file = objFunctionFilePtr_();
541  ios_base::fmtflags flags = file.flags();
542  flags |= std::cout.left;
543  file.flags(flags);
544  if (target_.valid())
545  {
546  file<< setw(width_) << "#target" << " "
547  << setw(width_) << target_() << endl;
548  }
549  if (targetLeft_.valid())
550  {
551  file<< setw(width_) << "#targetLeft" << " "
552  << setw(width_) << targetLeft_() << endl;
553  }
554  if (normalize_)
555  {
556  file<< setw(width_) << "#normFactor " << " "
557  << setw(width_) << normFactor_() << endl;
558  }
559  addHeaderInfo();
560  file<< setw(4) << "#" << " ";
561  file<< setw(width_) << "J" << " ";
562  file<< setw(width_) << "JCycle" << " ";
563  if (targetLeft_)
564  {
565  file<< setw(width_) << "JCycleLeft" << " ";
566  }
568  file<< endl;
569  }
570 
571  OFstream& file = objFunctionFilePtr_();
572  file<< setw(4) << mesh_.time().value() << " ";
573  file<< setw(width_) << J_ << " ";
574  file<< setw(width_) << JCycle() << " ";
575  if (targetLeft_)
576  {
577  file<< setw(width_) << JCycle(true) << " ";
578  }
579  addColumnValues();
580  file<< endl;
581  }
582 
583  return true;
584 }
585 
586 
588 {
589  if (Pstream::master())
590  {
591  // File is opened only upon invocation of the write function
592  // in order to avoid various instantiations of the same objective
593  // opening the same file
594  unsigned int width = IOstream::defaultPrecision() + 6;
596  {
598  }
601  << setw(width) << mesh_.time().value() << tab << J_ << endl;
602  }
603 }
604 
605 
607 {
608  if (Pstream::master())
609  {
611  {
613  }
614  }
615 }
616 
617 
618 void objective::writeMeanValue() const
619 {
620  if (Pstream::master())
621  {
622  // Write mean value if necessary
623  // Covers both steady and unsteady runs
624  if
625  (
628  )
629  {
630  // File is opened only upon invocation of the write function
631  // in order to avoid various instantiations of the same objective
632  // opening the same file
633  if (!meanValueFilePtr_)
634  {
636  }
637 
639  << mesh_.time().value() << tab << JMean_ << endl;
640  }
641  }
642 }
643 
644 
645 bool objective::writeData(Ostream& os) const
646 {
647  os.writeEntry("JMean", JMean_);
648  if (normFactor_)
649  {
650  os.writeEntry("normFactor", normFactor_());
651  }
652  os.writeEntry("weight", weight_);
653  return os.good();
654 }
655 
657 void objective::addHeaderInfo() const
658 {
659  // Does nothing in base
660 }
661 
663 void objective::addHeaderColumns() const
664 {
665  // Does nothing in base
666 }
667 
668 
669 void objective::addColumnValues() const
670 {
671  // Does nothing in base
672 }
673 
674 
675 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
676 
677 } // End namespace Foam
678 
679 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const Type & value() const noexcept
Return const reference to value.
void setInstantValueFilePtr() const
Set the output file ptr for the instantaneous value.
Definition: objective.C:65
dictionary dict
virtual void update_dSdbMultiplier()
Update d (normal dS) / db multiplier. Surface-based sensitivity term.
Definition: objective.H:503
virtual void updateNormalizationFactor()
Update normalization factors, for objectives in which the factor is not known a priori.
Definition: objective.C:294
autoPtr< boundaryVectorField > bdxdbMultPtr_
Term multiplying delta(x)/delta b at the boundary.
Definition: objective.H:171
rDeltaTY field()
autoPtr< scalar > target_
Target value, in case the objective is used as a constraint.
Definition: objective.H:110
A class for handling file names.
Definition: fileName.H:72
virtual void incrementIntegrationTimes(const scalar timeSpan)
Increment integration times.
Definition: objective.C:431
virtual void update_dJdb()
Definition: objective.H:488
const word adjointSolverName_
Definition: objective.H:65
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:91
virtual void writeMeanValue() const
Write mean objective function history.
Definition: objective.C:611
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:49
dictionary dict_
Definition: objective.H:64
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
bool hasdxdbDirectMult() const noexcept
Definition: objectiveI.H:199
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:209
bool hasBoundarydJdb() const noexcept
Definition: objectiveI.H:175
bool hasdJdbField() const noexcept
Definition: objectiveI.H:169
autoPtr< scalar > normFactor_
Normalization factor.
Definition: objective.H:101
virtual void update_dxdbMultiplier()
Update d (x) / db multiplier. Surface-based sensitivity term.
Definition: objective.H:515
Output to file stream, using an OSstream.
Definition: OFstream.H:49
autoPtr< scalar > targetLeft_
Target on the left hand-side of a double inequality, for double sided constraints.
Definition: objective.H:116
bool hasdndbMult() const noexcept
Definition: objectiveI.H:187
virtual void doNormalization()
Normalize all fields allocated by the objective.
Definition: objective.C:359
const word dictName("faMeshDefinition")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual scalar J()=0
Return the instantaneous objective function value.
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.
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:423
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
virtual void accumulateJMean()
Accumulate contribution for the mean objective value.
Definition: objective.C:325
autoPtr< vectorField3 > bEdgeContribution_
Contribution located in specific parts of a patch. Mainly intended for patch boundary edges contribut...
Definition: objective.H:185
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
autoPtr< boundaryVectorField > bdJdbPtr_
Term from material derivative.
Definition: objective.H:156
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
bool hasBoundaryEdgeContribution() const noexcept
Definition: objectiveI.H:205
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
const word objectiveName_
Definition: objective.H:67
virtual void update_gradDxDbMultiplier()
Update grad( dx/db multiplier). Volume-based sensitivity term.
Definition: objective.H:540
localIOdictionary is derived from IOdictionary but excludes parallel master reading.
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:37
bool hasGradDxDbMult() const noexcept
Definition: objectiveI.H:217
virtual void update_divDxDbMultiplier()
Update div( dx/db multiplier). Volume-based sensitivity term.
Definition: objective.H:534
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual bool write(const bool valid=true) const
Write objective function history.
Definition: objective.C:523
bool hasdSdbMult() const noexcept
Definition: objectiveI.H:181
static autoPtr< objective > New(const fvMesh &mesh, const dictionary &dict, const word &objectiveType, const word &adjointSolverName, const word &primalSolverName)
Return a reference to the selected turbulence model.
Definition: objective.C:215
autoPtr< scalar > integrationEndTimePtr_
Definition: objective.H:122
dynamicFvMesh & mesh
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:614
bool computeMeanFields_
Definition: objective.H:68
void setObjectiveFilePtr() const
Set the output file ptr.
Definition: objective.C:56
fileName globalPath() const
Return global path for the case = rootPath/globalCaseName. Same as TimePaths::globalPath() ...
Definition: Time.H:512
A class for handling words, derived from Foam::string.
Definition: word.H:63
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:405
void setMeanValueFilePtr() const
Set the output file ptr for the mean value.
Definition: objective.C:77
virtual void addHeaderColumns() const
Write headers for additional columns.
Definition: objective.C:656
virtual void update_dndbMultiplier()
Update d (normal) / db multiplier. Surface-based sensitivity term.
Definition: objective.H:509
bool hasdJdb() const noexcept
Definition: objectiveI.H:163
virtual void update_boundarydJdb()
Update objective function derivative term.
Definition: objective.H:497
const fvMesh & mesh_
Definition: objective.H:63
void negate(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1)
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
scalar JMean_
Average objective value.
Definition: objective.H:81
scalar J_
Objective function value and weight.
Definition: objective.H:76
#define DebugInfo
Report an information message using Foam::Info.
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
bool valid() const noexcept
Identical to good(), or bool operator.
Definition: autoPtr.H:415
virtual scalar weight() const
Return the objective function weight.
Definition: objective.C:347
virtual void update_dJdbField()
Definition: objective.H:491
virtual scalar JCycle(bool negate=false) const
Return the objective function of the optimisation cycle.
Definition: objective.C:258
Istream and Ostream manipulators taking arguments.
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
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
autoPtr< boundaryVectorField > bdndbMultPtr_
Term multiplying delta(n)/delta b.
Definition: objective.H:166
virtual bool writeData(Ostream &os) const
Write averaged objective for continuation.
Definition: objective.C:638
virtual void writeInstantaneousValue() const
Write objective function history at each primal solver iteration.
Definition: objective.C:580
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
virtual bool isWithinIntegrationTime() const
Check whether this is an objective integration time.
Definition: objective.C:409
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition: objective.H:194
virtual void addColumnValues() const
Write information to additional columns.
Definition: objective.C:662
autoPtr< OFstream > meanValueFilePtr_
File to keep the average objective values after the end of the primal solver.
Definition: objective.H:221
unsigned int width_
Default width of entries when writing in the objective files.
Definition: objective.H:226
virtual void update_boundaryEdgeContribution()
Update boundary edge contributions.
Definition: objective.H:528
An IOstream is an abstract base class for all input/output systems; be they streams, files, token lists etc.
Definition: IOstream.H:82
bool hasIntegrationStartTime() const noexcept
Definition: objectiveI.H:223
virtual void nullify()
Nullify adjoint contributions.
Definition: objective.C:468
autoPtr< scalarField > dJdbFieldPtr_
Contribution to field sensitivity derivatives.
Definition: objective.H:148
virtual void writeInstantaneousSeparator() const
Append a blank line after the end of optimisation cycle to the file holding the instantaneous objecti...
Definition: objective.C:599
const word & name() const
Return reference to name.
Definition: fvMesh.H:387
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool normalize() const
Is the objective normalized.
Definition: objective.C:353
scalar weight_
Objective weight.
Definition: objective.H:86
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
virtual void addHeaderInfo() const
Write any information that needs to go the header of the file.
Definition: objective.C:650
virtual bool readDict(const dictionary &dict)
Definition: objective.C:251
bool hasDivDxDbMult() const noexcept
Definition: objectiveI.H:211
bool hasdxdbMult() const noexcept
Definition: objectiveI.H:193
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
Definition: objective.H:199
autoPtr< scalar > integrationStartTimePtr_
Objective integration start and end times (for unsteady flows)
Definition: objective.H:121
virtual void update()=0
Update objective function derivatives.
Definition: objective.C:447
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
autoPtr< OFstream > instantValueFilePtr_
File to keep the objective values at each iteration of the primal solver.
Definition: objective.H:215
autoPtr< boundaryVectorField > bdxdbDirectMultPtr_
Term multiplying delta(x)/delta b at the boundary for objectives that directly depend on x...
Definition: objective.H:178
virtual void update_dxdbDirectMultiplier()
Update d (x) / db multiplier. Surface and volume-based sensitivity term.
Definition: objective.H:522
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
Definition: objective.H:161
fileName objFunctionFolder_
Output file variables.
Definition: objective.H:204
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:635
bool hasIntegrationEndTime() const noexcept
Definition: objectiveI.H:229
bool found
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127