33template<
class CloudType>
40 const word& patchName = patchData_[patchi].patchName();
42 forAll(nEscape_[patchi], injectori)
54template<
class CloudType>
62 patchData_(cloud.
mesh(), this->coeffDict()),
63 nEscape_(patchData_.size()),
64 massEscape_(nEscape_.size()),
65 nStick_(nEscape_.size()),
66 massStick_(nEscape_.size()),
67 writeFields_(this->coeffDict().getOrDefault(
"writeFields", false)),
69 massEscapePtr_(nullptr),
70 massStickPtr_(nullptr)
72 const bool outputByInjectorId
77 Info<<
" Interaction fields will be written to "
87 Info<<
" Interaction fields will not be written" <<
endl;
92 if (outputByInjectorId)
94 for (
const auto& inj : cloud.injectors())
96 injIdToIndex_.insert(inj.injectorID(), nInjectors++);
107 forAll(patchData_, patchi)
109 const word& interactionTypeName =
110 patchData_[patchi].interactionTypeName();
116 const word& patchName = patchData_[patchi].patchName();
118 <<
"Unknown patch interaction type "
119 << interactionTypeName <<
" for patch " << patchName
120 <<
". Valid selections are:"
125 nEscape_[patchi].setSize(nInjectors,
Zero);
126 massEscape_[patchi].setSize(nInjectors,
Zero);
127 nStick_[patchi].setSize(nInjectors,
Zero);
128 massStick_[patchi].setSize(nInjectors,
Zero);
133template<
class CloudType>
136 const LocalInteraction<CloudType>& pim
140 patchData_(pim.patchData_),
141 nEscape_(pim.nEscape_),
142 massEscape_(pim.massEscape_),
143 nStick_(pim.nStick_),
144 massStick_(pim.massStick_),
145 writeFields_(pim.writeFields_),
146 injIdToIndex_(pim.injIdToIndex_),
147 massEscapePtr_(nullptr),
148 massStickPtr_(nullptr)
154template<
class CloudType>
168 mesh.time().timeName(),
180 return *massEscapePtr_;
184template<
class CloudType>
189 const fvMesh&
mesh = this->owner().mesh();
198 mesh.time().timeName(),
210 return *massStickPtr_;
214template<
class CloudType>
222 const label patchi = patchData_.applyToPatch(
pp.index());
232 ? injIdToIndex_.lookup(
p.typeId(), 0)
237 this->wordToInteractionType
239 patchData_[patchi].interactionTypeName()
250 keepParticle =
false;
254 const scalar dm =
p.mass()*
p.nParticle();
256 nEscape_[patchi][idx]++;
257 massEscape_[patchi][idx] += dm;
261 const label pI =
pp.index();
262 const label fI =
pp.whichFace(
p.face());
263 massEscape().boundaryFieldRef()[pI][fI] += dm;
273 const scalar dm =
p.mass()*
p.nParticle();
275 nStick_[patchi][idx]++;
276 massStick_[patchi][idx] += dm;
280 const label pI =
pp.index();
281 const label fI =
pp.whichFace(
p.face());
282 massStick().boundaryFieldRef()[pI][fI] += dm;
294 this->owner().patchData(
p,
pp, nw, Up);
299 if (
mag(Up) > 0 &&
mag(
U) < this->Urmax())
302 <<
"Particle U the same as patch "
303 <<
" The particle has been removed" <<
nl <<
endl;
305 keepParticle =
false;
316 U -= (1.0 + patchData_[patchi].e())*Un*nw;
319 U -= patchData_[patchi].mu()*Ut;
329 <<
"Unknown interaction type "
330 << patchData_[patchi].interactionTypeName()
331 <<
"(" << it <<
") for patch "
332 << patchData_[patchi].patchName()
333 <<
". Valid selections are:" << this->interactionTypeNames_
345template<
class CloudType>
356 forAll(patchData_, patchi)
358 label lsd = nEscape_[patchi].size();
359 npe0[patchi].setSize(lsd,
Zero);
360 mpe0[patchi].setSize(lsd,
Zero);
361 nps0[patchi].setSize(lsd,
Zero);
362 mps0[patchi].setSize(lsd,
Zero);
366 this->getModelProperty(
"nEscape", npe0);
367 this->getModelProperty(
"massEscape", mpe0);
368 this->getModelProperty(
"nStick", nps0);
369 this->getModelProperty(
"massStick", mps0);
376 npe[i] = npe[i] + npe0[i];
383 mpe[i] = mpe[i] + mpe0[i];
390 nps[i] = nps[i] + nps0[i];
397 mps[i] = mps[i] + mps0[i];
400 if (injIdToIndex_.size())
404 labelList indexToInjector(injIdToIndex_.size());
407 indexToInjector[iter.val()] = iter.key();
410 forAll(patchData_, patchi)
412 forAll(mpe[patchi], indexi)
414 const word& patchName = patchData_[patchi].patchName();
416 Log_<<
" Parcel fate: patch " << patchName
417 <<
" (number, mass)" <<
nl
418 <<
" - escape (injector " << indexToInjector[indexi]
419 <<
" ) = " << npe[patchi][indexi]
420 <<
", " << mpe[patchi][indexi] <<
nl
421 <<
" - stick (injector " << indexToInjector[indexi]
422 <<
" ) = " << nps[patchi][indexi]
423 <<
", " << mps[patchi][indexi] <<
nl;
429 forAll(patchData_, patchi)
431 const word& patchName = patchData_[patchi].patchName();
433 Log_<<
" Parcel fate: patch " << patchName
434 <<
" (number, mass)" <<
nl
436 << npe[patchi][0] <<
", " << mpe[patchi][0] <<
nl
438 << nps[patchi][0] <<
", " << mps[patchi][0] <<
nl;
444 forAll(npe[patchi], injectori)
447 <<
tab << npe[patchi][injectori]
448 <<
tab << mpe[patchi][injectori]
449 <<
tab << nps[patchi][injectori]
450 <<
tab << mps[patchi][injectori];
454 this->file() <<
endl;
456 if (this->writeTime())
458 this->setModelProperty(
"nEscape", npe);
459 this->setModelProperty(
"massEscape", mpe);
460 this->setModelProperty(
"nStick", nps);
461 this->setModelProperty(
"massStick", mps);
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const CloudType & owner() const
Return const access to the owner cloud.
virtual bool writeTime() const
Flag to indicate when to write a property.
@ REGISTER
Request registration (bool: true).
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
void setSize(label n)
Alias for resize().
LocalInteraction(const dictionary &dict, CloudType &owner)
Construct from dictionary.
volScalarField & massStick()
Return access to the massStick field.
volScalarField & massEscape()
Return access to the massEscape field.
virtual bool correct(typename CloudType::parcelType &p, const polyPatch &pp, bool &keepParticle)
Apply velocity correction.
virtual void info()
Write patch interaction info.
virtual void writeFileHeader(Ostream &os)
Output file header information.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
const scalar & Urmax() const
Return Urmax.
static wordList interactionTypeNames_
PatchInteractionModel(CloudType &owner)
Construct null from owner.
virtual void info()
Write patch interaction info.
static interactionType wordToInteractionType(const word &itWord)
Convert word to interaction result.
virtual void writeFileHeader(Ostream &os)
Output file header information.
static void listGather(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather (reduce) list elements, applying bop to each list element.
A cloud is a registry collection of lagrangian particles.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
virtual OFstream & file()
Return access to the file (if only 1).
Mesh data needed to do the Finite Volume discretisation.
A patch is a list of labels that address the faces in the global face list.
bool getModelProperty(const word &entryName, Type &value) const
Retrieve generic property from the sub-model.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
#define Log_
Report write to Foam::Info if the class log switch is true.
#define WarningInFunction
Report a warning using Foam::Warning.
List< scalarList > scalarListList
List of scalarList.
DSMCCloud< dsmcParcel > CloudType
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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).
constexpr char tab
The tab '\t' character(0x09).
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.