streamlines.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-2025 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 \*---------------------------------------------------------------------------*/
25 
26 #include "streamlines.H"
27 #include "streamlinesCloud.H"
28 #include "meshSearch.H"
29 #include "sampledSet.H"
30 #include "globalIndex.H"
31 #include "interpolationCellPoint.H"
32 #include "PatchTools.H"
33 #include "writeFile.H"
34 #include "polyTopoChangeMap.H"
35 #include "polyMeshMap.H"
36 #include "polyDistributionMap.H"
37 #include "OSspecific.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 template<class Type>
47 {
48  List<List<Type>> gatheredField(Pstream::nProcs());
49  gatheredField[Pstream::myProcNo()] = field;
50  Pstream::gatherList(gatheredField);
51 
52  field =
53  ListListOps::combine<List<Type>>
54  (
55  gatheredField,
57  );
58 }
59 
60 }
61 
62 
63 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
64 
65 namespace Foam
66 {
67  namespace functionObjects
68  {
71  }
72 }
73 
76 {
77  "forward",
78  "backward",
79  "both"
80 };
81 
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
86 (
87  const word& name,
88  const Time& runTime,
89  const dictionary& dict
90 )
91 :
92  fvMeshFunctionObject(name, runTime, dict),
93  nSubCycle_(0),
94  searchEngine_(mesh_)
95 {
96  read(dict);
97 }
98 
99 
100 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
101 
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
109 {
110  Info<< type() << " " << name() << ":" << nl;
111 
112  dict.lookup("fields") >> fields_;
113 
114  UName_ = dict.lookupOrDefault("U", word("U"));
115 
116  writeAge_ = dict.lookupOrDefault<Switch>("writeAge", true);
117 
118  trackDirection_ = trackDirectionNames_[word(dict.lookup("direction"))];
119 
120  trackOutside_ = dict.lookupOrDefault<Switch>("outside", false);
121 
122  dict.lookup("lifeTime") >> lifeTime_;
123 
124  if (lifeTime_ < 1)
125  {
127  << "Illegal value " << lifeTime_ << " for lifeTime"
128  << exit(FatalError);
129  }
130 
131  bool subCycling = dict.found("nSubCycle");
132  bool fixedLength = dict.found("trackLength");
133  if (subCycling && fixedLength)
134  {
136  << "Cannot both specify automatic time stepping (through '"
137  << "nSubCycle' specification) and fixed track length (through '"
138  << "trackLength')"
139  << exit(FatalIOError);
140  }
141  if (subCycling)
142  {
143  nSubCycle_ = max(dict.lookup<scalar>("nSubCycle"), 1);
144  trackLength_ = vGreat;
145  Info<< " automatic track length specified through"
146  << " number of sub cycles : " << nSubCycle_ << nl << endl;
147  }
148  else
149  {
150  nSubCycle_ = 1;
151  dict.lookup("trackLength") >> trackLength_;
152  Info<< " fixed track length specified : "
153  << trackLength_ << nl << endl;
154  }
155 
156  interpolationScheme_ =
157  dict.lookupOrDefault
158  (
159  "interpolationScheme",
161  );
162 
163  cloudName_ = dict.lookupOrDefault<word>("cloudName", "streamlines");
164 
165  sampledSetPtr_ = sampledSet::New
166  (
167  "seedSampleSet",
168  mesh_,
169  searchEngine_,
170  dict.subDict("seedSampleSet")
171  );
172 
173  formatterPtr_ = setWriter::New(dict.lookup("setFormat"), dict);
174 
175  return true;
176 }
177 
178 
180 {
181  wordList allFields(fields_);
182  allFields.append(UName_);
183 
184  return allFields;
185 }
186 
187 
189 {
190  return true;
191 }
192 
193 
195 {
196  Info<< type() << " " << name() << " write:" << nl;
197 
198  // Create list of available fields
200  forAll(fields_, fieldi)
201  {
202  if
203  (
204  false
205  #define FoundTypeField(Type, nullArg) \
206  || foundObject<VolField<Type>>(fields_[fieldi])
208  #undef FoundTypeField
209  )
210  {
211  fieldNames.append(fields_[fieldi]);
212  }
213  else
214  {
215  cannotFindObject(fields_[fieldi]);
216  }
217  }
218 
219  // Lookup fields and construct interpolators
220  #define DeclareTypeInterpolator(Type, nullArg) \
221  PtrList<interpolation<Type>> Type##Interp(fieldNames.size());
223  #undef DeclareTypeInterpolator
224  forAll(fieldNames, fieldi)
225  {
226  #define ConstructTypeInterpolator(Type, nullArg) \
227  if (mesh_.foundObject<VolField<Type>>(fieldNames[fieldi])) \
228  { \
229  Type##Interp.set \
230  ( \
231  fieldi, \
232  interpolation<Type>::New \
233  ( \
234  interpolationScheme_, \
235  mesh_.lookupObject<VolField<Type>>(fieldNames[fieldi]) \
236  ) \
237  ); \
238  }
240  #undef ConstructTypeInterpolator
241  }
242 
243  // Create a velocity interpolator if it is not already available
244  const label UIndex = findIndex(fieldNames, UName_);
245  tmpNrc<interpolation<vector>> UInterp(nullptr);
246  if (UIndex == -1)
247  {
248  UInterp =
250  (
252  (
253  interpolationScheme_,
254  mesh_.lookupObject<volVectorField>(UName_)
255  ).ptr()
256  );
257  }
258 
259  // Do tracking to create sampled data
260  DynamicField<point> allPositions;
261  DynamicField<label> allTracks;
262  DynamicField<label> allTrackParts;
263  DynamicField<scalar> allAges;
264  #define DeclareAllTypes(Type, nullArg) \
265  List<DynamicField<Type>> all##Type##s(fieldNames.size());
267  #undef DeclareAllTypes
268  {
269  // Create a cloud and initialise with points from the sampled set
270  globalIndex gi(sampledSetPtr_().size());
272  (
273  mesh_,
274  cloudName_,
276  );
277  label nLocateBoundaryHits;
278  forAll(sampledSetPtr_(), i)
279  {
280  particles.addParticle
281  (
283  (
284  mesh_,
285  sampledSetPtr_().positions()[i],
286  sampledSetPtr_().cells()[i],
287  nLocateBoundaryHits,
288  lifeTime_,
289  gi.toGlobal(i)
290  )
291  );
292  }
293 
294  // Report the number of successful seeds
295  const label nSeeds = returnReduce(particles.size(), sumOp<label>());
296  Info << " Seeded " << nSeeds << " particles" << endl;
297 
298  // Create tracking data
300  (
301  particles,
302  #define TypeInterpolatorParameter(Type, nullArg) \
303  Type##Interp,
306  UIndex != -1 ? vectorInterp[UIndex] : UInterp(),
307  trackDirection_ == trackDirection::forward,
308  trackOutside_,
309  nSubCycle_,
310  trackLength_,
311  allPositions,
312  allTracks,
313  allTrackParts,
314  allAges
315  #define AllTypesParameter(Type, nullArg) \
316  , all##Type##s
318  #undef AllTypesParameter
319  );
320 
321  // Track
322  IDLList<streamlinesParticle> initialParticles;
323  if (trackDirection_ == trackDirection::both)
324  {
325  initialParticles = particles;
326  }
327 
328  particles.move(particles, td);
329 
330  if (trackDirection_ == trackDirection::both)
331  {
332  particles.IDLList<streamlinesParticle>::operator=(initialParticles);
333  td.trackForward_ = !td.trackForward_;
334  particles.move(particles, td);
335  }
336  }
337 
338  // Gather data on the master
339  if (Pstream::parRun())
340  {
341  gatherAndFlatten(allPositions);
342  gatherAndFlatten(allTracks);
343  gatherAndFlatten(allTrackParts);
344  gatherAndFlatten(allAges);
345  forAll(fieldNames, fieldi)
346  {
347  #define GatherAndFlattenAllTypes(Type, nullArg) \
348  if (Type##Interp.set(fieldi)) \
349  { \
350  gatherAndFlatten(all##Type##s[fieldi]); \
351  }
353  #undef GatherAndFlattenAllTypes
354  }
355  }
356 
357  // Report the total number of samples
358  Info<< " Sampled " << allPositions.size() << " locations" << endl;
359 
360  // Bin-sort by track and trackPart to build an ordering
361  labelList order(allPositions.size());
362  if (Pstream::master() && allPositions.size())
363  {
364  const label nTracks = max(allTracks) + 1;
365  const label trackParti0 = min(allTrackParts);
366  const label trackParti1 = max(allTrackParts) + 1;
367 
368  labelListList trackPartCounts
369  (
370  nTracks,
371  labelList(trackParti1 - trackParti0, 0)
372  );
373  forAll(allPositions, samplei)
374  {
375  const label tracki = allTracks[samplei];
376  const label trackParti = -trackParti0 + allTrackParts[samplei];
377  trackPartCounts[tracki][trackParti] ++;
378  }
379 
380  label offset = 0;
381  labelListList trackPartOffsets
382  (
383  nTracks,
384  labelList(trackParti1 - trackParti0, 0)
385  );
386  forAll(trackPartOffsets, tracki)
387  {
388  forAll(trackPartOffsets[tracki], trackParti)
389  {
390  trackPartOffsets[tracki][trackParti] += offset;
391  offset += trackPartCounts[tracki][trackParti];
392  }
393  }
394 
395  forAll(trackPartCounts, tracki)
396  {
397  trackPartCounts[tracki] = 0;
398  }
399 
400  forAll(allPositions, samplei)
401  {
402  const label tracki = allTracks[samplei];
403  const label trackParti = -trackParti0 + allTrackParts[samplei];
404 
405  order[samplei] =
406  trackPartOffsets[tracki][trackParti]
407  + trackPartCounts[tracki][trackParti];
408 
409  trackPartCounts[tracki][trackParti] ++;
410  }
411  }
412 
413  //auto reportTrackParts = [&]()
414  //{
415  // Info<< nl;
416  // forAll(allPositions, samplei)
417  // {
418  // if
419  // (
420  // samplei == 0
421  // || allTracks[samplei] != allTracks[samplei - 1]
422  // || allTrackParts[samplei] != allTrackParts[samplei - 1]
423  // )
424  // {
425  // Info<< "track #" << allTracks[samplei]
426  // << " part #" << allTrackParts[samplei]
427  // << " from i=" << samplei << " to ";
428  // }
429  // if
430  // (
431  // samplei == allPositions.size() - 1
432  // || allTracks[samplei + 1] != allTracks[samplei]
433  // || allTrackParts[samplei + 1] != allTrackParts[samplei]
434  // )
435  // {
436  // Info<< "i=" << samplei << nl;
437  // }
438  // }
439  //};
440 
441  //reportTrackParts();
442 
443  // Reorder
444  if (Pstream::master())
445  {
446  allPositions.rmap(allPositions, order);
447  allTracks.rmap(allTracks, order);
448  allTrackParts.rmap(allTrackParts, order);
449  allAges.rmap(allAges, order);
450  forAll(fieldNames, fieldi)
451  {
452  #define RMapAllTypes(Type, nullArg) \
453  if (Type##Interp.set(fieldi)) \
454  { \
455  all##Type##s[fieldi].rmap(all##Type##s[fieldi], order); \
456  }
458  #undef RMapAllTypes
459  }
460  }
461 
462  //reportTrackParts();
463 
464  // Relabel tracks and track parts into track labels only, and join the
465  // forward and backward track parts that are connected to the seed
466  if (Pstream::master())
467  {
468  label samplei = 0, tracki = 0;
469  forAll(allPositions, samplej)
470  {
471  const label trackj = allTracks[samplej];
472  const label trackPartj = allTrackParts[samplej];
473 
474  allPositions[samplei] = allPositions[samplej];
475  allTracks[samplei] = tracki;
476  allTrackParts[samplei] = 0;
477  allAges[samplei] = allAges[samplej];
478  forAll(fieldNames, fieldi)
479  {
480  #define ShuffleUpAllTypes(Type, nullArg) \
481  if (Type##Interp.set(fieldi)) \
482  { \
483  all##Type##s[fieldi][samplei] = \
484  all##Type##s[fieldi][samplej]; \
485  }
487  #undef ShuffleUpAllTypes
488  }
489 
490  const bool joinNewParts =
491  samplej != allPositions.size() - 1
492  && trackPartj == -1
493  && allTrackParts[samplej + 1] == 0;
494 
495  if (!joinNewParts) samplei ++;
496 
497  const bool newPart =
498  samplej == allPositions.size() - 1
499  || trackj != allTracks[samplej + 1]
500  || trackPartj != allTrackParts[samplej + 1];
501 
502  if (!joinNewParts && newPart) tracki ++;
503  }
504 
505  allPositions.resize(samplei);
506  allTracks.resize(samplei);
507  allTrackParts.resize(samplei);
508  allAges.resize(samplei);
509  forAll(fieldNames, fieldi)
510  {
511  #define ResizeAllTypes(Type, nullArg) \
512  if (Type##Interp.set(fieldi)) \
513  { \
514  all##Type##s[fieldi].resize(samplei); \
515  }
517  #undef ResizeAllTypes
518  }
519  }
520 
521  //reportTrackParts();
522 
523  // Write
524  if (Pstream::master() && allPositions.size())
525  {
526  // Make output directory
527  const fileName outputPath =
528  time_.globalPath()
530  /(mesh_.name() != polyMesh::defaultRegion ? mesh_.name() : word())
531  /name()
532  /time_.name();
533  mkDir(outputPath);
534 
535  // Pass data to the formatter to write
536  const label nValueSets = fieldNames.size() + writeAge_;
537  wordList valueSetNames(nValueSets);
538  #define DeclareTypeValueSets(Type, nullArg) \
539  UPtrList<const Field<Type>> Type##ValueSets(nValueSets);
541  #undef DeclareTypeValueSets
542  if (writeAge_)
543  {
544  valueSetNames[0] = "age";
545  scalarValueSets.set(0, &allAges);
546  }
547  forAll(fieldNames, fieldi)
548  {
549  valueSetNames[fieldi + writeAge_] = fieldNames[fieldi];
550 
551  #define SetTypeValueSetPtr(Type, nullArg) \
552  if (Type##Interp.set(fieldi)) \
553  { \
554  Type##ValueSets.set \
555  ( \
556  fieldi + writeAge_, \
557  &all##Type##s[fieldi] \
558  ); \
559  }
561  #undef SetTypeValueSetPtr
562  }
563  formatterPtr_->write
564  (
565  outputPath,
566  "tracks",
567  coordSet(allTracks, word::null, allPositions),
568  valueSetNames
569  #define TypeValueSetsParameter(Type, nullArg) , Type##ValueSets
572  );
573  }
574 
575  Info<< endl;
576 
577  return true;
578 }
579 
580 
582 {
583  if (&mesh == &mesh_)
584  {
585  searchEngine_.correct();
586  sampledSetPtr_->movePoints();
587  }
588 }
589 
590 
592 (
593  const polyTopoChangeMap& map
594 )
595 {
596  if (&map.mesh() == &mesh_)
597  {
598  searchEngine_.correct();
599 
600  sampledSetPtr_->topoChange(map);
601  }
602 }
603 
604 
606 (
607  const polyMeshMap& map
608 )
609 {
610  if (&map.mesh() == &mesh_)
611  {
612  searchEngine_.correct();
613 
614  sampledSetPtr_->mapMesh(map);
615  }
616 }
617 
618 
620 (
621  const polyDistributionMap& map
622 )
623 {
624  if (&map.mesh() == &mesh_)
625  {
626  searchEngine_.correct();
627 
628  sampledSetPtr_->distribute(map);
629  }
630 }
631 
632 
633 // ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
Dynamically sized Field.
Definition: DynamicField.H:72
void resize(const label)
Alter the addressed list size.
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:386
Generic GeometricField class.
Template class for intrusive linked lists.
Definition: ILList.H:67
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
Holds list of sampling positions.
Definition: coordSet.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A class for handling file names.
Definition: fileName.H:82
Abstract base-class for Time/database functionObjects.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
This functionObject tracks a particle cloud in the specified velocity field of an incompressible flow...
Definition: particles.H:111
Generates streamline data by sampling a set of user-specified fields along a particle track,...
Definition: streamlines.H:209
streamlines(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: streamlines.C:86
virtual wordList fields() const
Return the list of fields required.
Definition: streamlines.C:179
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: streamlines.C:592
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: streamlines.C:620
virtual ~streamlines()
Destructor.
Definition: streamlines.C:102
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: streamlines.C:606
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: streamlines.C:581
virtual bool execute()
Do nothing.
Definition: streamlines.C:188
virtual bool write()
Calculate and write the streamlines.
Definition: streamlines.C:194
virtual bool read(const dictionary &)
Read the field average data.
Definition: streamlines.C:108
static const NamedEnum< trackDirection, 3 > trackDirectionNames_
Track direction enumeration names.
Definition: streamlines.H:223
static const word outputPrefix
Directory prefix.
Definition: writeFile.H:72
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
static autoPtr< interpolation< Type > > New(const word &interpolationType, const VolField< Type > &psi)
Return a reference to the specified interpolation scheme.
Definition: interpolation.C:56
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
const polyMesh & mesh() const
Return polyMesh.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
const polyMesh & mesh() const
Return polyMesh.
Definition: polyMeshMap.H:75
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:270
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const polyMesh & mesh() const
Return polyMesh.
static autoPtr< sampledSet > New(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Return a reference to the selected sampledSet.
Definition: sampledSet.C:190
static autoPtr< setWriter > New(const word &writeType, const IOstream::streamFormat writeFormat=IOstream::ASCII, const IOstream::compressionType writeCompression=IOstream::UNCOMPRESSED)
Select given write options.
Definition: setWriter.C:286
A Cloud of streamlines particles.
Particle class that samples fields as it passes through. Used in streamlines calculation.
A class for managing temporary objects without reference counting.
Definition: tmpNrc.H:53
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
static List< word > fieldNames
Definition: globalFoam.H:46
const cellShapeList & cells
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
int order(const scalar s)
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
void offset(label &lst, const label o)
void gatherAndFlatten(DynamicField< Type > &field)
Definition: streamlines.C:46
static const char nl
Definition: Ostream.H:267
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
#define FoundTypeField(Type, nullArg)
#define DeclareTypeValueSets(Type, nullArg)
#define SetTypeValueSetPtr(Type, nullArg)
#define DeclareAllTypes(Type, nullArg)
#define AllTypesParameter(Type, nullArg)
#define ResizeAllTypes(Type, nullArg)
#define TypeValueSetsParameter(Type, nullArg)
#define DeclareTypeInterpolator(Type, nullArg)
#define GatherAndFlattenAllTypes(Type, nullArg)
#define ConstructTypeInterpolator(Type, nullArg)
#define RMapAllTypes(Type, nullArg)
#define TypeInterpolatorParameter(Type, nullArg)
#define ShuffleUpAllTypes(Type, nullArg)