Go to the documentation of this file.
33 template<
class CloudType>
37 return *cloudCopyPtr_;
41 template<
class CloudType>
48 template<
class CloudType>
52 return particleProperties_;
56 template<
class CloudType>
60 return outputProperties_;
64 template<
class CloudType>
67 return outputProperties_;
71 template<
class CloudType>
79 template<
class CloudType>
86 template<
class CloudType>
87 inline const typename CloudType::particleType::constantProperties&
94 template<
class CloudType>
95 inline typename CloudType::particleType::constantProperties&
102 template<
class CloudType>
106 return subModelProperties_;
110 template<
class CloudType>
117 template<
class CloudType>
124 template<
class CloudType>
131 template<
class CloudType>
138 template<
class CloudType>
145 template<
class CloudType>
152 template<
class CloudType>
161 template<
class CloudType>
169 template<
class CloudType>
177 template<
class CloudType>
185 template<
class CloudType>
193 template<
class CloudType>
197 return *dispersionModel_;
201 template<
class CloudType>
205 return *dispersionModel_;
209 template<
class CloudType>
213 return *patchInteractionModel_;
217 template<
class CloudType>
221 return *patchInteractionModel_;
225 template<
class CloudType>
229 return *stochasticCollisionModel_;
233 template<
class CloudType>
237 return *stochasticCollisionModel_;
241 template<
class CloudType>
245 return *surfaceFilmModel_;
249 template<
class CloudType>
253 return *surfaceFilmModel_;
257 template<
class CloudType>
261 return *packingModel_;
265 template<
class CloudType>
269 return *packingModel_;
273 template<
class CloudType>
277 return *dampingModel_;
281 template<
class CloudType>
285 return *dampingModel_;
289 template<
class CloudType>
293 return *isotropyModel_;
297 template<
class CloudType>
301 return *isotropyModel_;
305 template<
class CloudType>
309 return *UIntegrator_;
313 template<
class CloudType>
316 scalar sysMass = 0.0;
317 for (
const parcelType&
p : *
this)
319 sysMass +=
p.nParticle()*
p.mass();
326 template<
class CloudType>
332 for (
const parcelType&
p : *
this)
334 linearMomentum +=
p.nParticle()*
p.mass()*
p.U();
337 return linearMomentum;
341 template<
class CloudType>
345 scalar parPerParcel = 0;
347 for (
const parcelType&
p : *
this)
349 parPerParcel +=
p.nParticle();
356 template<
class CloudType>
360 scalar linearKineticEnergy = 0;
362 for (
const parcelType&
p : *
this)
364 linearKineticEnergy +=
p.nParticle()*0.5*
p.mass()*(
p.U() &
p.U());
367 return linearKineticEnergy;
371 template<
class CloudType>
380 for (
const parcelType&
p : *
this)
382 si +=
p.nParticle()*
pow(
p.d(), i);
383 sj +=
p.nParticle()*
pow(
p.d(), j);
386 reduce(si, sumOp<scalar>());
387 reduce(sj, sumOp<scalar>());
394 template<
class CloudType>
398 for (
const parcelType&
p : *
this)
409 template<
class CloudType>
416 template<
class CloudType>
420 if (!cellOccupancyPtr_)
422 buildCellOccupancy();
425 return *cellOccupancyPtr_;
429 template<
class CloudType>
433 return cellLengthScale_;
437 template<
class CloudType>
445 template<
class CloudType>
453 template<
class CloudType>
461 template<
class CloudType>
469 template<
class CloudType>
476 Pout<<
"UTrans min/max = " <<
min(UTrans()).value() <<
", "
477 <<
max(UTrans()).value() <<
nl
478 <<
"UCoeff min/max = " <<
min(UCoeff()).value() <<
", "
479 <<
max(UCoeff()).value() <<
endl;
488 if (solution_.coupled())
490 if (solution_.semiImplicit(
"U"))
492 volScalarField::Internal
493 Vdt(mesh_.V()*this->db().time().deltaT());
500 return UTrans()/Vdt -
fvm::Sp(UCoeff()/Vdt,
U) + UCoeff()/Vdt*
U;
507 fvm.source() = -UTrans()/(this->db().time().deltaT());
517 template<
class CloudType>
527 this->
name() +
":vDotSweep",
536 extrapolatedCalculatedFvPatchScalarField::typeName
541 for (
const parcelType&
p : *
this)
543 const label celli =
p.cell();
545 vDotSweep[celli] +=
p.nParticle()*
p.areaP()*
mag(
p.U() - U_[celli]);
555 template<
class CloudType>
565 this->
name() +
":theta",
574 extrapolatedCalculatedFvPatchScalarField::typeName
579 for (
const parcelType&
p : *
this)
581 const label celli =
p.cell();
583 theta[celli] +=
p.nParticle()*
p.volume();
593 template<
class CloudType>
603 this->
name() +
":alpha",
616 for (
const parcelType&
p : *
this)
618 const label celli =
p.cell();
620 alpha[celli] +=
p.nParticle()*
p.mass();
623 alpha /= (mesh_.V()*rho_);
629 template<
class CloudType>
639 this->
name() +
":rhoEff",
652 for (
const parcelType&
p : *
this)
654 const label celli =
p.cell();
656 rhoEff[celli] +=
p.nParticle()*
p.mass();
volScalarField::Internal & UCoeff()
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
tmp< fvVectorMatrix > SU(volVectorField &U, bool incompressible=false) const
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const integrationScheme & UIntegrator() const
functionType & functions()
A class for managing temporary objects.
Base class for dispersion modelling.
static constexpr const zero Zero
const dimensionSet dimDensity
const parcelType::constantProperties & constProps() const
scalar totalParticlePerParcel() const
const fvMesh & mesh() const
const dimensionedScalar alpha
const volScalarField & mu() const
const DampingModel< KinematicCloud< CloudType > > & dampingModel() const
scalar massInSystem() const
const tmp< volScalarField > rhoEff() const
const forceType & forces() const
const SurfaceFilmModel< KinematicCloud< CloudType > > & surfaceFilm() const
const KinematicCloud & cloudCopy() const
const InjectionModelList< KinematicCloud< CloudType > > & injectors() const
Ostream & endl(Ostream &os)
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Base class for collisional damping models.
const dimensionSet dimForce
Base for a set of schemes which integrate simple ODEs which arise from semi-implcit rate expressions.
label min(const labelHashSet &set, label minValue=labelMax)
const tmp< volScalarField > vDotSweep() const
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Templated base class for kinematic cloud.
Templated patch interaction model class.
fvMatrix< vector > fvVectorMatrix
Templated stochastic collision model class.
Generic templated field type.
List< DynamicList< parcelType * > > & cellOccupancy()
Templated wall surface film model class.
volVectorField::Internal & UTrans()
Stores all relevant solution info for cloud.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
Base class for packing models.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
reduce(hasMovingMesh, orOp< bool >())
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
const StochasticCollisionModel< KinematicCloud< CloudType > > & stochasticCollision() const
const dimensionSet & dimensions() const
const IOdictionary & particleProperties() const
const cloudSolution & solution() const
Calculate the matrix for implicit and explicit sources.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
void correctBoundaryConditions()
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
const IsotropyModel< KinematicCloud< CloudType > > & isotropyModel() const
const DispersionModel< KinematicCloud< CloudType > > & dispersion() const
Internal & ref(const bool updateAccessTime=true)
const PatchInteractionModel< KinematicCloud< CloudType > > & patchInteraction() const
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
const dimensionedVector & g() const
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const tmp< volScalarField > theta() const
A special matrix type and solver, designed for finite volume solutions of scalar equations....
const PackingModel< KinematicCloud< CloudType > > & packingModel() const
scalar Dij(const label i, const label j) const
word name(const expressions::valueTypeCode typeCode)
const dictionary & subModelProperties() const
Base class for collisional return-to-isotropy models.
const volVectorField & U() const
scalar linearKineticEnergyOfSystem() const
Generic GeometricField class.
vector linearMomentumOfSystem() const
const IOdictionary & outputProperties() const
void reset(const dimensionSet &ds)
List of injection models.
const tmp< volScalarField > alpha() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimless
const volScalarField & rho() const
const scalarField & cellLengthScale() const