30 template<
class CloudType>
38 patchData_(cloud.mesh(), this->coeffDict()),
39 nEscape_(patchData_.size(), 0),
40 massEscape_(patchData_.size(), 0.0),
41 nStick_(patchData_.size(), 0),
42 massStick_(patchData_.size(), 0.0),
43 writeFields_(this->coeffDict().lookupOrDefault(
"writeFields",
false)),
44 massEscapePtr_(
nullptr),
45 massStickPtr_(
nullptr)
49 word massEscapeName(this->owner().
name() +
":massEscape");
50 word massStickName(this->owner().
name() +
":massStick");
51 Info<<
" Interaction fields will be written to " << massEscapeName
52 <<
" and " << massStickName <<
endl;
59 Info<<
" Interaction fields will not be written" <<
endl;
65 const word& interactionTypeName =
66 patchData_[
patchi].interactionTypeName();
68 this->wordToInteractionType(interactionTypeName);
72 const word& patchName = patchData_[
patchi].patchName();
74 <<
"Unknown patch interaction type " 75 << interactionTypeName <<
" for patch " << patchName
76 <<
". Valid selections are:" 84 template<
class CloudType>
91 patchData_(pim.patchData_),
92 nEscape_(pim.nEscape_),
93 massEscape_(pim.massEscape_),
95 massStick_(pim.massStick_),
96 writeFields_(pim.writeFields_),
97 massEscapePtr_(
nullptr),
98 massStickPtr_(
nullptr)
104 template<
class CloudType>
111 template<
class CloudType>
114 if (!massEscapePtr_.valid())
124 this->owner().
name() +
":massEscape",
127 IOobject::READ_IF_PRESENT,
136 return massEscapePtr_();
140 template<
class CloudType>
143 if (!massStickPtr_.valid())
153 this->owner().
name() +
":massStick",
156 IOobject::READ_IF_PRESENT,
165 return massStickPtr_();
169 template<
class CloudType>
182 bool& moving = p.moving();
185 this->wordToInteractionType
187 patchData_[patchi].interactionTypeName()
198 scalar dm = p.mass()*p.nParticle();
200 keepParticle =
false;
204 massEscape_[
patchi] += dm;
209 massEscape().boundaryFieldRef()[pI][fI] += dm;
215 scalar dm = p.mass()*p.nParticle();
226 massStick().boundaryFieldRef()[pI][fI] += dm;
238 this->owner().patchData(p, pp, nw, Up);
248 U -= (1.0 + patchData_[
patchi].e())*Un*nw;
251 U -= patchData_[
patchi].mu()*Ut;
261 <<
"Unknown interaction type " 262 << patchData_[
patchi].interactionTypeName()
263 <<
"(" << it <<
") for patch " 264 << patchData_[
patchi].patchName()
265 <<
". Valid selections are:" << this->interactionTypeNames_
277 template<
class CloudType>
282 this->getModelProperty(
"nEscape", npe0);
285 this->getModelProperty(
"massEscape", mpe0);
288 this->getModelProperty(
"nStick", nps0);
291 this->getModelProperty(
"massStick", mps0);
313 os <<
" Parcel fate (number, mass) : patch " 314 << patchData_[i].patchName() <<
nl 315 <<
" - escape = " << npe[i]
316 <<
", " << mpe[i] <<
nl 317 <<
" - stick = " << nps[i]
318 <<
", " << mps[i] <<
nl;
321 if (this->writeTime())
323 this->setModelProperty(
"nEscape", npe);
326 this->setModelProperty(
"massEscape", mpe);
329 this->setModelProperty(
"nStick", nps);
332 this->setModelProperty(
"massStick", mps);
Patch interaction specified on a patch-by-patch basis.
volScalarField & massEscape()
Return access to the massEscape field.
#define forAll(list, i)
Loop across all elements in list.
errorManipArg< error, int > exit(error &err, const int errNo=1)
A list of keyword definitions, which are a keyword followed by any number of values (e...
virtual void info(Ostream &os)
Write patch interaction info to stream.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual ~LocalInteraction()
Destructor.
const Time & time() const
Return the top-level database.
Templated patch interaction model class.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
A class for handling words, derived from string.
LocalInteraction(const dictionary &dict, CloudType &owner)
Construct from dictionary.
errorManip< error > abort(error &err)
ParcelType parcelType
Type of parcel the cloud was instantiated for.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
volScalarField & massStick()
Return access to the massStick field.
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
label index() const
Return the index of this patch in the boundaryMesh.
A patch is a list of labels that address the faces in the global face list.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Templated base class for dsmc cloud.
virtual bool correct(typename CloudType::parcelType &p, const polyPatch &pp, bool &keepParticle)
Apply velocity correction.
label whichFace(const label l) const
Return label of face in patch from global face label.