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>
128 inline Foam::volScalarField::Boundary&
131 return q_.boundaryFieldRef();
135 template<
class ParcelType>
136 inline Foam::volVectorField::Boundary&
139 return fD_.boundaryFieldRef();
143 template<
class ParcelType>
144 inline Foam::volScalarField::Boundary&
147 return rhoN_.boundaryFieldRef();
151 template<
class ParcelType>
152 inline Foam::volScalarField::Boundary&
155 return rhoM_.boundaryFieldRef();
159 template<
class ParcelType>
160 inline Foam::volScalarField::Boundary&
163 return linearKE_.boundaryFieldRef();
167 template<
class ParcelType>
168 inline Foam::volScalarField::Boundary&
171 return internalE_.boundaryFieldRef();
175 template<
class ParcelType>
176 inline Foam::volScalarField::Boundary&
179 return iDof_.boundaryFieldRef();
183 template<
class ParcelType>
184 inline Foam::volVectorField::Boundary&
187 return momentum_.boundaryFieldRef();
191 template<
class ParcelType>
199 template<
class ParcelType>
207 template<
class ParcelType>
211 return binaryCollisionModel_;
215 template<
class ParcelType>
219 return binaryCollisionModel_();
223 template<
class ParcelType>
227 return wallInteractionModel_;
231 template<
class ParcelType>
235 return wallInteractionModel_();
239 template<
class ParcelType>
243 return inflowBoundaryModel_;
247 template<
class ParcelType>
251 return inflowBoundaryModel_();
255 template<
class ParcelType>
258 scalar sysMass = 0.0;
262 const ParcelType&
p = iter();
264 const typename ParcelType::constantProperties& cP = constProps
269 sysMass += cP.mass();
272 return nParticle_*sysMass;
276 template<
class ParcelType>
283 const ParcelType&
p = iter();
285 const typename ParcelType::constantProperties& cP = constProps
290 linearMomentum += cP.mass()*p.U();
293 return nParticle_*linearMomentum;
297 template<
class ParcelType>
301 scalar linearKineticEnergy = 0.0;
305 const ParcelType&
p = iter();
307 const typename ParcelType::constantProperties& cP = constProps
312 linearKineticEnergy += 0.5*cP.mass()*(p.U() & p.U());
315 return nParticle_*linearKineticEnergy;
319 template<
class ParcelType>
323 scalar internalEnergy = 0.0;
327 const ParcelType&
p = iter();
329 internalEnergy += p.Ei();
332 return nParticle_*internalEnergy;
336 template<
class ParcelType>
348 template<
class ParcelType>
361 template<
class ParcelType>
372 template<
class ParcelType>
385 template<
class ParcelType>
397 template<
class ParcelType>
411 template<
class ParcelType>
418 template<
class ParcelType>
425 template<
class ParcelType>
433 template<
class ParcelType>
440 template<
class ParcelType>
448 template<
class ParcelType>
456 template<
class ParcelType>
464 template<
class ParcelType>
472 template<
class ParcelType>
479 template<
class ParcelType>
const volVectorField & fD() const
Return force density at surface field.
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 volVectorField & momentum() const
Return the momentum density field.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Templated inflow boundary model class.
volScalarField::Boundary & rhoMBF()
Return non-const mass density boundary field reference.
Templated DSMC particle collision class.
const List< word > & typeIdList() const
Return the idList.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const IOdictionary & particleProperties() const
Return particle properties dictionary.
volScalarField::Boundary & iDofBF()
Return non-const internal degree of freedom density boundary.
const volScalarField & boundaryT() const
Return macroscopic temperature.
const InflowBoundaryModel< DSMCCloud< ParcelType > > & inflowBoundary() const
Return reference to wall interaction model.
volScalarField::Boundary & qBF()
Return non-const heat flux boundary field reference.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
scalar maxwellianMostProbableSpeed(scalar temperature, scalar mass) const
Most probable speed.
const volScalarField & dsmcRhoN() const
Return the field of number of DSMC particles.
scalar massInSystem() const
Total mass in system.
const volScalarField & rhoN() const
Return the real particle number density field.
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
const volScalarField & internalE() const
Return the internal energy density field.
A class for handling words, derived from string.
const volScalarField & q() const
Return heat flux at surface field.
Random & rndGen()
Return refernce to the random object.
void clear()
Clear the contents of the list.
volScalarField::Boundary & linearKEBF()
Return non-const linear kinetic energy density boundary.
const volScalarField & rhoM() const
Return the particle mass density field.
volVectorField::Boundary & fDBF()
Return non-const force density at boundary field reference.
scalar nParticle() const
Return the number of real particles represented by one.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
errorManip< error > abort(error &err)
volScalarField::Boundary & internalEBF()
Return non-const internal energy density boundary field.
Simple random number generator.
scalar internalEnergyOfSystem() const
Total internal energy in the system.
const volScalarField & iDof() const
Return the average internal degrees of freedom field.
void clear()
Clear the Cloud.
scalar maxwellianRMSSpeed(scalar temperature, scalar mass) const
RMS particle speed.
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
const volScalarField & linearKE() const
Return the total linear kinetic energy (translational and.
const dimensionedScalar k
Boltzmann constant.
const fvMesh & mesh() const
Return refernce to the mesh.
const List< DynamicList< ParcelType * > > & cellOccupancy() const
Return the cell occupancy addressing.
Mesh data needed to do the Finite Volume discretisation.
const volVectorField & boundaryU() const
Return macroscopic velocity.
const word & cloudName() const
Return the cloud type.
volScalarField::Boundary & rhoNBF()
Return non-const number density boundary field reference.
volScalarField & sigmaTcRMax()
Return the sigmaTcRMax field. non-const access to allow.
scalar maxwellianAverageSpeed(scalar temperature, scalar mass) const
Average particle speed.
A class for managing temporary objects.
Templated wall interaction model class.
vector linearMomentumOfSystem() const
Total linear momentum of the system.
const WallInteractionModel< DSMCCloud< ParcelType > > & wallInteraction() const
Return reference to wall interaction model.
scalarField & collisionSelectionRemainder()
Return the collision selection remainder field. non-const.
Templated base class for dsmc cloud.
const BinaryCollisionModel< DSMCCloud< ParcelType > > & binaryCollision() const
Return reference to binary elastic collision model.
volVectorField::Boundary & momentumBF()
Return non-const momentum density boundary field reference.