54class fieldInterpolator
58 const IOobjectList& objects_;
59 const wordRes& selectedFields_;
62 const interpolationWeights& interpolator_;
72 const IOobjectList& objects,
73 const wordRes& selectedFields,
76 const interpolationWeights& interpolator,
77 const wordList& timeNames,
84 selectedFields_(selectedFields),
87 interpolator_(interpolator),
88 timeNames_(timeNames),
92 template<
class GeoFieldType>
97template<
class GeoFieldType>
98void fieldInterpolator::interpolate()
105 selectedFields_.empty()
106 ? objects_.csorted<GeoFieldType>()
107 : objects_.csorted<GeoFieldType>(selectedFields_)
111 const word& fieldName =
io.name();
115 Info<<
" " << GeoFieldType::typeName <<
's';
118 Info<<
' ' << fieldName <<
'(';
120 const scalar deltaT = (ti1_.value() - ti_.value())/(divisions_ + 1);
122 for (
int j=0; j<divisions_; j++)
130 if (j < divisions_-1)
138 interpolator_.valueWeights
195int main(
int argc,
char *argv[])
199 "Interpolate fields between time-steps. Eg, for animation."
208 "The fields (or field) to be interpolated."
209 " Eg, '(U T p \"Y.*\")' or a single field 'U'"
215 "Specify number of temporal sub-divisions to create (default = 1)."
221 "The type of interpolation (linear or spline)"
232 if (selectedFields.empty())
234 Info<<
"Interpolating all fields" <<
nl <<
endl;
243 const int divisions =
args.getOrDefault<
int>(
"divisions", 1);
244 Info<<
"Using " << divisions <<
" per time interval" <<
nl <<
endl;
247 const word interpolationType =
248 args.getOrDefault<
word>(
"interpolationType",
"linear");
250 Info<<
"Using interpolation " << interpolationType <<
nl <<
endl;
259 timeVals[i] = timeDirs[i].value();
260 timeNames[i] = timeDirs[i].name();
274 Info<<
"Interpolating fields for times:" <<
endl;
276 for (label timei = 0; timei < timeDirs.
size() - 1; timei++)
278 runTime.setTime(timeDirs[timei], timei);
283 fieldInterpolator interpolator
Field reading functions for post-processing utilities.
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,...
const T & name() const noexcept
The name/key (const access).
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
static void addOption(const word &optName, const string ¶m="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
static void addNote(const string ¬e)
Add extra notes for the usage information.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name.
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.
A List of wordRe with additional matching capabilities.
A class for handling words, derived from Foam::string.
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< instant > instantList
List of instants.
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.
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.
GeometricField< sphericalTensor, fvsPatchField, surfaceMesh > surfaceSphericalTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
GeometricField< symmTensor, fvsPatchField, surfaceMesh > surfaceSymmTensorField
constexpr char nl
The newline '\n' character (0x0a).
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.