temporalInterpolate.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25  Interpolate fields between time-steps e.g. for animation.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "argList.H"
30 #include "timeSelector.H"
31 
32 #include "fvMesh.H"
33 #include "Time.H"
34 #include "volMesh.H"
35 #include "surfaceMesh.H"
36 #include "volFields.H"
37 #include "surfaceFields.H"
38 #include "pointFields.H"
39 #include "ReadFields.H"
40 #include "interpolationWeights.H"
41 #include "uniformInterpolate.H"
42 
43 using namespace Foam;
44 
45 class fieldInterpolator
46 {
47  Time& runTime_;
48  const fvMesh& mesh_;
49  const IOobjectList& objects_;
50  const HashSet<word>& selectedFields_;
51  instant ti_;
52  instant ti1_;
53  const interpolationWeights& interpolator_;
54  const wordList& timeNames_;
55  int divisions_;
56 
57 public:
58 
59  fieldInterpolator
60  (
61  Time& runTime,
62  const fvMesh& mesh,
63  const IOobjectList& objects,
64  const HashSet<word>& selectedFields,
65  const instant& ti,
66  const instant& ti1,
67  const interpolationWeights& interpolator,
68  const wordList& timeNames,
69  int divisions
70  )
71  :
72  runTime_(runTime),
73  mesh_(mesh),
74  objects_(objects),
75  selectedFields_(selectedFields),
76  ti_(ti),
77  ti1_(ti1),
78  interpolator_(interpolator),
79  timeNames_(timeNames),
80  divisions_(divisions)
81  {}
82 
83  template<class GeoFieldType>
84  void interpolate();
85 };
86 
87 
88 template<class GeoFieldType>
90 {
91  const word& fieldClassName = GeoFieldType::typeName;
92 
93  IOobjectList fields = objects_.lookupClass(fieldClassName);
94 
95  if (fields.size())
96  {
97  Info<< " " << fieldClassName << "s:";
98 
100  {
101  if
102  (
103  selectedFields_.empty()
104  || selectedFields_.found(fieldIter()->name())
105  )
106  {
107  Info<< " " << fieldIter()->name() << '(';
108 
109  scalar deltaT = (ti1_.value() - ti_.value())/(divisions_ + 1);
110 
111  for (int j=0; j<divisions_; j++)
112  {
113  instant timej = instant(ti_.value() + (j + 1)*deltaT);
114 
115  runTime_.setTime(instant(timej.name()), 0);
116 
117  Info<< timej.name();
118 
119  if (j < divisions_-1)
120  {
121  Info<< " ";
122  }
123 
124  // Calculate times to read and weights
125  labelList indices;
126  scalarField weights;
127  interpolator_.valueWeights
128  (
129  runTime_.value(),
130  indices,
131  weights
132  );
133 
134  const wordList selectedTimeNames
135  (
136  UIndirectList<word>(timeNames_, indices)()
137  );
138 
139  // Read on the objectRegistry all the required fields
140  ReadFields<GeoFieldType>
141  (
142  fieldIter()->name(),
143  mesh_,
144  selectedTimeNames
145  );
146 
147  GeoFieldType fieldj
148  (
149  uniformInterpolate<GeoFieldType>
150  (
151  IOobject
152  (
153  fieldIter()->name(),
154  runTime_.name(),
155  fieldIter()->db(),
158  false
159  ),
160  fieldIter()->name(),
161  selectedTimeNames,
162  weights
163  )
164  );
165 
166  fieldj.write();
167  }
168 
169  Info<< ')';
170  }
171  }
172 
173  Info<< endl;
174  }
175 }
176 
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 int main(int argc, char *argv[])
181 {
183  #include "addRegionOption.H"
185  (
186  "fields",
187  "list",
188  "specify a list of fields to be interpolated. Eg, '(U T p)' - "
189  "regular expressions not currently supported"
190  );
192  (
193  "divisions",
194  "integer",
195  "specify number of temporal sub-divisions to create (default = 1)."
196  );
198  (
199  "interpolationType",
200  "word",
201  "specify type of interpolation (linear or spline)"
202  );
203 
204  #include "setRootCase.H"
206 
207  HashSet<word> selectedFields;
208  if (args.optionFound("fields"))
209  {
210  args.optionLookup("fields")() >> selectedFields;
211  }
212  if (selectedFields.size())
213  {
214  Info<< "Interpolating fields " << selectedFields << nl << endl;
215  }
216  else
217  {
218  Info<< "Interpolating all fields" << nl << endl;
219  }
220 
221 
222  int divisions = 1;
223  if (args.optionFound("divisions"))
224  {
225  args.optionLookup("divisions")() >> divisions;
226  }
227  Info<< "Using " << divisions << " per time interval" << nl << endl;
228 
229 
230  const word interpolationType = args.optionLookupOrDefault<word>
231  (
232  "interpolationType",
233  "linear"
234  );
235  Info<< "Using interpolation " << interpolationType << nl << endl;
236 
237 
239 
240  scalarField timeVals(timeDirs.size());
241  wordList timeNames(timeDirs.size());
242  forAll(timeDirs, i)
243  {
244  timeVals[i] = timeDirs[i].value();
245  timeNames[i] = timeDirs[i].name();
246  }
247  autoPtr<interpolationWeights> interpolatorPtr
248  (
250  (
251  interpolationType,
252  timeVals
253  )
254  );
255 
256 
257  #include "createNamedMesh.H"
258 
259  Info<< "Interpolating fields for times:" << endl;
260 
261  for (label timei = 0; timei < timeDirs.size() - 1; timei++)
262  {
263  runTime.setTime(timeDirs[timei], timei);
264 
265  // Read objects in time directory
266  IOobjectList objects(mesh, runTime.name());
267 
268  fieldInterpolator interpolator
269  (
270  runTime,
271  mesh,
272  objects,
273  selectedFields,
274  timeDirs[timei],
275  timeDirs[timei+1],
276  interpolatorPtr(),
277  timeNames,
278  divisions
279  );
280 
281  // Interpolate vol fields
282  interpolator.interpolate<volScalarField>();
283  interpolator.interpolate<volVectorField>();
284  interpolator.interpolate<volSphericalTensorField>();
285  interpolator.interpolate<volSymmTensorField>();
286  interpolator.interpolate<volTensorField>();
287 
288  // Interpolate surface fields
289  interpolator.interpolate<surfaceScalarField>();
290  interpolator.interpolate<surfaceVectorField>();
291  interpolator.interpolate<surfaceSphericalTensorField>();
292  interpolator.interpolate<surfaceSymmTensorField>();
293  interpolator.interpolate<surfaceTensorField>();
294  }
295 
296  Info<< "End\n" << endl;
297 
298  return 0;
299 }
300 
301 
302 // ************************************************************************* //
Field reading functions for post-processing utilities.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Generic GeometricField class.
A HashTable with keys but without contents.
Definition: HashSet.H:62
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
A List with indirect addressing.
Definition: UIndirectList.H:60
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:120
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:243
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
An instant of time. Contains the time value and name.
Definition: instant.H:67
const word & name() const
Name (const access)
Definition: instant.H:126
Abstract base class for interpolating in 1D.
static autoPtr< interpolationWeights > New(const word &type, const scalarField &samples)
Return a reference to the selected interpolationWeights.
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:252
A class for handling words, derived from string.
Definition: word.H:62
int main(int argc, char *argv[])
Definition: financialFoam.C:44
static instantList timeDirs
Definition: globalFoam.H:44
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
static const char nl
Definition: Ostream.H:260
objects
Foam::argList args(argc, argv)
Foam::surfaceFields.