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-------------------------------------------------------------------------------
12License
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
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
43
44// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
45
46void objective::makeFolder()
47{
48 if (Pstream::master())
49 {
50 const Time& time = mesh_.time();
52 time.globalPath()/"optimisation"/type()/time.timeName();
53
55 }
56}
57
58
60{
62 (
64 );
65}
66
67
69{
71 (
72 new OFstream
73 (
75 )
76 );
77}
78
79
81{
83 (
84 new OFstream
85 (
87 )
88 );
89}
90
91
92// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
93
95{
96 return dict_;
97}
98
99
100// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
101
103(
104 const fvMesh& mesh,
105 const dictionary& dict,
106 const word& adjointSolverName,
107 const word& primalSolverName
108)
109:
111 (
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 divDxDbMultPtr_(nullptr),
160 gradDxDbMultPtr_(nullptr),
161
162 objFunctionFolder_("word"),
163 objFunctionFilePtr_(nullptr),
164 instantValueFilePtr_(nullptr),
165 meanValueFilePtr_(nullptr),
166 width_(IOstream::defaultPrecision() + 5)
167{
168 makeFolder();
169 // Read integration start and end times, if present.
170 // For unsteady runs only
171 if (dict.found("integrationStartTime"))
172 {
174 (
175 new scalar(dict.get<scalar>("integrationStartTime"))
176 );
177 }
178 if (dict.found("integrationEndTime"))
179 {
181 (
182 new scalar(dict.get<scalar>("integrationEndTime"))
183 );
184 }
185
186 // Set normalization factor, if present
187 if (normalize_)
188 {
189 scalar normFactor(Zero);
190 if (dict.readIfPresent("normFactor", normFactor))
191 {
192 normFactor_.reset(new scalar(normFactor));
193 }
194 else if (this->readIfPresent("normFactor", normFactor))
195 {
196 normFactor_.reset(new scalar(normFactor));
197 }
198 }
199}
200
201
202// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
203
205(
206 const fvMesh& mesh,
207 const dictionary& dict,
208 const word& objectiveType,
209 const word& adjointSolverName,
210 const word& primalSolverName
211)
212{
213 auto* ctorPtr = objectiveConstructorTable(objectiveType);
214
215 if (!ctorPtr)
216 {
218 (
219 dict,
220 "objective",
221 objectiveType,
222 *objectiveConstructorTablePtr_
223 ) << exit(FatalIOError);
224 }
225
226 return autoPtr<objective>
227 (
228 ctorPtr
229 (
230 mesh,
231 dict,
232 adjointSolverName,
233 primalSolverName
234 )
235 );
236}
237
238
239// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
240
241bool objective::readDict(const dictionary& dict)
242{
243 dict_ = dict;
244 return true;
245}
246
247
248scalar objective::JCycle() const
249{
250 scalar J(J_);
251 if
252 (
255 )
256 {
257 J = JMean_;
258 }
259
260 // Subtract target, in case the objective is used as a constraint
261 if (target_.valid())
262 {
263 J -= target_();
264 }
265
266 // Normalize here, in order to get the correct value for line search
267 if (normalize_ && normFactor_)
268 {
269 J /= normFactor_();
270 }
271
272 return J;
273}
274
275
277{
278 if (normalize_ && !normFactor_)
279 {
280 normFactor_.reset(new scalar(JCycle()));
281 }
282}
283
284
286{
288 {
289 const label iAverageIter = solverControl.averageIter();
290 if (iAverageIter == 0)
291 {
292 JMean_ = Zero;
293 }
294 scalar avIter(iAverageIter);
295 scalar oneOverItP1 = 1./(avIter + 1);
296 scalar mult = avIter*oneOverItP1;
297 JMean_ = JMean_*mult + J_*oneOverItP1;
298 }
299}
300
301
303{
305 {
306 const scalar time = mesh_.time().value();
308 {
309 const scalar dt = mesh_.time().deltaT().value();
310 const scalar elapsedTime = time - integrationStartTimePtr_();
311 const scalar denom = elapsedTime + dt;
312 JMean_ = (JMean_*elapsedTime + J_*dt)/denom;
313 }
314 }
315 else
316 {
318 << "Unallocated integration start or end time"
319 << exit(FatalError);
320 }
321}
322
323
324scalar objective::weight() const
325{
326 return weight_;
327}
328
329
331{
332 return normalize_;
333}
334
335
337{
338 if (normalize_ && normFactor_)
339 {
340 const scalar oneOverNorm(1./normFactor_());
341
342 if (hasdJdb())
343 {
344 dJdbPtr_().primitiveFieldRef() *= oneOverNorm;
345 }
346 if (hasBoundarydJdb())
347 {
348 bdJdbPtr_() *= oneOverNorm;
349 }
350 if (hasdSdbMult())
351 {
352 bdSdbMultPtr_() *= oneOverNorm;
353 }
354 if (hasdndbMult())
355 {
356 bdndbMultPtr_() *= oneOverNorm;
357 }
358 if (hasdxdbMult())
359 {
360 bdxdbMultPtr_() *= oneOverNorm;
361 }
362 if (hasdxdbDirectMult())
363 {
364 bdxdbDirectMultPtr_() *= oneOverNorm;
365 }
367 {
368 bEdgeContribution_() *= oneOverNorm;
369 }
370 if (hasDivDxDbMult())
371 {
372 divDxDbMultPtr_() *= oneOverNorm;
373 }
374 if (hasGradDxDbMult())
375 {
376 gradDxDbMultPtr_() *= oneOverNorm;
377 }
378 }
379}
380
381
383{
385 {
386 const scalar time = mesh_.time().value();
387 return
388 (
391 );
392 }
393 else
394 {
396 << "Unallocated integration start or end time"
397 << exit(FatalError);
398 }
399 return false;
400}
401
402
403void objective::incrementIntegrationTimes(const scalar timeSpan)
404{
406 {
407 integrationStartTimePtr_() += timeSpan;
408 integrationEndTimePtr_() += timeSpan;
409 }
410 else
411 {
413 << "Unallocated integration start or end time"
414 << exit(FatalError);
415 }
416}
417
418
420{
421 if (!dJdbPtr_)
422 {
423 // If pointer is not set, set it to a zero field
424 dJdbPtr_.reset
425 (
426 createZeroFieldPtr<scalar>
427 (
428 mesh_,
429 ("dJdb_" + objectiveName_),
430 dimensionSet(0, 5, -2, 0, 0, 0, 0)
431 )
432 );
433 }
434
435 return *dJdbPtr_;
436}
437
438
440{
441 if (!bdJdbPtr_)
442 {
443 bdJdbPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
444 }
445 return bdJdbPtr_()[patchI];
446}
447
448
450{
451 if (!bdSdbMultPtr_)
452 {
453 bdSdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
454 }
455 return bdSdbMultPtr_()[patchI];
456}
457
458
460{
461 if (!bdndbMultPtr_)
462 {
463 bdndbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
464 }
465 return bdndbMultPtr_()[patchI];
466}
467
468
470{
471 if (!bdxdbMultPtr_)
472 {
473 bdxdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
474 }
475 return bdxdbMultPtr_()[patchI];
476}
477
478
480{
482 {
483 bdxdbDirectMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
484 }
485 return bdxdbDirectMultPtr_()[patchI];
486}
487
488
490(
491 const label patchI,
492 const label edgeI
493)
494{
496 {
498 << "Unallocated boundaryEdgeMultiplier field"
499 << exit(FatalError);
500 }
501 return bEdgeContribution_()[patchI][edgeI];
502}
503
504
506{
507 if (!bdJdbPtr_)
508 {
509 bdJdbPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
510 }
511 return *bdJdbPtr_;
512}
513
514
516{
517 if (!bdSdbMultPtr_)
518 {
519 bdSdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
520 }
521 return *bdSdbMultPtr_;
522}
523
524
526{
527 if (!bdndbMultPtr_)
528 {
529 bdndbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
530 }
531 return *bdndbMultPtr_;
532}
533
534
536{
537 if (!bdxdbMultPtr_)
538 {
539 bdxdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
540 }
541 return *bdxdbMultPtr_;
542}
543
544
546{
548 {
549 bdxdbDirectMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
550 }
551 return *bdxdbDirectMultPtr_;
552}
553
554
556{
558 {
560 << "Unallocated boundaryEdgeMultiplier field"
561 << endl << endl
562 << exit(FatalError);
563 }
564 return *bEdgeContribution_;
565}
566
567
569{
570 if (!divDxDbMultPtr_)
571 {
572 // If pointer is not set, set it to a zero field
573 divDxDbMultPtr_.reset
574 (
575 createZeroFieldPtr<scalar>
576 (
577 mesh_,
578 ("divDxDbMult"+objectiveName_),
579 // Variable dimensions!!
580 // Dummy dimensionless. Only the internalField will be used
581 dimless
582 )
583 );
584 }
585 return *divDxDbMultPtr_;
586}
587
588
590{
591 if (!gradDxDbMultPtr_)
592 {
593 // If pointer is not set, set it to a zero field
594 gradDxDbMultPtr_.reset
595 (
596 createZeroFieldPtr<tensor>
597 (
598 mesh_,
599 ("gradDxDbMult"+objectiveName_),
600 // Variable dimensions!!
602 )
603 );
604 }
605 return *gradDxDbMultPtr_;
606}
607
608
610{
611 if (!nullified_)
612 {
613 if (hasdJdb())
614 {
615 dJdbPtr_() == dimensionedScalar(dJdbPtr_().dimensions(), Zero);
616 }
617 if (hasBoundarydJdb())
618 {
620 }
621 if (hasdSdbMult())
622 {
624 }
625 if (hasdndbMult())
626 {
628 }
629 if (hasdxdbMult())
630 {
632 }
633 if (hasdxdbDirectMult())
634 {
636 }
638 {
640 {
642 }
643 }
644 if (hasDivDxDbMult())
645 {
646 divDxDbMultPtr_() ==
647 dimensionedScalar(divDxDbMultPtr_().dimensions(), Zero);
648 }
649 if (hasGradDxDbMult())
650 {
652 dimensionedTensor(gradDxDbMultPtr_().dimensions(), Zero);
653 }
654
655 nullified_ = true;
656 }
657}
658
659
660bool objective::write(const bool valid) const
661{
662 if (Pstream::master())
663 {
664 // File is opened only upon invocation of the write function
665 // in order to avoid various instantiations of the same objective
666 // opening the same file
668 {
671 ios_base::fmtflags flags = file.flags();
672 flags |= std::cout.left;
673 file.flags(flags);
674 if (target_.valid())
675 {
676 file<< setw(width_) << "#target" << " "
677 << setw(width_) << target_() << endl;
678 }
679 if (normalize_)
680 {
681 file<< setw(width_) << "#normFactor " << " "
682 << setw(width_) << normFactor_() << endl;
683 }
685 file<< setw(4) << "#" << " ";
686 file<< setw(width_) << "J" << " ";
687 file<< setw(width_) << "JCycle" << " ";
689 file<< endl;
690 }
691
693 file<< setw(4) << mesh_.time().value() << " ";
694 file<< setw(width_) << J_ << " ";
695 file<< setw(width_) << JCycle() << " ";
697 file<< endl;
698 }
699
700 return true;
701}
702
703
705{
706 if (Pstream::master())
707 {
708 // File is opened only upon invocation of the write function
709 // in order to avoid various instantiations of the same objective
710 // opening the same file
712 {
714 }
715
716 instantValueFilePtr_() << mesh_.time().value() << tab << J_ << endl;
717 }
718}
719
720
722{
723 if (Pstream::master())
724 {
726 {
728 }
729 }
730}
731
732
734{
735 if (Pstream::master())
736 {
737 // Write mean value if necessary
738 // Covers both steady and unsteady runs
739 if
740 (
743 )
744 {
745 // File is opened only upon invocation of the write function
746 // in order to avoid various instantiations of the same objective
747 // opening the same file
749 {
751 }
752
754 << mesh_.time().value() << tab << JMean_ << endl;
755 }
756 }
757}
758
759
760bool objective::writeData(Ostream& os) const
761{
762 os.writeEntry("JMean", JMean_);
763 if (normFactor_)
764 {
765 os.writeEntry("normFactor", normFactor_());
766 }
767 return os.good();
768}
769
770
772{
773 // Does nothing in base
774}
775
776
778{
779 // Does nothing in base
780}
781
782
784{
785 // Does nothing in base
786}
787
788
789// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
790
791} // End namespace Foam
792
793// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
bool found
Generic GeometricBoundaryField class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:506
An IOstream is an abstract base class for all input/output systems; be they streams,...
Definition: IOstream.H:82
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:233
Output to file stream, using an OSstream.
Definition: OFstream.H:57
virtual ios_base::fmtflags flags() const
Get stream flags.
Definition: OSstream.C:288
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
fileName globalPath() const
Return global path for the case.
Definition: TimePathsI.H:80
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:55
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
bool valid() const noexcept
Identical to good(), or bool operator.
Definition: autoPtr.H:272
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
const Type & value() const
Return const reference to value.
A class for handling file names.
Definition: fileName.H:76
virtual bool write()
Write the output fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
localIOdictionary is derived from IOdictionary but excludes parallel master reading.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:62
virtual void writeInstantaneousSeparator() const
Definition: objective.C:721
virtual void writeInstantaneousValue() const
Write objective function history at each primal solver iteration.
Definition: objective.C:704
scalar JCycle() const
Definition: objective.C:248
autoPtr< OFstream > instantValueFilePtr_
Definition: objective.H:145
autoPtr< boundaryVectorField > bdJdbPtr_
Term from material derivative.
Definition: objective.H:106
const boundaryVectorField & dSdbMultiplier()
Multiplier of delta(n dS)/delta b for all patches.
Definition: objective.C:515
const volTensorField & gradDxDbMultiplier()
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
Definition: objective.C:589
void accumulateJMean()
Accumulate contribution for the mean objective value.
Definition: objective.C:302
const fvMesh & mesh_
Definition: objective.H:67
scalar JMean_
Average objective value.
Definition: objective.H:80
bool hasdSdbMult() const
Definition: objectiveI.H:51
bool hasdndbMult() const
Definition: objectiveI.H:57
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
Definition: objective.H:109
const word objectiveName_
Definition: objective.H:71
bool hasBoundaryEdgeContribution() const
Definition: objectiveI.H:75
bool hasGradDxDbMult() const
Definition: objectiveI.H:87
void setMeanValueFilePtr() const
Set the output file ptr for the mean value.
Definition: objective.C:80
virtual void writeMeanValue() const
Write mean objective function history.
Definition: objective.C:733
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:141
autoPtr< OFstream > meanValueFilePtr_
Definition: objective.H:149
virtual void addHeaderColumns() const
Write headers for additional columns.
Definition: objective.C:777
scalar weight() const
Return the objective function weight.
Definition: objective.C:324
autoPtr< boundaryVectorField > bdxdbDirectMultPtr_
Definition: objective.H:120
bool hasIntegrationStartTime() const
Definition: objectiveI.H:93
virtual void nullify()
Nullify adjoint contributions.
Definition: objective.C:609
dictionary dict_
Definition: objective.H:68
autoPtr< scalar > integrationStartTimePtr_
Objective integration start and end times (for unsteady flows)
Definition: objective.H:94
const boundaryVectorField & boundarydJdb()
Contribution to surface sensitivities for all patches.
Definition: objective.C:505
void setObjectiveFilePtr() const
Set the output file ptr.
Definition: objective.C:59
autoPtr< volScalarField > dJdbPtr_
Contribution to field sensitivity derivatives.
Definition: objective.H:100
autoPtr< scalar > normFactor_
Normalization factor.
Definition: objective.H:86
scalar weight_
Objective weight.
Definition: objective.H:83
bool computeMeanFields_
Definition: objective.H:72
const volScalarField & divDxDbMultiplier()
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
Definition: objective.C:568
const boundaryVectorField & dxdbDirectMultiplier()
Multiplier of delta(x)/delta b for all patches.
Definition: objective.C:545
unsigned int width_
Default width of entries when writing in the objective files.
Definition: objective.H:152
virtual scalar J()=0
Return the instantaneous objective function value.
autoPtr< scalar > integrationEndTimePtr_
Definition: objective.H:95
virtual void addHeaderInfo() const
Write any information that needs to go the header of the file.
Definition: objective.C:771
bool hasDivDxDbMult() const
Definition: objectiveI.H:81
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition: objective.H:132
const boundaryVectorField & dndbMultiplier()
Multiplier of delta(n dS)/delta b for all patches.
Definition: objective.C:525
scalar J_
Objective function value and weight.
Definition: objective.H:77
bool hasdxdbDirectMult() const
Definition: objectiveI.H:69
bool hasdxdbMult() const
Definition: objectiveI.H:63
void setInstantValueFilePtr() const
Set the output file ptr for the instantaneous value.
Definition: objective.C:68
void incrementIntegrationTimes(const scalar timeSpan)
Increment integration times.
Definition: objective.C:403
virtual void doNormalization()
Normalize all fields allocated by the objective.
Definition: objective.C:336
bool normalize() const
Is the objective normalized.
Definition: objective.C:330
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:94
const volScalarField & dJdb()
Contribution to field sensitivities.
Definition: objective.C:419
bool hasBoundarydJdb() const
Definition: objectiveI.H:45
bool hasIntegrationEndTime() const
Definition: objectiveI.H:99
virtual void updateNormalizationFactor()
Definition: objective.C:276
fileName objFunctionFolder_
Output file variables.
Definition: objective.H:138
bool hasdJdb() const
Definition: objectiveI.H:39
autoPtr< scalar > target_
Target value, in case the objective is used as a constraint.
Definition: objective.H:91
autoPtr< boundaryVectorField > bdndbMultPtr_
Term multiplying delta(n)/delta b.
Definition: objective.H:112
autoPtr< boundaryVectorField > bdxdbMultPtr_
Term multiplying delta(x)/delta b at the boundary.
Definition: objective.H:115
const word adjointSolverName_
Definition: objective.H:69
const boundaryVectorField & dxdbMultiplier()
Multiplier of delta(x)/delta b for all patches.
Definition: objective.C:535
virtual void addColumnValues() const
Write information to additional columns.
Definition: objective.C:783
const vectorField3 & boundaryEdgeMultiplier()
Multiplier located at patch boundary edges.
Definition: objective.C:555
bool isWithinIntegrationTime() const
Check whether this is an objective integration time.
Definition: objective.C:382
autoPtr< vectorField3 > bEdgeContribution_
Definition: objective.H:125
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
Definition: objective.H:135
Base class for solver control classes.
Definition: solverControl.H:52
bool doAverageIter() const
label & averageIter()
Return average iteration index reference.
splitCell * master() const
Definition: splitCell.H:113
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
rDeltaTY field()
bool
Definition: EEqn.H:20
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const word dictName("faMeshDefinition")
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
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:515
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensionSet pow2(const dimensionSet &ds)
Definition: dimensionSet.C:369
constexpr char tab
The tab '\t' character(0x09)
Definition: Ostream.H:52
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
A non-counting (dummy) refCount.
Definition: refCount.H:59