Go to the documentation of this file.
35 template<
class CloudType>
46 particleFluxAccumulators_()
56 if (isType<polyPatch>(patch))
66 this->coeffDict().subDict(
"numberDensities")
72 particleFluxAccumulators_.setSize(patches_.size());
87 numberDensities_.setSize(molecules.
size());
93 numberDensitiesDict.
lookup(molecules[i])
98 if (moleculeTypeIds_[i] == -1)
101 <<
"typeId " << molecules[i] <<
"not defined in cloud." <<
nl
106 numberDensities_ /=
cloud.nParticle();
113 template<
class CloudType>
120 template<
class CloudType>
135 pFA[facei].
setSize(patch.size(), 0);
141 template<
class CloudType>
148 const scalar deltaT =
mesh.time().deltaTValue();
154 label particlesInserted = 0;
158 cloud.boundaryT().boundaryField()
163 cloud.boundaryU().boundaryField()
181 label typeId = moleculeTypeIds_[i];
183 scalar mass =
cloud.constProps(typeId).mass();
188 <<
"Zero boundary temperature detected, check boundaryT "
189 <<
"condition." <<
nl
195 cloud.maxwellianMostProbableSpeed
219 exp(-
sqr(sCosTheta)) + sqrtPi*sCosTheta*(1 +
erf(sCosTheta))
229 const face&
f = patch[pFI];
233 label cellI =
mesh.faceOwner()[globalFaceIndex];
249 scalar previousCummulativeSum = 0.0;
257 + previousCummulativeSum;
259 previousCummulativeSum = cTriAFracs[triI];
264 cTriAFracs.last() = 1.0;
281 scalar faceTemperature = boundaryT[
patchi][pFI];
287 scalar& faceAccumulator = pFA[i][pFI];
294 if ((faceAccumulator - nI) >
rndGen.scalar01())
299 faceAccumulator -= nI;
301 label typeId = moleculeTypeIds_[i];
303 scalar mass =
cloud.constProps(typeId).mass();
305 for (
label i = 0; i < nI; i++)
310 scalar triSelection =
rndGen.scalar01();
313 label selectedTriI = -1;
319 if (cTriAFracs[triI] >= triSelection)
327 const tetIndices& faceTetIs = faceTets[selectedTriI];
333 scalar mostProbableSpeed
335 cloud.maxwellianMostProbableSpeed
342 scalar sCosTheta = (faceVelocity &
n)/mostProbableSpeed;
345 scalar uNormProbCoeffA =
346 sCosTheta +
sqrt(
sqr(sCosTheta) + 2.0);
348 scalar uNormProbCoeffB =
352 + sCosTheta*(sCosTheta -
sqrt(
sqr(sCosTheta) + 2.0))
356 scalar randomScaling = 3.0;
360 randomScaling =
mag(sCosTheta) + 1;
368 scalar uNormalThermal;
374 randomScaling*(2.0*
rndGen.scalar01() - 1);
376 uNormal = uNormalThermal + sCosTheta;
384 P = 2.0*uNormal/uNormProbCoeffA
385 *
exp(uNormProbCoeffB -
sqr(uNormalThermal));
388 }
while (P <
rndGen.scalar01());
396 + (t1 & faceVelocity)*t1
397 + (t2 & faceVelocity)*t2
398 + mostProbableSpeed*uNormal*
n;
400 scalar Ei =
cloud.equipartitionInternalEnergy
403 cloud.constProps(typeId).internalDegreesOfFreedom()
425 Info<<
" Particles inserted = "
426 << particlesInserted <<
endl;
Simple random number generator.
Point randomPoint(Random &rndGen) const
Return a random point on the triangle from a uniform.
const dimensionedScalar k
Boltzmann constant.
#define forAll(list, i)
Loop across all elements in list.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar exp(const dimensionedScalar &ds)
virtual ~FreeStream()
Destructor.
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar erf(const dimensionedScalar &ds)
label tetPt() const
Return the characterising tetPtI.
Mesh consisting of general polyhedral cells.
virtual void autoMap(const mapPolyMesh &)
Remap the particles to the correct cells following mesh change.
virtual void inflow()
Introduce particles.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the.
A patch is a list of labels that address the faces in the global face list.
FreeStream(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
Templated base class for dsmc cloud.
A list of keyword definitions, which are a keyword followed by any number of values (e....
errorManip< error > abort(error &err)
label start() const
Return start label of this patch in the polyMesh face list.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
scalar mag() const
Return scalar magnitude.
void setSize(const label)
Reset size of List.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
A cloud is a collection of lagrangian particles.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
const vectorField::subField faceCentres() const
Return face centres.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const vectorField::subField faceAreas() const
Return face normals.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A face is a list of labels corresponding to mesh vertices.
wordList toc() const
Return the table of contents.
void size(const label)
Override size to be inconsistent with allocated storage.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
cachedRandom rndGen(label(0), -1)
static List< tetIndices > faceTetIndices(const polyMesh &mesh, label fI, label cI)
Return the tet decomposition of the given face, with.