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-2019 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();
108 
109  evolveCloud(cloud, td);
110 
111  if (solution_.coupled())
112  {
113  cloud.relaxSources(cloud.cloudCopy());
114  }
115  }
116  else
117  {
118  cloud.preEvolve();
119 
120  evolveCloud(cloud, td);
121 
122  if (solution_.coupled())
123  {
124  cloud.scaleSources();
125  }
126  }
127 
128  cloud.info();
129 
130  cloud.postEvolve();
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  Info<< endl;
236 
237  if (debug)
238  {
239  this->writePositions();
240  }
241 
242  this->dispersion().cacheFields(false);
243 
244  forces_.cacheFields(false);
245 
246  functions_.postEvolve();
247 
248  solution_.nextIter();
249 
250  if (this->db().time().writeTime())
251  {
252  outputProperties_.writeObject
253  (
254  IOstream::ASCII,
255  IOstream::currentVersion,
256  this->db().time().writeCompression(),
257  true
258  );
259  }
260 }
261 
262 
263 template<class CloudType>
265 {
266  CloudType::cloudReset(c);
267 
268  rndGen_ = c.rndGen_;
269 
270  forces_.transfer(c.forces_);
271 
272  functions_.transfer(c.functions_);
273 
274  injectors_.transfer(c.injectors_);
275 
276  dispersionModel_.reset(c.dispersionModel_.ptr());
277  patchInteractionModel_.reset(c.patchInteractionModel_.ptr());
278  stochasticCollisionModel_.reset(c.stochasticCollisionModel_.ptr());
279  surfaceFilmModel_.reset(c.surfaceFilmModel_.ptr());
280 
281  UIntegrator_.reset(c.UIntegrator_.ptr());
282 }
283 
284 
285 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
286 
287 template<class CloudType>
289 (
290  const word& cloudName,
291  const volScalarField& rho,
292  const volVectorField& U,
293  const volScalarField& mu,
294  const dimensionedVector& g,
295  bool readFields
296 )
297 :
298  CloudType(rho.mesh(), cloudName, false),
299  kinematicCloud(),
300  cloudCopyPtr_(nullptr),
301  mesh_(rho.mesh()),
302  particleProperties_
303  (
304  IOobject
305  (
306  cloudName + "Properties",
307  rho.mesh().time().constant(),
308  rho.mesh(),
309  IOobject::MUST_READ_IF_MODIFIED,
310  IOobject::NO_WRITE
311  )
312  ),
313  outputProperties_
314  (
315  IOobject
316  (
317  cloudName + "OutputProperties",
318  mesh_.time().timeName(),
319  "uniform"/cloud::prefix/cloudName,
320  mesh_,
321  IOobject::READ_IF_PRESENT,
322  IOobject::NO_WRITE
323  )
324  ),
325  solution_(mesh_, particleProperties_.subDict("solution")),
326  constProps_(particleProperties_),
327  subModelProperties_
328  (
329  particleProperties_.subOrEmptyDict
330  (
331  "subModels",
332  keyType::REGEX,
333  solution_.active()
334  )
335  ),
336  rndGen_(Pstream::myProcNo()),
337  cellOccupancyPtr_(),
338  cellLengthScale_(mag(cbrt(mesh_.V()))),
339  rho_(rho),
340  U_(U),
341  mu_(mu),
342  g_(g),
343  pAmbient_(0.0),
344  forces_
345  (
346  *this,
347  mesh_,
348  subModelProperties_.subOrEmptyDict
349  (
350  "particleForces",
351  keyType::REGEX,
352  solution_.active()
353  ),
354  solution_.active()
355  ),
356  functions_
357  (
358  *this,
359  particleProperties_.subOrEmptyDict("cloudFunctions"),
360  solution_.active()
361  ),
362  injectors_
363  (
364  subModelProperties_.subOrEmptyDict("injectionModels"),
365  *this
366  ),
367  dispersionModel_(nullptr),
368  patchInteractionModel_(nullptr),
369  stochasticCollisionModel_(nullptr),
370  surfaceFilmModel_(nullptr),
371  UIntegrator_(nullptr),
372  UTrans_
373  (
375  (
376  IOobject
377  (
378  this->name() + ":UTrans",
379  this->db().time().timeName(),
380  this->db(),
381  IOobject::READ_IF_PRESENT,
382  IOobject::AUTO_WRITE
383  ),
384  mesh_,
386  )
387  ),
388  UCoeff_
389  (
391  (
392  IOobject
393  (
394  this->name() + ":UCoeff",
395  this->db().time().timeName(),
396  this->db(),
397  IOobject::READ_IF_PRESENT,
398  IOobject::AUTO_WRITE
399  ),
400  mesh_,
402  )
403  )
404 {
405  if (solution_.active())
406  {
407  setModels();
408 
409  if (readFields)
410  {
411  parcelType::readFields(*this);
412  this->deleteLostParticles();
413  }
414  }
415 
416  if (solution_.resetSourcesOnStartup())
417  {
418  resetSourceTerms();
419  }
420 }
421 
422 
423 template<class CloudType>
425 (
427  const word& name
428 )
429 :
430  CloudType(c.mesh_, name, c),
431  kinematicCloud(),
432  cloudCopyPtr_(nullptr),
433  mesh_(c.mesh_),
434  particleProperties_(c.particleProperties_),
435  outputProperties_(c.outputProperties_),
436  solution_(c.solution_),
437  constProps_(c.constProps_),
438  subModelProperties_(c.subModelProperties_),
439  rndGen_(c.rndGen_, true),
440  cellOccupancyPtr_(nullptr),
441  cellLengthScale_(c.cellLengthScale_),
442  rho_(c.rho_),
443  U_(c.U_),
444  mu_(c.mu_),
445  g_(c.g_),
446  pAmbient_(c.pAmbient_),
447  forces_(c.forces_),
448  functions_(c.functions_),
449  injectors_(c.injectors_),
450  dispersionModel_(c.dispersionModel_->clone()),
451  patchInteractionModel_(c.patchInteractionModel_->clone()),
452  stochasticCollisionModel_(c.stochasticCollisionModel_->clone()),
453  surfaceFilmModel_(c.surfaceFilmModel_->clone()),
454  UIntegrator_(c.UIntegrator_->clone()),
455  UTrans_
456  (
458  (
459  IOobject
460  (
461  this->name() + ":UTrans",
462  this->db().time().timeName(),
463  this->db(),
464  IOobject::NO_READ,
465  IOobject::NO_WRITE,
466  false
467  ),
468  c.UTrans_()
469  )
470  ),
471  UCoeff_
472  (
474  (
475  IOobject
476  (
477  name + ":UCoeff",
478  this->db().time().timeName(),
479  this->db(),
480  IOobject::NO_READ,
481  IOobject::NO_WRITE,
482  false
483  ),
484  c.UCoeff_()
485  )
486  )
487 {}
488 
489 
490 template<class CloudType>
492 (
493  const fvMesh& mesh,
494  const word& name,
496 )
497 :
499  kinematicCloud(),
500  cloudCopyPtr_(nullptr),
501  mesh_(mesh),
502  particleProperties_
503  (
504  IOobject
505  (
506  name + "Properties",
507  mesh.time().constant(),
508  mesh,
509  IOobject::NO_READ,
510  IOobject::NO_WRITE,
511  false
512  )
513  ),
514  outputProperties_
515  (
516  IOobject
517  (
518  name + "OutputProperties",
519  mesh_.time().timeName(),
520  "uniform"/cloud::prefix/name,
521  mesh_,
522  IOobject::NO_READ,
523  IOobject::NO_WRITE,
524  false
525  )
526  ),
527  solution_(mesh),
528  constProps_(),
529  subModelProperties_(dictionary::null),
530  rndGen_(),
531  cellOccupancyPtr_(nullptr),
532  cellLengthScale_(c.cellLengthScale_),
533  rho_(c.rho_),
534  U_(c.U_),
535  mu_(c.mu_),
536  g_(c.g_),
537  pAmbient_(c.pAmbient_),
538  forces_(*this, mesh),
539  functions_(*this),
540  injectors_(*this),
541  dispersionModel_(nullptr),
542  patchInteractionModel_(nullptr),
543  stochasticCollisionModel_(nullptr),
544  surfaceFilmModel_(nullptr),
545  UIntegrator_(nullptr),
546  UTrans_(nullptr),
547  UCoeff_(nullptr)
548 {}
549 
550 
551 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
552 
553 template<class CloudType>
555 (
556  parcelType& parcel,
557  const scalar lagrangianDt
558 )
559 {
560  parcel.rho() = constProps_.rho0();
561 }
562 
563 
564 template<class CloudType>
566 (
567  parcelType& parcel,
568  const scalar lagrangianDt,
569  const bool fullyDescribed
570 )
571 {
572  const scalar carrierDt = mesh_.time().deltaTValue();
573  parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
574 
575  if (parcel.typeId() == -1)
576  {
577  parcel.typeId() = constProps_.parcelTypeId();
578  }
579 }
580 
581 
582 template<class CloudType>
584 {
585  cloudCopyPtr_.reset
586  (
587  static_cast<KinematicCloud<CloudType>*>
588  (
589  clone(this->name() + "Copy").ptr()
590  )
591  );
592 }
593 
594 
595 template<class CloudType>
597 {
598  cloudReset(cloudCopyPtr_());
599  cloudCopyPtr_.clear();
600 }
601 
602 
603 template<class CloudType>
605 {
606  UTrans().field() = Zero;
607  UCoeff().field() = 0.0;
608 }
609 
610 
611 template<class CloudType>
612 template<class Type>
614 (
616  const DimensionedField<Type, volMesh>& field0,
617  const word& name
618 ) const
619 {
620  const scalar coeff = solution_.relaxCoeff(name);
621  field = field0 + coeff*(field - field0);
622 }
623 
624 
625 template<class CloudType>
626 template<class Type>
628 (
630  const word& name
631 ) const
632 {
633  const scalar coeff = solution_.relaxCoeff(name);
634  field *= coeff;
635 }
636 
637 
638 template<class CloudType>
640 (
641  const KinematicCloud<CloudType>& cloudOldTime
642 )
643 {
644  this->relax(UTrans_(), cloudOldTime.UTrans(), "U");
645  this->relax(UCoeff_(), cloudOldTime.UCoeff(), "U");
646 }
647 
648 
649 template<class CloudType>
651 {
652  this->scale(UTrans_(), "U");
653  this->scale(UCoeff_(), "U");
654 }
655 
656 
657 template<class CloudType>
659 {
660  // force calculation of mesh dimensions - needed for parallel runs
661  // with topology change due to lazy evaluation of valid mesh dimensions
662  label nGeometricD = mesh_.nGeometricD();
663 
664  Info<< "\nSolving " << nGeometricD << "-D cloud " << this->name() << endl;
665 
666  this->dispersion().cacheFields(true);
667  forces_.cacheFields(true);
668  updateCellOccupancy();
669 
670  pAmbient_ = constProps_.dict().template
671  lookupOrDefault<scalar>("pAmbient", pAmbient_);
672 
673  functions_.preEvolve();
674 }
675 
676 
677 template<class CloudType>
679 {
680  if (solution_.canEvolve())
681  {
682  typename parcelType::trackingData td(*this);
683 
684  solve(*this, td);
685  }
686 }
687 
688 
689 template<class CloudType>
690 template<class TrackCloudType>
692 (
693  TrackCloudType& cloud,
694  typename parcelType::trackingData& td
695 )
696 {
697  td.part() = parcelType::trackingData::tpLinearTrack;
698  CloudType::move(cloud, td, solution_.trackTime());
699 
700  updateCellOccupancy();
701 }
702 
703 
704 template<class CloudType>
706 (
707  const parcelType& p,
708  const polyPatch& pp,
709  vector& nw,
710  vector& Up
711 ) const
712 {
713  p.patchData(nw, Up);
714 
715  // If this is a wall patch, then there may be a non-zero tangential velocity
716  // component; the lid velocity in a lid-driven cavity case, for example. We
717  // want the particle to interact with this velocity, so we look it up in the
718  // velocity field and use it to set the wall-tangential component.
719  if (isA<wallPolyPatch>(pp))
720  {
721  const label patchi = pp.index();
722  const label patchFacei = pp.whichFace(p.face());
723 
724  // We only want to use the boundary condition value only if it is set
725  // by the boundary condition. If the boundary values are extrapolated
726  // (e.g., slip conditions) then they represent the motion of the fluid
727  // just inside the domain rather than that of the wall itself.
728  if (U_.boundaryField()[patchi].fixesValue())
729  {
730  const vector Uw1 = U_.boundaryField()[patchi][patchFacei];
731  const vector& Uw0 =
732  U_.oldTime().boundaryField()[patchi][patchFacei];
733 
734  const scalar f = p.currentTimeFraction();
735 
736  const vector Uw = Uw0 + f*(Uw1 - Uw0);
737 
738  const tensor nnw = nw*nw;
739 
740  Up = (nnw & Up) + Uw - (nnw & Uw);
741  }
742  }
743 }
744 
745 
746 template<class CloudType>
748 {
749  updateCellOccupancy();
750  injectors_.updateMesh();
751  cellLengthScale_ = mag(cbrt(mesh_.V()));
752 }
753 
754 
755 template<class CloudType>
757 {
759 
760  updateMesh();
761 }
762 
763 
764 template<class CloudType>
766 {
767  const vector linearMomentum =
768  returnReduce(linearMomentumOfSystem(), sumOp<vector>());
769 
770  const scalar linearKineticEnergy =
771  returnReduce(linearKineticEnergyOfSystem(), sumOp<scalar>());
772 
773  const label nTotParcel = returnReduce(this->size(), sumOp<label>());
774 
775  const scalar particlePerParcel =
776  (
777  nTotParcel
778  ? (returnReduce(totalParticlePerParcel(), sumOp<scalar>()) / nTotParcel)
779  : 0
780  );
781 
782  Info<< "Cloud: " << this->name() << nl
783  << " Current number of parcels = " << nTotParcel << nl
784  << " Current mass in system = "
785  << returnReduce(massInSystem(), sumOp<scalar>()) << nl
786  << " Linear momentum = " << linearMomentum << nl
787  << " |Linear momentum| = " << mag(linearMomentum) << nl
788  << " Linear kinetic energy = " << linearKineticEnergy << nl
789  << " Average particle per parcel = " << particlePerParcel << nl;
790 
791  injectors_.info(Info);
792  this->surfaceFilm().info(Info);
793  this->patchInteraction().info(Info);
794 }
795 
796 
797 template<class CloudType>
799 {
800  parcelType::readObjects(*this, obr);
801 }
802 
803 
804 template<class CloudType>
806 {
807  parcelType::writeObjects(*this, obr);
808 }
809 
810 
811 // ************************************************************************* //
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:706
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:650
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:805
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:692
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.
Definition: zero.H:128
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:57
Foam::KinematicCloud::relaxSources
void relaxSources(const KinematicCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
Definition: KinematicCloud.C:640
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:583
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::KinematicCloud::UTrans
volVectorField::Internal & UTrans()
Return reference to momentum source.
Definition: KinematicCloudI.H:392
interpolation.H
rho
rho
Definition: readInitialConditions.H:96
Foam::IOobject::info
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:519
solve
CEqn solve()
Foam::sumOp
Definition: ops.H:213
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::KinematicCloud::readObjects
virtual void readObjects(const objectRegistry &obr)
Read particle fields from objects in the obr registry.
Definition: KinematicCloud.C:798
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.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:678
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
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:62
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:566
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:66
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
subCycleTime.H
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
field
rDeltaTY field()
Foam::CloudType
DSMCCloud< dsmcParcel > CloudType
Definition: makeDSMCParcelBinaryCollisionModels.C:38
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
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:756
Foam::KinematicCloud::setParcelThermoProperties
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
Definition: KinematicCloud.C:555
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:604
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:398
Foam::nl
constexpr char nl
Definition: Ostream.H:372
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::KinematicCloud::cloudReset
void cloudReset(KinematicCloud< CloudType > &c)
Reset state of cloud.
Definition: KinematicCloud.C:264
Foam::KinematicCloud::postEvolve
void postEvolve()
Post-evolve.
Definition: KinematicCloud.C:233
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:747
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:614
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:596
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:765
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
Foam::patchIdentifier::index
label index() const
The index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:133
Foam::KinematicCloud::preEvolve
void preEvolve()
Pre-evolve.
Definition: KinematicCloud.C:658
Foam::KinematicCloud::scale
void scale(DimensionedField< Type, volMesh > &field, const word &name) const
Scale field.
Definition: KinematicCloud.C:628
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()