30 template<
class CloudType>
33 const label globalPatchi
38 if (patchIDs_[i] == globalPatchi)
48 template<
class CloudType>
65 template<
class CloudType>
76 p_(this->coeffDict().
template lookup<scalar>(
"p")),
77 psi_(this->coeffDict().
template lookupOrDefault<scalar>(
"psi", 2.0)),
78 K_(this->coeffDict().
template lookupOrDefault<scalar>(
"K", 2.0))
80 const wordList allPatchNames = owner.mesh().boundaryMesh().names();
91 <<
"Cannot find any patch names matching " << patchName[i]
95 uniquePatchIDs.
insert(patchIDs);
98 patchIDs_ = uniquePatchIDs.
toc();
105 template<
class CloudType>
113 patchIDs_(pe.patchIDs_),
122 template<
class CloudType>
129 template<
class CloudType>
134 QPtr_->primitiveFieldRef() = 0.0;
146 this->owner().
name() +
":Q",
149 IOobject::READ_IF_PRESENT,
160 template<
class CloudType>
169 const label localPatchi = applyToPatch(patchi);
171 if (localPatchi != -1)
175 this->owner().patchData(p, pp, nw, Up);
186 const scalar magU =
mag(U);
187 const vector UHat = U/magU;
194 scalar& Q = QPtr_->boundaryFieldRef()[
patchi][patchFacei];
197 const scalar coeff = p.nParticle()*p.mass()*
sqr(magU)/(p_*psi_*K_);
198 if (
tan(alpha) < K_/6)
200 Q += coeff*(
sin(2*alpha) - 6/K_*
sqr(
sin(alpha)));
204 Q += coeff*(K_*
sqr(
cos(alpha))/6);
dimensionedScalar acos(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
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...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
ParticleErosion(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Function object to create a field of eroded volume, Q, on a specified list of patches. The volume is calculated by the model of Finnie et al. The implementation follows the description given by the review of Yadav et al.
Ostream & endl(Ostream &os)
Add newline and flush stream.
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Holds information (coordinate and normal) regarding nearest wall point.
virtual void preEvolve()
Pre-evolve hook.
const Time & time() const
Return the top-level database.
bool insert(const Key &key)
Insert a new entry.
virtual ~ParticleErosion()
Destructor.
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar cos(const dimensionedScalar &ds)
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
A class for handling words, derived from string.
virtual void write()
Write post-processing info.
errorManip< error > abort(error &err)
dimensionedScalar sin(const dimensionedScalar &ds)
word name(const complex &)
Return a string representation of a complex.
#define WarningInFunction
Report a warning using Foam::Warning.
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.
const dimensionSet dimVolume
dimensioned< scalar > mag(const dimensioned< Type > &)
List< Key > toc() const
Return the table of contents.
label applyToPatch(const label globalPatchi) const
Returns local patchi if patch is in patchIds_ list.
virtual void postPatch(const parcelType &p, const polyPatch &pp, bool &keepParticle)
Post-patch hook.
A patch is a list of labels that address the faces in the global face list.
dimensionedScalar tan(const dimensionedScalar &ds)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Templated base class for dsmc cloud.
Templated cloud function object base class.
label whichFace(const label l) const
Return label of face in patch from global face label.