41template<
class ParcelType>
42void Foam::DSMCCloud<ParcelType>::buildConstProps()
44 Info<<
nl <<
"Constructing constant properties for" <<
endl;
45 constProps_.setSize(typeIdList_.size());
49 particleProperties_.subDict(
"moleculeProperties")
54 const word&
id = typeIdList_[i];
58 const dictionary& molDict(moleculeProperties.subDict(
id));
60 constProps_[i] =
typename ParcelType::constantProperties(molDict);
65template<
class ParcelType>
66void Foam::DSMCCloud<ParcelType>::buildCellOccupancy()
68 for (
auto& list : cellOccupancy_)
73 for (ParcelType&
p : *
this)
75 cellOccupancy_[
p.cell()].append(&
p);
80template<
class ParcelType>
81void Foam::DSMCCloud<ParcelType>::initialise
88 const scalar temperature
90 dsmcInitialiseDict.
get<scalar>(
"temperature")
97 dsmcInitialiseDict.
subDict(
"numberDensities")
100 List<word> molecules(numberDensitiesDict.toc());
106 numberDensities[i] = numberDensitiesDict.get<scalar>(molecules[i]);
109 numberDensities /= nParticle_;
111 forAll(mesh_.cells(), celli)
123 scalar tetVolume = tet.
mag();
127 const word& moleculeName(molecules[i]);
129 label typeId = typeIdList_.find(moleculeName);
134 <<
"typeId " << moleculeName <<
"not defined." <<
nl
138 const typename ParcelType::constantProperties& cP =
141 scalar numberDensity = numberDensities[i];
144 scalar particlesRequired = numberDensity*tetVolume;
147 label nParticlesToInsert = label(particlesRequired);
153 (particlesRequired - nParticlesToInsert)
154 > rndGen_.sample01<scalar>()
157 nParticlesToInsert++;
160 for (label pI = 0; pI < nParticlesToInsert; pI++)
164 vector U = equipartitionLinearVelocity
170 scalar Ei = equipartitionInternalEnergy
173 cP.internalDegreesOfFreedom()
178 addNewParcel(
p, celli,
U, Ei, typeId);
188 label mostAbundantType(
findMax(numberDensities));
190 const typename ParcelType::constantProperties& cP = constProps
195 sigmaTcRMax_.primitiveFieldRef() = cP.sigmaT()*maxwellianMostProbableSpeed
201 sigmaTcRMax_.correctBoundaryConditions();
205template<
class ParcelType>
206void Foam::DSMCCloud<ParcelType>::collisions()
208 if (!binaryCollision().active())
214 List<DynamicList<label>> subCells(8);
216 scalar deltaT =
mesh().time().deltaTValue();
218 label collisionCandidates = 0;
220 label collisions = 0;
222 forAll(cellOccupancy_, celli)
226 label nC(cellParcels.size());
240 List<label> whichSubCell(cellParcels.size());
242 const point& cC = mesh_.cellCentres()[celli];
246 const ParcelType&
p = *cellParcels[i];
247 vector relPos =
p.position() - cC;
250 pos0(relPos.x()) + 2*
pos0(relPos.y()) + 4*
pos0(relPos.z());
252 subCells[subCell].append(i);
253 whichSubCell[i] = subCell;
258 scalar sigmaTcRMax = sigmaTcRMax_[celli];
260 scalar selectedPairs =
261 collisionSelectionRemainder_[celli]
262 + 0.5*nC*(nC - 1)*nParticle_*sigmaTcRMax*deltaT
263 /mesh_.cellVolumes()[celli];
265 label nCandidates(selectedPairs);
266 collisionSelectionRemainder_[celli] = selectedPairs - nCandidates;
267 collisionCandidates += nCandidates;
269 for (label c = 0;
c < nCandidates;
c++)
275 label candidateP = rndGen_.position<label>(0, nC - 1);
278 label candidateQ = -1;
280 List<label> subCellPs = subCells[whichSubCell[candidateP]];
281 label nSC = subCellPs.
size();
291 label i = rndGen_.position<label>(0, nSC - 1);
292 candidateQ = subCellPs[i];
293 }
while (candidateP == candidateQ);
303 candidateQ = rndGen_.position<label>(0, nC - 1);
304 }
while (candidateP == candidateQ);
324 ParcelType& parcelP = *cellParcels[candidateP];
325 ParcelType& parcelQ = *cellParcels[candidateQ];
327 scalar sigmaTcR = binaryCollision().sigmaTcR
337 if (sigmaTcR > sigmaTcRMax_[celli])
339 sigmaTcRMax_[celli] = sigmaTcR;
342 if ((sigmaTcR/sigmaTcRMax) > rndGen_.sample01<scalar>())
344 binaryCollision().collide
360 sigmaTcRMax_.correctBoundaryConditions();
362 if (collisionCandidates)
364 Info<<
" Collisions = "
366 <<
" Acceptance rate = "
367 << scalar(collisions)/scalar(collisionCandidates) <<
nl
377template<
class ParcelType>
378void Foam::DSMCCloud<ParcelType>::resetFields()
393template<
class ParcelType>
394void Foam::DSMCCloud<ParcelType>::calculateFields()
398 scalarField& dsmcRhoN = dsmcRhoN_.primitiveFieldRef();
399 scalarField& linearKE = linearKE_.primitiveFieldRef();
400 scalarField& internalE = internalE_.primitiveFieldRef();
402 vectorField& momentum = momentum_.primitiveFieldRef();
404 for (
const ParcelType&
p : *
this)
406 const label celli =
p.cell();
409 rhoM[celli] += constProps(
p.typeId()).mass();
411 linearKE[celli] += 0.5*constProps(
p.typeId()).mass()*(
p.U() &
p.U());
412 internalE[celli] +=
p.Ei();
413 iDof[celli] += constProps(
p.typeId()).internalDegreesOfFreedom();
414 momentum[celli] += constProps(
p.typeId()).mass()*
p.U();
417 rhoN *= nParticle_/
mesh().cellVolumes();
418 rhoN_.correctBoundaryConditions();
420 rhoM *= nParticle_/
mesh().cellVolumes();
421 rhoM_.correctBoundaryConditions();
423 dsmcRhoN_.correctBoundaryConditions();
425 linearKE *= nParticle_/
mesh().cellVolumes();
426 linearKE_.correctBoundaryConditions();
428 internalE *= nParticle_/
mesh().cellVolumes();
429 internalE_.correctBoundaryConditions();
431 iDof *= nParticle_/
mesh().cellVolumes();
432 iDof_.correctBoundaryConditions();
434 momentum *= nParticle_/
mesh().cellVolumes();
435 momentum_.correctBoundaryConditions();
441template<
class ParcelType>
451 this->
addParticle(
new ParcelType(mesh_, position, celli,
U, Ei, typeId));
457template<
class ParcelType>
458Foam::DSMCCloud<ParcelType>::DSMCCloud
481 typeIdList_(particleProperties_.
lookup(
"typeIdList")),
482 nParticle_(particleProperties_.get<scalar>(
"nEquivalentParticles")),
483 cellOccupancy_(mesh_.nCells()),
488 this->
name() +
"SigmaTcRMax",
497 collisionSelectionRemainder_
501 IOobject::scopedName(this->
name(),
"collisionSelectionRemainder"),
653 binaryCollisionModel_
661 wallInteractionModel_
679 buildCellOccupancy();
683 forAll(collisionSelectionRemainder_, i)
685 collisionSelectionRemainder_[i] = rndGen_.
sample01<scalar>();
690 ParcelType::readFields(*
this);
695template<
class ParcelType>
696Foam::DSMCCloud<ParcelType>::DSMCCloud
719 typeIdList_(particleProperties_.
lookup(
"typeIdList")),
720 nParticle_(particleProperties_.get<scalar>(
"nEquivalentParticles")),
726 this->
name() +
"SigmaTcRMax",
737 collisionSelectionRemainder_
741 IOobject::scopedName(this->
name(),
"collisionSelectionRemainder"),
765 this->
name() +
"fD_",
778 this->
name() +
"rhoN_",
791 this->
name() +
"rhoM_",
804 this->
name() +
"dsmcRhoN_",
817 this->
name() +
"linearKE_",
830 this->
name() +
"internalE_",
843 this->
name() +
"iDof_",
856 this->
name() +
"momentum_",
899 binaryCollisionModel_(),
900 wallInteractionModel_(),
901 inflowBoundaryModel_()
905 initialise(dsmcInitialiseDict);
911template<
class ParcelType>
918template<
class ParcelType>
921 typename ParcelType::trackingData
td(*
this);
928 this->dumpParticlePositions();
932 this->inflowBoundary().inflow();
938 buildCellOccupancy();
948template<
class ParcelType>
951 label nDSMCParticles = this->size();
954 scalar nMol = nDSMCParticles*nParticle_;
956 vector linearMomentum = linearMomentumOfSystem();
959 scalar linearKineticEnergy = linearKineticEnergyOfSystem();
962 scalar internalEnergy = internalEnergyOfSystem();
966 <<
" Number of dsmc particles = "
972 Info<<
" Number of molecules = "
974 <<
" Mass in system = "
976 <<
" Average linear momentum = "
977 << linearMomentum/nMol <<
nl
978 <<
" |Average linear momentum| = "
979 <<
mag(linearMomentum)/nMol <<
nl
980 <<
" Average linear kinetic energy = "
981 << linearKineticEnergy/nMol <<
nl
982 <<
" Average internal energy = "
983 << internalEnergy/nMol <<
nl
984 <<
" Average total energy = "
985 << (internalEnergy + linearKineticEnergy)/nMol
991template<
class ParcelType>
1000 *rndGen_.GaussNormal<
vector>();
1004template<
class ParcelType>
1020 -
log(rndGen_.sample01<scalar>())
1026 const scalar a = 0.5*iDof - 1;
1027 scalar energyRatio = 0;
1032 energyRatio = 10*rndGen_.sample01<scalar>();
1033 P =
pow((energyRatio/a), a)*
exp(a - energyRatio);
1034 }
while (P < rndGen_.sample01<scalar>());
1040template<
class ParcelType>
1045 this->db().time().
path()/
"parcelPositions_"
1046 + this->
name() +
"_"
1047 + this->db().time().
timeName() +
".obj"
1050 for (
const ParcelType&
p : *
this)
1052 pObj<<
"v " <<
p.position().x()
1053 <<
' ' <<
p.position().y()
1054 <<
' ' <<
p.position().z()
1062template<
class ParcelType>
1068 cellOccupancy_.setSize(mesh_.nCells());
1069 buildCellOccupancy();
1072 this->inflowBoundary().autoMap(mapper);
const word cloudName(propsDict.get< word >("cloud"))
Templated DSMC particle collision class.
Base cloud calls templated on particle type.
Cloud(const polyMesh &mesh, const Foam::zero, const word &cloudName)
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
void addParticle(ParcelType *pPtr)
label size() const noexcept
The number of elements in list.
Virtual abstract base class for templated DSMCCloud.
DSMCBaseCloud()=default
Null constructor.
Templated base class for dsmc cloud.
const word & cloudName() const
Return the cloud type.
scalar massInSystem() const
Total mass in system.
virtual ~DSMCCloud()
Destructor.
vector linearMomentumOfSystem() const
Total linear momentum of the system.
const InflowBoundaryModel< DSMCCloud< ParcelType > > & inflowBoundary() const
Return reference to wall interaction model.
vector equipartitionLinearVelocity(scalar temperature, scalar mass)
Generate a random velocity sampled from the Maxwellian speed.
const volScalarField & iDof() const
Return the average internal degrees of freedom field.
virtual void autoMap(const mapPolyMesh &)
Remap the particles to the correct cells following mesh change.
scalar internalEnergyOfSystem() const
Total internal energy in the system.
void evolve()
Evolve the cloud (move, collide).
scalar equipartitionInternalEnergy(scalar temperature, direction internalDegreesOfFreedom)
Generate a random internal energy, sampled from the.
void clear()
Clear the Cloud.
const fvMesh & mesh() const
Return reference to the mesh.
void addNewParcel(const vector &position, const label celli, const vector &U, const scalar Ei, const label typeId)
Add new parcel.
void dumpParticlePositions() const
Dump particle positions to .obj file.
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
void info() const
Print cloud information.
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...
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const objectRegistry & db() const noexcept
Return the local objectRegistry.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Templated inflow boundary model class.
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
virtual void flush() override
Flush stream.
Inter-processor communications stream.
Type sample01()
Return a sample whose components lie in the range [0,1].
void size(const label n)
Older name for setAddressableSize.
Templated wall interaction model class.
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.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
Template invariant parts for fvPatchField.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const Time & time() const noexcept
Return time registry.
constant condensation/saturation model.
static List< tetIndices > cellTetIndices(const polyMesh &mesh, label cI)
Return the tet decomposition of the given cell, see.
Lookup type of boundary radiation properties.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
tetPointRef tet(const polyMesh &mesh) const
The tet geometry for this tet, where point0 is the cell centre.
scalar mag() const
Return volume.
Point randomPoint(Random &rndGen) const
Return a random point in the tetrahedron from a uniform distribution.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
const dimensionedScalar k
Boltzmann constant.
const dimensionedScalar c
Speed of light in a vacuum.
Different types of constants.
Namespace for handling debugging switches.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
dimensionedScalar pos0(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label findMax(const ListType &input, label start=0)
Linear search for the index of the max element, similar to std::max_element but for lists and returns...
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.