55 template<
class CloudType>
58 template<
class CloudType>
61 template<
class CloudType>
68 template<
class ParcelType>
71 public Cloud<ParcelType>,
78 const word cloudName_;
153 binaryCollisionModel_;
157 wallInteractionModel_;
161 inflowBoundaryModel_;
167 void buildConstProps();
170 void buildCellOccupancy();
173 void initialise(
const IOdictionary& dsmcInitialiseDict);
182 void calculateFields();
207 const word& cloudName,
263 inline const typename ParcelType::constantProperties&
273 inline volScalarField::GeometricBoundaryField&
qBF();
276 inline volVectorField::GeometricBoundaryField&
fDBF();
279 inline volScalarField::GeometricBoundaryField&
rhoNBF();
282 inline volScalarField::GeometricBoundaryField&
rhoMBF();
286 inline volScalarField::GeometricBoundaryField&
linearKEBF();
290 inline volScalarField::GeometricBoundaryField&
internalEBF();
294 inline volScalarField::GeometricBoundaryField&
iDofBF();
297 inline volVectorField::GeometricBoundaryField&
momentumBF();
457 const label tetFaceI,
volScalarField::GeometricBoundaryField & rhoNBF()
Return non-const number density boundary field reference.
const IOdictionary & particleProperties() const
Return particle properties dictionary.
Simple random number generator.
const WallInteractionModel< DSMCCloud< ParcelType > > & wallInteraction() const
Return reference to wall interaction model.
virtual ~DSMCCloud()
Destructor.
Mesh data needed to do the Finite Volume discretisation.
scalar maxwellianMostProbableSpeed(scalar temperature, scalar mass) const
Most probable speed.
volScalarField::GeometricBoundaryField & qBF()
Return non-const heat flux boundary field reference.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
volScalarField::GeometricBoundaryField & rhoMBF()
Return non-const mass density boundary field reference.
void clear()
Clear the Cloud.
const volVectorField & momentum() const
Return the momentum density field.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
volScalarField::GeometricBoundaryField & linearKEBF()
Return non-const linear kinetic energy density boundary.
scalarField & collisionSelectionRemainder()
Return the collision selection remainder field. non-const.
Virtual abstract base class for templated DSMCCloud.
const volScalarField & rhoN() const
Return the real particle number density field.
A class for handling words, derived from string.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const volScalarField & boundaryT() const
Return macroscopic temperature.
const InflowBoundaryModel< DSMCCloud< ParcelType > > & inflowBoundary() const
Return reference to wall interaction model.
volScalarField::GeometricBoundaryField & iDofBF()
Return non-const internal degree of freedom density boundary.
const volScalarField & linearKE() const
Return the total linear kinetic energy (translational and.
const volScalarField & iDof() const
Return the average internal degrees of freedom field.
const volScalarField & internalE() const
Return the internal energy density field.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
void addNewParcel(const vector &position, const vector &U, const scalar Ei, const label cellI, const label tetFaceI, const label tetPtI, const label typeId)
Add new parcel.
Templated wall interaction model class.
scalar equipartitionInternalEnergy(scalar temperature, direction internalDegreesOfFreedom)
Generate a random internal energy, sampled from the.
void dumpParticlePositions() const
Dump particle positions to .obj file.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
void evolve()
Evolve the cloud (move, collide)
vector equipartitionLinearVelocity(scalar temperature, scalar mass)
Generate a random velocity sampled from the Maxwellian speed.
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
const List< word > & typeIdList() const
Return the idList.
virtual void autoMap(const mapPolyMesh &)
Remap the particles to the correct cells following mesh change.
const volVectorField & boundaryU() const
Return macroscopic velocity.
vector linearMomentumOfSystem() const
Total linear momentum of the system.
scalar massInSystem() const
Total mass in system.
const volScalarField & dsmcRhoN() const
Return the field of number of DSMC particles.
const word & cloudName() const
Return the cloud type.
scalar nParticle() const
Return the number of real particles represented by one.
void info() const
Print cloud information.
volScalarField & sigmaTcRMax()
Return the sigmaTcRMax field. non-const access to allow.
volVectorField::GeometricBoundaryField & momentumBF()
Return non-const momentum density boundary field reference.
const volScalarField & rhoM() const
Return the particle mass density field.
Random & rndGen()
Return refernce to the random object.
Templated DSMC particle collision class.
volVectorField::GeometricBoundaryField & fDBF()
Return non-const force density at boundary field reference.
scalar maxwellianAverageSpeed(scalar temperature, scalar mass) const
Average particle speed.
const volScalarField & q() const
Return heat flux at surface field.
Base cloud calls templated on particle type.
const fvMesh & mesh() const
Return refernce to the mesh.
Templated base class for dsmc cloud.
Templated inflow boundary model class.
const BinaryCollisionModel< DSMCCloud< ParcelType > > & binaryCollision() const
Return reference to binary elastic collision model.
This function object reads fields from the time directories and adds them to the mesh database for fu...
const volVectorField & fD() const
Return force density at surface field.
scalar internalEnergyOfSystem() const
Total internal energy in the system.
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
const List< DynamicList< ParcelType * > > & cellOccupancy() const
Return the cell occupancy addressing.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
volScalarField::GeometricBoundaryField & internalEBF()
Return non-const internal energy density boundary field.
scalar maxwellianRMSSpeed(scalar temperature, scalar mass) const
RMS particle speed.