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-2020 PCOpt/NTUA
9  Copyright (C) 2013-2020 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();
51  objFunctionFolder_ =
52  time.globalPath()/"optimisation"/type()/time.timeName();
53 
54  mkDir(objFunctionFolder_);
55  }
56 }
57 
58 
59 void objective::setObjectiveFilePtr() const
60 {
61  objFunctionFilePtr_.reset
62  (
63  new OFstream(objFunctionFolder_/objectiveName_ + adjointSolverName_)
64  );
65 }
66 
67 
68 void objective::setInstantValueFilePtr() const
69 {
70  instantValueFilePtr_.reset
71  (
72  new OFstream
73  (
74  objFunctionFolder_/objectiveName_ + "Instant" + adjointSolverName_
75  )
76  );
77 }
78 
79 
80 void objective::setMeanValueFilePtr() const
81 {
82  meanValueFilePtr_.reset
83  (
84  new OFstream
85  (
86  objFunctionFolder_/objectiveName_ + "Mean" + adjointSolverName_
87  )
88  );
89 }
90 
91 
92 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
93 
94 const dictionary& objective::dict() const
95 {
96  return dict_;
97 }
98 
99 
100 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
101 
102 objective::objective
103 (
104  const fvMesh& mesh,
105  const dictionary& dict,
106  const word& adjointSolverName,
107  const word& primalSolverName
108 )
109 :
110  localIOdictionary
111  (
112  IOobject
113  (
114  dict.dictName(),
115  mesh.time().timeName(),
116  fileName("uniform")/fileName("objectives")/adjointSolverName,
117  mesh,
118  IOobject::READ_IF_PRESENT,
119  IOobject::AUTO_WRITE
120  ),
121  // avoid type checking since dictionary is read using the
122  // derived type name and type() will result in "objective"
123  // here
124  word::null
125  ),
126  mesh_(mesh),
127  dict_(dict),
128  adjointSolverName_(adjointSolverName),
129  primalSolverName_(primalSolverName),
130  objectiveName_(dict.dictName()),
131  computeMeanFields_(false), // is reset in derived classes
132  nullified_(false),
133  normalize_(dict.getOrDefault<bool>("normalize", false)),
134 
135  J_(Zero),
136  JMean_(this->getOrDefault<scalar>("JMean", Zero)),
137  weight_(Zero),
138  normFactor_(nullptr),
139  target_
140  (
141  dict.found("target") ?
142  autoPtr<scalar>::New(dict.get<scalar>("target")) :
143  nullptr
144  ),
145  integrationStartTimePtr_(nullptr),
146  integrationEndTimePtr_(nullptr),
147 
148  // Initialize pointers to nullptr.
149  // Not all of them are required for each objective function.
150  // Each child should allocate whatever is needed.
151 
152  dJdbPtr_(nullptr),
153  bdJdbPtr_(nullptr),
154  bdSdbMultPtr_(nullptr),
155  bdndbMultPtr_(nullptr),
156  bdxdbMultPtr_(nullptr),
157  bdxdbDirectMultPtr_(nullptr),
158  bEdgeContribution_(nullptr),
159  bdJdStressPtr_(nullptr),
160  divDxDbMultPtr_(nullptr),
161  gradDxDbMultPtr_(nullptr),
162 
163  objFunctionFolder_("word"),
164  objFunctionFilePtr_(nullptr),
165  instantValueFilePtr_(nullptr),
166  meanValueFilePtr_(nullptr),
167  width_(IOstream::defaultPrecision() + 5)
168 {
169  makeFolder();
170  // Read integration start and end times, if present.
171  // For unsteady runs only
172  if (dict.found("integrationStartTime"))
173  {
175  (
176  new scalar(dict.get<scalar>("integrationStartTime"))
177  );
178  }
179  if (dict.found("integrationEndTime"))
180  {
182  (
183  new scalar(dict.get<scalar>("integrationEndTime"))
184  );
185  }
186 
187  // Set normalization factor, if present
188  if (normalize_)
189  {
190  scalar normFactor(Zero);
191  if (dict.readIfPresent("normFactor", normFactor))
192  {
193  normFactor_.reset(new scalar(normFactor));
194  }
195  else if (this->readIfPresent("normFactor", normFactor))
196  {
197  normFactor_.reset(new scalar(normFactor));
198  }
199  }
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
204 
205 autoPtr<objective> objective::New
206 (
207  const fvMesh& mesh,
208  const dictionary& dict,
209  const word& objectiveType,
210  const word& adjointSolverName,
211  const word& primalSolverName
212 )
213 {
214  auto* ctorPtr = objectiveConstructorTable(objectiveType);
215 
216  if (!ctorPtr)
217  {
219  (
220  dict,
221  "objective",
222  objectiveType,
223  *objectiveConstructorTablePtr_
224  ) << exit(FatalIOError);
225  }
226 
227  return autoPtr<objective>
228  (
229  ctorPtr
230  (
231  mesh,
232  dict,
233  adjointSolverName,
234  primalSolverName
235  )
236  );
237 }
238 
239 
240 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
241 
242 bool objective::readDict(const dictionary& dict)
243 {
244  dict_ = dict;
245  return true;
246 }
247 
248 
249 scalar objective::JCycle() const
250 {
251  scalar J(J_);
252  if
253  (
254  computeMeanFields_
255  || (hasIntegrationStartTime() && hasIntegrationEndTime())
256  )
257  {
258  J = JMean_;
259  }
260 
261  // Subtract target, in case the objective is used as a constraint
262  if (target_.valid())
263  {
264  J -= target_();
265  }
266 
267  // Normalize here, in order to get the correct value for line search
268  if (normalize_ && normFactor_)
269  {
270  J /= normFactor_();
271  }
272 
273  return J;
274 }
275 
276 
277 void objective::updateNormalizationFactor()
278 {
279  if (normalize_ && !normFactor_)
280  {
281  normFactor_.reset(new scalar(JCycle()));
282  }
283 }
284 
285 
286 void objective::accumulateJMean(solverControl& solverControl)
287 {
288  if (solverControl.doAverageIter())
289  {
290  const label iAverageIter = solverControl.averageIter();
291  if (iAverageIter == 0)
292  {
293  JMean_ = Zero;
294  }
295  scalar avIter(iAverageIter);
296  scalar oneOverItP1 = 1./(avIter + 1);
297  scalar mult = avIter*oneOverItP1;
298  JMean_ = JMean_*mult + J_*oneOverItP1;
299  }
300 }
301 
302 
303 void objective::accumulateJMean()
304 {
305  if (hasIntegrationStartTime() && hasIntegrationEndTime())
306  {
307  const scalar time = mesh_.time().value();
308  if (isWithinIntegrationTime())
309  {
310  const scalar dt = mesh_.time().deltaT().value();
311  const scalar elapsedTime = time - integrationStartTimePtr_();
312  const scalar denom = elapsedTime + dt;
313  JMean_ = (JMean_*elapsedTime + J_*dt)/denom;
314  }
315  }
316  else
317  {
319  << "Unallocated integration start or end time"
320  << exit(FatalError);
321  }
322 }
323 
324 
325 scalar objective::weight() const
326 {
327  return weight_;
328 }
329 
330 
331 bool objective::normalize() const
332 {
333  return normalize_;
334 }
335 
336 
337 void objective::doNormalization()
338 {
339  if (normalize_ && normFactor_)
340  {
341  const scalar oneOverNorm(1./normFactor_());
342 
343  if (hasdJdb())
344  {
345  dJdbPtr_().primitiveFieldRef() *= oneOverNorm;
346  }
347  if (hasBoundarydJdb())
348  {
349  bdJdbPtr_() *= oneOverNorm;
350  }
351  if (hasdSdbMult())
352  {
353  bdSdbMultPtr_() *= oneOverNorm;
354  }
355  if (hasdndbMult())
356  {
357  bdndbMultPtr_() *= oneOverNorm;
358  }
359  if (hasdxdbMult())
360  {
361  bdxdbMultPtr_() *= oneOverNorm;
362  }
363  if (hasdxdbDirectMult())
364  {
365  bdxdbDirectMultPtr_() *= oneOverNorm;
366  }
367  if (hasBoundaryEdgeContribution())
368  {
369  bEdgeContribution_() *= oneOverNorm;
370  }
371  if (hasDivDxDbMult())
372  {
373  divDxDbMultPtr_() *= oneOverNorm;
374  }
375  if (hasGradDxDbMult())
376  {
377  gradDxDbMultPtr_() *= oneOverNorm;
378  }
379  if (hasBoundarydJdStress())
380  {
381  bdJdStressPtr_() *= oneOverNorm;
382  }
383  }
384 }
385 
386 
387 bool objective::isWithinIntegrationTime() const
388 {
389  if (hasIntegrationStartTime() && hasIntegrationEndTime())
390  {
391  const scalar time = mesh_.time().value();
392  return
393  (
394  time >= integrationStartTimePtr_()
395  && time <= integrationEndTimePtr_()
396  );
397  }
398  else
399  {
401  << "Unallocated integration start or end time"
402  << exit(FatalError);
403  }
404  return false;
405 }
406 
407 
408 void objective::incrementIntegrationTimes(const scalar timeSpan)
409 {
410  if (hasIntegrationStartTime() && hasIntegrationEndTime())
411  {
412  integrationStartTimePtr_() += timeSpan;
413  integrationEndTimePtr_() += timeSpan;
414  }
415  else
416  {
418  << "Unallocated integration start or end time"
419  << exit(FatalError);
420  }
421 }
422 
423 
424 const volScalarField& objective::dJdb()
425 {
426  if (!dJdbPtr_)
427  {
428  // If pointer is not set, set it to a zero field
429  dJdbPtr_.reset
430  (
431  createZeroFieldPtr<scalar>
432  (
433  mesh_,
434  ("dJdb_" + objectiveName_),
435  dimensionSet(0, 5, -2, 0, 0, 0, 0)
436  )
437  );
438  }
439 
440  return *dJdbPtr_;
441 }
442 
443 
444 const fvPatchVectorField& objective::boundarydJdb(const label patchI)
445 {
446  if (!bdJdbPtr_)
447  {
448  bdJdbPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
449  }
450  return bdJdbPtr_()[patchI];
451 }
452 
453 
454 const fvPatchVectorField& objective::dSdbMultiplier(const label patchI)
455 {
456  if (!bdSdbMultPtr_)
457  {
458  bdSdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
459  }
460  return bdSdbMultPtr_()[patchI];
461 }
462 
463 
464 const fvPatchVectorField& objective::dndbMultiplier(const label patchI)
465 {
466  if (!bdndbMultPtr_)
467  {
468  bdndbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
469  }
470  return bdndbMultPtr_()[patchI];
471 }
472 
473 
474 const fvPatchVectorField& objective::dxdbMultiplier(const label patchI)
475 {
476  if (!bdxdbMultPtr_)
477  {
478  bdxdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
479  }
480  return bdxdbMultPtr_()[patchI];
481 }
482 
483 
484 const fvPatchVectorField& objective::dxdbDirectMultiplier(const label patchI)
485 {
486  if (!bdxdbDirectMultPtr_)
487  {
488  bdxdbDirectMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
489  }
490  return bdxdbDirectMultPtr_()[patchI];
491 }
492 
493 
494 const vectorField& objective::boundaryEdgeMultiplier
495 (
496  const label patchI,
497  const label edgeI
498 )
499 {
500  if (!bdxdbDirectMultPtr_)
501  {
503  << "Unallocated boundaryEdgeMultiplier field"
504  << exit(FatalError);
505  }
506  return bEdgeContribution_()[patchI][edgeI];
507 }
508 
509 
510 const fvPatchTensorField& objective::boundarydJdStress(const label patchI)
511 {
512  if (!bdJdStressPtr_)
513  {
514  bdJdStressPtr_.reset(createZeroBoundaryPtr<tensor>(mesh_));
515  }
516  return bdJdStressPtr_()[patchI];
517 }
518 
519 
520 const boundaryVectorField& objective::boundarydJdb()
521 {
522  if (!bdJdbPtr_)
523  {
524  bdJdbPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
525  }
526  return *bdJdbPtr_;
527 }
528 
529 
530 const boundaryVectorField& objective::dSdbMultiplier()
531 {
532  if (!bdSdbMultPtr_)
533  {
534  bdSdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
535  }
536  return *bdSdbMultPtr_;
537 }
538 
539 
540 const boundaryVectorField& objective::dndbMultiplier()
541 {
542  if (!bdndbMultPtr_)
543  {
544  bdndbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
545  }
546  return *bdndbMultPtr_;
547 }
548 
549 
550 const boundaryVectorField& objective::dxdbMultiplier()
551 {
552  if (!bdxdbMultPtr_)
553  {
554  bdxdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
555  }
556  return *bdxdbMultPtr_;
557 }
558 
559 
560 const boundaryVectorField& objective::dxdbDirectMultiplier()
561 {
562  if (!bdxdbDirectMultPtr_)
563  {
564  bdxdbDirectMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
565  }
566  return *bdxdbDirectMultPtr_;
567 }
568 
569 
570 const vectorField3& objective::boundaryEdgeMultiplier()
571 {
572  if (!bdxdbDirectMultPtr_)
573  {
575  << "Unallocated boundaryEdgeMultiplier field"
576  << endl << endl
577  << exit(FatalError);
578  }
579  return *bEdgeContribution_;
580 }
581 
582 
583 const boundaryTensorField& objective::boundarydJdStress()
584 {
585  if (!bdJdStressPtr_)
586  {
587  bdJdStressPtr_.reset(createZeroBoundaryPtr<tensor>(mesh_));
588  }
589  return *bdJdStressPtr_;
590 }
591 
592 
593 const volScalarField& objective::divDxDbMultiplier()
594 {
595  if (!divDxDbMultPtr_)
596  {
597  // If pointer is not set, set it to a zero field
598  divDxDbMultPtr_.reset
599  (
600  createZeroFieldPtr<scalar>
601  (
602  mesh_,
603  ("divDxDbMult"+objectiveName_),
604  // Variable dimensions!!
605  // Dummy dimensionless. Only the internalField will be used
606  dimless
607  )
608  );
609  }
610  return *divDxDbMultPtr_;
611 }
612 
613 
614 const volTensorField& objective::gradDxDbMultiplier()
615 {
616  if (!gradDxDbMultPtr_)
617  {
618  // If pointer is not set, set it to a zero field
619  gradDxDbMultPtr_.reset
620  (
621  createZeroFieldPtr<tensor>
622  (
623  mesh_,
624  ("gradDxDbMult"+objectiveName_),
625  // Variable dimensions!!
626  dimensionSet(pow2(dimLength)/pow3(dimTime))
627  )
628  );
629  }
630  return *gradDxDbMultPtr_;
631 }
632 
633 
634 void objective::nullify()
635 {
636  if (!nullified_)
637  {
638  if (hasdJdb())
639  {
640  dJdbPtr_() == dimensionedScalar(dJdbPtr_().dimensions(), Zero);
641  }
642  if (hasBoundarydJdb())
643  {
644  bdJdbPtr_() == vector::zero;
645  }
646  if (hasdSdbMult())
647  {
648  bdSdbMultPtr_() == vector::zero;
649  }
650  if (hasdndbMult())
651  {
652  bdndbMultPtr_() == vector::zero;
653  }
654  if (hasdxdbMult())
655  {
656  bdxdbMultPtr_() == vector::zero;
657  }
658  if (hasdxdbDirectMult())
659  {
660  bdxdbDirectMultPtr_() == vector::zero;
661  }
662  if (hasBoundaryEdgeContribution())
663  {
664  for (Field<vectorField>& field : bEdgeContribution_())
665  {
667  }
668  }
669  if (hasBoundarydJdStress())
670  {
671  bdJdStressPtr_() == tensor::zero;
672  }
673  if (hasDivDxDbMult())
674  {
675  divDxDbMultPtr_() ==
676  dimensionedScalar(divDxDbMultPtr_().dimensions(), Zero);
677  }
678  if (hasGradDxDbMult())
679  {
680  gradDxDbMultPtr_() ==
681  dimensionedTensor(gradDxDbMultPtr_().dimensions(), Zero);
682  }
683 
684  nullified_ = true;
685  }
686 }
687 
688 
689 bool objective::write(const bool valid) const
690 {
691  if (Pstream::master())
692  {
693  // File is opened only upon invocation of the write function
694  // in order to avoid various instantiations of the same objective
695  // opening the same file
696  if (!objFunctionFilePtr_)
697  {
698  setObjectiveFilePtr();
699  OFstream& file = objFunctionFilePtr_();
700  ios_base::fmtflags flags = file.flags();
701  flags |= std::cout.left;
702  file.flags(flags);
703  if (target_.valid())
704  {
705  file<< setw(width_) << "#target" << " "
706  << setw(width_) << target_() << endl;
707  }
708  if (normalize_)
709  {
710  file<< setw(width_) << "#normFactor " << " "
711  << setw(width_) << normFactor_() << endl;
712  }
713  addHeaderInfo();
714  file<< setw(4) << "#" << " ";
715  file<< setw(width_) << "J" << " ";
716  file<< setw(width_) << "JCycle" << " ";
717  addHeaderColumns();
718  file<< endl;
719  }
720 
721  OFstream& file = objFunctionFilePtr_();
722  file<< setw(4) << mesh_.time().value() << " ";
723  file<< setw(width_) << J_ << " ";
724  file<< setw(width_) << JCycle() << " ";
725  addColumnValues();
726  file<< endl;
727  }
728 
729  return true;
730 }
731 
732 
733 void objective::writeInstantaneousValue() const
734 {
735  if (Pstream::master())
736  {
737  // File is opened only upon invocation of the write function
738  // in order to avoid various instantiations of the same objective
739  // opening the same file
740  if (!instantValueFilePtr_)
741  {
742  setInstantValueFilePtr();
743  }
744 
745  instantValueFilePtr_() << mesh_.time().value() << tab << J_ << endl;
746  }
747 }
748 
749 
750 void objective::writeInstantaneousSeparator() const
751 {
752  if (Pstream::master())
753  {
754  if (instantValueFilePtr_)
755  {
756  instantValueFilePtr_() << endl;
757  }
758  }
759 }
760 
761 
762 void objective::writeMeanValue() const
763 {
764  if (Pstream::master())
765  {
766  // Write mean value if necessary
767  // Covers both steady and unsteady runs
768  if
769  (
770  computeMeanFields_
771  || (hasIntegrationStartTime() && hasIntegrationEndTime())
772  )
773  {
774  // File is opened only upon invocation of the write function
775  // in order to avoid various instantiations of the same objective
776  // opening the same file
777  if (!meanValueFilePtr_)
778  {
779  setMeanValueFilePtr();
780  }
781 
782  meanValueFilePtr_()
783  << mesh_.time().value() << tab << JMean_ << endl;
784  }
785  }
786 }
787 
788 
789 bool objective::writeData(Ostream& os) const
790 {
791  os.writeEntry("JMean", JMean_);
792  if (normFactor_)
793  {
794  os.writeEntry("normFactor", normFactor_());
795  }
796  return os.good();
797 }
798 
799 
800 void objective::addHeaderInfo() const
801 {
802  // Does nothing in base
803 }
804 
805 
806 void objective::addHeaderColumns() const
807 {
808  // Does nothing in base
809 }
810 
811 
812 void objective::addColumnValues() const
813 {
814  // Does nothing in base
815 }
816 
817 
818 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
819 
820 } // End namespace Foam
821 
822 // ************************************************************************* //
Foam::objective::normFactor_
autoPtr< scalar > normFactor_
Normalization factor.
Definition: objective.H:86
Foam::volTensorField
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:66
Foam::objective::integrationStartTimePtr_
autoPtr< scalar > integrationStartTimePtr_
Objective integration start and end times (for unsteady flows)
Definition: objective.H:94
Foam::objective::normalize_
bool normalize_
Definition: objective.H:74
Foam::baseIOdictionary::writeData
virtual bool writeData(Ostream &) const
The writeData function required by regIOobject write operation.
Definition: baseIOdictionary.C:109
Foam::objective::objectiveName_
const word objectiveName_
Definition: objective.H:71
Foam::dictionary::New
static autoPtr< dictionary > New(Istream &is)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:69
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::pow2
dimensionSet pow2(const dimensionSet &ds)
Definition: dimensionSet.C:369
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
dictName
const word dictName("faMeshDefinition")
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::IOobject::time
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:493
Foam::objective::adjointSolverName_
const word adjointSolverName_
Definition: objective.H:69
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Foam::vectorField3
Field< Field< vectorField > > vectorField3
Definition: objectiveFwd.H:39
createZeroField.H
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::boundaryTensorField
volTensorField::Boundary boundaryTensorField
Definition: boundaryFieldsFwd.H:56
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
field
rDeltaTY field()
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::objective::integrationEndTimePtr_
autoPtr< scalar > integrationEndTimePtr_
Definition: objective.H:95
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::autoPtr::reset
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
Foam::fvPatchVectorField
fvPatchField< vector > fvPatchVectorField
Definition: fvPatchFieldsFwd.H:43
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::tab
constexpr char tab
Definition: Ostream.H:403
Foam::dictionary::write
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:206
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::objective::mesh_
const fvMesh & mesh_
Definition: objective.H:67
bool
bool
Definition: EEqn.H:20
Foam::TimePaths::globalPath
fileName globalPath() const
Return global path for the case.
Definition: TimePathsI.H:80
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::objective::dict_
dictionary dict_
Definition: objective.H:68
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::fvPatchTensorField
fvPatchField< tensor > fvPatchTensorField
Definition: fvPatchFieldsFwd.H:46
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
objective.H
Foam::boundaryVectorField
volVectorField::Boundary boundaryVectorField
Definition: boundaryFieldsFwd.H:55
Foam::dimensionedTensor
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
Definition: dimensionedTensor.H:52
Foam::normalize
quaternion normalize(const quaternion &q)
Return the normalized (unit) quaternion of the given quaternion.
Definition: quaternionI.H:661
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189