Loading...
Searching...
No Matches
KinematicCloudI.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
29#include "fvmSup.H"
30
31// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32
33template<class CloudType>
36{
37 return *cloudCopyPtr_;
38}
39
40
41template<class CloudType>
44 return mesh_;
45}
46
47
48template<class CloudType>
49inline const Foam::IOdictionary&
53}
54
55
56template<class CloudType>
57inline const Foam::IOdictionary&
62
63
64template<class CloudType>
67 return outputProperties_;
68}
69
70
71template<class CloudType>
72inline const Foam::cloudSolution&
77
78
79template<class CloudType>
82 return solution_;
83}
84
85
86template<class CloudType>
87inline const typename CloudType::particleType::constantProperties&
90 return constProps_;
91}
92
93
94template<class CloudType>
95inline typename CloudType::particleType::constantProperties&
98 return constProps_;
99}
100
101
102template<class CloudType>
103inline const Foam::dictionary&
108
109
110template<class CloudType>
112{
113 return rho_;
114}
115
116
117template<class CloudType>
119{
120 return U_;
121}
122
123
124template<class CloudType>
126{
127 return mu_;
128}
129
130
131template<class CloudType>
133{
134 return g_;
135}
136
137
138template<class CloudType>
140{
141 return pAmbient_;
142}
143
144
145template<class CloudType>
147{
148 return pAmbient_;
149}
150
151
152template<class CloudType>
153//inline const typename CloudType::parcelType::forceType&
157 return forces_;
158}
159
160
161template<class CloudType>
165 return forces_;
166}
167
168
169template<class CloudType>
173 return functions_;
174}
175
176
177template<class CloudType>
181 return injectors_;
182}
183
184
185template<class CloudType>
189 return injectors_;
190}
191
192
193template<class CloudType>
197 return *dispersionModel_;
198}
199
200
201template<class CloudType>
205 return *dispersionModel_;
206}
207
208
209template<class CloudType>
214}
215
216
217template<class CloudType>
222}
223
224
225template<class CloudType>
230}
231
232
233template<class CloudType>
238}
239
240
241template<class CloudType>
245 return *surfaceFilmModel_;
246}
247
248
249template<class CloudType>
253 return *surfaceFilmModel_;
254}
255
256
257template<class CloudType>
261 return *packingModel_;
262}
263
264
265template<class CloudType>
269 return *packingModel_;
270}
271
272
273template<class CloudType>
277 return *dampingModel_;
278}
279
280
281template<class CloudType>
285 return *dampingModel_;
286}
287
288
289template<class CloudType>
293 return *isotropyModel_;
294}
295
296
297template<class CloudType>
301 return *isotropyModel_;
302}
303
304
305template<class CloudType>
306inline const Foam::integrationScheme&
308{
309 return *UIntegrator_;
310}
311
312
313template<class CloudType>
314inline Foam::scalar Foam::KinematicCloud<CloudType>::massInSystem() const
315{
316 scalar sysMass = 0.0;
317 for (const parcelType& p : *this)
318 {
319 sysMass += p.nParticle()*p.mass();
320 }
322 return sysMass;
323}
324
325
326template<class CloudType>
327inline Foam::vector
329{
330 vector linearMomentum(Zero);
331
332 for (const parcelType& p : *this)
333 {
334 linearMomentum += p.nParticle()*p.mass()*p.U();
335 }
337 return linearMomentum;
338}
339
340
341template<class CloudType>
342inline Foam::scalar
344{
345 scalar parPerParcel = 0;
346
347 for (const parcelType& p : *this)
348 {
349 parPerParcel += p.nParticle();
350 }
352 return parPerParcel;
353}
354
355
356template<class CloudType>
357inline Foam::scalar
359{
360 scalar linearKineticEnergy = 0;
361
362 for (const parcelType& p : *this)
363 {
364 linearKineticEnergy += p.nParticle()*0.5*p.mass()*(p.U() & p.U());
366
367 return linearKineticEnergy;
368}
369
370
371template<class CloudType>
372inline Foam::scalar Foam::KinematicCloud<CloudType>::Dij
373(
374 const label i,
375 const label j
376) const
377{
378 scalar si = 0.0;
379 scalar sj = 0.0;
380 for (const parcelType& p : *this)
381 {
382 si += p.nParticle()*pow(p.d(), i);
383 sj += p.nParticle()*pow(p.d(), j);
384 }
385
386 reduce(si, sumOp<scalar>());
387 reduce(sj, sumOp<scalar>());
388 sj = max(sj, VSMALL);
389
390 return si/sj;
391}
392
393
394template<class CloudType>
395inline Foam::scalar Foam::KinematicCloud<CloudType>::Dmax() const
396{
397 scalar d = -GREAT;
398 for (const parcelType& p : *this)
399 {
400 d = max(d, p.d());
401 }
402
404
405 return max(0.0, d);
406}
407
408
409template<class CloudType>
412 return rndGen_;
413}
414
415
416template<class CloudType>
419{
420 if (!cellOccupancyPtr_)
421 {
422 buildCellOccupancy();
423 }
425 return *cellOccupancyPtr_;
426}
427
428
429template<class CloudType>
430inline const Foam::scalarField&
432{
433 return cellLengthScale_;
434}
435
436
437template<class CloudType>
439(
440 const parcelType& p,
441 const typename parcelType::trackingData& td
442)
443{
444 const scalar m = p.nParticle()*p.mass();
445
446 rhokTrans()[p.cell()] += m;
448 UTrans()[p.cell()] += m*p.U();
449}
450
451
452template<class CloudType>
456 return *rhokTrans_;
457}
458
459
460template<class CloudType>
464 return *rhokTrans_;
465}
466
467
468template<class CloudType>
472 return *UTrans_;
473}
474
475
476template<class CloudType>
480 return *UTrans_;
482
483
484template<class CloudType>
488 return *UCoeff_;
489}
490
491
492template<class CloudType>
496 return *UCoeff_;
497}
498
499
500template<class CloudType>
503{
504 if (debug)
505 {
506 Pout<< "rhokTrans min/max = " << min(rhokTrans()).value() << ", "
507 << max(rhokTrans()).value() << endl;
508 }
509
510 if (this->solution().coupled())
511 {
512 return rhokTrans()/this->db().time().deltaT()/this->mesh().V();
513 }
514
515 return tmp<volScalarField::Internal>::New
516 (
517 this->newIOobject(IOobject::scopedName(this->name(), "rhokTrans")),
518 this->mesh(),
520 (
521 rhokTrans().dimensions()/dimTime/dimVolume, Zero
523 );
524}
526
527template<class CloudType>
530const
531{
532 if (debug)
533 {
534 Pout<< "UTrans min/max = " << min(UTrans()).value() << ", "
535 << max(UTrans()).value() << nl
536 << "UCoeff min/max = " << min(UCoeff()).value() << ", "
537 << max(UCoeff()).value() << endl;
538 }
541 if (incompressible)
542 {
543 dim.reset(dimForce/dimDensity);
544 }
545
546 if (solution_.coupled())
548 if (solution_.semiImplicit("U"))
549 {
551 Vdt(mesh_.V()*this->db().time().deltaT());
553 if (incompressible)
554 {
555 Vdt.dimensions() *= dimDensity;
556 }
558 return UTrans()/Vdt - fvm::Sp(UCoeff()/Vdt, U) + UCoeff()/Vdt*U;
559 }
560 else
561 {
562 auto tfvm = tmp<fvVectorMatrix>::New(U, dim);
563 auto& fvm = tfvm.ref();
564
565 fvm.source() = -UTrans()/(this->db().time().deltaT());
566
567 return tfvm;
568 }
569 }
571 return tmp<fvVectorMatrix>::New(U, dim);
572}
573
574
575template<class CloudType>
578{
579 auto tvDotSweep = tmp<volScalarField>::New
580 (
581 this->newIOobject(IOobject::scopedName(this->name(), "vDotSweep")),
582 mesh_,
585 );
586 auto& vDotSweep = tvDotSweep.ref();
587
588 for (const parcelType& p : *this)
589 {
590 const label celli = p.cell();
592 vDotSweep[celli] += p.nParticle()*p.areaP()*mag(p.U() - U_[celli]);
593 }
594
595 vDotSweep.primitiveFieldRef() /= mesh_.V();
596 vDotSweep.correctBoundaryConditions();
598 return tvDotSweep;
599}
601
602template<class CloudType>
605{
606 auto ttheta = tmp<volScalarField>::New
607 (
608 this->newIOobject(IOobject::scopedName(this->name(), "theta")),
609 mesh_,
612 );
613 auto& theta = ttheta.ref();
614
615 for (const parcelType& p : *this)
616 {
617 const label celli = p.cell();
618
619 theta[celli] += p.nParticle()*p.volume();
620 }
621
622 theta.primitiveFieldRef() /= mesh_.V();
623 theta.correctBoundaryConditions();
625 return ttheta;
626}
627
628
629template<class CloudType>
632{
633 auto talpha = tmp<volScalarField>::New
634 (
635 this->newIOobject(IOobject::scopedName(this->name(), "alpha")),
636 mesh_,
638 );
639
640 scalarField& alpha = talpha.ref().primitiveFieldRef();
641 for (const parcelType& p : *this)
642 {
643 const label celli = p.cell();
644
645 alpha[celli] += p.nParticle()*p.mass();
646 }
647
648 alpha /= (mesh_.V()*rho_);
650 return talpha;
651}
652
653
654template<class CloudType>
657{
658 auto trhoEff = tmp<volScalarField>::New
659 (
660 this->newIOobject(IOobject::scopedName(this->name(), "rhoEff")),
661 mesh_,
663 );
664
665 scalarField& rhoEff = trhoEff.ref().primitiveFieldRef();
666 for (const parcelType& p : *this)
667 {
668 const label celli = p.cell();
669
670 rhoEff[celli] += p.nParticle()*p.mass();
671 }
672
673 rhoEff /= mesh_.V();
675 return trhoEff;
676}
677
678
679// ************************************************************************* //
Base class for collisional damping models.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const noexcept
Return dimensions.
Base class for dispersion modelling.
DimensionedField< scalar, volMesh > Internal
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
List of injection models.
Base class for collisional return-to-isotropy models.
Templated base class for kinematic cloud.
const parcelType::constantProperties & constProps() const
Return the constant properties.
const dictionary & subModelProperties() const
Return reference to the sub-models dictionary.
autoPtr< volScalarField::Internal > rhokTrans_
Mass for kinematic cloud.
cloudSolution solution_
Solution properties.
autoPtr< DampingModel< KinematicCloud< CloudType > > > dampingModel_
Damping model.
autoPtr< volVectorField::Internal > UTrans_
Momentum.
scalar massInSystem() const
Total mass in system.
scalar totalParticlePerParcel() const
Average particle per parcel.
forceType forces_
Optional particle forces.
tmp< volScalarField::Internal > Srhok() const
Return tmp mass source for kinematic.
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
const fvMesh & mesh_
References to the mesh and time databases.
autoPtr< SurfaceFilmModel< KinematicCloud< CloudType > > > surfaceFilmModel_
Surface film model.
const tmp< volScalarField > alpha() const
Return the particle mass fraction field.
const scalarField & cellLengthScale() const
Return the cell length scale.
scalar pAmbient_
Averaged ambient domain pressure.
const SurfaceFilmModel< KinematicCloud< CloudType > > & surfaceFilm() const
Return const-access to the surface film model.
scalarField cellLengthScale_
Cell length scale.
const PatchInteractionModel< KinematicCloud< CloudType > > & patchInteraction() const
Return const-access to the patch interaction model.
autoPtr< integrationScheme > UIntegrator_
Velocity integration.
vector linearMomentumOfSystem() const
Total linear momentum of the system.
const volVectorField & U() const
Return carrier gas velocity.
scalar pAmbient() const
Return const-access to the ambient pressure.
const KinematicCloud & cloudCopy() const
Return a reference to the cloud copy.
InjectionModelList< KinematicCloud< CloudType > > injectors_
Injector models.
const volScalarField & rho() const
Return carrier gas density.
autoPtr< IsotropyModel< KinematicCloud< CloudType > > > isotropyModel_
Exchange model.
autoPtr< DispersionModel< KinematicCloud< CloudType > > > dispersionModel_
Dispersion model.
autoPtr< PatchInteractionModel< KinematicCloud< CloudType > > > patchInteractionModel_
Patch interaction model.
const volVectorField & U_
Velocity [m/s].
const integrationScheme & UIntegrator() const
Return reference to velocity integration.
List< DynamicList< parcelType * > > & cellOccupancy()
Return the cell occupancy information for each.
const volScalarField & mu_
Dynamic viscosity [Pa.s].
const dimensionedVector & g_
Gravity.
functionType & functions()
Optional cloud function objects.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
const tmp< volScalarField > rhoEff() const
Return the particle effective density field.
const cloudSolution & solution() const
Return const access to the solution properties.
const dictionary subModelProperties_
Sub-models dictionary.
const tmp< volScalarField > vDotSweep() const
Volume swept rate of parcels per cell.
const PackingModel< KinematicCloud< CloudType > > & packingModel() const
Return const access to the packing model.
const DispersionModel< KinematicCloud< CloudType > > & dispersion() const
Return const-access to the dispersion model.
volScalarField::Internal & UCoeff()
Return coefficient for carrier phase U equation.
functionType functions_
Optional cloud function objects.
const StochasticCollisionModel< KinematicCloud< CloudType > > & stochasticCollision() const
Return const-access to the stochastic collision model.
parcelType::constantProperties constProps_
Parcel constant properties.
scalar Dmax() const
Max diameter.
autoPtr< List< DynamicList< parcelType * > > > cellOccupancyPtr_
Cell occupancy information for each parcel, (demand driven).
const DampingModel< KinematicCloud< CloudType > > & dampingModel() const
Return const access to the damping model.
const volScalarField & rho_
Density [kg/m3].
CloudFunctionObjectList< KinematicCloud< CloudType > > functionType
Function object type.
scalar Dij(const label i, const label j) const
Mean diameter Dij.
IOdictionary particleProperties_
Dictionary of particle properties.
const IOdictionary & particleProperties() const
Return particle properties dictionary.
volVectorField::Internal & UTrans()
Return reference to momentum source.
tmp< fvVectorMatrix > SU(volVectorField &U, bool incompressible=false) const
Return tmp momentum source term (compressible).
const dimensionedVector & g() const
Gravity.
void transferToCarrier(const parcelType &p, const typename parcelType::trackingData &td)
Transfer the effect of parcel to the carrier phase.
const fvMesh & mesh() const
Return reference to the mesh.
Random rndGen_
Random number generator - used by some injection routines.
const InjectionModelList< KinematicCloud< CloudType > > & injectors() const
Return const access to the injection model.
Random & rndGen() const
Return reference to the random object.
volScalarField::Internal & rhokTrans()
Return reference to mass for kinematic source.
const forceType & forces() const
Optional particle forces.
ParticleForceList< KinematicCloud< CloudType > > forceType
Force models type.
autoPtr< PackingModel< KinematicCloud< CloudType > > > packingModel_
Packing model.
const IOdictionary & outputProperties() const
Return output properties dictionary.
const IsotropyModel< KinematicCloud< CloudType > > & isotropyModel() const
Return const access to the isotropy model.
void buildCellOccupancy()
Build the cellOccupancy.
IOdictionary outputProperties_
Dictionary of output properties.
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
const volScalarField & mu() const
Return carrier gas dynamic viscosity.
autoPtr< StochasticCollisionModel< KinematicCloud< CloudType > > > stochasticCollisionModel_
Stochastic collision model.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
Base class for packing models.
Templated patch interaction model class.
Random number generator.
Definition Random.H:56
Templated stochastic collision model class.
Templated wall surface film model class.
Stores all relevant solution info for cloud.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Base for a set of schemes which integrate simple ODEs which arise from semi-implcit rate expressions.
IOobject newIOobject(const word &name, IOobjectOption ioOpt) const
Create an IOobject at the current time instance (timeName) with the specified options.
Selector class for relaxation factors, solver type and solution.
Definition solution.H:95
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
U
Definition pEqn.H:72
volScalarField & p
bool coupled
dynamicFvMesh & mesh
Calculate the finiteVolume matrix for implicit and explicit sources.
auto & name
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for handling debugging switches.
Definition debug.C:45
Namespace of functions to calculate implicit derivatives returning a matrix.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
const dimensionSet dimForce
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimVolume(pow3(dimLength))
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
const dimensionSet dimDensity
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha