temporalInterpolate.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-2015 OpenFOAM Foundation
9  Copyright (C) 2018-2020 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 Application
28  temporalInterpolate
29 
30 Group
31  grpPostProcessingUtilities
32 
33 Description
34  Interpolate fields between time-steps e.g. for animation.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "argList.H"
39 #include "timeSelector.H"
40 
41 #include "fvMesh.H"
42 #include "Time.H"
43 #include "volMesh.H"
44 #include "surfaceMesh.H"
45 #include "volFields.H"
46 #include "surfaceFields.H"
47 #include "pointFields.H"
48 #include "ReadFields.H"
49 #include "interpolationWeights.H"
50 #include "uniformInterpolate.H"
51 
52 using namespace Foam;
53 
54 class fieldInterpolator
55 {
56  Time& runTime_;
57  const fvMesh& mesh_;
58  const IOobjectList& objects_;
59  const wordRes& selectedFields_;
60  instant ti_;
61  instant ti1_;
62  const interpolationWeights& interpolator_;
63  const wordList& timeNames_;
64  int divisions_;
65 
66 public:
67 
68  fieldInterpolator
69  (
70  Time& runTime,
71  const fvMesh& mesh,
72  const IOobjectList& objects,
73  const wordRes& selectedFields,
74  const instant& ti,
75  const instant& ti1,
76  const interpolationWeights& interpolator,
77  const wordList& timeNames,
78  int divisions
79  )
80  :
81  runTime_(runTime),
82  mesh_(mesh),
83  objects_(objects),
84  selectedFields_(selectedFields),
85  ti_(ti),
86  ti1_(ti1),
87  interpolator_(interpolator),
88  timeNames_(timeNames),
89  divisions_(divisions)
90  {}
91 
92  template<class GeoFieldType>
93  void interpolate();
94 };
95 
96 
97 template<class GeoFieldType>
99 {
100  const word& clsName = GeoFieldType::typeName;
101 
102  const wordList fieldNames =
103  (
104  selectedFields_.empty()
105  ? objects_.sortedNames(clsName)
106  : objects_.sortedNames(clsName, selectedFields_)
107  );
108 
109  if (fieldNames.size())
110  {
111  Info<< " " << clsName << 's';
112  }
113 
114  for (const word& fieldName : fieldNames)
115  {
116  Info<< ' ' << fieldName << '(';
117 
118  const scalar deltaT = (ti1_.value() - ti_.value())/(divisions_ + 1);
119 
120  for (int j=0; j<divisions_; j++)
121  {
122  instant timej = instant(ti_.value() + (j + 1)*deltaT);
123 
124  runTime_.setTime(instant(timej.name()), 0);
125 
126  Info<< timej.name();
127 
128  if (j < divisions_-1)
129  {
130  Info<< " ";
131  }
132 
133  // Calculate times to read and weights
134  labelList indices;
135  scalarField weights;
136  interpolator_.valueWeights
137  (
138  runTime_.value(),
139  indices,
140  weights
141  );
142 
143  const wordList selectedTimeNames
144  (
145  UIndirectList<word>(timeNames_, indices)()
146  );
147 
148  //Info<< "For time " << runTime_.value()
149  // << " need times " << selectedTimeNames
150  // << " need weights " << weights << endl;
151 
152 
153  // Read on the objectRegistry all the required fields
154  ReadFields<GeoFieldType>
155  (
156  fieldName,
157  mesh_,
158  selectedTimeNames
159  );
160 
161  GeoFieldType fieldj
162  (
163  uniformInterpolate<GeoFieldType>
164  (
165  IOobject
166  (
167  fieldName,
168  runTime_.timeName(),
169  objects_[fieldName]->db(),
172  false
173  ),
174  fieldName,
175  selectedTimeNames,
176  weights
177  )
178  );
179 
180  fieldj.write();
181  }
182 
183  Info<< ')';
184  }
185 
186  if (fieldNames.size()) Info<< endl;
187 }
188 
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 int main(int argc, char *argv[])
193 {
195  (
196  "Interpolate fields between time-steps. Eg, for animation."
197  );
198 
200  #include "addRegionOption.H"
202  (
203  "fields",
204  "wordRes",
205  "The fields (or field) to be interpolated."
206  " Eg, '(U T p \"Y.*\")' or a single field 'U'"
207  );
209  (
210  "divisions",
211  "integer",
212  "Specify number of temporal sub-divisions to create (default = 1)."
213  );
215  (
216  "interpolationType",
217  "word",
218  "The type of interpolation (linear or spline)"
219  );
220 
221  argList::noFunctionObjects(); // Never use function objects
222 
223  #include "setRootCase.H"
224  #include "createTime.H"
225 
226  // Non-mandatory
227  const wordRes selectedFields(args.getList<wordRe>("fields", false));
228 
229  if (selectedFields.empty())
230  {
231  Info<< "Interpolating all fields" << nl << endl;
232  }
233  else
234  {
235  Info<< "Interpolating fields " << flatOutput(selectedFields)
236  << nl << endl;
237  }
238 
239 
240  const int divisions = args.getOrDefault<int>("divisions", 1);
241  Info<< "Using " << divisions << " per time interval" << nl << endl;
242 
243 
244  const word interpolationType =
245  args.getOrDefault<word>("interpolationType", "linear");
246 
247  Info<< "Using interpolation " << interpolationType << nl << endl;
248 
249 
251 
252  scalarField timeVals(timeDirs.size());
253  wordList timeNames(timeDirs.size());
254  forAll(timeDirs, i)
255  {
256  timeVals[i] = timeDirs[i].value();
257  timeNames[i] = timeDirs[i].name();
258  }
259  autoPtr<interpolationWeights> interpolatorPtr
260  (
262  (
263  interpolationType,
264  timeVals
265  )
266  );
267 
268 
269  #include "createNamedMesh.H"
270 
271  Info<< "Interpolating fields for times:" << endl;
272 
273  for (label timei = 0; timei < timeDirs.size() - 1; timei++)
274  {
275  runTime.setTime(timeDirs[timei], timei);
276 
277  // Read objects in time directory
278  IOobjectList objects(mesh, runTime.timeName());
279 
280  fieldInterpolator interpolator
281  (
282  runTime,
283  mesh,
284  objects,
285  selectedFields,
286  timeDirs[timei],
287  timeDirs[timei+1],
288  interpolatorPtr(),
289  timeNames,
290  divisions
291  );
292 
293  // Interpolate vol fields
294  interpolator.interpolate<volScalarField>();
295  interpolator.interpolate<volVectorField>();
296  interpolator.interpolate<volSphericalTensorField>();
297  interpolator.interpolate<volSymmTensorField>();
298  interpolator.interpolate<volTensorField>();
299 
300  // Interpolate surface fields
301  interpolator.interpolate<surfaceScalarField>();
302  interpolator.interpolate<surfaceVectorField>();
303  interpolator.interpolate<surfaceSphericalTensorField>();
304  interpolator.interpolate<surfaceSymmTensorField>();
305  interpolator.interpolate<surfaceTensorField>();
306  }
307 
308  Info<< "End\n" << endl;
309 
310  return 0;
311 }
312 
313 
314 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
volFields.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::argList::getOrDefault
T getOrDefault(const word &optName, const T &deflt) const
Get a value from the named option if present, or return default.
Definition: argListI.H:307
volMesh.H
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::argList::addNote
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:412
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
surfaceFields.H
Foam::surfaceFields.
Foam::wordRe
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:80
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::argList::noFunctionObjects
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
Definition: argList.C:473
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::Field< scalar >
uniformInterpolate.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
argList.H
addRegionOption.H
interpolationWeights.H
Foam::interpolationWeights
Abstract base class for interpolating in 1D.
Definition: interpolationWeights.H:58
surfaceMesh.H
createNamedMesh.H
Required Variables.
fieldNames
const wordRes fieldNames(propsDict.getOrDefault< wordRes >("fields", wordRes()))
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:216
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:55
Foam::Instant::name
const T & name() const
The name/key (const access)
Definition: Instant.H:111
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
setRootCase.H
ReadFields.H
Field reading functions for post-processing utilities.
Foam::interpolationWeights::New
static autoPtr< interpolationWeights > New(const word &type, const scalarField &samples)
Return a reference to the selected interpolationWeights.
Definition: interpolationWeights.C:53
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::timeSelector::addOptions
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:102
Foam::List< word >
Foam::Time::setTime
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:1003
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
timeSelector.H
createTime.H
Foam::argList::getList
List< T > getList(const label index) const
Get a List of values from the argument at index.
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::instant
An instant of time. Contains the time value and name.
Definition: instant.H:52
Foam::timeSelector::select0
static instantList select0(Time &runTime, const argList &args)
Definition: timeSelector.C:235
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::argList::addOption
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:335
Foam::IOobject::NO_READ
Definition: IOobject.H:188
args
Foam::argList args(argc, argv)
pointFields.H
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.