Time.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "Time.H"
30#include "PstreamReduceOps.H"
31#include "argList.H"
32#include "HashSet.H"
33#include "profiling.H"
34#include "IOdictionary.H"
35#include "registerSwitch.H"
36#include <sstream>
37
38// * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
39
40namespace Foam
41{
43}
44
45const Foam::Enum
46<
48>
50({
51 { stopAtControls::saEndTime, "endTime" },
52 { stopAtControls::saNoWriteNow, "noWriteNow" },
53 { stopAtControls::saWriteNow, "writeNow" },
54 { stopAtControls::saNextWrite, "nextWrite" },
55 // Leave saUnknown untabulated - fallback to flag unknown settings
56});
57
58
59const Foam::Enum
60<
62>
64({
65 { writeControls::wcNone, "none" },
66 { writeControls::wcTimeStep, "timeStep" },
67 { writeControls::wcRunTime, "runTime" },
68 { writeControls::wcAdjustableRunTime, "adjustable" },
69 { writeControls::wcAdjustableRunTime, "adjustableRunTime" },
70 { writeControls::wcClockTime, "clockTime" },
71 { writeControls::wcCpuTime, "cpuTime" },
72 // Leave wcUnknown untabulated - fallback to flag unknown settings
73});
74
75
77
79
80const int Foam::Time::maxPrecision_(3 - log10(SMALL));
81
83
85(
86 Foam::debug::infoSwitch("printExecutionFormat", 0)
87);
88
90(
91 "printExecutionFormat",
92 int,
94);
95
96
97// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
98
100{
101 bool adjustTime = false;
102 scalar timeToNextWrite = VGREAT;
103
105 {
106 adjustTime = true;
107 timeToNextWrite = max
108 (
109 0.0,
111 );
112 }
113
114 if (adjustTime)
115 {
116 scalar nSteps = timeToNextWrite/deltaT_;
117
118 // For tiny deltaT the label can overflow!
119 if (nSteps < scalar(labelMax))
120 {
121 // nSteps can be < 1 so make sure at least 1
122 label nStepsToNextWrite = max(1, round(nSteps));
123
124 scalar newDeltaT = timeToNextWrite/nStepsToNextWrite;
125
126 // Control the increase of the time step to within a factor of 2
127 // and the decrease within a factor of 5.
128 if (newDeltaT >= deltaT_)
129 {
130 deltaT_ = min(newDeltaT, 2.0*deltaT_);
131 }
132 else
133 {
134 deltaT_ = max(newDeltaT, 0.2*deltaT_);
135 }
136 }
137 }
138
139 functionObjects_.adjustTimeStep();
140}
141
142
144{
145 // default is to resume calculation from "latestTime"
146 const word startFrom = controlDict_.getOrDefault<word>
147 (
148 "startFrom",
149 "latestTime"
150 );
151
152 if (startFrom == "startTime")
153 {
154 controlDict_.readEntry("startTime", startTime_);
155 }
156 else
157 {
158 // Search directory for valid time directories
159 instantList timeDirs = findTimes(path(), constant());
160
161 const label nTimes = timeDirs.size();
162
163 if (startFrom == "firstTime")
164 {
165 if (nTimes > 1 && timeDirs.first().name() == constant())
166 {
167 startTime_ = timeDirs[1].value();
168 }
169 else if (nTimes)
170 {
171 startTime_ = timeDirs.first().value();
172 }
173 }
174 else if (startFrom == "latestTime")
175 {
176 if (nTimes)
177 {
178 startTime_ = timeDirs.last().value();
179 }
180 }
181 else
182 {
183 FatalIOErrorInFunction(controlDict_)
184 << "expected startTime, firstTime or latestTime"
185 << " found '" << startFrom << "'"
186 << exit(FatalIOError);
187 }
188 }
189
190 setTime(startTime_, 0);
191
192 readDict();
193 deltaTSave_ = deltaT_;
194 deltaT0_ = deltaT_;
195
196 // Check if time directory exists
197 // If not increase time precision to see if it is formatted differently.
198 // Note: do not use raw fileHandler().exists(...) since does not check
199 // alternative processorsDDD directories naming
200 if (fileHandler().filePath(timePath()).empty())
201 {
202 int oldPrecision = precision_;
203 int requiredPrecision = -1;
204 bool found = false;
206 for
207 (
208 precision_ = maxPrecision_;
209 precision_ > oldPrecision;
210 --precision_
211 )
212 {
213 // Update the time formatting
214 setTime(startTime_, 0);
215
216 word newTime(timeName());
217 if (newTime == oldTime)
218 {
219 break;
220 }
221 oldTime = newTime;
222
223 // Check the existence of the time directory with the new format
224 //found = fileHandler().exists(timePath(), false);
225 const fileName dirName(fileHandler().filePath(timePath()));
226 found = !dirName.empty();
227
228 if (found)
229 {
230 requiredPrecision = precision_;
231 }
232 }
233
234 if (requiredPrecision > 0)
235 {
236 // Update the time precision
237 precision_ = requiredPrecision;
238
239 // Update the time formatting
240 setTime(startTime_, 0);
241
243 << "Increasing the timePrecision from " << oldPrecision
244 << " to " << precision_
245 << " to support the formatting of the current time directory "
246 << timeName() << nl << endl;
247 }
248 else
249 {
250 // Could not find time directory so assume it is not present
251 precision_ = oldPrecision;
252
253 // Revert the time formatting
254 setTime(startTime_, 0);
255 }
256 }
257
258 if (Pstream::parRun())
259 {
260 scalar sumStartTime = startTime_;
261 reduce(sumStartTime, sumOp<scalar>());
262 if
263 (
264 mag(Pstream::nProcs()*startTime_ - sumStartTime)
265 > Pstream::nProcs()*deltaT_/10.0
266 )
267 {
268 FatalIOErrorInFunction(controlDict_)
269 << "Start time is not the same for all processors" << nl
270 << "processor " << Pstream::myProcNo() << " has startTime "
271 << startTime_ << exit(FatalIOError);
272 }
273 }
274
275 IOdictionary timeDict
276 (
278 (
279 "time",
280 timeName(),
281 "uniform",
282 *this,
285 false
286 )
287 );
288
289 // Read and set the deltaT only if time-step adjustment is active
290 // otherwise use the deltaT from the controlDict
291 if (controlDict_.getOrDefault("adjustTimeStep", false))
292 {
293 if (timeDict.readIfPresent("deltaT", deltaT_))
294 {
295 deltaTSave_ = deltaT_;
296 deltaT0_ = deltaT_;
297 }
298 }
299
300 timeDict.readIfPresent("deltaT0", deltaT0_);
301
302 if (timeDict.readIfPresent("index", startTimeIndex_))
303 {
304 timeIndex_ = startTimeIndex_;
305 }
306
307
308 // Check if values stored in time dictionary are consistent
309
310 // 1. Based on time name
311 bool checkValue = true;
312
313 string storedTimeName;
314 if (timeDict.readIfPresent("name", storedTimeName))
315 {
316 if (storedTimeName == timeName())
317 {
318 // Same time. No need to check stored value
319 checkValue = false;
320 }
321 }
322
323 // 2. Based on time value
324 // (consistent up to the current time writing precision so it won't
325 // trigger if we just change the write precision)
326 if (checkValue)
327 {
328 scalar storedTimeValue;
329 if (timeDict.readIfPresent("value", storedTimeValue))
330 {
331 word storedTimeName(timeName(storedTimeValue));
332
333 if (storedTimeName != timeName())
334 {
335 IOWarningInFunction(timeDict)
336 << "Time read from time dictionary " << storedTimeName
337 << " differs from actual time " << timeName() << '.' << nl
338 << " This may cause unexpected database behaviour."
339 << " If you are not interested" << nl
340 << " in preserving time state delete"
341 << " the time dictionary."
342 << endl;
343 }
344 }
345 }
346}
347
348
349void Foam::Time::setMonitoring(const bool forceProfiling)
350{
351 const dictionary* profilingDict = controlDict_.findDict("profiling");
352 if (!profilingDict)
353 {
354 // ... or from etc/controlDict
355 profilingDict = debug::controlDict().findDict("profiling");
356 }
357
358 // initialize profiling on request
359 // otherwise rely on profiling entry within controlDict
360 // and skip if 'active' keyword is explicitly set to false
361 if (forceProfiling)
362 {
363 profiling::initialize
364 (
366 (
367 "profiling",
368 timeName(),
369 "uniform",
370 *this,
373 ),
374 *this
375 );
376 }
377 else if
378 (
379 profilingDict
380 && profilingDict->getOrDefault("active", true)
381 )
382 {
383 profiling::initialize
384 (
385 *profilingDict,
387 (
388 "profiling",
389 timeName(),
390 "uniform",
391 *this,
394 ),
395 *this
396 );
397 }
398
399 // Time objects not registered so do like objectRegistry::checkIn ourselves.
400 if (runTimeModifiable_)
401 {
402 // Monitor all files that controlDict depends on
403 fileHandler().addWatches(controlDict_, controlDict_.files());
404 }
405
406 // Clear dependent files - not needed now
407 controlDict_.files().clear();
408}
409
410
411// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
412
414(
415 const word& ctrlDictName,
416 const fileName& rootPath,
417 const fileName& caseName,
418 const word& systemName,
419 const word& constantName,
420 const bool enableFunctionObjects,
421 const bool enableLibs
422)
423:
425 (
426 rootPath,
427 caseName,
428 systemName,
429 constantName
430 ),
431
432 objectRegistry(*this),
433
434 loopProfiling_(nullptr),
435 libs_(),
436
437 controlDict_
438 (
440 (
441 ctrlDictName,
442 system(),
443 *this,
444 IOobject::MUST_READ_IF_MODIFIED,
445 IOobject::NO_WRITE,
446 false
447 )
448 ),
449
450 startTimeIndex_(0),
451 startTime_(0),
452 endTime_(0),
453
454 stopAt_(saEndTime),
455 writeControl_(wcTimeStep),
456 writeInterval_(GREAT),
457 purgeWrite_(0),
458 subCycling_(0),
459 writeOnce_(false),
460 sigWriteNow_(*this, true),
461 sigStopAtWriteNow_(*this, true),
462 writeStreamOption_(IOstream::ASCII),
463 graphFormat_("raw"),
464 runTimeModifiable_(false),
465 functionObjects_(*this, false)
466{
467 if (enableFunctionObjects)
468 {
469 functionObjects_.on();
470 }
471
472 if (enableLibs)
473 {
474 libs_.open("libs", controlDict_);
475 }
476
477 // Explicitly set read flags on objectRegistry so anything constructed
478 // from it reads as well (e.g. fvSolution).
480
481 setControls();
483}
484
485
487(
488 const word& ctrlDictName,
489 const argList& args,
490 const word& systemName,
491 const word& constantName,
492 const bool enableFunctionObjects,
493 const bool enableLibs
494)
495:
496 TimePaths(args, systemName, constantName),
497
498 objectRegistry(*this),
499
500 loopProfiling_(nullptr),
501 libs_(),
502
503 controlDict_
504 (
506 (
507 ctrlDictName,
508 system(),
509 *this,
510 IOobject::MUST_READ_IF_MODIFIED,
511 IOobject::NO_WRITE,
512 false
513 )
514 ),
515
516 startTimeIndex_(0),
517 startTime_(0),
518 endTime_(0),
519
520 stopAt_(saEndTime),
521 writeControl_(wcTimeStep),
522 writeInterval_(GREAT),
523 purgeWrite_(0),
524 subCycling_(0),
525 writeOnce_(false),
526 sigWriteNow_(*this, true),
527 sigStopAtWriteNow_(*this, true),
528 writeStreamOption_(IOstream::ASCII),
529 graphFormat_("raw"),
530 runTimeModifiable_(false),
531 functionObjects_(*this, false)
532{
533 // Functions
534 //
535 // * '-withFunctionObjects' exists and used = enable
536 // * '-noFunctionObjects' exists and used = disable
537 // * default: no functions if there is no way to enable/disable them
538 if
539 (
540 argList::validOptions.found("withFunctionObjects")
541 ? args.found("withFunctionObjects")
542 : argList::validOptions.found("noFunctionObjects")
543 ? !args.found("noFunctionObjects")
544 : false
545 )
546 {
547 if (enableFunctionObjects)
548 {
549 functionObjects_.on();
550 }
551 }
552
553 // Libraries
554 //
555 // * enable by default unless '-no-libs' option was used
556 if (enableLibs && !args.found("no-libs"))
557 {
558 libs_.open("libs", controlDict_);
559 }
560
561 // Explicitly set read flags on objectRegistry so anything constructed
562 // from it reads as well (e.g. fvSolution).
564
565 setControls();
566
567 // '-profiling' = force profiling, ignore controlDict entry
568 setMonitoring(args.found("profiling"));
569}
570
571
573(
574 const dictionary& dict,
575 const fileName& rootPath,
576 const fileName& caseName,
577 const word& systemName,
578 const word& constantName,
579 const bool enableFunctionObjects,
580 const bool enableLibs
581)
582:
584 (
585 rootPath,
586 caseName,
587 systemName,
588 constantName
589 ),
590
591 objectRegistry(*this),
592
593 loopProfiling_(nullptr),
594 libs_(),
595
596 controlDict_
597 (
599 (
600 controlDictName,
601 system(),
602 *this,
603 IOobject::NO_READ,
604 IOobject::NO_WRITE,
605 false
606 ),
607 dict
608 ),
609
610 startTimeIndex_(0),
611 startTime_(0),
612 endTime_(0),
613
614 stopAt_(saEndTime),
615 writeControl_(wcTimeStep),
616 writeInterval_(GREAT),
617 purgeWrite_(0),
618 subCycling_(0),
619 writeOnce_(false),
620 sigWriteNow_(*this, true),
621 sigStopAtWriteNow_(*this, true),
622 writeStreamOption_(IOstream::ASCII),
623 graphFormat_("raw"),
624 runTimeModifiable_(false),
625 functionObjects_(*this, false)
626{
627 if (enableFunctionObjects)
628 {
629 functionObjects_.on();
630 }
631
632 if (enableLibs)
633 {
634 libs_.open("libs", controlDict_);
635 }
636
637
638 // Explicitly set read flags on objectRegistry so anything constructed
639 // from it reads as well (e.g. fvSolution).
641
642 // Since could not construct regIOobject with setting:
644
645 setControls();
647}
648
649
651(
652 const fileName& rootPath,
653 const fileName& caseName,
654 const word& systemName,
655 const word& constantName,
656 const bool enableFunctionObjects,
657 const bool enableLibs
658)
659:
661 (
662 rootPath,
663 caseName,
664 systemName,
665 constantName
666 ),
667
668 objectRegistry(*this),
669
670 loopProfiling_(nullptr),
671 libs_(),
672
673 controlDict_
674 (
676 (
677 controlDictName,
678 system(),
679 *this,
680 IOobject::NO_READ,
681 IOobject::NO_WRITE,
682 false
683 )
684 ),
685
686 startTimeIndex_(0),
687 startTime_(0),
688 endTime_(0),
689
690 stopAt_(saEndTime),
691 writeControl_(wcTimeStep),
692 writeInterval_(GREAT),
693 purgeWrite_(0),
694 subCycling_(0),
695 writeOnce_(false),
696 writeStreamOption_(IOstream::ASCII),
697 graphFormat_("raw"),
698 runTimeModifiable_(false),
699 functionObjects_(*this, false)
700{
701 if (enableFunctionObjects)
702 {
703 functionObjects_.on();
704 }
705
706 if (enableLibs)
707 {
708 libs_.open("libs", controlDict_);
709 }
710
711 setMonitoring(); // for profiling etc
712}
713
714
715// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
716
718{
719 return
721 (
722 fileName("."), // root-path
723 fileName("."), // case-name
724 false, // No enableFunctionObjects
725 false // No enableLibs
726 );
727}
728
729
731{
732 return
734 (
735 caseDir.path(), // root-path
736 caseDir.name(), // case-name
737 false, // No enableFunctionObjects
738 false // No enableLibs
739 );
740}
741
742
744{
745 return
747 (
749 args,
750 false, // No enableFunctionObjects
751 false // No enableLibs
752 );
753}
754
755
756// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
757
759{
760 loopProfiling_.reset(nullptr);
761
762 forAllReverse(controlDict_.watchIndices(), i)
763 {
764 fileHandler().removeWatch(controlDict_.watchIndices()[i]);
765 }
766
767 // Destroy function objects first
768 functionObjects_.clear();
769
770 // Clean up profiling
771 profiling::stop(*this);
772
773 // Ensure all owned objects are also cleaned up now
775}
776
777
778// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
779
780Foam::word Foam::Time::timeName(const scalar t, const int precision)
781{
782 std::ostringstream buf;
783 buf.setf(ios_base::fmtflags(format_), ios_base::floatfield);
784 buf.precision(precision);
785 buf << t;
786 return buf.str();
787}
788
789
791{
793}
794
795
797(
798 const fileName& dir,
799 const word& name,
800 const IOobject::readOption rOpt,
801 const word& stopInstance
802) const
803{
804 // Note: name might be empty!
805 IOobject startIO(name, timeName(), dir, *this, rOpt);
806
808 (
809 fileHandler().findInstance
810 (
811 startIO,
812 timeOutputValue(),
813 stopInstance
814 )
815 );
816 return io.instance();
817}
818
819
821(
822 const fileName& directory,
823 const instant& t
824) const
825{
826 // Simplified version: use findTimes (readDir + sort). The expensive
827 // bit is the readDir, not the sorting. Tbd: avoid calling findInstancePath
828 // from filePath.
829
830 instantList timeDirs = findTimes(path(), constant());
831 // Note:
832 // - times will include constant (with value 0) as first element.
833 // For backwards compatibility make sure to find 0 in preference
834 // to constant.
835 // - list is sorted so could use binary search
836
837 forAllReverse(timeDirs, i)
838 {
839 if (t.equal(timeDirs[i].value()))
840 {
841 return timeDirs[i].name();
842 }
843 }
844
845 return word::null;
846}
847
848
850{
851 return findInstancePath(path(), t);
852}
853
854
855Foam::label Foam::Time::startTimeIndex() const
856{
857 return startTimeIndex_;
858}
859
860
862{
863 return dimensionedScalar("startTime", dimTime, startTime_);
864}
865
866
868{
869 return dimensionedScalar("endTime", dimTime, endTime_);
870}
871
872
874{
875 return stopAt_;
876}
877
878
879bool Foam::Time::run() const
880{
881 loopProfiling_.reset(nullptr);
882
883 bool isRunning = value() < (endTime_ - 0.5*deltaT_);
884
885 if (!subCycling_)
886 {
887 // Only execute when the condition is no longer true
888 // ie, when exiting the control loop
889 if (!isRunning && timeIndex_ != startTimeIndex_)
890 {
891 // Ensure functionObjects execute on last time step
892 // (and hence write uptodate functionObjectProperties)
893 {
894 addProfiling(fo, "functionObjects.execute()");
895 functionObjects_.execute();
896 }
897 {
898 addProfiling(fo, "functionObjects.end()");
899 functionObjects_.end();
900 }
901 }
902 }
903
904 if (isRunning)
905 {
906 if (!subCycling_)
907 {
908 const_cast<Time&>(*this).readModifiedObjects();
909
910 if (timeIndex_ == startTimeIndex_)
911 {
912 addProfiling(functionObjects, "functionObjects.start()");
913 functionObjects_.start();
914 }
915 else
916 {
917 addProfiling(functionObjects, "functionObjects.execute()");
918 functionObjects_.execute();
919 }
920
921 // Check if the execution of functionObjects require re-reading
922 // any files. This moves effect of e.g. 'timeActivatedFileUpdate'
923 // one time step forward. Note that we cannot call
924 // readModifiedObjects from within timeActivatedFileUpdate since
925 // it might re-read the functionObjects themselves (and delete
926 // the timeActivatedFileUpdate one)
927 if (functionObjects_.filesModified())
928 {
929 const_cast<Time&>(*this).readModifiedObjects();
930 }
931 }
932
933 // Update the "is-running" status following the
934 // possible side-effects from functionObjects
935 isRunning = value() < (endTime_ - 0.5*deltaT_);
936
937 // (re)trigger profiling
938 if (profiling::active())
939 {
940 loopProfiling_.reset
941 (
942 new profilingTrigger("time.run() " + objectRegistry::name())
943 );
944 }
945 }
946
947 return isRunning;
948}
949
950
952{
953 const bool isRunning = run();
954
955 if (isRunning)
956 {
957 operator++();
958 }
959
960 return isRunning;
961}
962
963
964bool Foam::Time::end() const
965{
966 return value() > (endTime_ + 0.5*deltaT_);
967}
968
969
970bool Foam::Time::stopAt(const stopAtControls stopCtrl) const
971{
972 if (stopCtrl == stopAtControls::saUnknown)
973 {
974 return false;
975 }
976
977 const bool changed = (stopAt_ != stopCtrl);
978 stopAt_ = stopCtrl;
979 endTime_ = GREAT;
980
981 // Adjust endTime
982 if (stopCtrl == stopAtControls::saEndTime)
983 {
984 controlDict_.readEntry("endTime", endTime_);
985 }
986
987 return changed;
988}
989
990
992{
993 return controlDict_.getOrDefault("adjustTimeStep", false);
994}
995
996
998{
999 value() = t.value();
1000 dimensionedScalar::name() = t.dimensionedScalar::name();
1001 timeIndex_ = t.timeIndex_;
1002 fileHandler().setTime(*this);
1003}
1004
1005
1006void Foam::Time::setTime(const instant& inst, const label newIndex)
1007{
1008 value() = inst.value();
1009 dimensionedScalar::name() = inst.name();
1010 timeIndex_ = newIndex;
1011
1012 IOdictionary timeDict
1013 (
1014 IOobject
1015 (
1016 "time",
1017 timeName(),
1018 "uniform",
1019 *this,
1022 false
1023 )
1024 );
1025
1026 timeDict.readIfPresent("deltaT", deltaT_);
1027 timeDict.readIfPresent("deltaT0", deltaT0_);
1028 timeDict.readIfPresent("index", timeIndex_);
1029 fileHandler().setTime(*this);
1030}
1031
1032
1033void Foam::Time::setTime(const dimensionedScalar& newTime, const label newIndex)
1034{
1035 setTime(newTime.value(), newIndex);
1036}
1037
1038
1039void Foam::Time::setTime(const scalar newTime, const label newIndex)
1040{
1041 value() = newTime;
1042 dimensionedScalar::name() = timeName(timeToUserTime(newTime));
1043 timeIndex_ = newIndex;
1044 fileHandler().setTime(*this);
1045}
1046
1047
1049{
1050 setEndTime(endTime.value());
1051}
1052
1053
1054void Foam::Time::setEndTime(const scalar endTime)
1055{
1056 endTime_ = endTime;
1057}
1058
1059
1061(
1062 const dimensionedScalar& deltaT,
1063 const bool adjust
1064)
1065{
1066 setDeltaT(deltaT.value(), adjust);
1067}
1068
1069
1070void Foam::Time::setDeltaT(const scalar deltaT, const bool adjust)
1071{
1072 deltaT_ = deltaT;
1073 deltaTchanged_ = true;
1074
1075 if (adjust)
1076 {
1077 adjustDeltaT();
1078 }
1079}
1080
1081
1083{
1084 #ifdef FULLDEBUG
1085 if (prevTimeState_)
1086 {
1088 << "previous time state already set" << nl
1089 << exit(FatalError);
1090 }
1091 #endif
1092
1093 prevTimeState_.reset(new TimeState(*this));
1094
1095 setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
1096 deltaT_ /= nSubCycles;
1097 deltaT0_ /= nSubCycles;
1098 deltaTSave_ = deltaT0_;
1099
1100 subCycling_ = nSubCycles;
1101
1102 return prevTimeState();
1103}
1104
1105
1106void Foam::Time::subCycleIndex(const label index)
1107{
1108 // Only permit adjustment if sub-cycling was already active
1109 // and if the index is valid (positive, non-zero).
1110 // This avoids potential mixups for deleting.
1111
1112 if (subCycling_ && index > 0)
1113 {
1114 subCycling_ = index;
1115 }
1116}
1117
1118
1120{
1121 if (subCycling_)
1122 {
1123 TimeState::operator=(prevTimeState());
1124 prevTimeState_.reset(nullptr);
1125 }
1126
1127 subCycling_ = 0;
1128}
1129
1130
1131// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1132
1134{
1135 return operator+=(deltaT.value());
1136}
1137
1138
1140{
1141 setDeltaT(deltaT);
1142 return operator++();
1143}
1144
1145
1147{
1148 deltaT0_ = deltaTSave_;
1149 deltaTSave_ = deltaT_;
1150
1151 // Save old time value and name
1152 const scalar oldTimeValue = timeToUserTime(value());
1153 const word oldTimeName = dimensionedScalar::name();
1154
1155 // Increment time
1156 setTime(value() + deltaT_, timeIndex_ + 1);
1157
1158 if (!subCycling_)
1159 {
1160 // If the time is very close to zero reset to zero
1161 if (mag(value()) < 10*SMALL*deltaT_)
1162 {
1163 setTime(0.0, timeIndex_);
1164 }
1165
1166 if (sigStopAtWriteNow_.active() || sigWriteNow_.active())
1167 {
1168 // A signal might have been sent on one processor only
1169 // Reduce so all decide the same.
1170
1171 label flag = 0;
1172 if (sigStopAtWriteNow_.active() && stopAt_ == saWriteNow)
1173 {
1174 flag += 1;
1175 }
1176 if (sigWriteNow_.active() && writeOnce_)
1177 {
1178 flag += 2;
1179 }
1180 reduce(flag, maxOp<label>());
1181
1182 if (flag & 1)
1183 {
1184 stopAt_ = saWriteNow;
1185 }
1186 if (flag & 2)
1187 {
1188 writeOnce_ = true;
1189 }
1190 }
1191
1192 writeTime_ = false;
1193
1194 switch (writeControl_)
1195 {
1196 case wcNone:
1197 case wcUnknown:
1198 break;
1199
1200 case wcTimeStep:
1201 writeTime_ = !(timeIndex_ % label(writeInterval_));
1202 break;
1203
1204 case wcRunTime:
1205 case wcAdjustableRunTime:
1206 {
1207 const label writeIndex = label
1208 (
1209 ((value() - startTime_) + 0.5*deltaT_)
1210 / writeInterval_
1211 );
1212
1213 if (writeIndex > writeTimeIndex_)
1214 {
1215 writeTime_ = true;
1216 writeTimeIndex_ = writeIndex;
1217 }
1218 }
1219 break;
1220
1221 case wcCpuTime:
1222 {
1223 const label writeIndex = label
1224 (
1225 returnReduce(elapsedCpuTime(), maxOp<double>())
1226 / writeInterval_
1227 );
1228 if (writeIndex > writeTimeIndex_)
1229 {
1230 writeTime_ = true;
1231 writeTimeIndex_ = writeIndex;
1232 }
1233 }
1234 break;
1235
1236 case wcClockTime:
1237 {
1238 const label writeIndex = label
1239 (
1240 returnReduce(elapsedClockTime(), maxOp<double>())
1241 / writeInterval_
1242 );
1243 if (writeIndex > writeTimeIndex_)
1244 {
1245 writeTime_ = true;
1246 writeTimeIndex_ = writeIndex;
1247 }
1248 }
1249 break;
1250 }
1251
1252
1253 // Check if endTime needs adjustment to stop at the next run()/end()
1254 if (!end())
1255 {
1256 if (stopAt_ == saNoWriteNow)
1257 {
1258 endTime_ = value();
1259 }
1260 else if (stopAt_ == saWriteNow)
1261 {
1262 endTime_ = value();
1263 writeTime_ = true;
1264 }
1265 else if (stopAt_ == saNextWrite && writeTime_ == true)
1266 {
1267 endTime_ = value();
1268 }
1269 }
1270
1271 // Override writeTime if one-shot writing
1272 if (writeOnce_)
1273 {
1274 writeTime_ = true;
1275 writeOnce_ = false;
1276 }
1277
1278 // Adjust the precision of the time directory name if necessary
1279 if (writeTime_)
1280 {
1281 // User-time equivalent of deltaT
1282 const scalar userDeltaT =
1283 timeToUserTime(value()) - timeToUserTime(value() - deltaT_);
1284
1285 // Tolerance used when testing time equivalence
1286 const scalar timeTol =
1287 max(min(pow(scalar(10), -precision_), 0.1*userDeltaT), SMALL);
1288
1289 // Time value obtained by reading timeName
1290 scalar timeNameValue = -VGREAT;
1291
1292 // Check that new time representation differs from old one
1293 // reinterpretation of the word
1294 if
1295 (
1296 readScalar(dimensionedScalar::name(), timeNameValue)
1297 && (mag(timeNameValue - oldTimeValue - userDeltaT) > timeTol)
1298 )
1299 {
1300 int oldPrecision = precision_;
1301 while
1302 (
1303 precision_ < maxPrecision_
1304 && readScalar(dimensionedScalar::name(), timeNameValue)
1305 && (mag(timeNameValue - oldTimeValue - userDeltaT) > timeTol)
1306 )
1307 {
1308 precision_++;
1309 setTime(value(), timeIndex());
1310 }
1311
1312 if (precision_ != oldPrecision)
1313 {
1315 << "Increased the timePrecision from " << oldPrecision
1316 << " to " << precision_
1317 << " to distinguish between timeNames at time "
1319 << endl;
1320
1321 if (precision_ == maxPrecision_)
1322 {
1323 // Reached maxPrecision limit
1325 << "Current time name " << dimensionedScalar::name()
1326 << nl
1327 << " The maximum time precision has been reached"
1328 " which might result in overwriting previous"
1329 " results."
1330 << endl;
1331 }
1332
1333 // Check if round-off error caused time-reversal
1334 scalar oldTimeNameValue = -VGREAT;
1335 if
1336 (
1337 readScalar(oldTimeName, oldTimeNameValue)
1338 && (
1339 sign(timeNameValue - oldTimeNameValue)
1340 != sign(deltaT_)
1341 )
1342 )
1343 {
1345 << "Current time name " << dimensionedScalar::name()
1346 << " is set to an instance prior to the "
1347 "previous one "
1348 << oldTimeName << nl
1349 << " This might result in temporal "
1350 "discontinuities."
1351 << endl;
1352 }
1353 }
1354 }
1355 }
1356 }
1357
1358 return *this;
1359}
1360
1361
1363{
1364 return operator++();
1365}
1366
1367
1368// ************************************************************************* //
Inter-processor communication reduction functions.
bool found
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
readOption readOpt() const noexcept
The read option.
Definition: IOobjectI.H:164
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
readOption
Enumeration defining the read options.
Definition: IOobject.H:177
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:180
An IOstream is an abstract base class for all input/output systems; be they streams,...
Definition: IOstream.H:82
scalar value() const noexcept
The value (const access)
Definition: Instant.H:118
const T & name() const noexcept
The name/key (const access)
Definition: Instant.H:130
bool equal(scalar val) const noexcept
True if values are equal (includes SMALL for rounding)
Definition: Instant.H:142
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
Address the time paths without using the Time class.
Definition: TimePaths.H:57
The time value with time-stepping information, user-defined remapping, etc.
Definition: TimeState.H:54
scalar deltaT_
Definition: TimeState.H:60
label writeTimeIndex_
Definition: TimeState.H:58
label timeIndex_
Definition: TimeState.H:57
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
virtual bool isAdjustTimeStep() const
Return true if adjustTimeStep is true.
Definition: Time.C:991
virtual bool run() const
Return true if run should continue,.
Definition: Time.C:879
static int precision_
Time directory name precision.
Definition: Time.H:186
virtual void setDeltaT(const dimensionedScalar &deltaT, const bool adjust=true)
Reset time step, normally also calling adjustDeltaT()
Definition: Time.C:1061
static const Enum< stopAtControls > stopAtControlNames
Names for stopAtControls.
Definition: Time.H:119
virtual dimensionedScalar startTime() const
Return start time.
Definition: Time.C:861
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
writeControls
Write control options.
Definition: Time.H:85
@ wcAdjustableRunTime
"adjustable" / "adjustableRunTime"
Definition: Time.H:89
void setMonitoring(const bool forceProfiling=false)
Set file monitoring, profiling, etc.
Definition: Time.C:349
void adjustDeltaT()
Adjust the time step so that writing occurs at the specified time.
Definition: Time.C:99
scalar writeInterval_
Definition: Time.H:159
static int printExecutionFormat_
Style for "ExecutionTime = " output.
Definition: Time.H:127
static const int maxPrecision_
Maximum time directory name precision.
Definition: Time.H:189
word findInstancePath(const fileName &directory, const instant &t) const
Definition: Time.C:821
word findInstance(const fileName &dir, const word &name=word::null, const IOobject::readOption rOpt=IOobject::MUST_READ, const word &stopInstance=word::null) const
Definition: Time.C:797
fmtflags
Supported time directory name formats.
Definition: Time.H:108
static fmtflags format_
Time directory name format.
Definition: Time.H:183
stopAtControls
Definition: Time.H:98
writeControls writeControl_
Definition: Time.H:157
void setControls()
Set the controls from the current controlDict.
Definition: Time.C:143
virtual word timeName() const
Return current time name.
Definition: Time.C:790
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:226
virtual void setEndTime(const dimensionedScalar &endTime)
Reset end time.
Definition: Time.C:1048
static const Enum< writeControls > writeControlNames
Names for writeControls.
Definition: Time.H:116
virtual stopAtControls stopAt() const
Return the stop control information.
Definition: Time.C:873
virtual label startTimeIndex() const
Return start time index.
Definition: Time.C:855
virtual void endSubCycle()
Reset time after sub-cycling back to previous TimeState.
Definition: Time.C:1119
virtual void subCycleIndex(const label index)
Adjust the reported sub-cycle index.
Definition: Time.C:1106
virtual Time & operator++()
Prefix increment,.
Definition: Time.C:1146
scalar startTime_
Definition: Time.H:151
void readModifiedObjects()
Read the objects that have been modified.
Definition: TimeIO.C:467
virtual bool loop()
Return true if run should continue and if so increment time.
Definition: Time.C:951
virtual dimensionedScalar endTime() const
Return end time.
Definition: Time.C:867
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:997
virtual ~Time()
Destructor.
Definition: Time.C:758
virtual bool end() const
Return true if end of run,.
Definition: Time.C:964
virtual Time & operator+=(const dimensionedScalar &deltaT)
Set deltaT to that specified and increment time via operator++()
Definition: Time.C:1133
T & first()
Return the first element of the list.
Definition: UListI.H:202
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
T & last()
Return the last element of the list.
Definition: UListI.H:216
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:124
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:178
static HashTable< string > validOptions
A list of valid options.
Definition: argList.H:212
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
const scalar & value() const
Return const reference to value.
const word & name() const
Return const reference to name.
Particle-size distribution model wherein random samples are drawn from a given arbitrary probability ...
Definition: general.H:176
bool open(bool verbose=true)
void stop()
Stop parsing, freeing the allocated parser.
A class for handling file names.
Definition: fileName.H:76
static std::string path(const std::string &str)
Return directory path name (part before last /)
Definition: fileNameI.H:176
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:199
virtual void addWatches(regIOobject &, const fileNameList &) const
Helper: add watches for list of regIOobjects.
virtual bool removeWatch(const label) const
Remove watch on a file (using handle)
virtual void setTime(const Time &) const
Callback for time change.
bool adjustTimeStep()
Called at the end of Time::adjustDeltaT() if adjustTime is true.
void on()
Switch the function objects on.
void operator=(const ObukhovLength &)=delete
No copy assignment.
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name.
Definition: instant.H:56
Registry of regIOobjects.
void clear()
Clear all entries from the registry.
constant condensation/saturation model.
int myProcNo() const noexcept
Return processor number.
Triggers for starting/stopping code profiling.
static bool active()
True if profiling is allowed and is active.
Definition: profiling.C:111
Perform a subCycleTime on a field.
Definition: subCycle.H:127
A class for handling words, derived from Foam::string.
Definition: word.H:68
static const word null
An empty word.
Definition: word.H:80
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));if(U.nOldTimes()){ volVectorField *Uold=&U.oldTime();volScalarField *Kold=&K.oldTime();*Kold==0.5 *magSqr(*Uold);while(Uold->nOldTimes()) { Uold=&Uold-> oldTime()
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
runTimeSource setTime(sourceTimes[sourceTimeIndex], sourceTimeIndex)
word timeName
Definition: getTimeIndex.H:3
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
#define WarningInFunction
Report a warning using Foam::Warning.
int infoSwitch(const char *name, const int deflt=0)
Lookup info switch or add default value.
Definition: debug.C:231
dictionary & controlDict()
Definition: debug.C:143
Namespace for OpenFOAM.
const fileOperation & fileHandler()
Get current file handler.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar sign(const dimensionedScalar &ds)
int system(const std::string &command, const bool bg=false)
Execute the specified command via the shell.
Definition: MSwindows.C:1158
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
dimensionedScalar log10(const dimensionedScalar &ds)
constexpr label labelMax
Definition: label.H:61
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
IOerror FatalIOError
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
label timeIndex
Definition: getTimeIndex.H:30
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
#define registerInfoSwitch(Name, Type, SwitchVar)
dictionary dict
Foam::argList args(argc, argv)
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:346