36 template<
class CloudType>
47 particleFluxAccumulators_()
57 if (isType<polyPatch>(patch))
67 this->
coeffDict().subDict(
"numberDensities")
73 particleFluxAccumulators_.setSize(patches_.
size());
93 numberDensitiesDict.
lookup<scalar>(molecules[i]);
97 if (moleculeTypeIds_[i] == -1)
100 <<
"typeId " << molecules[i] <<
"not defined in cloud." <<
nl
105 numberDensities_ /=
cloud.nParticle();
112 template<
class CloudType>
119 template<
class CloudType>
140 template<
class CloudType>
154 label particlesInserted = 0;
158 cloud.boundaryT().boundaryField()
163 cloud.boundaryU().boundaryField()
166 label nLocateBoundaryHits = 0;
182 label typeId = moleculeTypeIds_[i];
184 scalar mass =
cloud.constProps(typeId).mass();
189 <<
"Zero boundary temperature detected, check boundaryT "
190 <<
"condition." <<
nl
196 cloud.maxwellianMostProbableSpeed
220 exp(-
sqr(sCosTheta)) + sqrtPi*sCosTheta*(1 +
erf(sCosTheta))
230 const face&
f = patch[pFI];
250 scalar previousCumulativeSum = 0.0;
258 + previousCumulativeSum;
260 previousCumulativeSum = cTriAFracs[triI];
265 cTriAFracs.
last() = 1.0;
282 scalar faceTemperature = boundaryT[
patchi][pFI];
288 scalar& faceAccumulator = pFA[i][pFI];
300 faceAccumulator -= nI;
302 label typeId = moleculeTypeIds_[i];
304 scalar mass =
cloud.constProps(typeId).mass();
306 for (
label i = 0; i < nI; i++)
314 label selectedTriI = -1;
320 if (cTriAFracs[triI] >= triSelection)
328 const tetIndices& faceTetIs = faceTets[selectedTriI];
334 scalar mostProbableSpeed
336 cloud.maxwellianMostProbableSpeed
343 scalar sCosTheta = (faceVelocity &
n)/mostProbableSpeed;
346 scalar uNormProbCoeffA =
347 sCosTheta +
sqrt(
sqr(sCosTheta) + 2.0);
349 scalar uNormProbCoeffB =
353 + sCosTheta*(sCosTheta -
sqrt(
sqr(sCosTheta) + 2.0))
357 scalar randomScaling = 3.0;
361 randomScaling =
mag(sCosTheta) + 1;
369 scalar uNormalThermal;
377 uNormal = uNormalThermal + sCosTheta;
385 P = 2.0*uNormal/uNormProbCoeffA
386 *
exp(uNormProbCoeffB -
sqr(uNormalThermal));
394 + (t1 & faceVelocity)*t1
395 + (t2 & faceVelocity)*t2
396 + mostProbableSpeed*uNormal*
n;
398 scalar Ei =
cloud.equipartitionInternalEnergy
401 cloud.constProps(typeId).internalDegreesOfFreedom()
421 if (nLocateBoundaryHits != 0)
424 <<
"Freestream inflow model for cloud " << this->owner().name()
425 <<
" did not accurately locate " << nLocateBoundaryHits
426 <<
" particles" <<
endl;
431 Info<<
" Particles inserted = "
432 << particlesInserted <<
endl;
#define forAll(list, i)
Loop across all elements in list.
Templated base class for dsmc cloud.
virtual void topoChange()
Update following mesh change.
FreeStream(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
virtual void inflow()
Introduce particles.
virtual ~FreeStream()
Destructor.
Generic GeometricBoundaryField class.
Templated inflow boundary model class.
const dictionary & coeffDict() const
Return the coefficients dictionary.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
void size(const label)
Override size to be inconsistent with allocated storage.
void setSize(const label)
Reset size of List.
scalar deltaTValue() const
Return time step value.
label size() const
Return the number of elements in the UList.
T & last()
Return the last element of the list.
A cloud is a collection of lagrangian particles.
A list of keyword definitions, which are a keyword followed by any number of values (e....
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
wordList toc() const
Return the table of contents.
Standard normal distribution. Not selectable.
virtual scalar sample() const
Sample the distribution.
A face is a list of labels corresponding to mesh vertices.
const Time & time() const
Return time.
static List< tetIndices > faceTetIndices(const polyMesh &mesh, label fI, label cI)
Return the tet decomposition of the given face, with.
Mesh consisting of general polyhedral cells.
virtual const labelList & faceOwner() const
Return face owner.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
virtual const pointField & points() const
Return raw points.
A patch is a list of labels that address the faces in the global face list.
const vectorField::subField faceAreas() const
Return face areas.
label start() const
Return start label of this patch in the polyMesh face list.
const vectorField::subField faceCentres() const
Return face centres.
const scalarField::subField magFaceAreas() const
Return face area magnitudes.
scalar scalar01()
Return a scalar uniformly distributed between zero and one.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
scalar mag() const
Return scalar magnitude.
Point randomPoint(randomGenerator &rndGen) const
Return a random point on the triangle from a uniform.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
dimensionedScalar exp(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.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
errorManip< error > abort(error &err)
dimensionedScalar erf(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
randomGenerator rndGen(653213)