33 template<
class ParcelType>
40 template<
class ParcelType>
47 template<
class ParcelType>
51 return particleProperties_;
55 template<
class ParcelType>
63 template<
class ParcelType>
70 template<
class ParcelType>
74 return cellOccupancy_;
78 template<
class ParcelType>
85 template<
class ParcelType>
89 return collisionSelectionRemainder_;
93 template<
class ParcelType>
101 template<
class ParcelType>
102 inline const typename ParcelType::constantProperties&
108 if (typeId < 0 || typeId >= constProps_.size())
111 <<
"constantProperties for requested typeId index "
112 << typeId <<
" do not exist" <<
nl
116 return constProps_[typeId];
120 template<
class ParcelType>
127 template<
class ParcelType>
135 template<
class ParcelType>
139 return q_.boundaryFieldRef();
143 template<
class ParcelType>
147 return fD_.boundaryFieldRef();
151 template<
class ParcelType>
155 return rhoN_.boundaryFieldRef();
159 template<
class ParcelType>
163 return rhoM_.boundaryFieldRef();
167 template<
class ParcelType>
171 return linearKE_.boundaryFieldRef();
175 template<
class ParcelType>
179 return internalE_.boundaryFieldRef();
183 template<
class ParcelType>
187 return iDof_.boundaryFieldRef();
191 template<
class ParcelType>
195 return momentum_.boundaryFieldRef();
199 template<
class ParcelType>
207 template<
class ParcelType>
215 template<
class ParcelType>
219 return binaryCollisionModel_;
223 template<
class ParcelType>
227 return binaryCollisionModel_();
231 template<
class ParcelType>
235 return wallInteractionModel_;
239 template<
class ParcelType>
243 return wallInteractionModel_();
247 template<
class ParcelType>
251 return inflowBoundaryModel_;
255 template<
class ParcelType>
259 return inflowBoundaryModel_();
263 template<
class ParcelType>
266 scalar sysMass = 0.0;
270 const ParcelType&
p = iter();
272 const typename ParcelType::constantProperties& cP = constProps
277 sysMass += cP.mass();
280 return nParticle_*sysMass;
284 template<
class ParcelType>
291 const ParcelType&
p = iter();
293 const typename ParcelType::constantProperties& cP = constProps
298 linearMomentum += cP.mass()*
p.U();
301 return nParticle_*linearMomentum;
305 template<
class ParcelType>
309 scalar linearKineticEnergy = 0.0;
313 const ParcelType&
p = iter();
315 const typename ParcelType::constantProperties& cP = constProps
320 linearKineticEnergy += 0.5*cP.mass()*(
p.U() &
p.U());
323 return nParticle_*linearKineticEnergy;
327 template<
class ParcelType>
331 scalar internalEnergy = 0.0;
335 const ParcelType&
p = iter();
337 internalEnergy +=
p.Ei();
340 return nParticle_*internalEnergy;
344 template<
class ParcelType>
356 template<
class ParcelType>
369 template<
class ParcelType>
380 template<
class ParcelType>
393 template<
class ParcelType>
405 template<
class ParcelType>
419 template<
class ParcelType>
426 template<
class ParcelType>
433 template<
class ParcelType>
441 template<
class ParcelType>
448 template<
class ParcelType>
456 template<
class ParcelType>
464 template<
class ParcelType>
472 template<
class ParcelType>
480 template<
class ParcelType>
487 template<
class ParcelType>
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Templated DSMC particle collision class.
Templated base class for dsmc cloud.
volScalarField::Boundary & rhoMBF()
Return non-const mass density boundary field reference.
const volScalarField & boundaryT() const
Return macroscopic temperature.
const word & cloudName() const
Return the cloud type.
scalar massInSystem() const
Total mass in system.
volScalarField::Boundary & rhoNBF()
Return non-const number density boundary field reference.
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
const volScalarField & q() const
Return heat flux at surface field.
const List< DynamicList< ParcelType * > > & cellOccupancy() const
Return the cell occupancy addressing.
volScalarField::Boundary & internalEBF()
Return non-const internal energy density boundary field.
volScalarField & sigmaTcRMax()
Return the sigmaTcRMax field. non-const access to allow.
vector linearMomentumOfSystem() const
Total linear momentum of the system.
const InflowBoundaryModel< DSMCCloud< ParcelType > > & inflowBoundary() const
Return reference to wall interaction model.
const volScalarField & iDof() const
Return the average internal degrees of freedom field.
const BinaryCollisionModel< DSMCCloud< ParcelType > > & binaryCollision() const
Return reference to binary elastic collision model.
scalar nParticle() const
Return the number of real particles represented by one.
const volVectorField & boundaryU() const
Return macroscopic velocity.
scalar maxwellianRMSSpeed(scalar temperature, scalar mass) const
RMS particle speed.
distributions::standardNormal & stdNormal()
Return reference to the standard normal distribution.
scalar maxwellianMostProbableSpeed(scalar temperature, scalar mass) const
Most probable speed.
volScalarField::Boundary & linearKEBF()
Return non-const linear kinetic energy density boundary.
volVectorField::Boundary & momentumBF()
Return non-const momentum density boundary field reference.
volScalarField::Boundary & qBF()
Return non-const heat flux boundary field reference.
volVectorField::Boundary & fDBF()
Return non-const force density at boundary field reference.
scalar internalEnergyOfSystem() const
Total internal energy in the system.
randomGenerator & rndGen()
Return reference to the random generator.
scalarField & collisionSelectionRemainder()
Return the collision selection remainder field. non-const.
const volVectorField & momentum() const
Return the momentum density field.
const List< word > & typeIdList() const
Return the idList.
volScalarField::Boundary & iDofBF()
Return non-const internal degree of freedom density boundary.
const IOdictionary & particleProperties() const
Return particle properties dictionary.
const volScalarField & linearKE() const
Return the total linear kinetic energy (translational and.
void clear()
Clear the Cloud.
const volScalarField & dsmcRhoN() const
Return the field of number of DSMC particles.
const fvMesh & mesh() const
Return references to the mesh.
const volScalarField & internalE() const
Return the internal energy density field.
scalar maxwellianAverageSpeed(scalar temperature, scalar mass) const
Average particle speed.
const volScalarField & rhoM() const
Return the particle mass density field.
const volScalarField & rhoN() const
Return the real particle number density field.
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
const volVectorField & fD() const
Return force density at surface field.
const WallInteractionModel< DSMCCloud< ParcelType > > & wallInteraction() const
Return reference to wall interaction model.
Generic GeometricBoundaryField class.
Generic GeometricField class.
void clear()
Clear the contents of the list.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Templated inflow boundary model class.
Templated wall interaction model class.
Standard normal distribution. Not selectable.
Mesh data needed to do the Finite Volume discretisation.
A class for managing temporary objects.
A class for handling words, derived from string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const dimensionedScalar k
Boltzmann constant.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
errorManip< error > abort(error &err)
dimensionedScalar sqrt(const dimensionedScalar &ds)