Loading...
Searching...
No Matches
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-2023 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
27Application
28 temporalInterpolate
29
30Group
31 grpPostProcessingUtilities
32
33Description
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"
50#include "uniformInterpolate.H"
51
52using namespace Foam;
53
54class 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
66public:
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
97template<class GeoFieldType>
98void fieldInterpolator::interpolate()
99{
100 label nFields = 0;
101 for
102 (
103 const IOobject& io :
104 (
105 selectedFields_.empty()
106 ? objects_.csorted<GeoFieldType>()
107 : objects_.csorted<GeoFieldType>(selectedFields_)
108 )
109 )
110 {
111 const word& fieldName = io.name();
112
113 if (!nFields)
114 {
115 Info<< " " << GeoFieldType::typeName << 's';
116 }
117
118 Info<< ' ' << fieldName << '(';
119
120 const scalar deltaT = (ti1_.value() - ti_.value())/(divisions_ + 1);
121
122 for (int j=0; j<divisions_; j++)
123 {
124 instant timej = instant(ti_.value() + (j + 1)*deltaT);
125
126 runTime_.setTime(instant(timej.name()), 0);
127
128 Info<< timej.name();
129
130 if (j < divisions_-1)
131 {
132 Info<< " ";
133 }
134
135 // Calculate times to read and weights
136 labelList indices;
137 scalarField weights;
138 interpolator_.valueWeights
139 (
140 runTime_.value(),
141 indices,
142 weights
143 );
144
145 const wordList selectedTimeNames
146 (
147 UIndirectList<word>(timeNames_, indices)()
148 );
149
150 //Info<< "For time " << runTime_.value()
151 // << " need times " << selectedTimeNames
152 // << " need weights " << weights << endl;
153
154
155 // Read on the objectRegistry all the required fields
157 (
158 fieldName,
159 mesh_,
160 selectedTimeNames
161 );
162
163 GeoFieldType fieldj
164 (
166 (
168 (
169 fieldName,
170 runTime_.timeName(),
171 io.db(),
175 ),
176 fieldName,
177 selectedTimeNames,
178 weights
179 )
180 );
181
182 fieldj.write();
183 }
184
185 Info<< ')';
186 ++nFields;
187 }
188
189 if (nFields) Info<< endl;
190}
191
192
193// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194
195int main(int argc, char *argv[])
196{
198 (
199 "Interpolate fields between time-steps. Eg, for animation."
200 );
201
203 #include "addRegionOption.H"
205 (
206 "fields",
207 "wordRes",
208 "The fields (or field) to be interpolated."
209 " Eg, '(U T p \"Y.*\")' or a single field 'U'"
210 );
212 (
213 "divisions",
214 "integer",
215 "Specify number of temporal sub-divisions to create (default = 1)."
216 );
218 (
219 "interpolationType",
220 "word",
221 "The type of interpolation (linear or spline)"
222 );
223
224 argList::noFunctionObjects(); // Never use function objects
225
226 #include "setRootCase.H"
227 #include "createTime.H"
228
229 // Non-mandatory
230 const wordRes selectedFields(args.getList<wordRe>("fields", false));
231
232 if (selectedFields.empty())
233 {
234 Info<< "Interpolating all fields" << nl << endl;
235 }
236 else
237 {
238 Info<< "Interpolating fields " << flatOutput(selectedFields)
239 << nl << endl;
240 }
241
242
243 const int divisions = args.getOrDefault<int>("divisions", 1);
244 Info<< "Using " << divisions << " per time interval" << nl << endl;
245
246
247 const word interpolationType =
248 args.getOrDefault<word>("interpolationType", "linear");
249
250 Info<< "Using interpolation " << interpolationType << nl << endl;
251
252
254
255 scalarField timeVals(timeDirs.size());
256 wordList timeNames(timeDirs.size());
257 forAll(timeDirs, i)
258 {
259 timeVals[i] = timeDirs[i].value();
260 timeNames[i] = timeDirs[i].name();
261 }
262 autoPtr<interpolationWeights> interpolatorPtr
263 (
265 (
266 interpolationType,
267 timeVals
268 )
269 );
270
271
272 #include "createNamedMesh.H"
273
274 Info<< "Interpolating fields for times:" << endl;
275
276 for (label timei = 0; timei < timeDirs.size() - 1; timei++)
277 {
278 runTime.setTime(timeDirs[timei], timei);
279
280 // Read objects in time directory
281 IOobjectList objects(mesh, runTime.timeName());
282
283 fieldInterpolator interpolator
284 (
285 runTime,
286 mesh,
287 objects,
288 selectedFields,
289 timeDirs[timei],
290 timeDirs[timei+1],
291 interpolatorPtr(),
292 timeNames,
293 divisions
294 );
295
296 // Interpolate vol fields
297 interpolator.interpolate<volScalarField>();
298 interpolator.interpolate<volVectorField>();
299 interpolator.interpolate<volSphericalTensorField>();
300 interpolator.interpolate<volSymmTensorField>();
301 interpolator.interpolate<volTensorField>();
302
303 // Interpolate surface fields
304 interpolator.interpolate<surfaceScalarField>();
305 interpolator.interpolate<surfaceVectorField>();
306 interpolator.interpolate<surfaceSphericalTensorField>();
307 interpolator.interpolate<surfaceSymmTensorField>();
308 interpolator.interpolate<surfaceTensorField>();
309 }
310
311 Info<< "End\n" << endl;
312
313 return 0;
314}
315
316
317// ************************************************************************* //
Field reading functions for post-processing utilities.
Required Classes.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const T & name() const noexcept
The name/key (const access).
Definition Instant.H:149
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
Definition argList.C:562
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:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name.
Definition instant.H:56
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 timeSelector options to argList::validOptions.
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options and also set the runTime to the first i...
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition wordRe.H:81
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
mesh interpolate(rAU)
dynamicFvMesh & mesh
engineTime & runTime
Required Classes.
const auto & io
Namespace for OpenFOAM.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
List< word > wordList
List of word.
Definition fileName.H:60
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< instant > instantList
List of instants.
Definition instantList.H:41
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< GeoField > uniformInterpolate(const HashPtrTable< GeoField, label, Hash< label > > &fields, const labelList &indices, const scalarField &weights)
Interpolate selected fields (given by indices and corresponding weights).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
GeometricField< tensor, fvsPatchField, surfaceMesh > surfaceTensorField
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
GeometricField< sphericalTensor, fvsPatchField, surfaceMesh > surfaceSphericalTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
GeometricField< symmTensor, fvsPatchField, surfaceMesh > surfaceSymmTensorField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.