33#include "surfaceInterpolate.H"
72 const label nFaces = fz.
size();
73 const label allFaces =
returnReduce(nFaces, sumOp<label>());
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_];
114 const label nFaces = fz.
size();
115 label nLocations = nInjectorLocations_;
120 scalar fraction = scalar(nFaces)/scalar(nGlobalFaces);
121 nLocations = ceil(fraction*nInjectorLocations_);
124 Pout<<
"nFaces:" << nFaces
125 <<
", nGlobalFaces:" << nGlobalFaces
126 <<
", fraction:" << fraction
127 <<
", nLocations:" << nLocations
132 pairPatchAgglomeration ppa
145 label nCoarseFaces = 0;
148 fineToCoarseAddr_ = ppa.restrictTopBottomAddressing();
149 nCoarseFaces =
max(fineToCoarseAddr_) + 1;
152 globalCoarseFaces_ = globalIndex(nCoarseFaces);
155 <<
" coarse faces" <<
endl;
188 patchIDs_.setSize(fz.
size(), -1);
189 patchFaceIDs_.setSize(fz.
size(), -1);
191 label nBlockedFaces = 0;
194 const label facei = fz[localFacei];
196 if (mesh_.isInternalFace(facei))
198 if (alphaf[facei] > alphaThreshold_)
200 blockedFaces[localFacei] =
true;
205 label patchi = mesh_.boundaryMesh().whichPatch(facei);
206 label patchFacei = -1;
208 const polyPatch&
pp = mesh_.boundaryMesh()[patchi];
214 patchFacei = (cpp->owner() ?
pp.whichFace(facei) : -1);
218 patchFacei =
pp.whichFace(facei);
221 if (patchFacei == -1)
225 else if (alphafp[patchFacei] > alphaThreshold_)
227 blockedFaces[localFacei] =
true;
230 patchIDs_[localFacei] = patchi;
231 patchFaceIDs_[localFacei] = patchFacei;
247 const label particlei = regionToParticleMap_[regioni];
250 if (
p.faceIHit != -1 && nInjectorLocations_)
253 label coarseFacei = fineToCoarseAddr_[
p.faceIHit];
254 p.faceIHit = globalCoarseFaces_.toGlobal(coarseFacei);
261 if ((pDiameter > minDiameter_) && (pDiameter < maxDiameter_))
266 const point position =
p.VC/(
p.V + ROOTVSMALL);
267 const vector U =
p.VU/(
p.V + ROOTVSMALL);
269 if (nInjectorLocations_)
274 injectedParticle* ip =
new injectedParticle
285 cloud_.addParticle(ip);
287 collectedVolume_ +=
p.V;
290 ++nCollectedParticles_;
303 const label nNewRegions,
312 labelList oldToNewRegion(particles_.size(), -1);
315 forAll(regionFaceIDs, facei)
317 label newRegioni = regionFaceIDs[facei];
318 label oldRegioni = regions0_[facei];
320 if (newRegioni != -1 && oldRegioni != -1)
324 newToNewRegion[newRegioni] =
325 max(newRegioni, oldToNewRegion[oldRegioni]);
326 oldToNewRegion[oldRegioni] = newRegioni;
334 label nParticle = -1;
336 Map<label> newRegionToParticleMap;
337 forAll(newToNewRegion, newRegioni0)
339 label newRegioni = newToNewRegion[newRegioni0];
340 if (newRegions.insert(newRegioni))
346 newRegionToParticleMap.insert(newRegioni0, nParticle);
353 List<eulerianParticle> newParticles(newRegionToParticleMap.size());
354 forAll(oldToNewRegion, oldRegioni)
356 label newRegioni = oldToNewRegion[oldRegioni];
357 if (newRegioni == -1)
361 <<
"Collecting particle from oldRegion:" << oldRegioni
364 collectParticle(time, oldRegioni);
369 label newParticlei = newRegionToParticleMap[newRegioni];
370 label oldParticlei = regionToParticleMap_[oldRegioni];
373 <<
"Combining newRegioni: " << newRegioni
374 <<
"(p:" << newParticlei <<
") and "
375 <<
"oldRegioni: " << oldRegioni
376 <<
"(p:" << oldParticlei <<
")"
379 newParticles[newParticlei] =
382 newParticles[newParticlei],
383 particles_[oldParticlei]
390 regionToParticleMap_ = newRegionToParticleMap;
411 const scalar deltaT = mesh_.time().deltaTValue();
412 const pointField& faceCentres = mesh_.faceCentres();
414 forAll(regionFaceIDs, localFacei)
416 const label newRegioni = regionFaceIDs[localFacei];
417 if (newRegioni != -1)
419 const label particlei = regionToParticleMap_[newRegioni];
420 const label meshFacei = fz[localFacei];
426 p.faceIHit = localFacei;
433 scalar magPhii =
mag(faceValue(
phi, localFacei, meshFacei));
434 vector Ufi = faceValue(
Uf, localFacei, meshFacei);
435 scalar dV = magPhii*deltaT;
437 p.VC += dV*faceCentres[meshFacei];
455 cloud_(mesh_,
"eulerianParticleCloud"),
461 alphaThreshold_(0.1),
465 nInjectorLocations_(0),
467 globalCoarseFaces_(),
470 regionToParticleMap_(),
471 minDiameter_(ROOTVSMALL),
473 nCollectedParticles_(getProperty<label>(
"nCollectedParticles", 0)),
474 collectedVolume_(getProperty<scalar>(
"collectedVolume", 0)),
475 nDiscardedParticles_(getProperty<label>(
"nDiscardedParticles", 0)),
476 discardedVolume_(getProperty<scalar>(
"discardedVolume", 0))
481 << name <<
" function object only applicable to 3-D cases"
500 dict.readEntry(
"faceZone", faceZoneName_);
501 dict.readEntry(
"alpha", alphaName_);
503 dict.readIfPresent(
"alphaThreshold", alphaThreshold_);
504 dict.readIfPresent(
"U", UName_);
505 dict.readIfPresent(
"rho", rhoName_);
506 dict.readIfPresent(
"phi", phiName_);
507 dict.readIfPresent(
"nLocations", nInjectorLocations_);
508 dict.readIfPresent(
"minDiameter", minDiameter_);
509 dict.readIfPresent(
"maxDiameter", maxDiameter_);
513 if (nInjectorLocations_)
541 const faceZone& fz = mesh_.faceZones()[zoneID_];
550 setBlockedFaces(talphaf(), fz, blockedFaces);
557 const label nRegionsNew = regionFaceIDs.
nRegions();
565 mesh_.time().value(),
571 accumulateParticleInfo(talphaf(), tphi(), regionFaceIDs, fz);
573 Log <<
" Collected particles : " << nCollectedParticles_ <<
nl
574 <<
" Collected volume : " << collectedVolume_ <<
nl
575 <<
" Discarded particles : " << nDiscardedParticles_ <<
nl
576 <<
" Discarded volume : " << discardedVolume_ <<
nl
590 setProperty(
"nCollectedParticles", nCollectedParticles_);
591 setProperty(
"collectedVolume", collectedVolume_);
592 setProperty(
"nDiscardedParticles", nDiscardedParticles_);
593 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.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvsPatchField< scalar >::calculatedType())
const Boundary & boundaryField() const noexcept
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.
@ NO_REGISTER
Do not request registration (bool: false).
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
A List with indirect addressing.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
A HashTable to objects of type <T> with a label key.
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
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(const bool on) noexcept
Set as parallel run on/off.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
label size() const noexcept
The number of entries in the list.
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,...
Lightweight class to store particle data derived from VOF calculations, with special handling for inp...
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 bool read(const dictionary &dict)
Read and set the function object if its data have changed.
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 nDiscardedParticles_
Total number of discarded particles.
label nInjectorLocations_
Number of sample locations to generate.
scalar alphaThreshold_
Value of phase fraction used to identify particle boundaries.
scalar minDiameter_
Minimum diameter (optional).
scalar collectedVolume_
Total collected volume [m3].
virtual void setBlockedFaces(const surfaceScalarField &alphaf, const faceZone &fz, boolList &blockedFaces)
Set the blocked faces, i.e. where alpha > alpha threshold value.
word UName_
Name of the velocity field, default = U.
globalIndex globalCoarseFaces_
Global coarse face addressing.
label nCollectedParticles_
Total number of collected particles.
Type faceValue(const GeometricField< Type, fvsPatchField, surfaceMesh > &field, const label localFaceI, const label globalFaceI) const
virtual void collectParticle(const scalar time, const label regioni)
Collect particles that have passed through the faceZone.
labelList patchFaceIDs_
Patch face indices where faceZone face intersect patch.
extractEulerianParticles(const word &name, const Time &runTime, const dictionary &dict)
Construct from components.
scalar maxDiameter_
Maximum diameter (optional).
injectedParticleCloud cloud_
Storage for collected particles.
Map< label > regionToParticleMap_
Map from region to index in particles_ list.
labelList regions0_
Region indices in faceZone faces from last iteration.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
word phiName_
Name of the flux field, default ="rho".
word rhoName_
Name of the density field, default = rho.
word faceZoneName_
Name of faceZone to sample.
List< eulerianParticle > particles_
Particle properties (partial, being accumulated).
virtual void initialiseBins()
Initialise the particle collection bins.
scalar discardedVolume_
Total discarded volume [m3].
word alphaName_
Name of phase fraction field.
virtual void checkFaceZone()
Check that the faceZone is valid.
label zoneID_
Index of the faceZone.
labelList fineToCoarseAddr_
Agglomeration addressing from fine to coarse (local proc only).
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
labelList patchIDs_
Patch indices where faceZone face intersect patch.
virtual void calculateAddressing(const label nRegionsNew, const scalar time, labelList ®ionFaceIDs)
Calculate the addressing between regions between iterations Returns the number of active regions (par...
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
void setProperty(const word &entryName, const Type &value)
Add generic property.
Type getProperty(const word &entryName, const Type &defaultValue=Type(Zero)) const
Retrieve generic property.
const Time & time() const
Return time database.
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
virtual bool read(const dictionary &dict)
Read.
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
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.
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.
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.
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
constexpr scalar pi(M_PI)
Namespace for handling debugging switches.
const std::string patch
OpenFOAM patch number as a std::string.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
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.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
vector point
Point is a vector.
List< bool > boolList
A List of bools.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
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.