KinematicCloud.C
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) 2016-2020 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 "KinematicCloud.H"
30 #include "integrationScheme.H"
31 #include "interpolation.H"
32 #include "subCycleTime.H"
33 
34 #include "InjectionModelList.H"
35 #include "DispersionModel.H"
36 #include "PatchInteractionModel.H"
38 #include "SurfaceFilmModel.H"
39 #include "profiling.H"
40 
41 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
42 
43 template<class CloudType>
45 {
46  dispersionModel_.reset
47  (
49  (
50  subModelProperties_,
51  *this
52  ).ptr()
53  );
54 
55  patchInteractionModel_.reset
56  (
58  (
59  subModelProperties_,
60  *this
61  ).ptr()
62  );
63 
64  stochasticCollisionModel_.reset
65  (
67  (
68  subModelProperties_,
69  *this
70  ).ptr()
71  );
72 
73  surfaceFilmModel_.reset
74  (
76  (
77  subModelProperties_,
78  *this
79  ).ptr()
80  );
81 
82  UIntegrator_.reset
83  (
85  (
86  "U",
87  solution_.integrationSchemes()
88  ).ptr()
89  );
90 }
91 
92 
93 template<class CloudType>
94 template<class TrackCloudType>
96 (
97  TrackCloudType& cloud,
98  typename parcelType::trackingData& td
99 )
100 {
101  addProfiling(prof, "cloud::solve");
102 
103  if (solution_.steadyState())
104  {
105  cloud.storeState();
106 
107  cloud.preEvolve(td);
108 
109  evolveCloud(cloud, td);
110 
111  if (solution_.coupled())
112  {
113  cloud.relaxSources(cloud.cloudCopy());
114  }
115  }
116  else
117  {
118  cloud.preEvolve(td);
119 
120  evolveCloud(cloud, td);
121 
122  if (solution_.coupled())
123  {
124  cloud.scaleSources();
125  }
126  }
127 
128  cloud.info();
129 
130  cloud.postEvolve(td);
131 
132  if (solution_.steadyState())
133  {
134  cloud.restoreState();
135  }
136 }
137 
138 
139 template<class CloudType>
141 {
142  if (cellOccupancyPtr_.empty())
143  {
144  cellOccupancyPtr_.reset
145  (
146  new List<DynamicList<parcelType*>>(mesh_.nCells())
147  );
148  }
149  else if (cellOccupancyPtr_().size() != mesh_.nCells())
150  {
151  // If the size of the mesh has changed, reset the
152  // cellOccupancy size
153 
154  cellOccupancyPtr_().setSize(mesh_.nCells());
155  }
156 
157  List<DynamicList<parcelType*>>& cellOccupancy = cellOccupancyPtr_();
158 
159  for (auto& list : cellOccupancy)
160  {
161  list.clear();
162  }
163 
164  for (parcelType& p : *this)
165  {
166  cellOccupancy[p.cell()].append(&p);
167  }
168 }
169 
170 
171 template<class CloudType>
173 {
174  // Only build the cellOccupancy if the pointer is set, i.e. it has
175  // been requested before.
176 
177  if (cellOccupancyPtr_.valid())
178  {
179  buildCellOccupancy();
180  }
181 }
182 
183 
184 template<class CloudType>
185 template<class TrackCloudType>
187 (
188  TrackCloudType& cloud,
189  typename parcelType::trackingData& td
190 )
191 {
192  if (solution_.coupled())
193  {
194  cloud.resetSourceTerms();
195  }
196 
197  if (solution_.transient())
198  {
199  label preInjectionSize = this->size();
200 
201  this->surfaceFilm().inject(cloud);
202 
203  // Update the cellOccupancy if the size of the cloud has changed
204  // during the injection.
205  if (preInjectionSize != this->size())
206  {
207  updateCellOccupancy();
208  preInjectionSize = this->size();
209  }
210 
211  injectors_.inject(cloud, td);
212 
213 
214  // Assume that motion will update the cellOccupancy as necessary
215  // before it is required.
216  cloud.motion(cloud, td);
217 
218  stochasticCollision().update(td, solution_.trackTime());
219  }
220  else
221  {
222 // this->surfaceFilm().injectSteadyState(cloud);
223 
224  injectors_.injectSteadyState(cloud, td, solution_.trackTime());
225 
226  td.part() = parcelType::trackingData::tpLinearTrack;
227  CloudType::move(cloud, td, solution_.trackTime());
228  }
229 }
230 
231 
232 template<class CloudType>
234 (
235  const typename parcelType::trackingData& td
236 )
237 {
238  Info<< endl;
239 
240  if (debug)
241  {
242  this->writePositions();
243  }
244 
245  this->dispersion().cacheFields(false);
246 
247  forces_.cacheFields(false);
248 
249  functions_.postEvolve(td);
250 
251  solution_.nextIter();
252 
253  if (this->db().time().writeTime())
254  {
255  outputProperties_.writeObject
256  (
258  (
259  IOstream::ASCII,
260  this->db().time().writeCompression()
261  ),
262  true
263  );
264  }
265 }
266 
267 
268 template<class CloudType>
270 {
271  CloudType::cloudReset(c);
272 
273  rndGen_ = c.rndGen_;
274 
275  forces_.transfer(c.forces_);
276 
277  functions_.transfer(c.functions_);
278 
279  injectors_.transfer(c.injectors_);
280 
281  dispersionModel_.reset(c.dispersionModel_.ptr());
282  patchInteractionModel_.reset(c.patchInteractionModel_.ptr());
283  stochasticCollisionModel_.reset(c.stochasticCollisionModel_.ptr());
284  surfaceFilmModel_.reset(c.surfaceFilmModel_.ptr());
285 
286  UIntegrator_.reset(c.UIntegrator_.ptr());
287 }
288 
289 
290 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
291 
292 template<class CloudType>
294 (
295  const word& cloudName,
296  const volScalarField& rho,
297  const volVectorField& U,
298  const volScalarField& mu,
299  const dimensionedVector& g,
300  bool readFields
301 )
302 :
303  CloudType(rho.mesh(), cloudName, false),
304  kinematicCloud(),
305  cloudCopyPtr_(nullptr),
306  mesh_(rho.mesh()),
307  particleProperties_
308  (
309  IOobject
310  (
311  cloudName + "Properties",
312  mesh_.time().constant(),
313  mesh_,
314  IOobject::MUST_READ_IF_MODIFIED,
315  IOobject::NO_WRITE
316  )
317  ),
318  outputProperties_
319  (
320  IOobject
321  (
322  cloudName + "OutputProperties",
323  mesh_.time().timeName(),
324  "uniform"/cloud::prefix/cloudName,
325  mesh_,
326  IOobject::READ_IF_PRESENT,
327  IOobject::NO_WRITE
328  )
329  ),
330  solution_(mesh_, particleProperties_.subDict("solution")),
331  constProps_(particleProperties_),
332  subModelProperties_
333  (
334  particleProperties_.subOrEmptyDict
335  (
336  "subModels",
337  keyType::REGEX,
338  solution_.active()
339  )
340  ),
341  rndGen_(Pstream::myProcNo()),
342  cellOccupancyPtr_(),
343  cellLengthScale_(mag(cbrt(mesh_.V()))),
344  rho_(rho),
345  U_(U),
346  mu_(mu),
347  g_(g),
348  pAmbient_(0.0),
349  forces_
350  (
351  *this,
352  mesh_,
353  subModelProperties_.subOrEmptyDict
354  (
355  "particleForces",
356  keyType::REGEX,
357  solution_.active()
358  ),
359  solution_.active()
360  ),
361  functions_
362  (
363  *this,
364  particleProperties_.subOrEmptyDict("cloudFunctions"),
365  solution_.active()
366  ),
367  injectors_
368  (
369  subModelProperties_.subOrEmptyDict("injectionModels"),
370  *this
371  ),
372  dispersionModel_(nullptr),
373  patchInteractionModel_(nullptr),
374  stochasticCollisionModel_(nullptr),
375  surfaceFilmModel_(nullptr),
376  UIntegrator_(nullptr),
377  UTrans_
378  (
379  new volVectorField::Internal
380  (
381  IOobject
382  (
383  this->name() + ":UTrans",
384  this->db().time().timeName(),
385  this->db(),
386  IOobject::READ_IF_PRESENT,
387  IOobject::AUTO_WRITE
388  ),
389  mesh_,
391  )
392  ),
393  UCoeff_
394  (
395  new volScalarField::Internal
396  (
397  IOobject
398  (
399  this->name() + ":UCoeff",
400  this->db().time().timeName(),
401  this->db(),
402  IOobject::READ_IF_PRESENT,
403  IOobject::AUTO_WRITE
404  ),
405  mesh_,
407  )
408  )
409 {
410  if (solution_.active())
411  {
412  setModels();
413 
414  if (readFields)
415  {
416  parcelType::readFields(*this);
417  this->deleteLostParticles();
418  }
419  }
420 
421  if (solution_.resetSourcesOnStartup())
422  {
423  resetSourceTerms();
424  }
425 }
426 
427 
428 template<class CloudType>
430 (
432  const word& name
433 )
434 :
435  CloudType(c.mesh_, name, c),
436  kinematicCloud(),
437  cloudCopyPtr_(nullptr),
438  mesh_(c.mesh_),
439  particleProperties_(c.particleProperties_),
440  outputProperties_(c.outputProperties_),
441  solution_(c.solution_),
442  constProps_(c.constProps_),
443  subModelProperties_(c.subModelProperties_),
444  rndGen_(c.rndGen_, true),
445  cellOccupancyPtr_(nullptr),
446  cellLengthScale_(c.cellLengthScale_),
447  rho_(c.rho_),
448  U_(c.U_),
449  mu_(c.mu_),
450  g_(c.g_),
451  pAmbient_(c.pAmbient_),
452  forces_(c.forces_),
453  functions_(c.functions_),
454  injectors_(c.injectors_),
455  dispersionModel_(c.dispersionModel_->clone()),
456  patchInteractionModel_(c.patchInteractionModel_->clone()),
457  stochasticCollisionModel_(c.stochasticCollisionModel_->clone()),
458  surfaceFilmModel_(c.surfaceFilmModel_->clone()),
459  UIntegrator_(c.UIntegrator_->clone()),
460  UTrans_
461  (
463  (
464  IOobject
465  (
466  this->name() + ":UTrans",
467  this->db().time().timeName(),
468  this->db(),
469  IOobject::NO_READ,
470  IOobject::NO_WRITE,
471  false
472  ),
473  c.UTrans_()
474  )
475  ),
476  UCoeff_
477  (
479  (
480  IOobject
481  (
482  name + ":UCoeff",
483  this->db().time().timeName(),
484  this->db(),
485  IOobject::NO_READ,
486  IOobject::NO_WRITE,
487  false
488  ),
489  c.UCoeff_()
490  )
491  )
492 {}
493 
494 
495 template<class CloudType>
497 (
498  const fvMesh& mesh,
499  const word& name,
501 )
502 :
504  kinematicCloud(),
505  cloudCopyPtr_(nullptr),
506  mesh_(mesh),
507  particleProperties_
508  (
509  IOobject
510  (
511  name + "Properties",
512  mesh_.time().constant(),
513  mesh_,
514  IOobject::NO_READ,
515  IOobject::NO_WRITE,
516  false
517  )
518  ),
519  outputProperties_
520  (
521  IOobject
522  (
523  name + "OutputProperties",
524  mesh_.time().timeName(),
525  "uniform"/cloud::prefix/name,
526  mesh_,
527  IOobject::NO_READ,
528  IOobject::NO_WRITE,
529  false
530  )
531  ),
532  solution_(mesh),
533  constProps_(),
534  subModelProperties_(dictionary::null),
535  rndGen_(),
536  cellOccupancyPtr_(nullptr),
537  cellLengthScale_(c.cellLengthScale_),
538  rho_(c.rho_),
539  U_(c.U_),
540  mu_(c.mu_),
541  g_(c.g_),
542  pAmbient_(c.pAmbient_),
543  forces_(*this, mesh),
544  functions_(*this),
545  injectors_(*this),
546  dispersionModel_(nullptr),
547  patchInteractionModel_(nullptr),
548  stochasticCollisionModel_(nullptr),
549  surfaceFilmModel_(nullptr),
550  UIntegrator_(nullptr),
551  UTrans_(nullptr),
552  UCoeff_(nullptr)
553 {}
554 
555 
556 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
557 
558 template<class CloudType>
560 (
561  parcelType& parcel,
562  const scalar lagrangianDt
563 )
564 {
565  parcel.rho() = constProps_.rho0();
566 }
567 
568 
569 template<class CloudType>
571 (
572  parcelType& parcel,
573  const scalar lagrangianDt,
574  const bool fullyDescribed
575 )
576 {
577  const scalar carrierDt = mesh_.time().deltaTValue();
578  parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
579 
580  if (parcel.typeId() == -1)
581  {
582  parcel.typeId() = constProps_.parcelTypeId();
583  }
584 }
585 
586 
587 template<class CloudType>
589 {
590  cloudCopyPtr_.reset
591  (
592  static_cast<KinematicCloud<CloudType>*>
593  (
594  clone(this->name() + "Copy").ptr()
595  )
596  );
597 }
598 
599 
600 template<class CloudType>
602 {
603  cloudReset(cloudCopyPtr_());
604  cloudCopyPtr_.clear();
605 }
606 
607 
608 template<class CloudType>
610 {
611  UTrans().field() = Zero;
612  UCoeff().field() = 0.0;
613 }
614 
615 
616 template<class CloudType>
617 template<class Type>
619 (
621  const DimensionedField<Type, volMesh>& field0,
622  const word& name
623 ) const
624 {
625  const scalar coeff = solution_.relaxCoeff(name);
626  field = field0 + coeff*(field - field0);
627 }
628 
629 
630 template<class CloudType>
631 template<class Type>
633 (
635  const word& name
636 ) const
637 {
638  const scalar coeff = solution_.relaxCoeff(name);
639  field *= coeff;
640 }
641 
642 
643 template<class CloudType>
645 (
646  const KinematicCloud<CloudType>& cloudOldTime
647 )
648 {
649  this->relax(UTrans_(), cloudOldTime.UTrans(), "U");
650  this->relax(UCoeff_(), cloudOldTime.UCoeff(), "U");
651 }
652 
653 
654 template<class CloudType>
656 {
657  this->scale(UTrans_(), "U");
658  this->scale(UCoeff_(), "U");
659 }
660 
661 
662 template<class CloudType>
664 (
665  const typename parcelType::trackingData& td
666 )
667 {
668  // force calculation of mesh dimensions - needed for parallel runs
669  // with topology change due to lazy evaluation of valid mesh dimensions
670  label nGeometricD = mesh_.nGeometricD();
671 
672  Info<< "\nSolving " << nGeometricD << "-D cloud " << this->name() << endl;
673 
674  this->dispersion().cacheFields(true);
675  forces_.cacheFields(true);
676  updateCellOccupancy();
677 
678  pAmbient_ = constProps_.dict().template
679  getOrDefault<scalar>("pAmbient", pAmbient_);
680 
681  functions_.preEvolve(td);
682 }
683 
684 
685 template<class CloudType>
687 {
688  if (solution_.canEvolve())
689  {
690  typename parcelType::trackingData td(*this);
691 
692  solve(*this, td);
693  }
694 }
695 
696 
697 template<class CloudType>
698 template<class TrackCloudType>
700 (
701  TrackCloudType& cloud,
702  typename parcelType::trackingData& td
703 )
704 {
705  td.part() = parcelType::trackingData::tpLinearTrack;
706  CloudType::move(cloud, td, solution_.trackTime());
707 
708  updateCellOccupancy();
709 }
710 
711 
712 template<class CloudType>
714 (
715  const parcelType& p,
716  const polyPatch& pp,
717  vector& nw,
718  vector& Up
719 ) const
720 {
721  p.patchData(nw, Up);
722 
723  // If this is a wall patch, then there may be a non-zero tangential velocity
724  // component; the lid velocity in a lid-driven cavity case, for example. We
725  // want the particle to interact with this velocity, so we look it up in the
726  // velocity field and use it to set the wall-tangential component.
727  if (isA<wallPolyPatch>(pp))
728  {
729  const label patchi = pp.index();
730  const label patchFacei = pp.whichFace(p.face());
731 
732  // We only want to use the boundary condition value only if it is set
733  // by the boundary condition. If the boundary values are extrapolated
734  // (e.g., slip conditions) then they represent the motion of the fluid
735  // just inside the domain rather than that of the wall itself.
736  if (U_.boundaryField()[patchi].fixesValue())
737  {
738  const vector Uw1 = U_.boundaryField()[patchi][patchFacei];
739  const vector& Uw0 =
740  U_.oldTime().boundaryField()[patchi][patchFacei];
741 
742  const scalar f = p.currentTimeFraction();
743 
744  const vector Uw = Uw0 + f*(Uw1 - Uw0);
745 
746  const tensor nnw = nw*nw;
747 
748  Up = (nnw & Up) + Uw - (nnw & Uw);
749  }
750  }
751 }
752 
753 
754 template<class CloudType>
756 {
757  updateCellOccupancy();
758  injectors_.updateMesh();
759  cellLengthScale_ = mag(cbrt(mesh_.V()));
760 }
761 
762 
763 template<class CloudType>
765 {
767 
768  updateMesh();
769 }
770 
771 
772 template<class CloudType>
774 {
775  const vector linearMomentum =
776  returnReduce(linearMomentumOfSystem(), sumOp<vector>());
777 
778  const scalar linearKineticEnergy =
779  returnReduce(linearKineticEnergyOfSystem(), sumOp<scalar>());
780 
781  const label nTotParcel = returnReduce(this->size(), sumOp<label>());
782 
783  const scalar particlePerParcel =
784  (
785  nTotParcel
786  ? (returnReduce(totalParticlePerParcel(), sumOp<scalar>()) / nTotParcel)
787  : 0
788  );
789 
790  Info<< "Cloud: " << this->name() << nl
791  << " Current number of parcels = " << nTotParcel << nl
792  << " Current mass in system = "
793  << returnReduce(massInSystem(), sumOp<scalar>()) << nl
794  << " Linear momentum = " << linearMomentum << nl
795  << " |Linear momentum| = " << mag(linearMomentum) << nl
796  << " Linear kinetic energy = " << linearKineticEnergy << nl
797  << " Average particle per parcel = " << particlePerParcel << nl;
798 
799  injectors_.info(Info);
800  this->surfaceFilm().info(Info);
801  this->patchInteraction().info(Info);
802 }
803 
804 
805 template<class CloudType>
807 {
808  parcelType::readObjects(*this, obr);
809 }
810 
811 
812 template<class CloudType>
814 {
815  parcelType::writeObjects(*this, obr);
816 }
817 
818 
819 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::KinematicCloud::patchData
void patchData(const parcelType &p, const polyPatch &pp, vector &normal, vector &Up) const
Calculate the patch normal and velocity to interact with,.
Definition: KinematicCloud.C:714
profiling.H
PatchInteractionModel.H
Foam::Tensor< scalar >
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::KinematicCloud::scaleSources
void scaleSources()
Apply scaling to (transient) cloud sources.
Definition: KinematicCloud.C:655
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::KinematicCloud::writeObjects
virtual void writeObjects(objectRegistry &obr) const
Write particle fields as objects into the obr registry.
Definition: KinematicCloud.C:813
cloudName
const word cloudName(propsDict.get< word >("cloud"))
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::KinematicCloud::setModels
void setModels()
Set cloud sub-models.
Definition: KinematicCloud.C:44
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::KinematicCloud::motion
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
Definition: KinematicCloud.C:700
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Foam::DispersionModel
Base class for dispersion modelling.
Definition: KinematicCloud.H:86
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
Foam::KinematicCloud::relaxSources
void relaxSources(const KinematicCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
Definition: KinematicCloud.C:645
Foam::KinematicCloud::UCoeff
volScalarField::Internal & UCoeff()
Return coefficient for carrier phase U equation.
Definition: KinematicCloudI.H:408
Foam::KinematicCloud::storeState
void storeState()
Store the current cloud state.
Definition: KinematicCloud.C:588
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::KinematicCloud::UTrans
volVectorField::Internal & UTrans()
Return reference to momentum source.
Definition: KinematicCloudI.H:392
interpolation.H
rho
rho
Definition: readInitialConditions.H:88
Foam::IOobject::info
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:532
solve
CEqn solve()
Foam::sumOp
Definition: ops.H:213
Foam::KinematicCloud::readObjects
virtual void readObjects(const objectRegistry &obr)
Read particle fields from objects in the obr registry.
Definition: KinematicCloud.C:806
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::keyType
A class for handling keywords in dictionaries.
Definition: keyType.H:60
Foam::kinematicCloud
Virtual abstract base class for templated KinematicCloud.
Definition: kinematicCloud.H:52
Foam::KinematicCloud::evolve
void evolve()
Evolve the cloud.
Definition: KinematicCloud.C:686
Foam::KinematicCloud
Templated base class for kinematic cloud.
Definition: KinematicCloud.H:103
Foam::PatchInteractionModel
Templated patch interaction model class.
Definition: KinematicCloud.H:89
Foam::StochasticCollisionModel
Templated stochastic collision model class.
Definition: KinematicCloud.H:95
cellOccupancy
const List< DynamicList< molecule * > > & cellOccupancy
Definition: calculateMDFields.H:1
Foam::KinematicCloud::checkParcelProperties
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Definition: KinematicCloud.C:571
Foam::SurfaceFilmModel
Templated wall surface film model class.
Definition: KinematicCloud.H:92
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:67
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
surfaceFilm
regionModels::surfaceFilmModel & surfaceFilm
Definition: createFieldRefs.H:3
Foam::CollidingParcel
Wrapper around kinematic parcel types to add collision modelling.
Definition: CollidingParcel.H:61
Foam::IOstreamOption
The IOstreamOption is a simple container for options an IOstream can normally have.
Definition: IOstreamOption.H:63
subCycleTime.H
field
rDeltaTY field()
Foam::CloudType
DSMCCloud< dsmcParcel > CloudType
Definition: makeDSMCParcelBinaryCollisionModels.C:38
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::ILList
Template class for intrusive linked lists.
Definition: ILList.H:52
addProfiling
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
Definition: profilingTrigger.H:114
Foam::KinematicCloud::preEvolve
void preEvolve(const typename parcelType::trackingData &td)
Pre-evolve.
Definition: KinematicCloud.C:664
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
KinematicCloud.H
Foam::dimensioned< vector >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::KinematicCloud::autoMap
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: KinematicCloud.C:764
Foam::KinematicCloud::setParcelThermoProperties
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
Definition: KinematicCloud.C:560
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
Foam::KinematicCloud::resetSourceTerms
void resetSourceTerms()
Reset the cloud source terms.
Definition: KinematicCloud.C:609
DispersionModel.H
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
U
U
Definition: pEqn.H:72
Foam::KinematicCloud::solve
void solve(TrackCloudType &cloud, typename parcelType::trackingData &td)
Solve the cloud - calls all evolution functions.
Definition: KinematicCloud.C:96
Foam::polyPatch::whichFace
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:399
Foam::nl
constexpr char nl
Definition: Ostream.H:385
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::KinematicCloud::cloudReset
void cloudReset(KinematicCloud< CloudType > &c)
Reset state of cloud.
Definition: KinematicCloud.C:269
Foam::Pstream
Inter-processor communications stream.
Definition: Pstream.H:56
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:102
Foam::KinematicCloud::updateMesh
void updateMesh()
Update mesh.
Definition: KinematicCloud.C:755
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::CollidingParcel::trackingData
ParcelType::trackingData trackingData
Use base tracking data.
Definition: CollidingParcel.H:127
Foam::KinematicCloud::relax
void relax(DimensionedField< Type, volMesh > &field, const DimensionedField< Type, volMesh > &field0, const word &name) const
Relax field.
Definition: KinematicCloud.C:619
Foam::Cloud
Base cloud calls templated on particle type.
Definition: Cloud.H:54
StochasticCollisionModel.H
integrationScheme.H
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::KinematicCloud::restoreState
void restoreState()
Reset the current cloud to the previously stored state.
Definition: KinematicCloud.C:601
SurfaceFilmModel.H
Foam::KinematicCloud::updateCellOccupancy
void updateCellOccupancy()
Update (i.e. build) the cellOccupancy if it has.
Definition: KinematicCloud.C:172
Foam::KinematicCloud::info
void info()
Print cloud information.
Definition: KinematicCloud.C:773
Foam::KinematicCloud::buildCellOccupancy
void buildCellOccupancy()
Build the cellOccupancy.
Definition: KinematicCloud.C:140
InjectionModelList.H
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::KinematicCloud::evolveCloud
void evolveCloud(TrackCloudType &cloud, typename parcelType::trackingData &td)
Evolve the cloud.
Definition: KinematicCloud.C:187
constant
constant condensation/saturation model.
Foam::patchIdentifier::index
label index() const
The index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:158
Foam::KinematicCloud::scale
void scale(DimensionedField< Type, volMesh > &field, const word &name) const
Scale field.
Definition: KinematicCloud.C:633
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
relax
UEqn relax()
Foam::KinematicCloud::postEvolve
void postEvolve(const typename parcelType::trackingData &td)
Post-evolve.
Definition: KinematicCloud.C:234