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 -------------------------------------------------------------------------------
11 License
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 
33 template<class CloudType>
36 {
37  return *cloudCopyPtr_;
38 }
39 
40 
41 template<class CloudType>
43 {
44  return mesh_;
45 }
46 
47 
48 template<class CloudType>
49 inline const Foam::IOdictionary&
51 {
52  return particleProperties_;
53 }
54 
55 
56 template<class CloudType>
57 inline const Foam::IOdictionary&
59 {
60  return outputProperties_;
61 }
62 
63 
64 template<class CloudType>
66 {
67  return outputProperties_;
68 }
69 
70 
71 template<class CloudType>
72 inline const Foam::cloudSolution&
74 {
75  return solution_;
76 }
77 
78 
79 template<class CloudType>
81 {
82  return solution_;
83 }
84 
85 
86 template<class CloudType>
87 inline const typename CloudType::particleType::constantProperties&
89 {
90  return constProps_;
91 }
92 
93 
94 template<class CloudType>
95 inline typename CloudType::particleType::constantProperties&
97 {
98  return constProps_;
99 }
100 
101 
102 template<class CloudType>
103 inline const Foam::dictionary&
105 {
106  return subModelProperties_;
107 }
108 
109 
110 template<class CloudType>
112 {
113  return rho_;
114 }
115 
116 
117 template<class CloudType>
119 {
120  return U_;
121 }
122 
123 
124 template<class CloudType>
126 {
127  return mu_;
128 }
129 
130 
131 template<class CloudType>
133 {
134  return g_;
135 }
136 
137 
138 template<class CloudType>
139 inline Foam::scalar Foam::KinematicCloud<CloudType>::pAmbient() const
140 {
141  return pAmbient_;
142 }
143 
144 
145 template<class CloudType>
146 inline Foam::scalar& Foam::KinematicCloud<CloudType>::pAmbient()
147 {
148  return pAmbient_;
149 }
150 
151 
152 template<class CloudType>
153 //inline const typename CloudType::parcelType::forceType&
154 inline const typename Foam::KinematicCloud<CloudType>::forceType&
156 {
157  return forces_;
158 }
159 
160 
161 template<class CloudType>
164 {
165  return forces_;
166 }
167 
168 
169 template<class CloudType>
172 {
173  return functions_;
174 }
175 
176 
177 template<class CloudType>
180 {
181  return injectors_;
182 }
183 
184 
185 template<class CloudType>
188 {
189  return injectors_;
190 }
191 
192 
193 template<class CloudType>
196 {
197  return *dispersionModel_;
198 }
199 
200 
201 template<class CloudType>
204 {
205  return *dispersionModel_;
206 }
207 
208 
209 template<class CloudType>
212 {
213  return *patchInteractionModel_;
214 }
215 
216 
217 template<class CloudType>
220 {
221  return *patchInteractionModel_;
222 }
223 
224 
225 template<class CloudType>
228 {
229  return *stochasticCollisionModel_;
230 }
231 
232 
233 template<class CloudType>
236 {
237  return *stochasticCollisionModel_;
238 }
239 
240 
241 template<class CloudType>
244 {
245  return *surfaceFilmModel_;
246 }
247 
248 
249 template<class CloudType>
252 {
253  return *surfaceFilmModel_;
254 }
255 
256 
257 template<class CloudType>
260 {
261  return *packingModel_;
262 }
263 
264 
265 template<class CloudType>
268 {
269  return *packingModel_;
270 }
271 
272 
273 template<class CloudType>
276 {
277  return *dampingModel_;
278 }
279 
280 
281 template<class CloudType>
284 {
285  return *dampingModel_;
286 }
287 
288 
289 template<class CloudType>
292 {
293  return *isotropyModel_;
294 }
295 
296 
297 template<class CloudType>
300 {
301  return *isotropyModel_;
302 }
303 
304 
305 template<class CloudType>
306 inline const Foam::integrationScheme&
308 {
309  return *UIntegrator_;
310 }
311 
312 
313 template<class CloudType>
314 inline 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 
326 template<class CloudType>
327 inline 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 
341 template<class CloudType>
342 inline 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 
356 template<class CloudType>
357 inline 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 }
369 
370 
371 template<class CloudType>
372 inline 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 
394 template<class CloudType>
395 inline 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 
409 template<class CloudType>
411 {
412  return rndGen_;
413 }
414 
415 
416 template<class CloudType>
419 {
420  if (!cellOccupancyPtr_)
421  {
422  buildCellOccupancy();
423  }
425  return *cellOccupancyPtr_;
426 }
427 
428 
429 template<class CloudType>
430 inline const Foam::scalarField&
432 {
433  return cellLengthScale_;
434 }
435 
436 
437 template<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 
452 template<class CloudType>
455 {
456  return *rhokTrans_;
457 }
458 
459 
460 template<class CloudType>
463 {
464  return *rhokTrans_;
465 }
466 
467 
468 template<class CloudType>
471 {
472  return *UTrans_;
473 }
474 
475 
476 template<class CloudType>
479 {
480  return *UTrans_;
481 }
482 
483 
484 template<class CloudType>
487 {
488  return *UCoeff_;
489 }
490 
491 
492 template<class CloudType>
495 {
496  return *UCoeff_;
497 }
498 
499 
500 template<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 
516  (
517  this->newIOobject(IOobject::scopedName(this->name(), "rhokTrans")),
518  this->mesh(),
520  (
521  rhokTrans().dimensions()/dimTime/dimVolume, Zero
522  )
523  );
524 }
525 
526 
527 template<class CloudType>
530 const
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  }
539 
540  dimensionSet dim(dimForce);
541  if (incompressible)
542  {
543  dim.reset(dimForce/dimDensity);
544  }
545 
546  if (solution_.coupled())
547  {
548  if (solution_.semiImplicit("U"))
549  {
550  volScalarField::Internal
551  Vdt(mesh_.V()*this->db().time().deltaT());
552 
553  if (incompressible)
554  {
555  Vdt.dimensions() *= dimDensity;
556  }
557 
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 
575 template<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();
591 
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 }
600 
601 
602 template<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 
629 template<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 
654 template<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();
674 
675  return trhoEff;
676 }
677 
678 
679 // ************************************************************************* //
tmp< fvVectorMatrix > SU(volVectorField &U, bool incompressible=false) const
Return tmp momentum source term (compressible)
const SurfaceFilmModel< KinematicCloud< CloudType > > & surfaceFilm() const
Return const-access to the surface film model.
const StochasticCollisionModel< KinematicCloud< CloudType > > & stochasticCollision() const
Return const-access to the stochastic collision model.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar totalParticlePerParcel() const
Average particle per parcel.
const scalarField & cellLengthScale() const
Return the cell length scale.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
scalar Dmax() const
Max diameter.
const InjectionModelList< KinematicCloud< CloudType > > & injectors() const
Return const access to the injection model.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
const KinematicCloud & cloudCopy() const
Return a reference to the cloud copy.
const volScalarField & rho() const
Return carrier gas density.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
const fvMesh & mesh() const
Return reference to the mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
volVectorField::Internal & UTrans()
Return reference to momentum source.
Base class for packing models.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const IsotropyModel< KinematicCloud< CloudType > > & isotropyModel() const
Return const access to the isotropy model.
const tmp< volScalarField > rhoEff() const
Return the particle effective density field.
const DampingModel< KinematicCloud< CloudType > > & dampingModel() const
Return const access to the damping model.
const dimensionSet dimless
Dimensionless.
List< DynamicList< parcelType * > > & cellOccupancy()
Return the cell occupancy information for each.
const parcelType::constantProperties & constProps() const
Return the constant properties.
Templated patch interaction model class.
scalar massInSystem() const
Total mass in system.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
const forceType & forces() const
Optional particle forces.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
volScalarField::Internal & rhokTrans()
Return reference to mass for kinematic source.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
volScalarField::Internal & UCoeff()
Return coefficient for carrier phase U equation.
const PackingModel< KinematicCloud< CloudType > > & packingModel() const
Return const access to the packing model.
tmp< volScalarField::Internal > Srhok() const
Return tmp mass source for kinematic.
List of injection models.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
Random number generator.
Definition: Random.H:55
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
const dimensionSet dimForce
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
const tmp< volScalarField > vDotSweep() const
Volume swept rate of parcels per cell.
const DispersionModel< KinematicCloud< CloudType > > & dispersion() const
Return const-access to the dispersion model.
scalar pAmbient() const
Return const-access to the ambient pressure.
Templated wall surface film model class.
const dimensionSet dimDensity
const cloudSolution & solution() const
Return const access to the solution properties.
Base for a set of schemes which integrate simple ODEs which arise from semi-implcit rate expressions...
void transferToCarrier(const parcelType &p, const typename parcelType::trackingData &td)
Transfer the effect of parcel to the carrier phase.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Base class for collisional return-to-isotropy models.
const IOdictionary & particleProperties() const
Return particle properties dictionary.
Templated base class for kinematic cloud.
Stores all relevant solution info for cloud.
Definition: cloudSolution.H:49
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const integrationScheme & UIntegrator() const
Return reference to velocity integration.
Templated stochastic collision model class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Base class for dispersion modelling.
vector linearMomentumOfSystem() const
Total linear momentum of the system.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionedVector & g() const
Gravity.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const IOdictionary & outputProperties() const
Return output properties dictionary.
const volVectorField & U() const
Return carrier gas velocity.
const PatchInteractionModel< KinematicCloud< CloudType > > & patchInteraction() const
Return const-access to the patch interaction model.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
bool coupled
scalar Dij(const label i, const label j) const
Mean diameter Dij.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
const volScalarField & mu() const
Return carrier gas dynamic viscosity.
const word extrapolatedCalculatedType
A combined zero-gradient and calculated patch field type.
const tmp< volScalarField > alpha() const
Return the particle mass fraction field.
functionType & functions()
Optional cloud function objects.
Base class for collisional damping models.
Random & rndGen() const
Return reference to the random object.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
Calculate the finiteVolume matrix for implicit and explicit sources.
const dictionary & subModelProperties() const
Return reference to the sub-models dictionary.
zeroField Sp
Definition: alphaSuSp.H:2
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127