53 ListListOps::combine<List<Type>>
67 namespace functionObjects
112 dict.lookup(
"fields") >> fields_;
114 UName_ =
dict.lookupOrDefault(
"U",
word(
"U"));
116 writeAge_ =
dict.lookupOrDefault<
Switch>(
"writeAge",
true);
118 trackDirection_ = trackDirectionNames_[
word(
dict.lookup(
"direction"))];
120 trackOutside_ =
dict.lookupOrDefault<
Switch>(
"outside",
false);
122 dict.lookup(
"lifeTime") >> lifeTime_;
127 <<
"Illegal value " << lifeTime_ <<
" for lifeTime"
131 bool subCycling =
dict.found(
"nSubCycle");
132 bool fixedLength =
dict.found(
"trackLength");
133 if (subCycling && fixedLength)
136 <<
"Cannot both specify automatic time stepping (through '"
137 <<
"nSubCycle' specification) and fixed track length (through '"
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;
151 dict.lookup(
"trackLength") >> trackLength_;
152 Info<<
" fixed track length specified : "
153 << trackLength_ <<
nl <<
endl;
156 interpolationScheme_ =
159 "interpolationScheme",
163 cloudName_ =
dict.lookupOrDefault<
word>(
"cloudName",
"streamlines");
170 dict.subDict(
"seedSampleSet")
215 cannotFindObject(fields_[fieldi]);
220 #define DeclareTypeInterpolator(Type, nullArg) \
221 PtrList<interpolation<Type>> Type##Interp(fieldNames.size());
223 #undef DeclareTypeInterpolator
226 #define ConstructTypeInterpolator(Type, nullArg) \
227 if (mesh_.foundObject<VolField<Type>>(fieldNames[fieldi])) \
232 interpolation<Type>::New \
234 interpolationScheme_, \
235 mesh_.lookupObject<VolField<Type>>(fieldNames[fieldi]) \
240 #undef ConstructTypeInterpolator
253 interpolationScheme_,
264 #define DeclareAllTypes(Type, nullArg) \
265 List<DynamicField<Type>> all##Type##s(fieldNames.size());
267 #undef DeclareAllTypes
277 label nLocateBoundaryHits;
278 forAll(sampledSetPtr_(), i)
285 sampledSetPtr_().positions()[i],
286 sampledSetPtr_().
cells()[i],
296 Info <<
" Seeded " << nSeeds <<
" particles" <<
endl;
306 UIndex != -1 ? vectorInterp[UIndex] : UInterp(),
307 trackDirection_ == trackDirection::forward,
323 if (trackDirection_ == trackDirection::both)
330 if (trackDirection_ == trackDirection::both)
347 #define GatherAndFlattenAllTypes(Type, nullArg) \
348 if (Type##Interp.set(fieldi)) \
350 gatherAndFlatten(all##Type##s[fieldi]); \
353 #undef GatherAndFlattenAllTypes
358 Info<<
" Sampled " << allPositions.
size() <<
" locations" <<
endl;
364 const label nTracks =
max(allTracks) + 1;
365 const label trackParti0 =
min(allTrackParts);
366 const label trackParti1 =
max(allTrackParts) + 1;
373 forAll(allPositions, samplei)
375 const label tracki = allTracks[samplei];
376 const label trackParti = -trackParti0 + allTrackParts[samplei];
377 trackPartCounts[tracki][trackParti] ++;
386 forAll(trackPartOffsets, tracki)
388 forAll(trackPartOffsets[tracki], trackParti)
390 trackPartOffsets[tracki][trackParti] +=
offset;
391 offset += trackPartCounts[tracki][trackParti];
395 forAll(trackPartCounts, tracki)
397 trackPartCounts[tracki] = 0;
400 forAll(allPositions, samplei)
402 const label tracki = allTracks[samplei];
403 const label trackParti = -trackParti0 + allTrackParts[samplei];
406 trackPartOffsets[tracki][trackParti]
407 + trackPartCounts[tracki][trackParti];
409 trackPartCounts[tracki][trackParti] ++;
448 allTrackParts.
rmap(allTrackParts,
order);
452 #define RMapAllTypes(Type, nullArg) \
453 if (Type##Interp.set(fieldi)) \
455 all##Type##s[fieldi].rmap(all##Type##s[fieldi], order); \
468 label samplei = 0, tracki = 0;
469 forAll(allPositions, samplej)
471 const label trackj = allTracks[samplej];
472 const label trackPartj = allTrackParts[samplej];
474 allPositions[samplei] = allPositions[samplej];
475 allTracks[samplei] = tracki;
476 allTrackParts[samplei] = 0;
477 allAges[samplei] = allAges[samplej];
480 #define ShuffleUpAllTypes(Type, nullArg) \
481 if (Type##Interp.set(fieldi)) \
483 all##Type##s[fieldi][samplei] = \
484 all##Type##s[fieldi][samplej]; \
487 #undef ShuffleUpAllTypes
490 const bool joinNewParts =
491 samplej != allPositions.
size() - 1
493 && allTrackParts[samplej + 1] == 0;
495 if (!joinNewParts) samplei ++;
498 samplej == allPositions.
size() - 1
499 || trackj != allTracks[samplej + 1]
500 || trackPartj != allTrackParts[samplej + 1];
502 if (!joinNewParts && newPart) tracki ++;
505 allPositions.
resize(samplei);
506 allTracks.
resize(samplei);
507 allTrackParts.
resize(samplei);
511 #define ResizeAllTypes(Type, nullArg) \
512 if (Type##Interp.set(fieldi)) \
514 all##Type##s[fieldi].resize(samplei); \
517 #undef ResizeAllTypes
538 #define DeclareTypeValueSets(Type, nullArg) \
539 UPtrList<const Field<Type>> Type##ValueSets(nValueSets);
541 #undef DeclareTypeValueSets
544 valueSetNames[0] =
"age";
545 scalarValueSets.set(0, &allAges);
549 valueSetNames[fieldi + writeAge_] =
fieldNames[fieldi];
551 #define SetTypeValueSetPtr(Type, nullArg) \
552 if (Type##Interp.set(fieldi)) \
554 Type##ValueSets.set \
556 fieldi + writeAge_, \
557 &all##Type##s[fieldi] \
561 #undef SetTypeValueSetPtr
585 searchEngine_.correct();
586 sampledSetPtr_->movePoints();
596 if (&map.
mesh() == &mesh_)
598 searchEngine_.correct();
600 sampledSetPtr_->topoChange(map);
610 if (&map.
mesh() == &mesh_)
612 searchEngine_.correct();
614 sampledSetPtr_->mapMesh(map);
624 if (&map.
mesh() == &mesh_)
626 searchEngine_.correct();
628 sampledSetPtr_->distribute(map);
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.
Macros for easy insertion into run-time selection tables.
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
Generic GeometricField class.
Template class for intrusive linked lists.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void append(const T &)
Append an element at the end of the list.
void size(const label)
Override size to be inconsistent with allocated storage.
Initialise the NamedEnum HashTable from the static list of names.
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,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static bool master(const label communicator=0)
Am I the master process.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
static bool & parRun()
Is this a parallel run?
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Holds list of sampling positions.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
A class for handling file names.
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...
Generates streamline data by sampling a set of user-specified fields along a particle track,...
streamlines(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
virtual wordList fields() const
Return the list of fields required.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual ~streamlines()
Destructor.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
virtual bool execute()
Do nothing.
virtual bool write()
Calculate and write the streamlines.
virtual bool read(const dictionary &)
Read the field average data.
static const NamedEnum< trackDirection, 3 > trackDirectionNames_
Track direction enumeration names.
static const word outputPrefix
Directory prefix.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
label toGlobal(const label i) const
From local to global.
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.
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.
const polyMesh & mesh() const
Return polyMesh.
Mesh consisting of general polyhedral cells.
static word defaultRegion
Return the default region name.
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.
static autoPtr< setWriter > New(const word &writeType, const IOstream::streamFormat writeFormat=IOstream::ASCII, const IOstream::compressionType writeCompression=IOstream::UNCOMPRESSED)
Select given write options.
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.
A class for handling words, derived from string.
static const word null
An empty word.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
static List< word > fieldNames
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)
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< label > labelList
A List of labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Ostream & endl(Ostream &os)
Add newline and flush stream.
int order(const scalar s)
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)
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,.
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
void offset(label &lst, const label o)
void gatherAndFlatten(DynamicField< Type > &field)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
#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)