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-2020 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\*---------------------------------------------------------------------------*/
28
29#include "fvmSup.H"
30
31// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32
33template<class CloudType>
36{
37 return *cloudCopyPtr_;
38}
39
40
41template<class CloudType>
43{
44 return mesh_;
45}
46
47
48template<class CloudType>
49inline const Foam::IOdictionary&
51{
52 return particleProperties_;
53}
54
55
56template<class CloudType>
57inline const Foam::IOdictionary&
59{
60 return outputProperties_;
61}
62
63
64template<class CloudType>
66{
67 return outputProperties_;
68}
69
70
71template<class CloudType>
72inline const Foam::cloudSolution&
74{
75 return solution_;
76}
77
78
79template<class CloudType>
81{
82 return solution_;
83}
84
85
86template<class CloudType>
89{
90 return constProps_;
91}
92
93
94template<class CloudType>
97{
98 return constProps_;
99}
100
101
102template<class CloudType>
103inline const Foam::dictionary&
105{
106 return subModelProperties_;
107}
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&
156{
157 return forces_;
158}
159
160
161template<class CloudType>
164{
165 return forces_;
166}
167
168
169template<class CloudType>
172{
173 return functions_;
174}
175
176
177template<class CloudType>
180{
181 return injectors_;
182}
183
184
185template<class CloudType>
188{
189 return injectors_;
190}
191
192
193template<class CloudType>
196{
197 return *dispersionModel_;
198}
199
200
201template<class CloudType>
204{
205 return *dispersionModel_;
206}
207
208
209template<class CloudType>
212{
213 return *patchInteractionModel_;
214}
215
216
217template<class CloudType>
220{
221 return *patchInteractionModel_;
222}
223
224
225template<class CloudType>
228{
229 return *stochasticCollisionModel_;
230}
231
232
233template<class CloudType>
236{
237 return *stochasticCollisionModel_;
238}
239
240
241template<class CloudType>
244{
245 return *surfaceFilmModel_;
246}
247
248
249template<class CloudType>
252{
253 return *surfaceFilmModel_;
254}
255
256
257template<class CloudType>
260{
261 return *packingModel_;
262}
263
264
265template<class CloudType>
268{
269 return *packingModel_;
270}
271
272
273template<class CloudType>
276{
277 return *dampingModel_;
278}
279
280
281template<class CloudType>
284{
285 return *dampingModel_;
286}
287
288
289template<class CloudType>
292{
293 return *isotropyModel_;
294}
295
296
297template<class CloudType>
300{
301 return *isotropyModel_;
302}
303
304
305template<class CloudType>
306inline const Foam::integrationScheme&
308{
309 return *UIntegrator_;
310}
311
312
313template<class CloudType>
315{
316 scalar sysMass = 0.0;
317 for (const parcelType& p : *this)
318 {
319 sysMass += p.nParticle()*p.mass();
320 }
321
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 }
336
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 }
351
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());
365 }
367 return linearKineticEnergy;
368}
370
371template<class CloudType>
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)
382 si += p.nParticle()*pow(p.d(), i);
383 sj += p.nParticle()*pow(p.d(), j);
384 }
386 reduce(si, sumOp<scalar>());
387 reduce(sj, sumOp<scalar>());
388 sj = max(sj, VSMALL);
389
390 return si/sj;
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);
407
408
409template<class CloudType>
411{
412 return rndGen_;
413}
414
416template<class CloudType>
419{
420 if (!cellOccupancyPtr_)
421 {
422 buildCellOccupancy();
423 }
425 return *cellOccupancyPtr_;
426}
428
429template<class CloudType>
430inline const Foam::scalarField&
432{
433 return cellLengthScale_;
435
436
437template<class CloudType>
441 return *UTrans_;
442}
443
444
445template<class CloudType>
448{
449 return *UTrans_;
450}
452
453template<class CloudType>
456{
457 return *UCoeff_;
458}
460
461template<class CloudType>
464{
465 return *UCoeff_;
466}
468
469template<class CloudType>
472const
473{
474 if (debug)
475 {
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;
481
482 dimensionSet dim(dimForce);
483 if (incompressible)
485 dim.reset(dimForce/dimDensity);
486 }
487
488 if (solution_.coupled())
490 if (solution_.semiImplicit("U"))
491 {
493 Vdt(mesh_.V()*this->db().time().deltaT());
494
495 if (incompressible)
496 {
497 Vdt.dimensions() *= dimDensity;
498 }
499
500 return UTrans()/Vdt - fvm::Sp(UCoeff()/Vdt, U) + UCoeff()/Vdt*U;
502 else
503 {
504 tmp<fvVectorMatrix> tfvm(new fvVectorMatrix(U, dim));
505 fvVectorMatrix& fvm = tfvm.ref();
506
507 fvm.source() = -UTrans()/(this->db().time().deltaT());
508
509 return tfvm;
510 }
511 }
512
513 return tmp<fvVectorMatrix>::New(U, dim);
514}
516
517template<class CloudType>
520{
521 tmp<volScalarField> tvDotSweep
522 (
524 (
526 (
527 this->name() + ":vDotSweep",
528 this->db().time().timeName(),
529 this->db(),
530 IOobject::NO_READ,
531 IOobject::NO_WRITE,
532 false
533 ),
534 mesh_,
535 dimensionedScalar(dimless/dimTime, Zero),
536 extrapolatedCalculatedFvPatchScalarField::typeName
538 );
539
540 volScalarField& vDotSweep = tvDotSweep.ref();
541 for (const parcelType& p : *this)
542 {
543 const label celli = p.cell();
544
545 vDotSweep[celli] += p.nParticle()*p.areaP()*mag(p.U() - U_[celli]);
546 }
547
548 vDotSweep.primitiveFieldRef() /= mesh_.V();
549 vDotSweep.correctBoundaryConditions();
550
551 return tvDotSweep;
552}
554
555template<class CloudType>
558{
560 (
564 (
565 this->name() + ":theta",
566 this->db().time().timeName(),
567 this->db(),
568 IOobject::NO_READ,
569 IOobject::NO_WRITE,
570 false
571 ),
572 mesh_,
573 dimensionedScalar(dimless, Zero),
574 extrapolatedCalculatedFvPatchScalarField::typeName
575 )
576 );
578 volScalarField& theta = ttheta.ref();
579 for (const parcelType& p : *this)
580 {
581 const label celli = p.cell();
582
583 theta[celli] += p.nParticle()*p.volume();
584 }
586 theta.primitiveFieldRef() /= mesh_.V();
588
589 return ttheta;
590}
591
592
593template<class CloudType>
596{
598 (
600 (
602 (
603 this->name() + ":alpha",
604 this->db().time().timeName(),
605 this->db(),
606 IOobject::NO_READ,
607 IOobject::NO_WRITE,
608 false
609 ),
610 mesh_,
612 )
613 );
614
615 scalarField& alpha = talpha.ref().primitiveFieldRef();
616 for (const parcelType& p : *this)
617 {
618 const label celli = p.cell();
619
620 alpha[celli] += p.nParticle()*p.mass();
621 }
622
623 alpha /= (mesh_.V()*rho_);
624
625 return talpha;
626}
627
628
629template<class CloudType>
632{
633 tmp<volScalarField> trhoEff
634 (
636 (
638 (
639 this->name() + ":rhoEff",
640 this->db().time().timeName(),
641 this->db(),
642 IOobject::NO_READ,
643 IOobject::NO_WRITE,
644 false
645 ),
646 mesh_,
648 )
649 );
650
651 scalarField& rhoEff = trhoEff.ref().primitiveFieldRef();
652 for (const parcelType& p : *this)
653 {
654 const label celli = p.cell();
655
656 rhoEff[celli] += p.nParticle()*p.mass();
657 }
658
659 rhoEff /= mesh_.V();
660
661 return trhoEff;
662}
663
664
665// ************************************************************************* //
Y[inertIndex] max(0.0)
reduce(hasMovingMesh, orOp< bool >())
Class to hold DSMC particle constant properties.
Definition: DSMCParcel.H:82
Base class for collisional damping models.
Definition: DampingModel.H:66
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Base class for dispersion modelling.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
void correctBoundaryConditions()
Correct boundary field.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
List of injection models.
Base class for collisional return-to-isotropy models.
Definition: IsotropyModel.H:66
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.
scalar massInSystem() const
Total mass in system.
scalar totalParticlePerParcel() const
Average particle per parcel.
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
const tmp< volScalarField > alpha() const
Return the particle mass fraction field.
const scalarField & cellLengthScale() const
Return the cell length scale.
const SurfaceFilmModel< KinematicCloud< CloudType > > & surfaceFilm() const
Return const-access to the surface film model.
const PatchInteractionModel< KinematicCloud< CloudType > > & patchInteraction() const
Return const-access to the patch interaction model.
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.
const volScalarField & rho() const
Return carrier gas density.
const integrationScheme & UIntegrator() const
Return reference to velocity integration.
List< DynamicList< parcelType * > > & cellOccupancy()
Return the cell occupancy information for each.
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 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.
const StochasticCollisionModel< KinematicCloud< CloudType > > & stochasticCollision() const
Return const-access to the stochastic collision model.
scalar Dmax() const
Max diameter.
const DampingModel< KinematicCloud< CloudType > > & dampingModel() const
Return const access to the damping model.
scalar Dij(const label i, const label j) const
Mean diameter Dij.
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.
const fvMesh & mesh() const
Return reference to the mesh.
const InjectionModelList< KinematicCloud< CloudType > > & injectors() const
Return const access to the injection model.
Random & rndGen() const
Return reference to the random object.
const forceType & forces() const
Optional particle forces.
const IOdictionary & outputProperties() const
Return output properties dictionary.
const IsotropyModel< KinematicCloud< CloudType > > & isotropyModel() const
Return const access to the isotropy model.
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
const volScalarField & mu() const
Return carrier gas dynamic viscosity.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
Base class for packing models.
Definition: PackingModel.H:71
Templated patch interaction model class.
Random number generator.
Definition: Random.H:60
Templated stochastic collision model class.
Templated wall surface film model class.
Stores all relevant solution info for cloud.
Definition: cloudSolution.H:54
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Field< Type > & source() noexcept
Definition: fvMatrix.H:458
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Base for a set of schemes which integrate simple ODEs which arise from semi-implcit rate expressions.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
U
Definition: pEqn.H:72
volScalarField & p
Calculate the finiteVolume matrix for implicit and explicit sources.
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimless
Dimensionless.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
const dimensionSet dimForce
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const dimensionSet dimDensity
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & alpha