37template<
class CloudType>
48 particleFluxAccumulators_()
68 this->
coeffDict().subDict(
"numberDensities")
74 particleFluxAccumulators_.setSize(patches_.size());
87 moleculeTypeIds_.setSize(molecules.
size());
89 numberDensities_.setSize(molecules.
size());
93 numberDensities_[i] = numberDensitiesDict.
get<scalar>(molecules[i]);
95 moleculeTypeIds_[i] =
cloud.typeIdList().
find(molecules[i]);
97 if (moleculeTypeIds_[i] == -1)
100 <<
"typeId " << molecules[i] <<
"not defined in cloud." <<
nl
105 numberDensities_ /=
cloud.nParticle();
112template<
class CloudType>
119template<
class CloudType>
127 label patchi = patches_[
p];
140template<
class CloudType>
153 label particlesInserted = 0;
157 cloud.boundaryT().boundaryField()
162 cloud.boundaryU().boundaryField()
168 label patchi = patches_[
p];
180 label typeId = moleculeTypeIds_[i];
182 scalar mass =
cloud.constProps(typeId).mass();
184 if (
min(boundaryT[patchi]) < SMALL)
187 <<
"Zero boundary temperature detected, check boundaryT "
188 <<
"condition." <<
nl
194 cloud.maxwellianMostProbableSpeed
208 (boundaryU[patchi] & -
patch.faceAreas()/
mag(
patch.faceAreas()))
215 mag(
patch.faceAreas())*numberDensities_[i]*deltaT
218 exp(-
sqr(sCosTheta)) + sqrtPi*sCosTheta*(1 +
erf(sCosTheta))
230 label globalFaceIndex = pFI +
patch.start();
236 scalar fA =
mag(
patch.faceAreas()[pFI]);
248 scalar previousCummulativeSum = 0.0;
256 + previousCummulativeSum;
258 previousCummulativeSum = cTriAFracs[triI];
263 cTriAFracs.last() = 1.0;
277 scalar faceTemperature = boundaryT[patchi][pFI];
279 const vector& faceVelocity = boundaryU[patchi][pFI];
283 scalar& faceAccumulator = pFA[i][pFI];
286 label nI =
max(label(faceAccumulator), 0);
295 faceAccumulator -= nI;
297 label typeId = moleculeTypeIds_[i];
299 scalar mass =
cloud.constProps(typeId).mass();
301 for (label i = 0; i < nI; i++)
309 label selectedTriI = -1;
315 if (cTriAFracs[triI] >= triSelection)
323 const tetIndices& faceTetIs = faceTets[selectedTriI];
329 scalar mostProbableSpeed
331 cloud.maxwellianMostProbableSpeed
338 scalar sCosTheta = (faceVelocity &
n)/mostProbableSpeed;
341 scalar uNormProbCoeffA =
342 sCosTheta +
sqrt(
sqr(sCosTheta) + 2.0);
344 scalar uNormProbCoeffB =
348 + sCosTheta*(sCosTheta -
sqrt(
sqr(sCosTheta) + 2.0))
352 scalar randomScaling = 3.0;
356 randomScaling =
mag(sCosTheta) + 1;
364 scalar uNormalThermal;
372 uNormal = uNormalThermal + sCosTheta;
380 P = 2.0*uNormal/uNormProbCoeffA
381 *
exp(uNormProbCoeffB -
sqr(uNormalThermal));
387 sqrt(physicoChemical::k.value()*faceTemperature/mass)
392 + (t1 & faceVelocity)*t1
393 + (t2 & faceVelocity)*t2
394 + mostProbableSpeed*uNormal*
n;
396 scalar Ei =
cloud.equipartitionInternalEnergy
399 cloud.constProps(typeId).internalDegreesOfFreedom()
402 cloud.addNewParcel(
p, celli,
U, Ei, typeId);
410 Info<<
" Particles inserted = "
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
FreeStream(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
virtual void inflow()
Introduce particles.
virtual void autoMap(const mapPolyMesh &)
Remap the particles to the correct cells following mesh change.
virtual ~FreeStream()
Destructor.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
void transfer(HashTable< T, Key, Hash > &rhs)
Transfer contents into this table.
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
InflowBoundaryModel(CloudType &owner)
Construct null from owner.
const dictionary & coeffDict() const
Return the coefficients dictionary.
const dictionary & dict() const
Return the owner cloud dictionary.
const CloudType & owner() const
Return const access the owner cloud object.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void setSize(label n)
Alias for resize().
Type GaussNormal()
Return a sample whose components are normally distributed with zero mean and unity variance N(0,...
Type sample01()
Return a sample whose components lie in the range [0,1].
scalar deltaTValue() const noexcept
Return time step value.
void size(const label n)
Older name for setAddressableSize.
T & last()
Access last element of the list, position [size()-1].
A cloud is a registry collection of lagrangian particles.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
wordList toc() const
Return the table of contents.
A face is a list of labels corresponding to mesh vertices.
const Time & time() const
Return the top-level database.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
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.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
virtual const labelList & faceOwner() const
Return face owner.
virtual const pointField & points() const
Return raw points.
A patch is a list of labels that address the faces in the global face list.
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
The triangle geometry for the face for this tet. The normal of the tri points out of the cell.
scalar mag() const
The magnitude of the triangle area.
Point randomPoint(Random &rndGen) const
Return a random point on the triangle from a uniform distribution.
const polyBoundaryMesh & patches
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
constexpr scalar pi(M_PI)
const dimensionedScalar k
Boltzmann constant.
const std::string patch
OpenFOAM patch number as a std::string.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
dimensionedScalar erf(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
errorManip< error > abort(error &err)
vector point
Point is a vector.
static constexpr const zero Zero
Global zero (0).
bool isType(const U &obj)
Check if typeid of the object and Type are identical.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.