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-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 
30 #include "objective.H"
31 #include "createZeroField.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 defineTypeNameAndDebug(objective, 0);
41 defineRunTimeSelectionTable(objective, objective);
42 
43 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
44 
45 void objective::makeFolder()
46 {
47  if (Pstream::master())
48  {
49  const Time& time = mesh_.time();
51  time.globalPath()/"optimisation"/type()/time.timeName();
52 
54  }
55 }
56 
57 
59 {
61  (
63  );
64 }
65 
66 
68 {
70  (
71  new OFstream
72  (
74  )
75  );
76 }
77 
78 
80 {
81  meanValueFilePtr_.reset
82  (
83  new OFstream
84  (
86  )
87  );
88 }
89 
90 
91 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
92 
94 {
95  return dict_;
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100 
101 objective::objective
102 (
103  const fvMesh& mesh,
104  const dictionary& dict,
105  const word& adjointSolverName,
106  const word& primalSolverName
107 )
108 :
109  mesh_(mesh),
110  dict_(dict),
111  adjointSolverName_(adjointSolverName),
112  primalSolverName_(primalSolverName),
113  objectiveName_(dict.dictName()),
114  computeMeanFields_(false), // is reset in derived classes
115  nullified_(false),
116 
117  J_(Zero),
118  JMean_(Zero),
119  weight_(Zero),
120 
121  integrationStartTimePtr_(nullptr),
122  integrationEndTimePtr_(nullptr),
123 
124  // Initialize pointers to nullptr.
125  // Not all of them are required for each objective function.
126  // Each child should allocate whatever is needed.
127 
128  dJdbPtr_(nullptr),
129  bdJdbPtr_(nullptr),
130  bdSdbMultPtr_(nullptr),
131  bdndbMultPtr_(nullptr),
132  bdxdbMultPtr_(nullptr),
133  bdxdbDirectMultPtr_(nullptr),
134  bEdgeContribution_(nullptr),
135  bdJdStressPtr_(nullptr),
136  divDxDbMultPtr_(nullptr),
137  gradDxDbMultPtr_(nullptr),
138 
139  objFunctionFolder_("word"),
140  objFunctionFilePtr_(nullptr),
141  instantValueFilePtr_(nullptr),
142  meanValueFilePtr_(nullptr)
143 {
144  makeFolder();
145  // Read integration start and end times, if present.
146  // For unsteady runs only
147  if (dict.found("integrationStartTime"))
148  {
149  integrationStartTimePtr_.reset
150  (
151  new scalar(dict.get<scalar>("integrationStartTime"))
152  );
153  }
154  if (dict.found("integrationEndTime"))
155  {
156  integrationEndTimePtr_.reset
157  (
158  new scalar(dict.get<scalar>("integrationEndTime"))
159  );
160  }
161 
162  // Read JMean from dictionary, if present
163  IOobject headObjectiveIODict
164  (
165  "objectiveDict" + objectiveName_,
166  mesh_.time().timeName(),
167  "uniform",
168  mesh_,
171  );
172  if (headObjectiveIODict.typeHeaderOk<IOdictionary>(false))
173  {
174  JMean_ = IOdictionary(headObjectiveIODict).get<scalar>("JMean");
175  }
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
180 
182 (
183  const fvMesh& mesh,
184  const dictionary& dict,
185  const word& objectiveType,
186  const word& adjointSolverName,
187  const word& primalSolverName
188 )
189 {
190  auto cstrIter = objectiveConstructorTablePtr_->cfind(objectiveType);
191 
192  if (!cstrIter.found())
193  {
195  (
196  dict,
197  "objective",
198  objectiveType,
199  *objectiveConstructorTablePtr_
200  ) << exit(FatalIOError);
201  }
202 
203  return autoPtr<objective>
204  (
205  cstrIter()
206  (
207  mesh,
208  dict,
209  adjointSolverName,
210  primalSolverName
211  )
212  );
213 }
214 
215 
216 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
217 
219 {
220  dict_ = dict;
221  return true;
222 }
223 
224 
225 scalar objective::JCycle() const
226 {
227  scalar J(J_);
228  if
229  (
232  )
233  {
234  J = JMean_;
235  }
236  return J;
237 }
238 
239 
241 {
242  // Does nothing in base
243 }
244 
245 
247 {
249  {
250  const label iAverageIter = solverControl.averageIter();
251  if (iAverageIter == 0)
252  {
253  JMean_ = Zero;
254  }
255  scalar avIter(iAverageIter);
256  scalar oneOverItP1 = 1./(avIter + 1);
257  scalar mult = avIter*oneOverItP1;
258  JMean_ = JMean_*mult + J_*oneOverItP1;
259  }
260 }
261 
262 
264 {
266  {
267  const scalar time = mesh_.time().value();
269  {
270  const scalar dt = mesh_.time().deltaT().value();
271  const scalar elapsedTime = time - integrationStartTimePtr_();
272  const scalar denom = elapsedTime + dt;
273  JMean_ = (JMean_*elapsedTime + J_*dt)/denom;
274  }
275  }
276  else
277  {
279  << "Unallocated integration start or end time"
280  << exit(FatalError);
281  }
282 }
283 
284 
285 scalar objective::weight() const
286 {
287  return weight_;
288 }
289 
290 
292 {
294  {
295  const scalar time = mesh_.time().value();
296  return
297  (
298  time >= integrationStartTimePtr_()
299  && time <= integrationEndTimePtr_()
300  );
301  }
302  else
303  {
305  << "Unallocated integration start or end time"
306  << exit(FatalError);
307  }
308  return false;
309 }
310 
311 
312 void objective::incrementIntegrationTimes(const scalar timeSpan)
313 {
315  {
316  integrationStartTimePtr_() += timeSpan;
317  integrationEndTimePtr_() += timeSpan;
318  }
319  else
320  {
322  << "Unallocated integration start or end time"
323  << exit(FatalError);
324  }
325 }
326 
327 
329 {
330  if (dJdbPtr_.empty())
331  {
332  // If pointer is not set, set it to a zero field
333  dJdbPtr_.reset
334  (
335  createZeroFieldPtr<scalar>
336  (
337  mesh_,
338  ("dJdb_" + objectiveName_),
339  dimensionSet(0, 5, -2, 0, 0, 0, 0)
340  )
341  );
342  }
343 
344  return dJdbPtr_();
345 }
346 
347 
349 {
350  if (bdJdbPtr_.empty())
351  {
352  bdJdbPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
353  }
354  return bdJdbPtr_()[patchI];
355 }
356 
357 
359 {
360  if (bdSdbMultPtr_.empty())
361  {
362  bdSdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
363  }
364  return bdSdbMultPtr_()[patchI];
365 }
366 
367 
369 {
370  if (bdndbMultPtr_.empty())
371  {
372  bdndbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
373  }
374  return bdndbMultPtr_()[patchI];
375 }
376 
377 
379 {
380  if (bdxdbMultPtr_.empty())
381  {
382  bdxdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
383  }
384  return bdxdbMultPtr_()[patchI];
385 }
386 
387 
389 {
390  if (bdxdbDirectMultPtr_.empty())
391  {
392  bdxdbDirectMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
393  }
394  return bdxdbDirectMultPtr_()[patchI];
395 }
396 
397 
399 (
400  const label patchI,
401  const label edgeI
402 )
403 {
404  if (bdxdbDirectMultPtr_.empty())
405  {
407  << "Unallocated boundaryEdgeMultiplier field"
408  << exit(FatalError);
409  }
410  return bEdgeContribution_()[patchI][edgeI];
411 }
412 
413 
415 {
416  if (bdJdStressPtr_.empty())
417  {
418  bdJdStressPtr_.reset(createZeroBoundaryPtr<tensor>(mesh_));
419  }
420  return bdJdStressPtr_()[patchI];
421 }
422 
423 
425 {
426  if (bdJdbPtr_.empty())
427  {
428  bdJdbPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
429  }
430  return bdJdbPtr_();
431 }
432 
433 
435 {
436  if (bdSdbMultPtr_.empty())
437  {
438  bdSdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
439  }
440  return bdSdbMultPtr_();
441 }
442 
443 
445 {
446  if (bdndbMultPtr_.empty())
447  {
448  bdndbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
449  }
450  return bdndbMultPtr_();
451 }
452 
453 
455 {
456  if (bdxdbMultPtr_.empty())
457  {
458  bdxdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
459  }
460  return bdxdbMultPtr_();
461 }
462 
463 
465 {
466  if (bdxdbDirectMultPtr_.empty())
467  {
468  bdxdbDirectMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
469  }
470  return bdxdbDirectMultPtr_();
471 }
472 
473 
475 {
476  if (bdxdbDirectMultPtr_.empty())
477  {
479  << "Unallocated boundaryEdgeMultiplier field"
480  << endl << endl
481  << exit(FatalError);
482  }
483  return bEdgeContribution_();
484 }
485 
486 
488 {
489  if (bdJdStressPtr_.empty())
490  {
491  bdJdStressPtr_.reset(createZeroBoundaryPtr<tensor>(mesh_));
492  }
493  return bdJdStressPtr_();
494 }
495 
496 
498 {
499  if (divDxDbMultPtr_.empty())
500  {
501  // If pointer is not set, set it to a zero field
502  divDxDbMultPtr_.reset
503  (
504  createZeroFieldPtr<scalar>
505  (
506  mesh_,
507  ("divDxDbMult"+objectiveName_),
508  // Variable dimensions!!
509  // Dummy dimensionless. Only the internalField will be used
510  dimless
511  )
512  );
513  }
514  return divDxDbMultPtr_();
515 }
516 
517 
519 {
520  if (gradDxDbMultPtr_.empty())
521  {
522  // If pointer is not set, set it to a zero field
523  gradDxDbMultPtr_.reset
524  (
525  createZeroFieldPtr<tensor>
526  (
527  mesh_,
528  ("gradDxDbMult"+objectiveName_),
529  // Variable dimensions!!
531  )
532  );
533  }
534  return gradDxDbMultPtr_();
535 }
536 
537 
539 {
540  if (!nullified_)
541  {
542  if (hasdJdb())
543  {
544  dJdbPtr_() == dimensionedScalar(dJdbPtr_().dimensions(), Zero);
545  }
546  if (hasBoundarydJdb())
547  {
548  bdJdbPtr_() == vector::zero;
549  }
550  if (hasdSdbMult())
551  {
553  }
554  if (hasdndbMult())
555  {
557  }
558  if (hasdxdbMult())
559  {
561  }
562  if (hasdxdbDirectMult())
563  {
565  }
567  {
569  {
571  }
572  }
573  if (hasBoundarydJdStress())
574  {
575  bdJdStressPtr_() == tensor::zero;
576  }
577  if (hasDivDxDbMult())
578  {
579  divDxDbMultPtr_() ==
580  dimensionedScalar(divDxDbMultPtr_().dimensions(), Zero);
581  }
582  if (hasGradDxDbMult())
583  {
584  gradDxDbMultPtr_() ==
585  dimensionedTensor(gradDxDbMultPtr_().dimensions(), Zero);
586  }
587 
588  nullified_ = true;
589  }
590 }
591 
592 
593 void objective::write() const
594 {
595  if (Pstream::master())
596  {
597  // File is opened only upon invocation of the write function
598  // in order to avoid various instantiations of the same objective
599  // opening the same file
600  if (objFunctionFilePtr_.empty())
601  {
603  }
604 
605  objFunctionFilePtr_() << mesh_.time().value() << tab << J_ << endl;
606  }
607 }
608 
609 
611 {
612  if (Pstream::master())
613  {
614  // File is opened only upon invocation of the write function
615  // in order to avoid various instantiations of the same objective
616  // opening the same file
617  if (instantValueFilePtr_.empty())
618  {
620  }
621 
622  instantValueFilePtr_() << mesh_.time().value() << tab << J_ << endl;
623  }
624 }
625 
626 
628 {
629  if (Pstream::master())
630  {
631  // Write mean value if necessary
632  // Covers both steady and unsteady runs
633  if
634  (
637  )
638  {
639  // File is opened only upon invocation of the write function
640  // in order to avoid various instantiations of the same objective
641  // opening the same file
642  if (meanValueFilePtr_.empty())
643  {
645  }
646 
648  << mesh_.time().value() << tab << JMean_ << endl;
649  }
650  }
651  // Write mean value under time/uniform, to allow for lineSearch to work
652  // appropriately in continuation runs, when field averaging is used
653  IOdictionary objectiveDict
654  (
655  IOobject
656  (
657  "objectiveDict" + objectiveName_,
658  mesh_.time().timeName(),
659  "uniform",
660  mesh_,
663  )
664  );
665  objectiveDict.add<scalar>("JMean", JMean_);
666  objectiveDict.regIOobject::write();
667 }
668 
669 
670 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
671 
672 } // End namespace Foam
673 
674 // ************************************************************************* //
Foam::fvPatchField< vector >
Foam::objective::nullify
virtual void nullify()
Nullify adjoint contributions.
Definition: objective.C:538
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::objective::bdxdbDirectMultPtr_
autoPtr< boundaryVectorField > bdxdbDirectMultPtr_
Definition: objective.H:105
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::objective::bdJdbPtr_
autoPtr< boundaryVectorField > bdJdbPtr_
Term from material derivative.
Definition: objective.H:91
Foam::objective::boundaryEdgeMultiplier
const vectorField3 & boundaryEdgeMultiplier()
Multiplier located at patch boundary edges.
Definition: objective.C:474
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::objective::hasBoundaryEdgeContribution
bool hasBoundaryEdgeContribution() const
Definition: objectiveI.H:75
Foam::solverControl::doAverageIter
bool doAverageIter() const
Definition: solverControlI.H:81
Foam::objective::integrationStartTimePtr_
autoPtr< scalar > integrationStartTimePtr_
Definition: objective.H:78
Foam::objective::gradDxDbMultPtr_
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
Definition: objective.H:123
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::objective::bdJdStressPtr_
autoPtr< boundaryTensorField > bdJdStressPtr_
For use with discrete-like sensitivities.
Definition: objective.H:113
Foam::TimeState::deltaT
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:54
Foam::objective::objectiveName_
const word objectiveName_
Definition: objective.H:68
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::pow2
dimensionSet pow2(const dimensionSet &ds)
Definition: dimensionSet.C:367
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionary.C:359
Foam::objective::dndbMultiplier
const boundaryVectorField & dndbMultiplier()
Multiplier of delta(n dS)/delta b for all patches.
Definition: objective.C:444
Foam::objective::hasdSdbMult
bool hasdSdbMult() const
Definition: objectiveI.H:51
Foam::objective::dxdbDirectMultiplier
const boundaryVectorField & dxdbDirectMultiplier()
Multiplier of delta(x)/delta b for all patches.
Definition: objective.C:464
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::objective::New
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:182
Foam::objective::setInstantValueFilePtr
void setInstantValueFilePtr() const
Set the output file ptr for the instantaneous value.
Definition: objective.C:67
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
Foam::objective::bdxdbMultPtr_
autoPtr< boundaryVectorField > bdxdbMultPtr_
Term multiplying delta(x)/delta b at the boundary.
Definition: objective.H:100
Foam::objective::setObjectiveFilePtr
void setObjectiveFilePtr() const
Set the output file ptr.
Definition: objective.C:58
Foam::FatalIOError
IOerror FatalIOError
Foam::solverControl
Base class for solver control classes.
Definition: solverControl.H:51
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:404
Foam::objective::JCycle
scalar JCycle() const
Definition: objective.C:225
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
Foam::objective::updateNormalizationFactor
virtual void updateNormalizationFactor()
Definition: objective.C:240
Foam::objective::hasdxdbMult
bool hasdxdbMult() const
Definition: objectiveI.H:63
Foam::objective::objFunctionFilePtr_
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:129
Foam::objective::hasdJdb
bool hasdJdb() const
Definition: objectiveI.H:39
Foam::objective::setMeanValueFilePtr
void setMeanValueFilePtr() const
Set the output file ptr for the mean value.
Definition: objective.C:79
Foam::objective::adjointSolverName_
const word adjointSolverName_
Definition: objective.H:66
Foam::objective::writeMeanValue
virtual void writeMeanValue() const
Write mean objective function history.
Definition: objective.C:627
Foam::objective::writeInstantaneousValue
virtual void writeInstantaneousValue() const
Write objective function history at each primal solver iteration.
Definition: objective.C:610
Foam::objective::dJdbPtr_
autoPtr< volScalarField > dJdbPtr_
Definition: objective.H:85
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::objective::hasIntegrationStartTime
bool hasIntegrationStartTime() const
Definition: objectiveI.H:99
Foam::objective::hasIntegrationEndTime
bool hasIntegrationEndTime() const
Definition: objectiveI.H:105
Foam::objective::hasDivDxDbMult
bool hasDivDxDbMult() const
Definition: objectiveI.H:81
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:380
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::objective::bEdgeContribution_
autoPtr< vectorField3 > bEdgeContribution_
Definition: objective.H:110
Foam::Field< vector >
Foam::objective::accumulateJMean
void accumulateJMean()
Accumulate contribution for the mean objective value.
Definition: objective.C:263
Foam::objective::J_
scalar J_
Definition: objective.H:73
Foam::objective::hasdndbMult
bool hasdndbMult() const
Definition: objectiveI.H:57
createZeroField.H
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::objective::divDxDbMultiplier
const volScalarField & divDxDbMultiplier()
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
Definition: objective.C:497
Foam::objective::gradDxDbMultiplier
const volTensorField & gradDxDbMultiplier()
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
Definition: objective.C:518
Foam::objective::bdndbMultPtr_
autoPtr< boundaryVectorField > bdndbMultPtr_
Term multiplying delta(n)/delta b.
Definition: objective.H:97
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
Foam::objective::dxdbMultiplier
const boundaryVectorField & dxdbMultiplier()
Multiplier of delta(x)/delta b for all patches.
Definition: objective.C:454
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
field
rDeltaTY field()
Foam::objective::weight_
scalar weight_
Definition: objective.H:75
Foam::objective::divDxDbMultPtr_
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition: objective.H:120
Foam::objective::meanValueFilePtr_
autoPtr< OFstream > meanValueFilePtr_
Definition: objective.H:137
Foam::solverControl::averageIter
label & averageIter()
Return average iteration index reference.
Definition: solverControlI.H:63
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionary.H:458
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
Foam::objective::integrationEndTimePtr_
autoPtr< scalar > integrationEndTimePtr_
Definition: objective.H:79
Foam::objective::boundarydJdStress
const boundaryTensorField & boundarydJdStress()
Objective partial deriv wrt stress tensor.
Definition: objective.C:487
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::objective::instantValueFilePtr_
autoPtr< OFstream > instantValueFilePtr_
Definition: objective.H:133
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::objective::hasBoundarydJdb
bool hasBoundarydJdb() const
Definition: objectiveI.H:45
Foam::objective::hasBoundarydJdStress
bool hasBoundarydJdStress() const
Definition: objectiveI.H:93
Foam::objective::bdSdbMultPtr_
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
Definition: objective.H:94
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:99
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::objective::nullified_
bool nullified_
Definition: objective.H:70
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::objective::computeMeanFields_
bool computeMeanFields_
Definition: objective.H:69
Foam::GeometricField::Boundary
Definition: GeometricField.H:114
Foam::tab
constexpr char tab
Definition: Ostream.H:371
Foam::objective::dSdbMultiplier
const boundaryVectorField & dSdbMultiplier()
Multiplier of delta(n dS)/delta b for all patches.
Definition: objective.C:434
Foam::objective::objFunctionFolder_
fileName objFunctionFolder_
Output file variables.
Definition: objective.H:126
Foam::objective::hasGradDxDbMult
bool hasGradDxDbMult() const
Definition: objectiveI.H:87
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:64
Foam::objective::boundarydJdb
const boundaryVectorField & boundarydJdb()
Contribution to surface sensitivities for all patches.
Definition: objective.C:424
Foam::objective::dict
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:93
Foam::objective::dict_
dictionary dict_
Definition: objective.H:65
Foam::objective::JMean_
scalar JMean_
Definition: objective.H:74
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::dictionary::add
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:703
Foam::objective::weight
scalar weight() const
Return the objective function weight.
Definition: objective.C:285
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
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::objective::hasdxdbDirectMult
bool hasdxdbDirectMult() const
Definition: objectiveI.H:69
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
objective.H
Foam::dimensionedTensor
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
Definition: dimensionedTensor.H:51
Foam::objective::incrementIntegrationTimes
void incrementIntegrationTimes(const scalar timeSpan)
Increment integration times.
Definition: objective.C:312
Foam::objective::readDict
virtual bool readDict(const dictionary &dict)
Definition: objective.C:218
Foam::objective::isWithinIntegrationTime
bool isWithinIntegrationTime() const
Check whether this is an objective integration time.
Definition: objective.C:291
Foam::objective::J
virtual scalar J()=0
Return the instantaneous objective function value.
Foam::objective::dJdb
const volScalarField & dJdb()
Contribution to field sensitivities.
Definition: objective.C:328
Foam::objective::write
virtual void write() const
Write objective function history.
Definition: objective.C:593