TimePaths.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-2013 OpenFOAM Foundation
9  Copyright (C) 2016-2024 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 
29 #include "TimePaths.H"
30 #include "argList.H"
31 #include "fileOperation.H"
32 #include "IOstreams.H"
33 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 bool Foam::TimePaths::detectProcessorCase()
37 {
38  if (processorCase_)
39  {
40  return processorCase_;
41  }
42 
43  // Look for "processor", but should really check for following digits too
44  const auto sep = globalCaseName_.rfind('/');
45  const auto pos = globalCaseName_.find
46  (
47  "processor",
48  (sep == string::npos ? 0 : sep)
49  );
50 
51  if (pos == 0)
52  {
53  globalCaseName_ = ".";
54  processorCase_ = true;
55  }
56  else if (pos != string::npos && sep != string::npos && sep == pos-1)
57  {
58  globalCaseName_.resize(sep);
59  processorCase_ = true;
60  }
61 
62  return processorCase_;
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67 
69 (
70  const bool processorCase,
71  const fileName& rootPath,
72  const bool distributed,
73  const fileName& globalCaseName,
74  const fileName& caseName,
75  const word& systemDirName,
76  const word& constantDirName
77 )
78 :
79  processorCase_(processorCase),
80  distributed_(distributed),
81  rootPath_(rootPath),
82  globalCaseName_(globalCaseName),
83  case_(caseName),
84  system_(systemDirName),
85  constant_(constantDirName)
86 {
87  // Find out from case name whether it is a processor directory
88  // and set processorCase flag so file searching goes up one level.
89  detectProcessorCase();
90 }
91 
92 
94 (
95  const argList& args,
96  const word& systemDirName,
97  const word& constantDirName
98 )
99 :
100  TimePaths
101  (
102  args.runControl().parRun(), // processorCase
103  args.rootPath(),
104  args.runControl().distributed(),
105  args.globalCaseName(),
106  args.caseName(),
107  systemDirName,
108  constantDirName
109  )
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 (
117  const fileName& directory,
118  const word& constantDirName
119 )
120 {
121  return fileHandler().findTimes(directory, constantDirName);
122 }
123 
124 
126 {
127  return findTimes(path(), constant());
128 }
129 
130 
132 (
133  const UList<instant>& timeDirs,
134  const instant& t
135 )
136 {
137  // Note:
138  // - timeDirs will include constant (with value 0) as first element.
139  // For backwards compatibility make sure to find 0 in preference
140  // to constant.
141  // - list is sorted so could use binary search
142 
143  forAllReverse(timeDirs, i)
144  {
145  if (t.equal(timeDirs[i].value()))
146  {
147  return timeDirs[i].name();
148  }
149  }
150 
151  return word();
152 }
153 
154 
155 // Foam::word Foam::Time::findInstancePath
156 // (
157 // const fileName& directory,
158 // const instant& t
159 // ) const
160 // {
161 // // Simplified version: use findTimes (readDir + sort).
162 // // The expensive bit is the readDir, not the sorting.
163 // // TBD: avoid calling findInstancePath from filePath.
164 //
165 // instantList timeDirs = findTimes(directory, constant());
166 //
167 // return findInstancePath(timeDirs, i);
168 // }
169 
170 
172 {
173  // Simplified version: use findTimes (readDir + sort).
174  // The expensive bit is the readDir, not the sorting.
175  // TBD: avoid calling findInstancePath from filePath.
177  instantList timeDirs = findTimes(path(), constant());
178  return findInstancePath(timeDirs, t);
179 }
180 
181 
183 (
184  const UList<instant>& timeDirs,
185  const scalar t,
186  const word& constantDirName
187 )
188 {
189  const label nTimes = timeDirs.size();
190 
191  label nearestIndex = -1;
192  scalar deltaT = GREAT;
193 
194  for (label timei=0; timei < nTimes; ++timei)
195  {
196  if (timeDirs[timei].name() == constantDirName) continue;
197 
198  const scalar diff = mag(timeDirs[timei].value() - t);
199  if (diff < deltaT)
200  {
201  deltaT = diff;
202  nearestIndex = timei;
203  }
204  }
205 
206  return nearestIndex;
207 }
208 
209 
211 {
212  instantList timeDirs = findTimes(path(), constant());
213 
214  const label nTimes = timeDirs.size();
215 
216  if (nTimes == 0)
217  {
218  // Cannot really fail at this point, but for some safety...
219  return instant(0, constant());
220  }
221  else if (nTimes == 1)
222  {
223  // Only one time (likely "constant") so return it
224  return timeDirs.front();
225  }
226  else if (t < timeDirs[1].value())
227  {
228  return timeDirs[1];
229  }
230  else if (t > timeDirs.back().value())
231  {
232  return timeDirs.back();
233  }
234 
235  label nearestIndex = 0; // Failsafe value
236  scalar deltaT = GREAT;
237 
238  for (label timei=1; timei < nTimes; ++timei)
239  {
240  const scalar diff = mag(timeDirs[timei].value() - t);
241  if (diff < deltaT)
242  {
243  deltaT = diff;
244  nearestIndex = timei;
245  }
246  }
247 
248  return timeDirs[nearestIndex];
249 }
250 
251 
252 // ************************************************************************* //
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:373
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
Address the time paths without using the Time class.
Definition: TimePaths.H:51
A class for handling file names.
Definition: fileName.H:72
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
static label findClosestTimeIndex(const UList< instant > &timeDirs, const scalar t, const word &constantDirName="constant")
Search instantList for the time index closest to the specified time.
Definition: TimePaths.C:176
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
dimensionedScalar pos(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:118
static word findInstancePath(const UList< instant > &timeDirs, const instant &t)
Search instantList for matching time value, return the instance name or word::null if nothing is equa...
Definition: TimePaths.C:125
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
static instantList findTimes(const fileName &directory, const word &constantDirName="constant")
Search a given directory for valid time directories.
Definition: TimePaths.C:109
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name...
Definition: instant.H:53
instantList times() const
Search the case for valid time directories.
Definition: TimePaths.C:118
TimePaths(const bool processorCase, const fileName &rootPath, const bool distributed, const fileName &globalCaseName, const fileName &caseName, const word &systemDirName="system", const word &constantDirName="constant")
Construct from all components.
Definition: TimePaths.C:62
instant findClosestTime(const scalar t) const
Search the case for the time closest to the given time.
Definition: TimePaths.C:203
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:437
List< instant > instantList
List of instants.
Definition: instantList.H:41
bool equal(scalar val) const noexcept
True if values are equal (includes SMALL for rounding)
Definition: Instant.H:155
Foam::argList args(argc, argv)