33#include "surfaceInterpolate.H"
44namespace functionObjects
79 <<
": Number of faceZone faces (" << allFaces
80 <<
") is less than the number of requested locations ("
87 <<
" faces : " << allFaces <<
nl
100 if (!nInjectorLocations_)
105 const faceZone& fz = mesh_.faceZones()[zoneID_];
115 label nLocations = nInjectorLocations_;
120 scalar fraction = scalar(
nFaces)/scalar(nGlobalFaces);
121 nLocations = ceil(fraction*nInjectorLocations_);
125 <<
", nGlobalFaces:" << nGlobalFaces
126 <<
", fraction:" << fraction
127 <<
", nLocations:" << nLocations
145 label nCoarseFaces = 0;
149 nCoarseFaces =
max(fineToCoarseAddr_) + 1;
155 <<
" coarse faces" <<
endl;
192 patchIDs_.setSize(fz.
size(), -1);
193 patchFaceIDs_.setSize(fz.
size(), -1);
195 label nBlockedFaces = 0;
198 const label facei = fz[localFacei];
200 if (mesh_.isInternalFace(facei))
202 if (alphaf[facei] > alphaThreshold_)
204 blockedFaces[localFacei] =
true;
209 label patchi = mesh_.boundaryMesh().whichPatch(facei);
210 label patchFacei = -1;
214 const auto* cpp = isA<coupledPolyPatch>(pp);
218 patchFacei = (cpp->owner() ? pp.
whichFace(facei) : -1);
220 else if (!isA<emptyPolyPatch>(pp))
225 if (patchFacei == -1)
229 else if (alphafp[patchFacei] > alphaThreshold_)
231 blockedFaces[localFacei] =
true;
234 patchIDs_[localFacei] = patchi;
235 patchFaceIDs_[localFacei] = patchFacei;
251 const label particlei = regionToParticleMap_[regioni];
254 if (
p.faceIHit != -1 && nInjectorLocations_)
257 label coarseFacei = fineToCoarseAddr_[
p.faceIHit];
258 p.faceIHit = globalCoarseFaces_.toGlobal(coarseFacei);
265 if ((pDiameter > minDiameter_) && (pDiameter < maxDiameter_))
270 const point position =
p.VC/(
p.V + ROOTVSMALL);
271 const vector U =
p.VU/(
p.V + ROOTVSMALL);
273 if (nInjectorLocations_)
289 cloud_.addParticle(ip);
291 collectedVolume_ +=
p.V;
294 ++nCollectedParticles_;
299 ++nDiscardedParticles_;
300 discardedVolume_ +=
p.V;
307 const label nNewRegions,
316 labelList oldToNewRegion(particles_.size(), -1);
319 forAll(regionFaceIDs, facei)
321 label newRegioni = regionFaceIDs[facei];
322 label oldRegioni = regions0_[facei];
324 if (newRegioni != -1 && oldRegioni != -1)
328 newToNewRegion[newRegioni] =
329 max(newRegioni, oldToNewRegion[oldRegioni]);
330 oldToNewRegion[oldRegioni] = newRegioni;
338 label nParticle = -1;
341 forAll(newToNewRegion, newRegioni0)
343 label newRegioni = newToNewRegion[newRegioni0];
344 if (newRegions.
insert(newRegioni))
350 newRegionToParticleMap.
insert(newRegioni0, nParticle);
358 forAll(oldToNewRegion, oldRegioni)
360 label newRegioni = oldToNewRegion[oldRegioni];
361 if (newRegioni == -1)
365 <<
"Collecting particle from oldRegion:" << oldRegioni
368 collectParticle(time, oldRegioni);
373 label newParticlei = newRegionToParticleMap[newRegioni];
374 label oldParticlei = regionToParticleMap_[oldRegioni];
377 <<
"Combining newRegioni: " << newRegioni
378 <<
"(p:" << newParticlei <<
") and "
379 <<
"oldRegioni: " << oldRegioni
380 <<
"(p:" << oldParticlei <<
")"
383 newParticles[newParticlei] =
386 newParticles[newParticlei],
387 particles_[oldParticlei]
393 particles_.transfer(newParticles);
394 regionToParticleMap_ = newRegionToParticleMap;
398 regions0_ = regionFaceIDs;
415 const scalar deltaT = mesh_.time().deltaTValue();
416 const pointField& faceCentres = mesh_.faceCentres();
418 forAll(regionFaceIDs, localFacei)
420 const label newRegioni = regionFaceIDs[localFacei];
421 if (newRegioni != -1)
423 const label particlei = regionToParticleMap_[newRegioni];
424 const label meshFacei = fz[localFacei];
430 p.faceIHit = localFacei;
437 scalar magPhii =
mag(faceValue(
phi, localFacei, meshFacei));
438 vector Ufi = faceValue(
Uf, localFacei, meshFacei);
439 scalar dV = magPhii*deltaT;
441 p.VC += dV*faceCentres[meshFacei];
459 cloud_(mesh_,
"eulerianParticleCloud"),
460 faceZoneName_(
word::null),
465 alphaThreshold_(0.1),
469 nInjectorLocations_(0),
471 globalCoarseFaces_(),
474 regionToParticleMap_(),
475 minDiameter_(ROOTVSMALL),
477 nCollectedParticles_(getProperty<label>(
"nCollectedParticles", 0)),
478 collectedVolume_(getProperty<scalar>(
"collectedVolume", 0)),
479 nDiscardedParticles_(getProperty<label>(
"nDiscardedParticles", 0)),
480 discardedVolume_(getProperty<scalar>(
"discardedVolume", 0))
485 <<
name <<
" function object only applicable to 3-D cases"
504 dict.readEntry(
"faceZone", faceZoneName_);
505 dict.readEntry(
"alpha", alphaName_);
507 dict.readIfPresent(
"alphaThreshold", alphaThreshold_);
508 dict.readIfPresent(
"U", UName_);
509 dict.readIfPresent(
"rho", rhoName_);
510 dict.readIfPresent(
"phi", phiName_);
511 dict.readIfPresent(
"nLocations", nInjectorLocations_);
512 dict.readIfPresent(
"minDiameter", minDiameter_);
513 dict.readIfPresent(
"maxDiameter", maxDiameter_);
517 if (nInjectorLocations_)
540 typeName +
":alphaf",
544 const faceZone& fz = mesh_.faceZones()[zoneID_];
553 setBlockedFaces(alphaf, fz, blockedFaces);
560 const label nRegionsNew = regionFaceIDs.
nRegions();
568 mesh_.time().value(),
574 accumulateParticleInfo(alphaf, tphi(), regionFaceIDs, fz);
576 Log <<
" Collected particles : " << nCollectedParticles_ <<
nl
577 <<
" Collected volume : " << collectedVolume_ <<
nl
578 <<
" Discarded particles : " << nDiscardedParticles_ <<
nl
579 <<
" Discarded volume : " << discardedVolume_ <<
nl
580 <<
" Particles in progress : " << particles_.size() <<
nl
593 setProperty(
"nCollectedParticles", nCollectedParticles_);
594 setProperty(
"collectedVolume", collectedVolume_);
595 setProperty(
"nDiscardedParticles", nDiscardedParticles_);
596 setProperty(
"discardedVolume", discardedVolume_);
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
label size() const noexcept
The number of elements in table.
A List with indirect addressing.
void setSize(const label n)
Alias for resize()
A HashTable to objects of type <T> with a label key.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
A list of faces which address into the list of points.
static void listCombineAllGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
void size(const label n)
Older name for setAddressableSize.
static bool & parRun() noexcept
Test if this a parallel run.
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
wordList names() const
A list of the zone names.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A subset of mesh faces organised as a primitive patch.
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
virtual const word & type() const =0
Runtime type information.
Generates particle size information from Eulerian calculations, e.g. VoF.
virtual void accumulateParticleInfo(const surfaceScalarField &alphaf, const surfaceScalarField &phi, const labelList ®ionFaceIDs, const faceZone &fz)
Process latest region information.
virtual tmp< surfaceScalarField > phiU() const
Return the volumetric flux.
label nInjectorLocations_
Number of sample locations to generate.
virtual void setBlockedFaces(const surfaceScalarField &alphaf, const faceZone &fz, boolList &blockedFaces)
Set the blocked faces, i.e. where alpha > alpha threshold value.
virtual void collectParticle(const scalar time, const label regioni)
Collect particles that have passed through the faceZone.
labelList regions0_
Region indices in faceZone faces from last iteration.
word faceZoneName_
Name of faceZone to sample.
virtual void initialiseBins()
Initialise the particle collection bins.
virtual void checkFaceZone()
Check that the faceZone is valid.
label zoneID_
Index of the faceZone.
virtual bool execute()
Execute.
virtual bool write()
Write.
virtual bool read(const dictionary &)
Read the field min/max data.
virtual void calculateAddressing(const label nRegionsNew, const scalar time, labelList ®ionFaceIDs)
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
Computes the magnitude of an input field.
Base class for writing single files from the function objects.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Primarily stores particle properties so that it can be injected at a later time. Note that this store...
Primitive patch pair agglomerate method.
const labelList & restrictTopBottomAddressing() const
Return restriction from top level to bottom level.
void agglomerate()
Agglomerate patch.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
A patch is a list of labels that address the faces in the global face list.
label whichFace(const label l) const
Return label of face in patch from global face label.
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Splits a patch into regions based on a mask field. Result is a globally consistent label list of regi...
label nRegions() const noexcept
Return the global number of regions.
splitCell * master() const
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
autoPtr< surfaceVectorField > Uf
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
constexpr scalar pi(M_PI)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Ostream & endl(Ostream &os)
Add newline and flush stream.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a)
#define forAll(list, i)
Loop across all elements in list.
A non-counting (dummy) refCount.