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-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "KinematicCloud.H"
30#include "integrationScheme.H"
31#include "interpolation.H"
32#include "subCycleTime.H"
33
34#include "InjectionModelList.H"
35#include "DispersionModel.H"
38#include "SurfaceFilmModel.H"
39#include "profiling.H"
40
41#include "PackingModel.H"
42#include "ParticleStressModel.H"
43#include "DampingModel.H"
44#include "IsotropyModel.H"
45#include "TimeScaleModel.H"
46
47// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
48
49template<class CloudType>
51{
52 dispersionModel_.reset
53 (
55 (
56 subModelProperties_,
57 *this
58 ).ptr()
59 );
60
61 patchInteractionModel_.reset
62 (
64 (
65 subModelProperties_,
66 *this
67 ).ptr()
68 );
69
70 stochasticCollisionModel_.reset
71 (
73 (
74 subModelProperties_,
75 *this
76 ).ptr()
77 );
78
79 surfaceFilmModel_.reset
80 (
82 (
83 subModelProperties_,
84 *this
85 ).ptr()
86 );
87
88 packingModel_.reset
89 (
91 (
92 subModelProperties_,
93 *this
94 ).ptr()
95 );
96
97 dampingModel_.reset
98 (
100 (
101 subModelProperties_,
102 *this
103 ).ptr()
104 );
105
106 isotropyModel_.reset
107 (
109 (
110 subModelProperties_,
111 *this
112 ).ptr()
113 );
114
115 UIntegrator_.reset
116 (
117 integrationScheme::New
118 (
119 "U",
120 solution_.integrationSchemes()
121 ).ptr()
122 );
123}
124
125
126template<class CloudType>
127template<class TrackCloudType>
129(
130 TrackCloudType& cloud,
131 typename parcelType::trackingData& td
132)
133{
134 addProfiling(prof, "cloud::solve");
135
136 if (solution_.steadyState())
137 {
138 cloud.storeState();
139
140 cloud.preEvolve(td);
141
142 evolveCloud(cloud, td);
143
144 if (solution_.coupled())
145 {
146 cloud.relaxSources(cloud.cloudCopy());
147 }
148 }
149 else
150 {
151 cloud.preEvolve(td);
152
153 evolveCloud(cloud, td);
154
155 if (solution_.coupled())
156 {
157 cloud.scaleSources();
158 }
159 }
160
161 cloud.info();
162
163 cloud.postEvolve(td);
164
165 if (solution_.steadyState())
166 {
167 cloud.restoreState();
168 }
169}
170
171
172template<class CloudType>
174{
175 if (!cellOccupancyPtr_)
176 {
177 cellOccupancyPtr_.reset
178 (
179 new List<DynamicList<parcelType*>>(mesh_.nCells())
180 );
181 }
182 else if (cellOccupancyPtr_().size() != mesh_.nCells())
183 {
184 // If the size of the mesh has changed, reset the
185 // cellOccupancy size
186
187 cellOccupancyPtr_().setSize(mesh_.nCells());
188 }
189
190 List<DynamicList<parcelType*>>& cellOccupancy = cellOccupancyPtr_();
191
192 for (auto& list : cellOccupancy)
193 {
194 list.clear();
195 }
196
197 for (parcelType& p : *this)
198 {
199 cellOccupancy[p.cell()].append(&p);
200 }
201}
202
203
204template<class CloudType>
206{
207 // Only build the cellOccupancy if the pointer is set, i.e. it has
208 // been requested before.
209
210 if (cellOccupancyPtr_)
211 {
212 buildCellOccupancy();
213 }
214}
215
216
217template<class CloudType>
218template<class TrackCloudType>
220(
221 TrackCloudType& cloud,
222 typename parcelType::trackingData& td
223)
224{
225 if (solution_.coupled())
226 {
227 cloud.resetSourceTerms();
228 }
229
230 if (solution_.transient())
231 {
232 label preInjectionSize = this->size();
233
234 this->surfaceFilm().inject(cloud);
235
236 // Update the cellOccupancy if the size of the cloud has changed
237 // during the injection.
238 if (preInjectionSize != this->size())
239 {
240 updateCellOccupancy();
241 preInjectionSize = this->size();
242 }
243
244 injectors_.inject(cloud, td);
245
246 // Assume that motion will update the cellOccupancy as necessary
247 // before it is required.
248 cloud.motion(cloud, td);
249
250 stochasticCollision().update(td, solution_.trackTime());
251 }
252 else
253 {
254// this->surfaceFilm().injectSteadyState(cloud);
255
256 injectors_.injectSteadyState(cloud, td, solution_.trackTime());
257
258 td.part() = parcelType::trackingData::tpLinearTrack;
259 CloudType::move(cloud, td, solution_.trackTime());
260 }
261}
262
263
264template<class CloudType>
266(
267 const typename parcelType::trackingData& td
268)
269{
270 Info<< endl;
271
272 if (debug)
273 {
274 this->writePositions();
275 }
276
277 this->dispersion().cacheFields(false);
278
279 this->patchInteraction().postEvolve();
280
281 forces_.cacheFields(false);
282
283 functions_.postEvolve(td);
284
285 solution_.nextIter();
286
287 if (this->db().time().writeTime())
288 {
289 outputProperties_.writeObject
290 (
292 (
293 IOstream::ASCII,
294 this->db().time().writeCompression()
295 ),
296 true
297 );
298 }
300 if (this->dampingModel().active())
301 {
302 this->dampingModel().cacheFields(false);
303 }
304 if (this->packingModel().active())
305 {
306 this->packingModel().cacheFields(false);
308}
309
310
311template<class CloudType>
313{
314 CloudType::cloudReset(c);
315
316 rndGen_ = c.rndGen_;
317
318 forces_.transfer(c.forces_);
319
320 functions_.transfer(c.functions_);
321
322 injectors_.transfer(c.injectors_);
323
324 dispersionModel_.reset(c.dispersionModel_.ptr());
325 patchInteractionModel_.reset(c.patchInteractionModel_.ptr());
326 stochasticCollisionModel_.reset(c.stochasticCollisionModel_.ptr());
327 surfaceFilmModel_.reset(c.surfaceFilmModel_.ptr());
328
329 packingModel_.reset(c.packingModel_.ptr());
330 dampingModel_.reset(c.dampingModel_.ptr());
331 isotropyModel_.reset(c.isotropyModel_.ptr());
332
333 UIntegrator_.reset(c.UIntegrator_.ptr());
334}
335
336
337// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
338
339template<class CloudType>
341(
342 const word& cloudName,
343 const volScalarField& rho,
344 const volVectorField& U,
345 const volScalarField& mu,
346 const dimensionedVector& g,
347 bool readFields
348)
349:
350 CloudType(rho.mesh(), cloudName, false),
352 cloudCopyPtr_(nullptr),
353 mesh_(rho.mesh()),
354 particleProperties_
355 (
357 (
358 cloudName + "Properties",
359 mesh_.time().constant(),
360 mesh_,
361 IOobject::MUST_READ_IF_MODIFIED,
362 IOobject::NO_WRITE
363 )
364 ),
365 outputProperties_
366 (
368 (
369 cloudName + "OutputProperties",
370 mesh_.time().timeName(),
371 "uniform"/cloud::prefix/cloudName,
372 mesh_,
373 IOobject::READ_IF_PRESENT,
374 IOobject::NO_WRITE
375 )
376 ),
377 solution_(mesh_, particleProperties_.subDict("solution")),
378 constProps_(particleProperties_),
379 subModelProperties_
380 (
381 particleProperties_.subOrEmptyDict
382 (
383 "subModels",
384 keyType::REGEX,
385 solution_.active()
386 )
387 ),
388 rndGen_(Pstream::myProcNo()),
389 cellOccupancyPtr_(),
390 cellLengthScale_(mag(cbrt(mesh_.V()))),
391 rho_(rho),
392 U_(U),
393 mu_(mu),
394 g_(g),
395 pAmbient_(0.0),
396 forces_
397 (
398 *this,
399 mesh_,
400 subModelProperties_.subOrEmptyDict
401 (
402 "particleForces",
403 keyType::REGEX,
404 solution_.active()
405 ),
406 solution_.active()
407 ),
408 functions_
409 (
410 *this,
411 particleProperties_.subOrEmptyDict("cloudFunctions"),
412 solution_.active()
413 ),
414 injectors_
415 (
416 subModelProperties_.subOrEmptyDict("injectionModels"),
417 *this
418 ),
419 dispersionModel_(nullptr),
420 patchInteractionModel_(nullptr),
421 stochasticCollisionModel_(nullptr),
422 surfaceFilmModel_(nullptr),
423
424 packingModel_(nullptr),
425 dampingModel_(nullptr),
426 isotropyModel_(nullptr),
427
428 UIntegrator_(nullptr),
429 UTrans_
430 (
431 new volVectorField::Internal
432 (
434 (
435 this->name() + ":UTrans",
436 this->db().time().timeName(),
437 this->db(),
438 IOobject::READ_IF_PRESENT,
439 IOobject::AUTO_WRITE
440 ),
441 mesh_,
443 )
444 ),
445 UCoeff_
446 (
447 new volScalarField::Internal
448 (
450 (
451 this->name() + ":UCoeff",
452 this->db().time().timeName(),
453 this->db(),
454 IOobject::READ_IF_PRESENT,
455 IOobject::AUTO_WRITE
456 ),
457 mesh_,
459 )
460 )
461{
462 if (solution_.active())
463 {
464 setModels();
465
466 if (readFields)
467 {
469 this->deleteLostParticles();
470 }
471 }
472
474 {
476 }
477}
478
479
480template<class CloudType>
482(
484 const word& name
485)
486:
487 CloudType(c.mesh_, name, c),
489 cloudCopyPtr_(nullptr),
490 mesh_(c.mesh_),
491 particleProperties_(c.particleProperties_),
492 outputProperties_(c.outputProperties_),
493 solution_(c.solution_),
494 constProps_(c.constProps_),
495 subModelProperties_(c.subModelProperties_),
496 rndGen_(c.rndGen_, true),
497 cellOccupancyPtr_(nullptr),
498 cellLengthScale_(c.cellLengthScale_),
499 rho_(c.rho_),
500 U_(c.U_),
501 mu_(c.mu_),
502 g_(c.g_),
503 pAmbient_(c.pAmbient_),
504 forces_(c.forces_),
505 functions_(c.functions_),
506 injectors_(c.injectors_),
507 dispersionModel_(c.dispersionModel_->clone()),
508 patchInteractionModel_(c.patchInteractionModel_->clone()),
509 stochasticCollisionModel_(c.stochasticCollisionModel_->clone()),
510 surfaceFilmModel_(c.surfaceFilmModel_->clone()),
511
512 packingModel_(c.packingModel_->clone()),
513 dampingModel_(c.dampingModel_->clone()),
514 isotropyModel_(c.isotropyModel_->clone()),
515
516 UIntegrator_(c.UIntegrator_->clone()),
517 UTrans_
518 (
519 new volVectorField::Internal
520 (
522 (
523 this->name() + ":UTrans",
524 this->db().time().timeName(),
525 this->db(),
526 IOobject::NO_READ,
527 IOobject::NO_WRITE,
528 false
529 ),
530 c.UTrans_()
531 )
532 ),
533 UCoeff_
534 (
535 new volScalarField::Internal
536 (
538 (
539 name + ":UCoeff",
540 this->db().time().timeName(),
541 this->db(),
542 IOobject::NO_READ,
543 IOobject::NO_WRITE,
544 false
545 ),
546 c.UCoeff_()
547 )
548 )
549{}
550
551
552template<class CloudType>
554(
555 const fvMesh& mesh,
556 const word& name,
558)
559:
562 cloudCopyPtr_(nullptr),
563 mesh_(mesh),
564 particleProperties_
565 (
567 (
568 name + "Properties",
569 mesh_.time().constant(),
570 mesh_,
571 IOobject::NO_READ,
572 IOobject::NO_WRITE,
573 false
574 )
575 ),
576 outputProperties_
577 (
579 (
580 name + "OutputProperties",
581 mesh_.time().timeName(),
582 "uniform"/cloud::prefix/name,
583 mesh_,
584 IOobject::NO_READ,
585 IOobject::NO_WRITE,
586 false
587 )
588 ),
589 solution_(mesh),
590 constProps_(),
591 subModelProperties_(),
592 rndGen_(),
593 cellOccupancyPtr_(nullptr),
594 cellLengthScale_(c.cellLengthScale_),
595 rho_(c.rho_),
596 U_(c.U_),
597 mu_(c.mu_),
598 g_(c.g_),
599 pAmbient_(c.pAmbient_),
600 forces_(*this, mesh),
601 functions_(*this),
602 injectors_(*this),
603 dispersionModel_(nullptr),
604 patchInteractionModel_(nullptr),
605 stochasticCollisionModel_(nullptr),
606 surfaceFilmModel_(nullptr),
607
608 packingModel_(nullptr),
609 dampingModel_(nullptr),
610 isotropyModel_(nullptr),
611
612 UIntegrator_(nullptr),
613 UTrans_(nullptr),
614 UCoeff_(nullptr)
615{}
617
618// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
619
620template<class CloudType>
622(
623 parcelType& parcel,
624 const scalar lagrangianDt
625)
626{
627 // If rho0 is given in the const properties
628 if (constProps_.rho0() != -1)
630 parcel.rho() = constProps_.rho0();
631 }
632}
633
634
635template<class CloudType>
637(
638 parcelType& parcel,
639 const scalar lagrangianDt,
640 const bool fullyDescribed
641)
643 const scalar carrierDt = mesh_.time().deltaTValue();
644 parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
645
646 if (parcel.typeId() == -1)
647 {
648 parcel.typeId() = constProps_.parcelTypeId();
649 }
650
651 if (parcel.rho() == -1)
654 << "The kinematic cloud needs rho0 in the constantProperties "
655 << " dictionary. " << nl
656 << abort(FatalError);
657 }
658}
659
661template<class CloudType>
663{
664 cloudCopyPtr_.reset
665 (
666 static_cast<KinematicCloud<CloudType>*>
667 (
668 clone(this->name() + "Copy").ptr()
669 )
670 );
671}
673
674template<class CloudType>
677 cloudReset(cloudCopyPtr_());
678 cloudCopyPtr_.clear();
679}
680
681
682template<class CloudType>
684{
685 UTrans().field() = Zero;
686 UCoeff().field() = 0.0;
687}
689
690template<class CloudType>
691template<class Type>
693(
696 const word& name
697) const
698{
699 const scalar coeff = solution_.relaxCoeff(name);
700 field = field0 + coeff*(field - field0);
701}
702
703
704template<class CloudType>
705template<class Type>
707(
709 const word& name
710) const
711{
712 const scalar coeff = solution_.relaxCoeff(name);
713 field *= coeff;
714}
715
716
717template<class CloudType>
719(
720 const KinematicCloud<CloudType>& cloudOldTime
721)
722{
723 this->relax(UTrans_(), cloudOldTime.UTrans(), "U");
724 this->relax(UCoeff_(), cloudOldTime.UCoeff(), "U");
725}
726
727
728template<class CloudType>
730{
731 this->scale(UTrans_(), "U");
732 this->scale(UCoeff_(), "U");
733}
734
735
736template<class CloudType>
738(
739 const typename parcelType::trackingData& td
740)
741{
742 // force calculation of mesh dimensions - needed for parallel runs
743 // with topology change due to lazy evaluation of valid mesh dimensions
744 label nGeometricD = mesh_.nGeometricD();
745
746 Info<< "\nSolving" << nGeometricD << "-D cloud " << this->name() << endl;
747
748 this->dispersion().cacheFields(true);
749 forces_.cacheFields(true);
750
751 pAmbient_ = constProps_.dict().template
752 getOrDefault<scalar>("pAmbient", pAmbient_);
753
754 if (this->dampingModel().active() || this->packingModel().active())
755 {
756 const_cast<typename parcelType::trackingData&>(td).updateAverages(*this);
757 }
758
759 if (this->dampingModel().active())
760 {
761 this->dampingModel().cacheFields(true);
762 }
763 if (this->packingModel().active())
764 {
765 this->packingModel().cacheFields(true);
766 }
767
768 updateCellOccupancy();
769
770 functions_.preEvolve(td);
771}
772
773
774template<class CloudType>
776{
777 if (solution_.canEvolve())
778 {
779 typename parcelType::trackingData td(*this);
780 solve(*this, td);
781 }
782}
783
784
785template<class CloudType>
786template<class TrackCloudType>
788(
789 TrackCloudType& cloud,
790 typename parcelType::trackingData& td
791)
792{
793 td.part() = parcelType::trackingData::tpLinearTrack;
794 CloudType::move(cloud, td, solution_.trackTime());
795
796 if (isotropyModel_->active())
797 {
798 td.updateAverages(cloud);
799 isotropyModel_->calculate();
800 }
801
802 updateCellOccupancy();
803}
804
805
806template<class CloudType>
808(
809 const parcelType& p,
810 const polyPatch& pp,
811 vector& nw,
812 vector& Up
813) const
814{
815 p.patchData(nw, Up);
816
817 // If this is a wall patch, then there may be a non-zero tangential velocity
818 // component; the lid velocity in a lid-driven cavity case, for example. We
819 // want the particle to interact with this velocity, so we look it up in the
820 // velocity field and use it to set the wall-tangential component.
821 if (isA<wallPolyPatch>(pp))
822 {
823 const label patchi = pp.index();
824 const label patchFacei = pp.whichFace(p.face());
825
826 // We only want to use the boundary condition value only if it is set
827 // by the boundary condition. If the boundary values are extrapolated
828 // (e.g., slip conditions) then they represent the motion of the fluid
829 // just inside the domain rather than that of the wall itself.
830 if (U_.boundaryField()[patchi].fixesValue())
831 {
832 const vector Uw1(U_.boundaryField()[patchi][patchFacei]);
833 const vector& Uw0 =
834 U_.oldTime().boundaryField()[patchi][patchFacei];
835
836 const scalar f = p.currentTimeFraction();
837
838 const vector Uw(Uw0 + f*(Uw1 - Uw0));
839
840 const tensor nnw(nw*nw);
841
842 Up = (nnw & Up) + Uw - (nnw & Uw);
843 }
844 }
845}
846
847
848template<class CloudType>
850{
851 updateCellOccupancy();
852 injectors_.updateMesh();
853 cellLengthScale_ = mag(cbrt(mesh_.V()));
854}
855
856
857template<class CloudType>
859{
861
862 updateMesh();
863}
864
865
866template<class CloudType>
868{
869 const vector linearMomentum =
870 returnReduce(linearMomentumOfSystem(), sumOp<vector>());
871
872 const scalar linearKineticEnergy =
873 returnReduce(linearKineticEnergyOfSystem(), sumOp<scalar>());
874
875 const label nTotParcel = returnReduce(this->size(), sumOp<label>());
876
877 const scalar particlePerParcel =
878 (
879 nTotParcel
880 ? (returnReduce(totalParticlePerParcel(), sumOp<scalar>()) / nTotParcel)
881 : 0
882 );
883
884 Info<< "Cloud: " << this->name() << nl
885 << " Current number of parcels = " << nTotParcel << nl
886 << " Current mass in system = "
887 << returnReduce(massInSystem(), sumOp<scalar>()) << nl
888 << " Linear momentum = " << linearMomentum << nl
889 << " |Linear momentum| = " << mag(linearMomentum) << nl
890 << " Linear kinetic energy = " << linearKineticEnergy << nl
891 << " Average particle per parcel = " << particlePerParcel << nl;
892
893 injectors_.info(Info);
894 this->surfaceFilm().info(Info);
895 this->patchInteraction().info(Info);
896
897 if (this->packingModel().active())
898 {
899 tmp<volScalarField> alpha = this->theta();
900
901 if (this->db().time().writeTime())
902 {
903 alpha().write();
904 }
905
906 const scalar alphaMin = gMin(alpha().primitiveField());
907 const scalar alphaMax = gMax(alpha().primitiveField());
908
909 Info<< " Min cell volume fraction = " << alphaMin << endl;
910 Info<< " Max cell volume fraction = " << alphaMax << endl;
911
912 if (alphaMax < SMALL)
913 {
914 return;
915 }
916
917 scalar nMin = GREAT;
918
919 forAll(this->mesh().cells(), celli)
920 {
921 const label n = this->cellOccupancy()[celli].size();
922
923 if (n > 0)
924 {
925 const scalar nPack = n*alphaMax/alpha()[celli];
926
927 if (nPack < nMin)
928 {
929 nMin = nPack;
930 }
931 }
932 }
933
934 reduce(nMin, minOp<scalar>());
935
936 Info<< " Min dense number of parcels = " << nMin << endl;
937 }
938}
939
940
941template<class CloudType>
943{
944 parcelType::readObjects(*this, obr);
945}
946
947
948template<class CloudType>
950{
951 parcelType::writeObjects(*this, obr);
952}
953
954
955// ************************************************************************* //
label n
const List< DynamicList< molecule * > > & cellOccupancy
const uniformDimensionedVectorField & g
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition: Cloud.C:118
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: Cloud.C:308
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
Base class for collisional damping models.
Definition: DampingModel.H:66
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Base class for dispersion modelling.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
Template class for intrusive linked lists.
Definition: ILList.H:69
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
InfoProxy< IOobject > info() const
Return info proxy, for printing information to a stream.
Definition: IOobject.H:690
The IOstreamOption is a simple container for options an IOstream can normally have.
Base class for collisional return-to-isotropy models.
Definition: IsotropyModel.H:66
Templated base class for kinematic cloud.
cloudSolution solution_
Solution properties.
void setModels()
Set cloud sub-models.
void relaxSources(const KinematicCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
void storeState()
Store the current cloud state.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
virtual void readObjects(const objectRegistry &obr)
Read particle fields from objects in the obr registry.
void cloudReset(KinematicCloud< CloudType > &c)
Reset state of cloud.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
void scaleSources()
Apply scaling to (transient) cloud sources.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
void updateCellOccupancy()
Update (i.e. build) the cellOccupancy if it has.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
void evolve()
Evolve the cloud.
volScalarField::Internal & UCoeff()
Return coefficient for carrier phase U equation.
void evolveCloud(TrackCloudType &cloud, typename parcelType::trackingData &td)
Evolve the cloud.
volVectorField::Internal & UTrans()
Return reference to momentum source.
void preEvolve(const typename parcelType::trackingData &td)
Pre-evolve.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
void resetSourceTerms()
Reset the cloud source terms.
void buildCellOccupancy()
Build the cellOccupancy.
void updateMesh()
Update mesh.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
virtual void postEvolve()
Post-evolve hook.
scalar scale
Overall scale factor.
Definition: PDRparams.H:140
Base class for packing models.
Definition: PackingModel.H:71
Templated patch interaction model class.
Inter-processor communications stream.
Definition: Pstream.H:63
Templated stochastic collision model class.
Templated wall surface film model class.
const Switch resetSourcesOnStartup() const
Return const access to the reset sources flag.
const Switch active() const
Return the active flag.
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
const motionSolver & motion() const
Return the motionSolver.
virtual void move()=0
Class used to pass data into container.
void relax()
Relax matrix (for steady-state solution).
Definition: faMatrix.C:581
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
Allows specification of different writing frequency of objects registered to the database.
Definition: writeObjects.H:142
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A class for handling keywords in dictionaries.
Definition: keyType.H:71
Virtual abstract base class for templated KinematicCloud.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Registry of regIOobjects.
const FieldField< Field, Type > & patchData() const
label index() const noexcept
The index of this patch in the boundaryMesh.
constant condensation/saturation model.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:451
virtual bool write(const bool valid=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
volScalarField & p
UEqn relax()
regionModels::surfaceFilmModel & surfaceFilm
rDeltaTY field()
const volScalarField & mu
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const cellShapeList & cells
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
dimensionedScalar cbrt(const dimensionedScalar &ds)
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
Type gMax(const FieldField< Field, Type > &f)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
labelList f(nPoints)
volScalarField & alpha
CEqn solve()
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
const word cloudName(propsDict.get< word >("cloud"))