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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvmSup.H"
27 #include "SortableList.H"
28 
29 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30 
31 template<class CloudType>
34 {
35  return cloudCopyPtr_();
36 }
37 
38 
39 template<class CloudType>
41 {
42  return mesh_;
43 }
44 
45 
46 template<class CloudType>
47 inline const Foam::IOdictionary&
49 {
50  return particleProperties_;
51 }
52 
53 
54 template<class CloudType>
55 inline const Foam::IOdictionary&
57 {
58  return outputProperties_;
59 }
60 
61 
62 template<class CloudType>
64 {
65  return outputProperties_;
66 }
67 
68 
69 template<class CloudType>
70 inline const Foam::cloudSolution&
72 {
73  return solution_;
74 }
75 
76 
77 template<class CloudType>
79 {
80  return solution_;
81 }
82 
83 
84 template<class CloudType>
85 inline const typename CloudType::particleType::constantProperties&
87 {
88  return constProps_;
89 }
90 
91 
92 template<class CloudType>
93 inline typename CloudType::particleType::constantProperties&
95 {
96  return constProps_;
97 }
98 
99 
100 template<class CloudType>
101 inline const Foam::dictionary&
103 {
104  return subModelProperties_;
105 }
106 
107 
108 template<class CloudType>
110 {
111  return rho_;
112 }
113 
114 
115 template<class CloudType>
117 {
118  return U_;
119 }
120 
121 
122 template<class CloudType>
124 {
125  return mu_;
126 }
127 
128 
129 template<class CloudType>
131 {
132  return g_;
133 }
134 
135 
136 template<class CloudType>
137 inline Foam::scalar Foam::KinematicCloud<CloudType>::pAmbient() const
138 {
139  return pAmbient_;
140 }
141 
142 
143 template<class CloudType>
144 inline Foam::scalar& Foam::KinematicCloud<CloudType>::pAmbient()
145 {
146  return pAmbient_;
147 }
148 
149 
150 template<class CloudType>
151 //inline const typename CloudType::parcelType::forceType&
152 inline const typename Foam::KinematicCloud<CloudType>::forceType&
154 {
155  return forces_;
156 }
157 
158 
159 template<class CloudType>
162 {
163  return forces_;
164 }
165 
166 
167 template<class CloudType>
170 {
171  return functions_;
172 }
173 
174 
175 template<class CloudType>
178 {
179  return injectors_;
180 }
181 
182 
183 template<class CloudType>
186 {
187  return injectors_;
188 }
189 
190 
191 template<class CloudType>
194 {
195  return dispersionModel_;
196 }
197 
198 
199 template<class CloudType>
202 {
203  return dispersionModel_();
204 }
205 
206 
207 template<class CloudType>
210 {
211  return patchInteractionModel_;
212 }
213 
214 
215 template<class CloudType>
218 {
219  return patchInteractionModel_();
220 }
221 
222 
223 template<class CloudType>
226 {
227  return stochasticCollisionModel_();
228 }
229 
230 
231 template<class CloudType>
234 {
235  return stochasticCollisionModel_();
236 }
237 
238 
239 template<class CloudType>
242 {
243  return surfaceFilmModel_();
244 }
245 
246 
247 template<class CloudType>
250 {
251  return surfaceFilmModel_();
252 }
253 
254 
255 template<class CloudType>
256 inline const Foam::vectorIntegrationScheme&
258 {
259  return UIntegrator_;
260 }
261 
262 
263 template<class CloudType>
265 {
266  return this->size();
267 }
268 
269 
270 template<class CloudType>
272 {
273  scalar sysMass = 0.0;
274  forAllConstIter(typename KinematicCloud<CloudType>, *this, iter)
275  {
276  const parcelType& p = iter();
277  sysMass += p.nParticle()*p.mass();
278  }
279 
280  return sysMass;
281 }
282 
283 
284 template<class CloudType>
285 inline Foam::vector
287 {
288  vector linearMomentum(vector::zero);
289 
290  forAllConstIter(typename KinematicCloud<CloudType>, *this, iter)
291  {
292  const parcelType& p = iter();
293 
294  linearMomentum += p.nParticle()*p.mass()*p.U();
295  }
296 
297  return linearMomentum;
298 }
299 
300 
301 template<class CloudType>
302 inline Foam::scalar
304 {
305  scalar linearKineticEnergy = 0.0;
306 
307  forAllConstIter(typename KinematicCloud<CloudType>, *this, iter)
308  {
309  const parcelType& p = iter();
310 
311  linearKineticEnergy += p.nParticle()*0.5*p.mass()*(p.U() & p.U());
312  }
313 
314  return linearKineticEnergy;
315 }
316 
317 
318 template<class CloudType>
319 inline Foam::scalar Foam::KinematicCloud<CloudType>::Dij
320 (
321  const label i,
322  const label j
323 ) const
324 {
325  scalar si = 0.0;
326  scalar sj = 0.0;
327  forAllConstIter(typename KinematicCloud<CloudType>, *this, iter)
328  {
329  const parcelType& p = iter();
330  si += p.nParticle()*pow(p.d(), i);
331  sj += p.nParticle()*pow(p.d(), j);
332  }
333 
334  reduce(si, sumOp<scalar>());
335  reduce(sj, sumOp<scalar>());
336  sj = max(sj, VSMALL);
337 
338  return si/sj;
339 }
340 
341 
342 template<class CloudType>
343 inline Foam::scalar Foam::KinematicCloud<CloudType>::Dmax() const
344 {
345  scalar d = -GREAT;
346  forAllConstIter(typename KinematicCloud<CloudType>, *this, iter)
347  {
348  const parcelType& p = iter();
349  d = max(d, p.d());
350  }
351 
352  reduce(d, maxOp<scalar>());
353 
354  return max(0.0, d);
355 }
356 
357 
358 template<class CloudType>
360 (
361  const scalar fraction
362 ) const
363 {
364  if ((fraction < 0) || (fraction > 1))
365  {
367  << "fraction should be in the range 0 < fraction < 1"
368  << exit(FatalError);
369  }
370 
371  scalar distance = 0.0;
372 
373  const label nParcel = this->size();
374  globalIndex globalParcels(nParcel);
375  const label nParcelSum = globalParcels.size();
376 
377  if (nParcelSum == 0)
378  {
379  return distance;
380  }
381 
382  // lists of parcels mass and distance from initial injection point
383  List<List<scalar> > procMass(Pstream::nProcs());
384  List<List<scalar> > procDist(Pstream::nProcs());
385 
386  List<scalar>& mass = procMass[Pstream::myProcNo()];
387  List<scalar>& dist = procDist[Pstream::myProcNo()];
388 
389  mass.setSize(nParcel);
390  dist.setSize(nParcel);
391 
392  label i = 0;
393  scalar mSum = 0.0;
394  forAllConstIter(typename KinematicCloud<CloudType>, *this, iter)
395  {
396  const parcelType& p = iter();
397  scalar m = p.nParticle()*p.mass();
398  scalar d = mag(p.position() - p.position0());
399  mSum += m;
400 
401  mass[i] = m;
402  dist[i] = d;
403 
404  i++;
405  }
406 
407  // calculate total mass across all processors
408  reduce(mSum, sumOp<scalar>());
409  Pstream::gatherList(procMass);
410  Pstream::gatherList(procDist);
411 
412  if (Pstream::master())
413  {
414  // flatten the mass lists
415  List<scalar> allMass(nParcelSum, 0.0);
416  SortableList<scalar> allDist(nParcelSum, 0.0);
417  for (label procI = 0; procI < Pstream::nProcs(); procI++)
418  {
420  (
421  allMass,
422  globalParcels.localSize(procI),
423  globalParcels.offset(procI)
424  ).assign(procMass[procI]);
425 
426  // flatten the distance list
428  (
429  allDist,
430  globalParcels.localSize(procI),
431  globalParcels.offset(procI)
432  ).assign(procDist[procI]);
433  }
434 
435  // sort allDist distances into ascending order
436  // note: allMass masses are left unsorted
437  allDist.sort();
438 
439  if (nParcelSum > 1)
440  {
441  const scalar mLimit = fraction*mSum;
442  const labelList& indices = allDist.indices();
443 
444  if (mLimit > (mSum - allMass[indices.last()]))
445  {
446  distance = allDist.last();
447  }
448  else
449  {
450  // assuming that 'fraction' is generally closer to 1 than 0,
451  // loop through in reverse distance order
452  const scalar mThreshold = (1.0 - fraction)*mSum;
453  scalar mCurrent = 0.0;
454  label i0 = 0;
455 
456  forAllReverse(indices, i)
457  {
458  label indI = indices[i];
459 
460  mCurrent += allMass[indI];
461 
462  if (mCurrent > mThreshold)
463  {
464  i0 = i;
465  break;
466  }
467  }
468 
469  if (i0 == indices.size() - 1)
470  {
471  distance = allDist.last();
472  }
473  else
474  {
475  // linearly interpolate to determine distance
476  scalar alpha = (mCurrent - mThreshold)/allMass[indices[i0]];
477  distance =
478  allDist[i0] + alpha*(allDist[i0+1] - allDist[i0]);
479  }
480  }
481  }
482  else
483  {
484  distance = allDist.first();
485  }
486  }
487 
488  Pstream::scatter(distance);
489 
490  return distance;
491 }
492 
493 
494 template<class CloudType>
496 {
497  return rndGen_;
498 }
499 
500 
501 template<class CloudType>
504 {
505  if (cellOccupancyPtr_.empty())
506  {
507  buildCellOccupancy();
508  }
509 
510  return cellOccupancyPtr_();
511 }
512 
513 
514 template<class CloudType>
515 inline const Foam::scalarField&
517 {
518  return cellLengthScale_;
519 }
520 
521 
522 template<class CloudType>
525 {
526  return UTrans_();
527 }
528 
529 
530 template<class CloudType>
533 {
534  return UTrans_();
535 }
536 
537 
538 template<class CloudType>
541 {
542  return UCoeff_();
543 }
544 
545 
546 template<class CloudType>
549 {
550  return UCoeff_();
551 }
552 
553 
554 template<class CloudType>
557 {
558  if (debug)
559  {
560  Info<< "UTrans min/max = " << min(UTrans()).value() << ", "
561  << max(UTrans()).value() << nl
562  << "UCoeff min/max = " << min(UCoeff()).value() << ", "
563  << max(UCoeff()).value() << endl;
564  }
565 
566  if (solution_.coupled())
567  {
568  if (solution_.semiImplicit("U"))
569  {
571  Vdt(mesh_.V()*this->db().time().deltaT());
572 
573  return UTrans()/Vdt - fvm::Sp(UCoeff()/Vdt, U) + UCoeff()/Vdt*U;
574  }
575  else
576  {
578  fvVectorMatrix& fvm = tfvm();
579 
580  fvm.source() = -UTrans()/(this->db().time().deltaT());
581 
582  return tfvm;
583  }
584  }
585 
587 }
588 
589 
590 template<class CloudType>
593 {
594  tmp<volScalarField> tvDotSweep
595  (
596  new volScalarField
597  (
598  IOobject
599  (
600  this->name() + ":vDotSweep",
601  this->db().time().timeName(),
602  this->db(),
603  IOobject::NO_READ,
604  IOobject::NO_WRITE,
605  false
606  ),
607  mesh_,
608  dimensionedScalar("zero", dimless/dimTime, 0.0),
609  zeroGradientFvPatchScalarField::typeName
610  )
611  );
612 
613  volScalarField& vDotSweep = tvDotSweep();
614  forAllConstIter(typename KinematicCloud<CloudType>, *this, iter)
615  {
616  const parcelType& p = iter();
617  const label cellI = p.cell();
618 
619  vDotSweep[cellI] += p.nParticle()*p.areaP()*mag(p.U() - U_[cellI]);
620  }
621 
622  vDotSweep.internalField() /= mesh_.V();
623  vDotSweep.correctBoundaryConditions();
624 
625  return tvDotSweep;
626 }
627 
628 
629 template<class CloudType>
632 {
633  tmp<volScalarField> ttheta
634  (
635  new volScalarField
636  (
637  IOobject
638  (
639  this->name() + ":theta",
640  this->db().time().timeName(),
641  this->db(),
642  IOobject::NO_READ,
643  IOobject::NO_WRITE,
644  false
645  ),
646  mesh_,
647  dimensionedScalar("zero", dimless, 0.0),
648  zeroGradientFvPatchScalarField::typeName
649  )
650  );
651 
652  volScalarField& theta = ttheta();
653  forAllConstIter(typename KinematicCloud<CloudType>, *this, iter)
654  {
655  const parcelType& p = iter();
656  const label cellI = p.cell();
657 
658  theta[cellI] += p.nParticle()*p.volume();
659  }
660 
661  theta.internalField() /= mesh_.V();
663 
664  return ttheta;
665 }
666 
667 
668 template<class CloudType>
671 {
672  tmp<volScalarField> talpha
673  (
674  new volScalarField
675  (
676  IOobject
677  (
678  this->name() + ":alpha",
679  this->db().time().timeName(),
680  this->db(),
681  IOobject::NO_READ,
682  IOobject::NO_WRITE,
683  false
684  ),
685  mesh_,
686  dimensionedScalar("zero", dimless, 0.0)
687  )
688  );
689 
690  scalarField& alpha = talpha().internalField();
691  forAllConstIter(typename KinematicCloud<CloudType>, *this, iter)
692  {
693  const parcelType& p = iter();
694  const label cellI = p.cell();
695 
696  alpha[cellI] += p.nParticle()*p.mass();
697  }
698 
699  alpha /= (mesh_.V()*rho_);
700 
701  return talpha;
702 }
703 
704 
705 template<class CloudType>
708 {
709  tmp<volScalarField> trhoEff
710  (
711  new volScalarField
712  (
713  IOobject
714  (
715  this->name() + ":rhoEff",
716  this->db().time().timeName(),
717  this->db(),
718  IOobject::NO_READ,
719  IOobject::NO_WRITE,
720  false
721  ),
722  mesh_,
723  dimensionedScalar("zero", dimDensity, 0.0)
724  )
725  );
726 
727  scalarField& rhoEff = trhoEff().internalField();
728  forAllConstIter(typename KinematicCloud<CloudType>, *this, iter)
729  {
730  const parcelType& p = iter();
731  const label cellI = p.cell();
732 
733  rhoEff[cellI] += p.nParticle()*p.mass();
734  }
735 
736  rhoEff /= mesh_.V();
737 
738  return trhoEff;
739 }
740 
741 
742 // ************************************************************************* //
Foam::KinematicCloud::UCoeff
DimensionedField< scalar, volMesh > & UCoeff()
Return coefficient for carrier phase U equation.
Definition: KinematicCloudI.H:540
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::maxOp
Definition: ops.H:172
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
p
p
Definition: pEqn.H:62
Foam::KinematicCloud::functions
functionType & functions()
Optional cloud function objects.
Definition: KinematicCloudI.H:169
Foam::KinematicCloud::penetration
scalar penetration(const scalar fraction) const
Penetration for fraction [0-1] of the current total mass.
Definition: KinematicCloudI.H:360
Foam::SortableList::sort
void sort()
(stable) sort the list (if changed after construction time)
Definition: SortableList.C:95
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::DispersionModel
Definition: KinematicCloud.H:79
Foam::dimDensity
const dimensionSet dimDensity
Foam::KinematicCloud::constProps
const parcelType::constantProperties & constProps() const
Return the constant properties.
Definition: KinematicCloudI.H:86
Foam::KinematicCloud::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: KinematicCloudI.H:40
Foam::KinematicCloud::SU
tmp< fvVectorMatrix > SU(volVectorField &U) const
Return tmp momentum source term.
Definition: KinematicCloudI.H:556
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::KinematicCloud::mu
const volScalarField & mu() const
Return carrier gas dynamic viscosity.
Definition: KinematicCloudI.H:123
Foam::KinematicCloud::massInSystem
scalar massInSystem() const
Total mass in system.
Definition: KinematicCloudI.H:271
Foam::KinematicCloud::rhoEff
const tmp< volScalarField > rhoEff() const
Return the particle effective density field.
Definition: KinematicCloudI.H:707
Foam::KinematicCloud::forces
const forceType & forces() const
Optional particle forces.
Definition: KinematicCloudI.H:153
Foam::KinematicCloud::surfaceFilm
const SurfaceFilmModel< KinematicCloud< CloudType > > & surfaceFilm() const
Return const-access to the surface film model.
Definition: KinematicCloudI.H:241
Foam::KinematicCloud::cloudCopy
const KinematicCloud & cloudCopy() const
Return a reference to the cloud copy.
Definition: KinematicCloudI.H:33
Foam::globalIndex::localSize
label localSize() const
My local size.
Definition: globalIndexI.H:60
Foam::KinematicCloud::injectors
const InjectionModelList< KinematicCloud< CloudType > > & injectors() const
Return const access to the injection model.
Definition: KinematicCloudI.H:177
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::dimForce
const dimensionSet dimForce
Foam::KinematicCloud::nParcels
label nParcels() const
Total number of parcels.
Definition: KinematicCloudI.H:264
U
U
Definition: pEqn.H:46
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::KinematicCloud::vDotSweep
const tmp< volScalarField > vDotSweep() const
Volume swept rate of parcels per cell.
Definition: KinematicCloudI.H:592
SortableList.H
Foam::IntegrationScheme
Top level model for Integration schemes.
Definition: IntegrationScheme.H:51
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::KinematicCloud
Templated base class for kinematic cloud.
Definition: KinematicCloud.H:96
Foam::PatchInteractionModel
Templated patch interaction model class.
Definition: KinematicCloud.H:82
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
Foam::StochasticCollisionModel
Templated stochastic collision model class.
Definition: KinematicCloud.H:88
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::KinematicCloud::cellOccupancy
List< DynamicList< parcelType * > > & cellOccupancy()
Return the cell occupancy information for each.
Definition: KinematicCloudI.H:503
Foam::cachedRandom
Random number generator.
Definition: cachedRandom.H:63
Foam::SurfaceFilmModel
Templated wall surface film model class.
Definition: KinematicCloud.H:85
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::KinematicCloud::UTrans
DimensionedField< vector, volMesh > & UTrans()
Return reference to momentum source.
Definition: KinematicCloudI.H:524
Foam::cloudSolution
Stores all relevant solution info for cloud.
Definition: cloudSolution.H:51
Foam::UList::assign
void assign(const UList< T > &)
Assign elements to those from UList.
Definition: UList.C:37
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::globalIndex::offset
label offset(const label procI) const
Start of procI data.
Definition: globalIndexI.H:48
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:418
Foam::SortableList
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: List.H:65
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::KinematicCloud::UIntegrator
const vectorIntegrationScheme & UIntegrator() const
Return reference to velocity integration.
Definition: KinematicCloudI.H:257
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::KinematicCloud::stochasticCollision
const StochasticCollisionModel< KinematicCloud< CloudType > > & stochasticCollision() const
Return const-access to the stochastic collision model.
Definition: KinematicCloudI.H:225
Foam::KinematicCloud::particleProperties
const IOdictionary & particleProperties() const
Return particle properties dictionary.
Definition: KinematicCloudI.H:48
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Foam::KinematicCloud::solution
const cloudSolution & solution() const
Return const access to the solution properties.
Definition: KinematicCloudI.H:71
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::CloudFunctionObjectList< KinematicCloud< CloudType > >
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::List::setSize
void setSize(const label)
Reset size of List.
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::KinematicCloud::dispersion
const DispersionModel< KinematicCloud< CloudType > > & dispersion() const
Return const-access to the dispersion model.
Definition: KinematicCloudI.H:193
Foam::sumOp
Definition: ops.H:162
Foam::KinematicCloud::pAmbient
scalar pAmbient() const
Return const-access to the ambient pressure.
Definition: KinematicCloudI.H:137
Foam::globalIndex::size
label size() const
Global sum of localSizes.
Definition: globalIndexI.H:66
Foam::KinematicCloud::patchInteraction
const PatchInteractionModel< KinematicCloud< CloudType > > & patchInteraction() const
Return const-access to the patch interaction model.
Definition: KinematicCloudI.H:209
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::KinematicCloud::g
const dimensionedVector & g() const
Gravity.
Definition: KinematicCloudI.H:130
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:291
Foam::KinematicCloud::theta
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
Definition: KinematicCloudI.H:631
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::SortableList::indices
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
Definition: SortableList.H:90
Foam::KinematicCloud::Dij
scalar Dij(const label i, const label j) const
Mean diameter Dij.
Definition: KinematicCloudI.H:320
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
Foam::KinematicCloud::subModelProperties
const dictionary & subModelProperties() const
Return reference to the sub-models dictionary.
Definition: KinematicCloudI.H:102
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::KinematicCloud::U
const volVectorField & U() const
Return carrier gas velocity.
Definition: KinematicCloudI.H:116
Foam::KinematicCloud::linearKineticEnergyOfSystem
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
Definition: KinematicCloudI.H:303
Foam::fvc::Sp
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::KinematicCloud::linearMomentumOfSystem
vector linearMomentumOfSystem() const
Total linear momentum of the system.
Definition: KinematicCloudI.H:286
Foam::KinematicCloud::Dmax
scalar Dmax() const
Max diameter.
Definition: KinematicCloudI.H:343
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::KinematicCloud::outputProperties
const IOdictionary & outputProperties() const
Return output properties dictionary.
Definition: KinematicCloudI.H:56
Foam::KinematicCloud::rndGen
cachedRandom & rndGen()
Return reference to the random object.
Definition: KinematicCloudI.H:495
Foam::ParticleForceList< KinematicCloud< CloudType > >
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::InjectionModelList
List of injection models.
Definition: KinematicCloud.H:76
Foam::KinematicCloud::alpha
const tmp< volScalarField > alpha() const
Return the particle mass fraction field.
Definition: KinematicCloudI.H:670
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::KinematicCloud::rho
const volScalarField & rho() const
Return carrier gas density.
Definition: KinematicCloudI.H:109
Foam::KinematicCloud::cellLengthScale
const scalarField & cellLengthScale() const
Return the cell length scale.
Definition: KinematicCloudI.H:516