timeControlFunctionObject.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) 2016 OpenFOAM Foundation
9  Copyright (C) 2017-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
30 #include "polyMesh.H"
31 #include "mapPolyMesh.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace functionObjects
39 {
40  defineTypeNameAndDebug(timeControl, 0);
41 }
42 }
43 
46 ({
47  { controlMode::TIME, "time" },
48  { controlMode::TRIGGER, "trigger" },
49  { controlMode::TIME_OR_TRIGGER, "timeOrTrigger" },
50  { controlMode::TIME_AND_TRIGGER, "timeAndTrigger" }
51 });
52 
53 
54 // * * * * * * * * * * * * * * * Private Members * * * * * * * * * * * * * * //
55 
56 void Foam::functionObjects::timeControl::readControls()
57 {
58  if (dict_.readIfPresent("timeStart", timeStart_))
59  {
60  timeStart_ = time_.userTimeToTime(timeStart_);
61  }
62  if (dict_.readIfPresent("timeEnd", timeEnd_))
63  {
64  timeEnd_ = time_.userTimeToTime(timeEnd_);
65  }
66 
67  controlMode_ =
68  controlModeNames_.getOrDefault("controlMode", dict_, controlMode::TIME);
69  dict_.readIfPresent("triggerStart", triggerStart_);
70  dict_.readIfPresent("triggerEnd", triggerEnd_);
71 
72  deltaTCoeff_ = GREAT;
73  if (dict_.readIfPresent("deltaTCoeff", deltaTCoeff_))
74  {
75  nStepsToStartTimeChange_ = labelMax;
76  }
77  else
78  {
79  nStepsToStartTimeChange_ = 3;
80  dict_.readIfPresent
81  (
82  "nStepsToStartTimeChange",
83  nStepsToStartTimeChange_
84  );
85  }
86 }
87 
88 
89 bool Foam::functionObjects::timeControl::active() const
90 {
91  label triggeri = time_.functionObjects().triggerIndex();
92 
93  bool inTime =
94  time_.value() >= (timeStart_ - 0.5*time_.deltaTValue())
95  && time_.value() <= (timeEnd_ + 0.5*time_.deltaTValue());
96 
97  bool inTrigger = triggeri >= triggerStart_ && triggeri <= triggerEnd_;
98 
100  << name() << " mode:" << controlModeNames_[controlMode_] << nl
101  << " - time:" << time_.value()
102  << " timeStart:" << timeStart_
103  << " timeEnd:" << timeEnd_
104  << " inTime:" << inTime << nl
105  << " - triggeri:" << triggeri
106  << " triggerStart:" << triggerStart_
107  << " triggerEnd:" << triggerEnd_
108  << " inTrigger:" << inTrigger
109  << endl;
110 
111  switch (controlMode_)
112  {
113  case controlMode::TIME:
114  {
115  return inTime;
116  }
117  case controlMode::TRIGGER:
118  {
119  return inTrigger;
120  }
121  case controlMode::TIME_OR_TRIGGER:
122  {
123  return inTime || inTrigger;
124  }
125  case controlMode::TIME_AND_TRIGGER:
126  {
127  return inTime && inTrigger;
128  }
129  default:
130  {
132  << "Unhandled enumeration: " << controlModeNames_[controlMode_]
133  << abort(FatalError);
134  }
135  }
136 
137  return false;
138 }
139 
140 
141 Foam::scalar Foam::functionObjects::timeControl::calcExpansion
142 (
143  const scalar startRatio,
144  const scalar y,
145  const label n
146 )
147 {
148  scalar ratio = startRatio;
149  // This function is used to calculate the 'expansion ratio' or
150  // 'time step growth factor'. Inputs:
151  // - ratio: wanted ratio
152  // - y: overall time required (divided by deltaT) to get from 1 to ratio
153  // - n: number of steps in which to get to wanted ratio
154 
155  // Let a0 = first term in geometric progression (GP), an = last term,
156  // n = number of terms, ratio = ratio of GP, Sn = sum of series for
157  // first n terms.
158  // We have an=a0*ratio^(n-1) Sn = a0(ratio^n -1)/(ratio -1)
159  // -> Sn=an/ratio^(n-1)*(ratio^n -1)/(ratio -1), assume y=Sn/an
160  // -> y*ratio^(n-1)=(ratio^n -1)/(ratio -1)
161  // -> y*ratio^n - y*ratio^(n-1) -ratio^n +1 = 0
162 
163  // Solve using Newton-Raphson iteration to determine the expansion
164  // ratio giving the desired overall timestep change.
165  // Note max iteration count to avoid hanging; function generally
166  // converges in 5 iterations or so.
167  for (label iter = 0; iter < 100; iter++)
168  {
169  // Dimensionless equation
170  scalar f = (y-1)*pow(ratio, n)+1-y*pow(ratio, n-1);
171  scalar dfdratio = (y-1)*n*pow(ratio, n-1)-y*(n-1)*pow(ratio, n-2);
172  scalar dratio = f/(dfdratio + SMALL);
173  ratio -= dratio;
174  // Exit if function is satisfied
175  if (mag(f) < 1e-6)
176  {
177  return ratio;
178  }
179  }
180 
181  if (debug)
182  {
184  << "Did not converge to find new timestep growth factor given "
185  << "overall factor " << y << " and " << n << " steps to do it in."
186  << endl << " Returning current best guess " << ratio << endl;
187  }
188  return ratio;
189 }
190 
191 
192 void Foam::functionObjects::timeControl::calcDeltaTCoeff
193 (
194  scalar& requiredDeltaTCoeff,
195  const scalar wantedDT,
196  const label nSteps,
197  const scalar presentTime,
198  const scalar timeToNextWrite,
199  const bool rampDirectionUp
200 )
201 {
202  const scalar writeInterval = writeControl_.interval();
203 
204 
205  // Set new deltaT approximated to last step value
206  // adjust to include geometric series sum of nSteps terms
207  scalar newDeltaT = deltaT0_;
208 
209  if (seriesDTCoeff_ != GREAT)
210  {
211  newDeltaT *= seriesDTCoeff_;
212  }
213 
214  // Recalculate new deltaTCoeff_ based on rounded steps
215  requiredDeltaTCoeff = Foam::exp(Foam::log(wantedDT/newDeltaT)/nSteps);
216 
217  // Calculate total time required with given dT increment
218  // to reach specified dT value
219  scalar requiredTimeInterval = newDeltaT;
220 
221  // For nSteps = 1 during we get coeff = 1.0. Adding SMALL in denominator
222  // might lead to SMALL requiredTimeInterval, and further unnecessary
223  // calculations
224  if (requiredDeltaTCoeff != 1.0)
225  {
226  requiredTimeInterval *=
227  (pow(requiredDeltaTCoeff, nSteps) - 1)
228  /(requiredDeltaTCoeff - 1);
229  }
230 
231  // Nearest time which is multiple of wantedDT, after
232  // requiredTimeInterval
233  scalar timeToNextMultiple = -presentTime;
234 
235  if (rampDirectionUp)
236  {
237  timeToNextMultiple +=
238  ceil
239  (
240  (presentTime + requiredTimeInterval)
241  /wantedDT
242  )
243  *wantedDT;
244  }
245  else
246  {
247  timeToNextMultiple +=
248  floor
249  (
250  (presentTime - requiredTimeInterval)
251  /wantedDT
252  )
253  *wantedDT;
254  }
255 
256  // Check nearest time and decide parameters
257  if (timeToNextWrite <= timeToNextMultiple)
258  {
259  // As the ramping can't be achieved, use current, unadapted deltaT
260  // (from e.g. maxCo calculation in solver) and adjust when
261  // writeInterval is closer than this deltaT
262  newDeltaT = deltaT0_;
263  if (timeToNextWrite < newDeltaT)
264  {
265  newDeltaT = timeToNextWrite;
266  }
267 
268  // Corner case when required time span is less than writeInterval
269  if (requiredTimeInterval > writeInterval)
270  {
272  << "With given ratio needed time span "
273  << requiredTimeInterval
274  << " exceeds available writeInterval "
275  << writeInterval << nl
276  << "Disabling all future time step ramping"
277  << endl;
278  deltaTCoeff_ = GREAT;
279  newDeltaT = wantedDT;
280  }
281  else
282  {
283  // Reset the series ratio, as this could change later
284  seriesDTCoeff_ = GREAT;
285  if (debug)
286  {
287  Info<< "Disabling ramping until time "
288  << presentTime+timeToNextWrite << endl;
289  }
290  }
291  requiredDeltaTCoeff = newDeltaT/deltaT0_;
292  }
293  else
294  {
295  // Find the minimum number of time steps (with constant deltaT
296  // variation) to get us to the writeInterval and multiples of wantedDT.
297  // Increases the number of nSteps to do the adjustment in until it
298  // has reached one that is possible. Returns the new expansionratio.
299 
300  scalar y = timeToNextMultiple/wantedDT;
301  label requiredSteps = nSteps;
302 
303  scalar ratioEstimate = deltaTCoeff_;
304  scalar ratioMax = deltaTCoeff_;
305 
306  if (seriesDTCoeff_ != GREAT)
307  {
308  ratioEstimate = seriesDTCoeff_;
309  }
310 
311  if (!rampDirectionUp)
312  {
313  ratioEstimate = 1/ratioEstimate;
314  ratioMax = 1/ratioMax;
315  requiredSteps *= -1;
316  }
317  // Provision for fall-back value if we can't satisfy requirements
318  bool searchConverged = false;
319  for (label iter = 0; iter < 100; iter++)
320  {
321  // Calculate the new expansion and the overall time step changes
322  const scalar newRatio = calcExpansion
323  (
324  ratioEstimate,
325  y,
326  requiredSteps
327  );
328  const scalar deltaTLoop =
329  wantedDT
330  / pow(newRatio, mag(requiredSteps)-1);
331  scalar firstDeltaRatio = deltaTLoop/deltaT0_;
332  // Avoid division by zero for ratio = 1.0
333  scalar Sn =
334  deltaTLoop
335  *(pow(newRatio, mag(requiredSteps))-1)
336  /(newRatio-1+SMALL);
337 
338  if (debug)
339  {
340  Info<< " nSteps " << requiredSteps
341  << " ratioEstimate " << ratioEstimate
342  << " a0 " << deltaTLoop
343  << " firstDeltaRatio " << firstDeltaRatio
344  << " Sn " << Sn << " Sn+Time " << Sn+presentTime
345  << " seriesRatio "<< newRatio << endl;
346  }
347 
348  // If series ratio becomes unity check next deltaT multiple
349  // Do the same if first ratio is decreasing when ramping up!
350  if
351  (
352  (firstDeltaRatio < 1.0 && rampDirectionUp)
353  || (firstDeltaRatio > 1.0 && !rampDirectionUp)
354  || newRatio == 1.0
355  )
356  {
357  y += 1;
358  requiredSteps = mag(nSteps);
359  if (debug)
360  {
361  Info<< "firstDeltaRatio " << firstDeltaRatio << " rampDir"
362  << rampDirectionUp << " newRatio " << newRatio
363  << " y " << y << " steps " << requiredSteps <<endl;
364  }
365  continue;
366  }
367 
368  if (firstDeltaRatio > ratioMax && rampDirectionUp)
369  {
370  requiredSteps++;
371  if (debug)
372  {
373  Info<< "First ratio " << firstDeltaRatio
374  << " exceeds threshold " << ratioMax << nl
375  << "Increasing required steps " << requiredSteps
376  << endl;
377  }
378  }
379  else if (firstDeltaRatio < ratioMax && !rampDirectionUp)
380  {
381  requiredSteps--;
382  if (debug)
383  {
384  Info<< "First ratio " << firstDeltaRatio
385  << " exceeds threshold " << ratioMax << nl
386  << "Decreasing required steps " << requiredSteps
387  << endl;
388  }
389  }
390  else if
391  (
392  (newRatio > ratioMax && rampDirectionUp)
393  || (newRatio < ratioMax && !rampDirectionUp)
394  )
395  {
396  y += 1;
397  requiredSteps = nSteps;
398  if (debug)
399  {
400  Info<< "Series ratio " << newRatio << "exceeds threshold "
401  << ratioMax << nl << "Consider next deltaT multiple "
402  << y*wantedDT + presentTime << endl;
403  }
404  }
405  else
406  {
407  seriesDTCoeff_ = newRatio;
408  // Reset future series expansion at last step
409  if (requiredSteps == 1)
410  {
411  Sn = deltaT0_*firstDeltaRatio;
412  seriesDTCoeff_ = GREAT;
413  }
414  // Situations where we achieve wantedDT but fail to achieve
415  // multiple of writeInterval
416  scalar jumpError =
417  remainder(Sn + presentTime, wantedDT) / wantedDT;
418 
419  if (mag(jumpError) > ROOTSMALL)
420  {
421  requiredSteps = label(timeToNextWrite/wantedDT);
422  firstDeltaRatio = timeToNextWrite/requiredSteps/deltaT0_;
423  }
424  if (debug)
425  {
426  Info<< "All conditions satisfied deltaT0_:" << deltaT0_
427  << " calculated deltaTCoeff:" << firstDeltaRatio
428  << " series ratio set to:" << seriesDTCoeff_ << endl;
429  }
430  searchConverged = true;
431  requiredDeltaTCoeff = firstDeltaRatio;
432  break;
433  }
434  }
435 
436  if (!searchConverged)
437  {
438  if (debug)
439  {
440  // Not converged. Return fall-back value
442  << "Did not converge to find new timestep growth factor"
443  << " given overall factor " << y
444  << " and " << requiredSteps << " steps to do it in."
445  << nl << " Falling back to non-adjusted deltatT "
446  << deltaT0_ << endl;
447  }
448  requiredDeltaTCoeff = 1.0;
449  }
450  }
451 }
452 
453 
454 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
455 
456 Foam::functionObjects::timeControl::timeControl
457 (
458  const word& name,
459  const Time& runTime,
460  const dictionary& dict
461 )
462 :
463  timeFunctionObject(name, runTime),
464  dict_(dict),
465  controlMode_(controlMode::TIME),
466  timeStart_(-VGREAT),
467  timeEnd_(VGREAT),
468  triggerStart_(labelMin),
469  triggerEnd_(labelMax),
470  nStepsToStartTimeChange_(labelMax),
471  executeControl_(runTime, dict, "execute"),
472  writeControl_(runTime, dict, "write"),
473  foPtr_(functionObject::New(name, runTime, dict_)),
474  executeTimeIndex_(-1),
475  deltaT0_(0),
476  seriesDTCoeff_(GREAT)
477 {
478  readControls();
479 }
480 
481 
482 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
483 
485 {
486  if
487  (
489  || Foam::timeControl::entriesPresent(dict, "output") // backwards compat
491  || dict.found("timeStart")
492  || dict.found("timeEnd")
493  || dict.found("triggerStart")
494  || dict.found("triggerEnd")
495  )
496  {
497  return true;
498  }
499 
500  return false;
501 }
502 
503 
505 {
506  // Store old deltaT for reference
507  deltaT0_ = time_.deltaTValue();
508 
509  if (active() && (postProcess || executeControl_.execute()))
510  {
511  executeTimeIndex_ = time_.timeIndex();
512  foPtr_->execute();
513  }
514 
515  return true;
516 }
517 
518 
519 bool Foam::functionObjects::timeControl::execute(const label subIndex)
520 {
521  if (active())
522  {
523  // Call underlying function object directly
524  foPtr_->execute(subIndex);
525  }
526 
527  return true;
528 }
529 
530 
532 {
533  if (active() && (postProcess || writeControl_.execute()))
534  {
535  // Ensure written results reflect the current state
536  if (executeTimeIndex_ != time_.timeIndex())
537  {
538  executeTimeIndex_ = time_.timeIndex();
539  foPtr_->execute();
540  }
541 
542  foPtr_->write();
543  }
544 
545  return true;
546 }
547 
548 
550 {
551  if (active() && (executeControl_.execute() || writeControl_.execute()))
552  {
553  foPtr_->end();
554  }
555 
556  return true;
557 }
558 
559 
561 {
562  if (active())
563  {
564  if
565  (
566  writeControl_.control()
568  )
569  {
570  const scalar presentTime = time_.value();
571 
572  //const scalar oldDeltaT = time_.deltaTValue();
573 
574  // Call underlying fo to do time step adjusting
575  foPtr_->adjustTimeStep();
576  const scalar wantedDT = time_.deltaTValue();
577 
578  //Pout<< "adjustTimeStep by " << foPtr_->name()
579  // << " from " << oldDeltaT << " to " << wantedDT << endl;
580 
581  const label writeTimeIndex = writeControl_.executionIndex();
582  const scalar writeInterval = writeControl_.interval();
583 
584  // Get time to next writeInterval
585  scalar timeToNextWrite =
586  (writeTimeIndex+1)*writeInterval
587  - (presentTime - time_.startTime().value());
588  if (timeToNextWrite <= 0.0)
589  {
590  timeToNextWrite = writeInterval;
591  }
592 
593  // The target DT (coming from fixed dT or maxCo etc) may not be
594  // an integral multiple of writeInterval. If so ignore any step
595  // calculation and give priority to writeInterval.
596  // Note looser tolerance on the relative intervalError - SMALL is
597  // too strict.
598  scalar intervalError =
599  remainder(writeInterval, wantedDT)/writeInterval;
600  if
601  (
602  mag(intervalError) > ROOTSMALL
603  || deltaTCoeff_ == GREAT
604  )
605  {
606  scalar deltaT = time_.deltaTValue();
607  scalar nSteps = timeToNextWrite/deltaT;
608 
609  // For tiny deltaT the label can overflow!
610  if
611  (
612  nSteps < scalar(labelMax)
613  && (
614  deltaTCoeff_ != GREAT
615  || nSteps < nStepsToStartTimeChange_
616  )
617  )
618  {
619  // nSteps can be < 1 so make sure at least 1
620 
621  label nStepsToNextWrite = max(1, round(nSteps));
622  scalar newDeltaT = timeToNextWrite/nStepsToNextWrite;
623 
624  // Backwards compatibility: clip deltaT to 0.2 .. 2
625  scalar clipThreshold = 2;
626  if (deltaTCoeff_ != GREAT)
627  {
628  clipThreshold = deltaTCoeff_;
629  }
630  // Adjust time step
631  if (newDeltaT >= deltaT)
632  {
633  deltaT = min(newDeltaT, clipThreshold*deltaT);
634  }
635  else
636  {
637  clipThreshold = 1/clipThreshold;
638  deltaT = max(newDeltaT, clipThreshold*deltaT);
639  }
640 
641  const_cast<Time&>(time_).setDeltaT(deltaT, false);
642  }
643  }
644  else
645  {
646  // Set initial approximation to user defined ratio.
647  scalar requiredDeltaTCoeff = deltaTCoeff_;
648 
649  // Use if we have series expansion ratio from last attempt.
650  // Starting from user defined coefficient, might get different
651  // steps and convergence time (could be recursive - worst case)
652  // This is more likely to happen when coefficient limit is high
653  if (seriesDTCoeff_ != GREAT)
654  {
655  requiredDeltaTCoeff = seriesDTCoeff_;
656  }
657  // Avoid divide by zero if we need ratio = 1
658  if (1/Foam::log(requiredDeltaTCoeff) > scalar(labelMax))
659  {
660  requiredDeltaTCoeff = deltaTCoeff_;
661  }
662 
663  // Try and observe the deltaTCoeff and the writeInterval
664  // and the wanted deltaT. Below code calculates the
665  // minimum number of time steps in which to end up on
666  // writeInterval+n*deltaT with the minimum possible deltaTCoeff.
667  // This is non-trivial!
668 
669  // Approx steps required from present dT to required dT
670  label nSteps = 0;
671  if (wantedDT < deltaT0_)
672  {
673  nSteps = label
674  (
675  floor
676  (
677  Foam::log(wantedDT/deltaT0_)
678  /Foam::log(requiredDeltaTCoeff + SMALL)
679  )
680  );
681  }
682  else
683  {
684  nSteps = label
685  (
686  ceil
687  (
688  Foam::log(wantedDT/deltaT0_)
689  /Foam::log(requiredDeltaTCoeff + SMALL)
690  )
691  );
692  }
693  // Negative steps -> ramp down needed,
694  // zero steps -> we are at writeInterval-multiple time!
695  if (nSteps < 1)
696  {
697  // Not enough steps to do gradual time step adjustment
698  // if earlier deltaT larger than wantedDT
699  if (wantedDT < deltaT0_)
700  {
701  requiredDeltaTCoeff = 1/requiredDeltaTCoeff;
702  calcDeltaTCoeff
703  (
704  requiredDeltaTCoeff,
705  wantedDT,
706  nSteps,
707  presentTime,
708  timeToNextWrite,
709  false
710  );
711  }
712  else
713  {
714  if (timeToNextWrite > wantedDT)
715  {
716  requiredDeltaTCoeff = wantedDT/deltaT0_;
717  }
718  else
719  {
720  requiredDeltaTCoeff = timeToNextWrite/deltaT0_;
721  }
722  }
723  // When allowed deltaTCoeff can't give required ramp
724  if (requiredDeltaTCoeff > deltaTCoeff_ && debug)
725  {
727  << "Required deltaTCoeff "<< requiredDeltaTCoeff
728  << " is larger than allowed value "<< deltaTCoeff_
729  << endl;
730  }
731  }
732  else
733  {
734  calcDeltaTCoeff
735  (
736  requiredDeltaTCoeff,
737  wantedDT,
738  nSteps,
739  presentTime,
740  timeToNextWrite,
741  true
742  );
743  }
744 
745  // Calculate deltaT to be used based on ramp
746  scalar newDeltaT = deltaT0_*requiredDeltaTCoeff;
747 
748  // Set time, deltaT adjustments for writeInterval purposes
749  // are already taken care. Hence disable auto-update
750  const_cast<Time&>(time_).setDeltaT(newDeltaT, false);
751 
752  if (seriesDTCoeff_ < 1.0)
753  {
754  requiredDeltaTCoeff = 1/requiredDeltaTCoeff;
755  seriesDTCoeff_ = 1/seriesDTCoeff_;
756  }
757  }
758  }
759  else
760  {
761  foPtr_->adjustTimeStep();
762  const scalar wantedDT = time_.deltaTValue();
763 
764  if (deltaTCoeff_ != GREAT)
765  {
766  // Clip time step change to deltaTCoeff
767 
768  scalar requiredDeltaTCoeff =
769  (
770  max
771  (
772  1.0/deltaTCoeff_,
773  min(deltaTCoeff_, wantedDT/deltaT0_)
774  )
775  );
776 
777  scalar newDeltaT = deltaT0_*requiredDeltaTCoeff;
778 
779  // Set time, deltaT adjustments for writeInterval purposes
780  // are already taken care. Hence disable auto-update
781  const_cast<Time&>( time_).setDeltaT(newDeltaT, false);
782  }
783  else
784  {
785  // Set the DT without any ramping
786  const_cast<Time&>( time_).setDeltaT(wantedDT, false);
787  }
788  }
789  }
790 
791  return true;
792 }
793 
794 
795 bool Foam::functionObjects::timeControl::read(const dictionary& dict)
796 {
797  if (dict != dict_)
798  {
799  dict_ = dict;
800 
801  writeControl_.read(dict);
802  executeControl_.read(dict);
803  readControls();
804 
805  // Forward to underlying function object
806  return foPtr_->read(dict);
807  }
808 
809  return false;
810 }
811 
814 {
815  return active() ? foPtr_->filesModified() : false;
816 }
817 
818 
820 {
821  if (active())
822  {
823  foPtr_->updateMesh(mpm);
824  }
825 }
826 
827 
829 {
830  if (active())
831  {
832  foPtr_->movePoints(mesh);
833  }
834 }
835 
836 
837 // ************************************************************************* //
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
static bool entriesPresent(const dictionary &dict, const word &prefix)
Identify if a timeControl object is present in the dictionary.
Definition: timeControl.C:80
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:49
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
static scalar Sn(const scalar a, const scalar x)
Definition: invIncGamma.C:64
static bool entriesPresent(const dictionary &dict)
Helper function to identify if a timeControl object is present.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static const Enum< controlMode > controlModeNames_
virtual bool execute()
Called at each ++ or += of the time-loop.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
constexpr label labelMin
Definition: label.H:54
Abstract base-class for Time/database function objects.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
EnumType getOrDefault(const word &key, const dictionary &dict, const EnumType deflt, const bool warnOnly=false) const
Find the key in the dictionary and return the corresponding enumeration element based on its name...
Definition: Enum.C:173
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
Macros for easy insertion into run-time selection tables.
virtual bool write()
Called at each ++ or += of the time-loop.
scalar y
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
dimensionedScalar exp(const dimensionedScalar &ds)
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: TimeState.C:42
A class for handling words, derived from Foam::string.
Definition: word.H:63
#define DebugInFunction
Report an information message using Foam::Info.
virtual void updateMesh(const mapPolyMesh &mpm)
Update for changes of mesh.
virtual bool filesModified() const
Did any file get changed during execution?
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
int debug
Static debugging option.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
labelList f(nPoints)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual void movePoints(const polyMesh &mesh)
Update for changes of mesh.
#define WarningInFunction
Report a warning using Foam::Warning.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr label labelMax
Definition: label.H:55
label n
Currently identical to "runTime".
Definition: timeControl.H:71
virtual bool adjustTimeStep()
Called at the end of Time::adjustDeltaT() if adjustTime is true.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
virtual bool read(const dictionary &)
Read and set the function object if its data have changed.
virtual bool end()
Called when Time::run() determines that the time-loop exits.
Namespace for OpenFOAM.
const Time & time_
Reference to the time database.