33 template<
class CloudType>
36 const label globalPatchi
41 if (patchIndices_[i] == globalPatchi)
51 template<
class CloudType>
68 template<
class CloudType>
79 p_(this->coeffDict().template lookup<scalar>(
"p")),
80 psi_(this->coeffDict().template lookupOrDefault<scalar>(
"psi", 2.0)),
81 K_(this->coeffDict().template lookupOrDefault<scalar>(
"K", 2.0))
94 <<
"Cannot find any patch names matching " << patchName[i]
98 uniquePatchIDs.
insert(patchIDs);
101 patchIndices_ = uniquePatchIDs.
toc();
108 template<
class CloudType>
116 patchIndices_(pe.patchIndices_),
125 template<
class CloudType>
132 template<
class CloudType>
137 QPtr_->primitiveFieldRef() = 0.0;
149 this->owner().
name() +
":Q",
163 template<
class CloudType>
167 if (!
p.onBoundaryFace(mesh))
return;
173 if (localPatchi != -1)
177 this->owner().patchData(
p, pp, nw, Up);
188 const scalar magU =
mag(
U);
196 scalar& Q = QPtr_->boundaryFieldRef()[
patchi][patchFacei];
199 const scalar coeff =
p.nParticle()*
p.mass()*
sqr(magU)/(p_*psi_*K_);
#define forAll(list, i)
Loop across all elements in list.
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Templated cloud function object base class.
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
const fvMesh & mesh() const
Return references to the mesh.
Generic GeometricField class.
bool insert(const Key &key)
Insert a new entry.
List< Key > toc() const
Return the table of contents.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Function object to create a field of eroded volume, Q, on a specified list of patches....
virtual ~ParticleErosion()
Destructor.
virtual void preFace(const parcelType &p)
Pre-face hook.
label applyToPatch(const label globalPatchi) const
Returns local patchi if patch is in patchIds_ list.
virtual void write()
Write post-processing info.
virtual void preEvolve()
Pre-evolve hook.
ParticleErosion(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
bool empty() const
Return true if the UList is empty (ie, size() is zero)
A list of keyword definitions, which are a keyword followed by any number of values (e....
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
const polyMesh & mesh() const
Return reference to polyMesh.
label index() const
Return the index of this patch in the boundaryMesh.
wordList names() const
Return the list of patch names.
const polyBoundaryMesh & boundaryMesh() const
Return boundary 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.
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
A class for handling words, derived from string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define WarningInFunction
Report a warning using Foam::Warning.
dimensionedScalar tan(const dimensionedScalar &ds)
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 findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
word name(const bool)
Return a word representation of a bool.
errorManip< error > abort(error &err)
dimensionedScalar sin(const dimensionedScalar &ds)
const dimensionSet dimVolume
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)