Go to the documentation of this file.
34 template<
class CloudType>
38 return *cloudCopyPtr_;
42 template<
class CloudType>
49 template<
class CloudType>
53 return particleProperties_;
57 template<
class CloudType>
61 return outputProperties_;
65 template<
class CloudType>
68 return outputProperties_;
72 template<
class CloudType>
80 template<
class CloudType>
87 template<
class CloudType>
88 inline const typename CloudType::particleType::constantProperties&
95 template<
class CloudType>
96 inline typename CloudType::particleType::constantProperties&
103 template<
class CloudType>
107 return subModelProperties_;
111 template<
class CloudType>
118 template<
class CloudType>
125 template<
class CloudType>
132 template<
class CloudType>
139 template<
class CloudType>
146 template<
class CloudType>
153 template<
class CloudType>
162 template<
class CloudType>
170 template<
class CloudType>
178 template<
class CloudType>
186 template<
class CloudType>
194 template<
class CloudType>
198 return *dispersionModel_;
202 template<
class CloudType>
206 return *dispersionModel_;
210 template<
class CloudType>
214 return *patchInteractionModel_;
218 template<
class CloudType>
222 return *patchInteractionModel_;
226 template<
class CloudType>
230 return *stochasticCollisionModel_;
234 template<
class CloudType>
238 return *stochasticCollisionModel_;
242 template<
class CloudType>
246 return *surfaceFilmModel_;
250 template<
class CloudType>
254 return *surfaceFilmModel_;
258 template<
class CloudType>
262 return *UIntegrator_;
266 template<
class CloudType>
269 scalar sysMass = 0.0;
272 sysMass +=
p.nParticle()*
p.mass();
279 template<
class CloudType>
287 linearMomentum +=
p.nParticle()*
p.mass()*
p.U();
290 return linearMomentum;
294 template<
class CloudType>
298 scalar parPerParcel = 0;
302 parPerParcel +=
p.nParticle();
309 template<
class CloudType>
313 scalar linearKineticEnergy = 0;
317 linearKineticEnergy +=
p.nParticle()*0.5*
p.mass()*(
p.U() &
p.U());
320 return linearKineticEnergy;
324 template<
class CloudType>
335 si +=
p.nParticle()*
pow(
p.d(), i);
336 sj +=
p.nParticle()*
pow(
p.d(), j);
341 sj =
max(sj, VSMALL);
347 template<
class CloudType>
362 template<
class CloudType>
369 template<
class CloudType>
373 if (cellOccupancyPtr_.empty())
375 buildCellOccupancy();
378 return *cellOccupancyPtr_;
382 template<
class CloudType>
386 return cellLengthScale_;
390 template<
class CloudType>
398 template<
class CloudType>
406 template<
class CloudType>
414 template<
class CloudType>
422 template<
class CloudType>
428 Pout<<
"UTrans min/max = " <<
min(UTrans()).value() <<
", "
429 <<
max(UTrans()).value() <<
nl
430 <<
"UCoeff min/max = " <<
min(UCoeff()).value() <<
", "
431 <<
max(UCoeff()).value() <<
endl;
434 if (solution_.coupled())
436 if (solution_.semiImplicit(
"U"))
439 Vdt(mesh_.V()*this->db().time().deltaT());
441 return UTrans()/Vdt -
fvm::Sp(UCoeff()/Vdt,
U) + UCoeff()/Vdt*
U;
448 fvm.
source() = -UTrans()/(this->db().time().deltaT());
458 template<
class CloudType>
468 this->
name() +
":vDotSweep",
477 extrapolatedCalculatedFvPatchScalarField::typeName
484 const label celli =
p.cell();
486 vDotSweep[celli] +=
p.nParticle()*
p.areaP()*
mag(
p.U() - U_[celli]);
496 template<
class CloudType>
506 this->
name() +
":theta",
515 extrapolatedCalculatedFvPatchScalarField::typeName
522 const label celli =
p.cell();
524 theta[celli] +=
p.nParticle()*
p.volume();
534 template<
class CloudType>
544 this->
name() +
":alpha",
559 const label celli =
p.cell();
561 alpha[celli] +=
p.nParticle()*
p.mass();
564 alpha /= (mesh_.V()*rho_);
570 template<
class CloudType>
580 this->
name() +
":rhoEff",
595 const label celli =
p.cell();
597 rhoEff[celli] +=
p.nParticle()*
p.mass();
int debug
Static debugging option.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
const DispersionModel< KinematicCloud< CloudType > > & dispersion() const
Return const-access to the dispersion model.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const integrationScheme & UIntegrator() const
Return reference to velocity integration.
const IOdictionary & outputProperties() const
Return output properties dictionary.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
functionType & functions()
Optional cloud function objects.
scalar pAmbient() const
Return const-access to the ambient pressure.
A class for managing temporary objects.
Base class for dispersion modelling.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimDensity
scalar totalParticlePerParcel() const
Average particle per parcel.
const fvMesh & mesh() const
Return reference to the mesh.
tmp< fvVectorMatrix > SU(volVectorField &U) const
Return tmp momentum source term.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
const volScalarField & mu() const
Return carrier gas dynamic viscosity.
scalar massInSystem() const
Total mass in system.
const tmp< volScalarField > rhoEff() const
Return the particle effective density field.
volScalarField::Internal & UCoeff()
Return coefficient for carrier phase U equation.
const KinematicCloud & cloudCopy() const
Return a reference to the cloud copy.
Ostream & endl(Ostream &os)
Add newline and flush stream.
volVectorField::Internal & UTrans()
Return reference to momentum source.
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
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)
Find the min value in labelHashSet, optionally limited by second argument.
const forceType & forces() const
Optional particle forces.
const StochasticCollisionModel< KinematicCloud< CloudType > > & stochasticCollision() const
Return const-access to the stochastic collision model.
const parcelType::constantProperties & constProps() const
Return the constant properties.
const tmp< volScalarField > vDotSweep() const
Volume swept rate of parcels per cell.
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.
List< DynamicList< parcelType * > > & cellOccupancy()
Return the cell occupancy information for each.
Templated wall surface film model class.
word name(const complex &c)
Return string representation of complex.
Stores all relevant solution info for cloud.
Wrapper around kinematic parcel types to add collision modelling.
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)
Find the max value in labelHashSet, optionally limited by second argument.
const cloudSolution & solution() const
Return const access to the solution properties.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
reduce(hasMovingMesh, orOp< bool >())
Mesh data needed to do the Finite Volume discretisation.
const IOdictionary & particleProperties() const
Return particle properties dictionary.
Calculate the matrix for implicit and explicit sources.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
void correctBoundaryConditions()
Correct boundary field.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Random & rndGen() const
Return reference to the random object.
const PatchInteractionModel< KinematicCloud< CloudType > > & patchInteraction() const
Return const-access to the patch interaction model.
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
Gravity.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
const InjectionModelList< KinematicCloud< CloudType > > & injectors() const
Return const access to the injection model.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
scalar Dij(const label i, const label j) const
Mean diameter Dij.
const dictionary & subModelProperties() const
Return reference to the sub-models dictionary.
const volVectorField & U() const
Return carrier gas velocity.
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
vector linearMomentumOfSystem() const
Total linear momentum of the system.
scalar Dmax() const
Max diameter.
const SurfaceFilmModel< KinematicCloud< CloudType > > & surfaceFilm() const
Return const-access to the surface film model.
List of injection models.
const tmp< volScalarField > alpha() const
Return the particle mass fraction field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const volScalarField & rho() const
Return carrier gas density.
const scalarField & cellLengthScale() const
Return the cell length scale.