Loading...
Searching...
No Matches
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-2023 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"
44#include "IsotropyModel.H"
45#include "TimeScaleModel.H"
46
47// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
48
49template<class CloudType>
51{
53 (
54 DispersionModel<KinematicCloud<CloudType>>::New
55 (
57 *this
58 ).ptr()
59 );
60
62 (
63 PatchInteractionModel<KinematicCloud<CloudType>>::New
64 (
66 *this
67 ).ptr()
68 );
69
71 (
72 StochasticCollisionModel<KinematicCloud<CloudType>>::New
73 (
75 *this
76 ).ptr()
77 );
78
80 (
81 SurfaceFilmModel<KinematicCloud<CloudType>>::New
82 (
84 *this
85 ).ptr()
86 );
87
88 packingModel_.reset
89 (
90 PackingModel<KinematicCloud<CloudType>>::New
91 (
93 *this
94 ).ptr()
95 );
96
97 dampingModel_.reset
98 (
99 DampingModel<KinematicCloud<CloudType>>::New
100 (
102 *this
103 ).ptr()
104 );
105
106 isotropyModel_.reset
107 (
108 IsotropyModel<KinematicCloud<CloudType>>::New
109 (
111 *this
112 ).ptr()
113 );
114
115 UIntegrator_.reset
116 (
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 log = solution_.log();
137
138 if (solution_.steadyState())
139 {
140 cloud.storeState();
141
142 cloud.preEvolve(td);
143
144 evolveCloud(cloud, td);
145
146 if (solution_.coupled())
147 {
148 cloud.relaxSources(cloud.cloudCopy());
149 }
150 }
151 else
152 {
153 cloud.preEvolve(td);
154
155 evolveCloud(cloud, td);
156
157 if (solution_.coupled())
158 {
159 cloud.scaleSources();
160 }
161 }
162
163 cloud.info();
164
165 cloud.postEvolve(td);
166
167 if (solution_.steadyState())
169 cloud.restoreState();
170 }
171}
172
173
174template<class CloudType>
176{
177 if (!cellOccupancyPtr_)
178 {
179 cellOccupancyPtr_.reset
180 (
181 new List<DynamicList<parcelType*>>(mesh_.nCells())
182 );
183 }
184 else if (cellOccupancyPtr_().size() != mesh_.nCells())
185 {
186 // If the size of the mesh has changed, reset the
187 // cellOccupancy size
188
189 cellOccupancyPtr_().setSize(mesh_.nCells());
190 }
191
192 List<DynamicList<parcelType*>>& cellOccupancy = cellOccupancyPtr_();
193
194 for (auto& list : cellOccupancy)
195 {
196 list.clear();
197 }
198
199 for (parcelType& p : *this)
201 cellOccupancy[p.cell()].append(&p);
202 }
203}
204
205
206template<class CloudType>
208{
209 // Only build the cellOccupancy if the pointer is set, i.e. it has
210 // been requested before.
211
212 if (cellOccupancyPtr_)
213 {
215 }
216}
217
218
219template<class CloudType>
220template<class TrackCloudType>
222(
223 TrackCloudType& cloud,
224 typename parcelType::trackingData& td
225)
226{
227 if (solution_.coupled())
228 {
229 cloud.resetSourceTerms();
230 }
231
232 if (solution_.transient())
233 {
234 label preInjectionSize = this->size();
235
236 this->surfaceFilm().inject(cloud);
237
238 // Update the cellOccupancy if the size of the cloud has changed
239 // during the injection.
240 if (preInjectionSize != this->size())
241 {
242 updateCellOccupancy();
243 preInjectionSize = this->size();
244 }
245
246 injectors_.inject(cloud, td);
247
248 // Assume that motion will update the cellOccupancy as necessary
249 // before it is required.
250 cloud.motion(cloud, td);
251
252 stochasticCollision().update(td, solution_.trackTime());
253 }
254 else
255 {
256// this->surfaceFilm().injectSteadyState(cloud);
257
258 injectors_.injectSteadyState(cloud, td, solution_.trackTime());
259
260 td.part() = parcelType::trackingData::tpLinearTrack;
261 CloudType::move(cloud, td, solution_.trackTime());
262 }
263}
264
265
266template<class CloudType>
268(
269 const typename parcelType::trackingData& td
270)
271{
272 Log_<< endl;
273
274 if (debug)
275 {
276 this->writePositions();
277 }
278
279 this->dispersion().cacheFields(false);
280
281 this->patchInteraction().postEvolve();
282
283 forces_.cacheFields(false);
284
285 functions_.postEvolve(td);
286
287 solution_.nextIter();
288
289 if (this->db().time().writeTime())
290 {
291 outputProperties_.writeObject
292 (
293 IOstreamOption
294 (
295 IOstreamOption::ASCII,
296 this->db().time().writeCompression()
297 ),
298 true
299 );
300 }
301
302 if (this->dampingModel().active())
303 {
304 this->dampingModel().cacheFields(false);
305 }
306 if (this->packingModel().active())
308 this->packingModel().cacheFields(false);
309 }
310}
311
312
313template<class CloudType>
314void Foam::KinematicCloud<CloudType>::cloudReset(KinematicCloud<CloudType>& c)
315{
316 CloudType::cloudReset(c);
317
318 rndGen_ = c.rndGen_;
319
320 forces_.transfer(c.forces_);
321
322 functions_.transfer(c.functions_);
323
324 injectors_.transfer(c.injectors_);
325
326 dispersionModel_.reset(c.dispersionModel_.ptr());
327 patchInteractionModel_.reset(c.patchInteractionModel_.ptr());
328 stochasticCollisionModel_.reset(c.stochasticCollisionModel_.ptr());
329 surfaceFilmModel_.reset(c.surfaceFilmModel_.ptr());
330
331 packingModel_.reset(c.packingModel_.ptr());
332 dampingModel_.reset(c.dampingModel_.ptr());
333 isotropyModel_.reset(c.isotropyModel_.ptr());
334
335 UIntegrator_.reset(c.UIntegrator_.ptr());
336}
337
338
339// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
341template<class CloudType>
342Foam::KinematicCloud<CloudType>::KinematicCloud
343(
344 const word& cloudName,
345 const volScalarField& rho,
346 const volVectorField& U,
347 const volScalarField& mu,
348 const dimensionedVector& g,
350)
351:
352 CloudType(rho.mesh(), cloudName, false),
354 cloudCopyPtr_(nullptr),
355 mesh_(rho.mesh()),
357 (
359 (
360 cloudName + "Properties",
361 mesh_.time().constant(),
362 mesh_,
366 )
367 ),
369 (
372 cloudName + "OutputProperties",
373 mesh_.time().timeName(),
374 "uniform"/cloud::prefix/cloudName,
375 mesh_,
379 )
380 ),
381 solution_(mesh_, particleProperties_.subDict("solution")),
384 (
385 particleProperties_.subOrEmptyDict
386 (
387 "subModels",
388 keyType::REGEX,
389 solution_.active()
390 )
391 ),
392 rndGen_(Pstream::myProcNo()),
395 rho_(rho),
396 U_(U),
397 mu_(mu),
398 g_(g),
399 pAmbient_(0.0),
400 forces_
401 (
402 *this,
403 mesh_,
404 subModelProperties_.subOrEmptyDict
405 (
406 "particleForces",
407 keyType::REGEX,
408 solution_.active()
409 ),
410 solution_.active()
411 ),
413 (
414 *this,
415 particleProperties_.subOrEmptyDict("cloudFunctions"),
416 solution_.active()
417 ),
419 (
420 subModelProperties_.subOrEmptyDict("injectionModels"),
421 *this
422 ),
423 dispersionModel_(nullptr),
424 patchInteractionModel_(nullptr),
426 surfaceFilmModel_(nullptr),
427
428 packingModel_(nullptr),
429 dampingModel_(nullptr),
430 isotropyModel_(nullptr),
431
432 UIntegrator_(nullptr),
434 (
435 new volScalarField::Internal
436 (
438 (
439 IOobject::scopedName(this->name(), "rhokTrans"),
440 this->db().time().timeName(),
441 this->db(),
445 ),
446 mesh_,
448 )
449 ),
450 UTrans_
451 (
452 new volVectorField::Internal
453 (
455 (
456 IOobject::scopedName(this->name(), "UTrans"),
457 this->db().time().timeName(),
458 this->db(),
462 ),
463 mesh_,
465 )
466 ),
467 UCoeff_
468 (
469 new volScalarField::Internal
470 (
472 (
473 IOobject::scopedName(this->name(), "UCoeff"),
474 this->db().time().timeName(),
475 this->db(),
479 ),
480 mesh_,
482 )
483 ),
484 log(true)
485{
486 if (solution_.active())
487 {
488 setModels();
489
490 if (readFields)
491 {
492 parcelType::readFields(*this);
493 this->deleteLostParticles();
494 }
495 }
496
497 if (solution_.resetSourcesOnStartup())
499 resetSourceTerms();
500 }
501}
502
503
504template<class CloudType>
505Foam::KinematicCloud<CloudType>::KinematicCloud
506(
508 const word& name
509)
510:
511 CloudType(c.mesh_, name, c),
513 cloudCopyPtr_(nullptr),
514 mesh_(c.mesh_),
515 particleProperties_(c.particleProperties_),
516 outputProperties_(c.outputProperties_),
517 solution_(c.solution_),
518 constProps_(c.constProps_),
519 subModelProperties_(c.subModelProperties_),
520 rndGen_(c.rndGen_, true),
521 cellOccupancyPtr_(nullptr),
522 cellLengthScale_(c.cellLengthScale_),
523 rho_(c.rho_),
524 U_(c.U_),
525 mu_(c.mu_),
526 g_(c.g_),
527 pAmbient_(c.pAmbient_),
528 forces_(c.forces_),
529 functions_(c.functions_),
530 injectors_(c.injectors_),
531 dispersionModel_(c.dispersionModel_->clone()),
532 patchInteractionModel_(c.patchInteractionModel_->clone()),
533 stochasticCollisionModel_(c.stochasticCollisionModel_->clone()),
534 surfaceFilmModel_(c.surfaceFilmModel_->clone()),
535
536 packingModel_(c.packingModel_->clone()),
537 dampingModel_(c.dampingModel_->clone()),
538 isotropyModel_(c.isotropyModel_->clone()),
539
540 UIntegrator_(c.UIntegrator_->clone()),
541 rhokTrans_
542 (
543 new volScalarField::Internal
544 (
546 (
547 IOobject::scopedName(this->name(), "rhokTrans"),
548 this->db().time().timeName(),
549 this->db(),
550 IOobject::NO_READ,
551 IOobject::NO_WRITE,
552 IOobject::NO_REGISTER
553 ),
554 c.rhokTrans_()
555 )
556 ),
557 UTrans_
558 (
559 new volVectorField::Internal
560 (
562 (
563 IOobject::scopedName(this->name(), "UTrans"),
564 this->db().time().timeName(),
565 this->db(),
566 IOobject::NO_READ,
567 IOobject::NO_WRITE,
568 IOobject::NO_REGISTER
569 ),
570 c.UTrans_()
571 )
572 ),
573 UCoeff_
574 (
575 new volScalarField::Internal
576 (
578 (
579 IOobject::scopedName(this->name(), "UCoeff"),
580 this->db().time().timeName(),
581 this->db(),
582 IOobject::NO_READ,
583 IOobject::NO_WRITE,
584 IOobject::NO_REGISTER
585 ),
586 c.UCoeff_()
588 ),
589 log(c.log)
590{}
591
592
593template<class CloudType>
594Foam::KinematicCloud<CloudType>::KinematicCloud
595(
596 const fvMesh& mesh,
597 const word& name,
599)
600:
601 CloudType(mesh, name, IDLList<parcelType>()),
603 cloudCopyPtr_(nullptr),
604 mesh_(mesh),
605 particleProperties_
606 (
608 (
609 name + "Properties",
610 mesh_.time().constant(),
611 mesh_,
612 IOobject::NO_READ,
613 IOobject::NO_WRITE,
614 IOobject::NO_REGISTER
615 )
616 ),
617 outputProperties_
618 (
620 (
621 name + "OutputProperties",
622 mesh_.time().timeName(),
623 "uniform"/cloud::prefix/name,
624 mesh_,
625 IOobject::NO_READ,
626 IOobject::NO_WRITE,
627 IOobject::NO_REGISTER
628 )
629 ),
630 solution_(mesh),
631 constProps_(),
632 subModelProperties_(),
633 rndGen_(),
634 cellOccupancyPtr_(nullptr),
635 cellLengthScale_(c.cellLengthScale_),
636 rho_(c.rho_),
637 U_(c.U_),
638 mu_(c.mu_),
639 g_(c.g_),
640 pAmbient_(c.pAmbient_),
641 forces_(*this, mesh),
642 functions_(*this),
643 injectors_(*this),
644 dispersionModel_(nullptr),
645 patchInteractionModel_(nullptr),
646 stochasticCollisionModel_(nullptr),
647 surfaceFilmModel_(nullptr),
648
649 packingModel_(nullptr),
650 dampingModel_(nullptr),
651 isotropyModel_(nullptr),
652
653 UIntegrator_(nullptr),
654 rhokTrans_(nullptr),
655 UTrans_(nullptr),
656 UCoeff_(nullptr),
658{}
659
660
661// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
662
663template<class CloudType>
665(
666 parcelType& parcel,
667 const scalar lagrangianDt
668)
669{
670 // If rho0 is given in the const properties
671 if (constProps_.rho0() != -1)
673 parcel.rho() = constProps_.rho0();
674 }
675}
676
677
678template<class CloudType>
680(
681 parcelType& parcel,
682 const scalar lagrangianDt,
683 const bool fullyDescribed
684)
685{
686 const scalar carrierDt = mesh_.time().deltaTValue();
687 parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
688
689 if (parcel.typeId() == -1)
690 {
691 parcel.typeId() = constProps_.parcelTypeId();
692 }
693
694 if (parcel.rho() == -1)
695 {
697 << "The kinematic cloud needs rho0 in the constantProperties "
698 << " dictionary. " << nl
699 << abort(FatalError);
700 }
701}
702
703
704template<class CloudType>
706{
707 cloudCopyPtr_.reset
708 (
709 static_cast<KinematicCloud<CloudType>*>
710 (
711 clone(this->name() + "Copy").ptr()
712 )
713 );
714}
715
716
717template<class CloudType>
720 cloudReset(cloudCopyPtr_());
721 cloudCopyPtr_.clear();
722}
723
724
725template<class CloudType>
727{
728 rhokTrans().field() = Zero;
729 UTrans().field() = Zero;
730 UCoeff().field() = Zero;
731}
732
733
734template<class CloudType>
735template<class Type>
737(
740 const word& name
741) const
742{
743 const scalar coeff = solution_.relaxCoeff(name);
744 field = field0 + coeff*(field - field0);
745}
746
747
748template<class CloudType>
749template<class Type>
751(
753 const word& name
754) const
756 const scalar coeff = solution_.relaxCoeff(name);
757 field *= coeff;
758}
759
760
761template<class CloudType>
763(
764 const KinematicCloud<CloudType>& cloudOldTime
765)
766{
767 this->relax(rhokTrans_(), cloudOldTime.rhokTrans(), "rhok");
768 this->relax(UTrans_(), cloudOldTime.UTrans(), "U");
769 this->relax(UCoeff_(), cloudOldTime.UCoeff(), "U");
770}
771
772
773template<class CloudType>
775{
776 this->scale(rhokTrans_(), "rhok");
777 this->scale(UTrans_(), "U");
778 this->scale(UCoeff_(), "U");
779}
780
781
782template<class CloudType>
784(
785 const typename parcelType::trackingData& td
786)
787{
788 // force calculation of mesh dimensions - needed for parallel runs
789 // with topology change due to lazy evaluation of valid mesh dimensions
790 label nGeometricD = mesh_.nGeometricD();
791
792 Log_<< "\nSolving" << nGeometricD << "-D cloud " << this->name() << endl;
793
794 this->dispersion().cacheFields(true);
795 forces_.cacheFields(true);
796
797 pAmbient_ = constProps_.dict().template
798 getOrDefault<scalar>("pAmbient", pAmbient_);
799
800 if (this->dampingModel().active() || this->packingModel().active())
801 {
802 const_cast<typename parcelType::trackingData&>(td).updateAverages(*this);
803 }
804
805 if (this->dampingModel().active())
806 {
807 this->dampingModel().cacheFields(true);
808 }
809 if (this->packingModel().active())
810 {
811 this->packingModel().cacheFields(true);
812 }
813
815
816 functions_.preEvolve(td);
817}
818
819
820template<class CloudType>
822{
823 if (solution_.canEvolve())
824 {
825 typename parcelType::trackingData td(*this);
826 solve(*this, td);
827 }
828}
829
830
831template<class CloudType>
832template<class TrackCloudType>
834(
835 TrackCloudType& cloud,
836 typename parcelType::trackingData& td
837)
838{
839 td.part() = parcelType::trackingData::tpLinearTrack;
840 CloudType::move(cloud, td, solution_.trackTime());
841
842 if (isotropyModel_->active())
843 {
844 td.updateAverages(cloud);
845 isotropyModel_->calculate();
847
848 updateCellOccupancy();
849}
850
852template<class CloudType>
854(
855 const parcelType& p,
856 const polyPatch& pp,
857 vector& nw,
858 vector& Up
859) const
860{
861 p.patchData(nw, Up);
862
863 // If this is a wall patch, then there may be a non-zero tangential velocity
864 // component; the lid velocity in a lid-driven cavity case, for example. We
865 // want the particle to interact with this velocity, so we look it up in the
866 // velocity field and use it to set the wall-tangential component.
868 {
869 const label patchi = pp.index();
870 const label patchFacei = pp.whichFace(p.face());
872 // We only want to use the boundary condition value only if it is set
873 // by the boundary condition. If the boundary values are extrapolated
874 // (e.g., slip conditions) then they represent the motion of the fluid
875 // just inside the domain rather than that of the wall itself.
876 if (U_.boundaryField()[patchi].fixesValue())
878 const vector Uw1(U_.boundaryField()[patchi][patchFacei]);
879 const vector& Uw0 =
880 U_.oldTime().boundaryField()[patchi][patchFacei];
881
882 const scalar f = p.currentTimeFraction();
883
884 const vector Uw(Uw0 + f*(Uw1 - Uw0));
885
886 const tensor nnw(nw*nw);
887
888 Up = (nnw & Up) + Uw - (nnw & Uw);
889 }
890 }
891}
892
893
894template<class CloudType>
896{
898 injectors_.updateMesh();
900}
901
903template<class CloudType>
905{
908 updateMesh();
909}
910
911
912template<class CloudType>
914{
915 const vector linearMomentum =
917
918 const scalar linearKineticEnergy =
920
921 const label nTotParcel = returnReduce(this->size(), sumOp<label>());
922
923 const scalar particlePerParcel =
924 (
925 nTotParcel
927 : 0
928 );
929
930 Log_<< "Cloud: " << this->name() << nl
931 << " Current number of parcels = " << nTotParcel << nl
932 << " Current mass in system = "
934 << " Linear momentum = " << linearMomentum << nl
935 << " |Linear momentum| = " << mag(linearMomentum) << nl
936 << " Linear kinetic energy = " << linearKineticEnergy << nl
937 << " Average particle per parcel = " << particlePerParcel << nl;
938
939
940 injectors_.info();
941 this->surfaceFilm().info();
942 this->patchInteraction().info();
943
944 if (this->packingModel().active())
945 {
947
948 if (this->db().time().writeTime())
949 {
950 alpha().write();
951 }
952
953 auto limits = gMinMax(alpha().primitiveField());
954
955 Log_<< " Min cell volume fraction = " << limits.min() << nl
956 << " Max cell volume fraction = " << limits.max() << endl;
957
958 if (limits.max() < SMALL)
959 {
960 return;
962
963 scalar nMin = GREAT;
964
965 forAll(this->mesh().cells(), celli)
967 const label n = this->cellOccupancy()[celli].size();
968
969 if (n > 0)
970 {
971 const scalar nPack = n*limits.max()/alpha()[celli];
972
973 if (nPack < nMin)
974 {
975 nMin = nPack;
976 }
977 }
978 }
979
980 reduce(nMin, minOp<scalar>());
982 Log_<< " Min dense number of parcels = " << nMin << endl;
983 }
984}
985
986
987template<class CloudType>
992
993
994template<class CloudType>
996{
997 parcelType::writeObjects(*this, obr);
998}
999
1000
1001// ************************************************************************* //
label n
const List< DynamicList< molecule * > > & cellOccupancy
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
void writePositions() const
Write positions to <cloudName>_positions.obj file.
Definition Cloud.C:329
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
Definition Cloud.C:137
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition Cloud.C:297
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
Definition Cloud.C:125
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
const word & cloudName() const
Base class for collisional damping models.
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:68
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
Definition IOobject.H:1041
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A simple container for options an IOstream can normally have.
@ ASCII
"ascii" (normal default)
Base class for collisional return-to-isotropy models.
Templated base class for kinematic cloud.
autoPtr< volScalarField::Internal > rhokTrans_
cloudSolution solution_
Solution properties.
autoPtr< DampingModel< KinematicCloud< CloudType > > > dampingModel_
Damping model.
void setModels()
Set cloud sub-models.
void relaxSources(const KinematicCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
forceType forces_
Optional particle forces.
void scale(DimensionedField< Type, volMesh > &field, const word &name) const
Scale field.
const fvMesh & mesh_
References to the mesh and time databases.
void storeState()
Store the current cloud state.
autoPtr< SurfaceFilmModel< KinematicCloud< CloudType > > > surfaceFilmModel_
Surface film model.
void patchData(const parcelType &p, const polyPatch &pp, vector &normal, vector &Up) const
Calculate the patch normal and velocity to interact with,.
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.
const PatchInteractionModel< KinematicCloud< CloudType > > & patchInteraction() const
Return const-access to the patch interaction model.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
autoPtr< integrationScheme > UIntegrator_
Velocity integration.
const volVectorField & U() const
Return carrier gas velocity.
void relax(DimensionedField< Type, volMesh > &field, const DimensionedField< Type, volMesh > &field0, const word &name) const
Relax field.
InjectionModelList< KinematicCloud< CloudType > > injectors_
Injector models.
const volScalarField & rho() const
Return carrier gas density.
autoPtr< IsotropyModel< KinematicCloud< CloudType > > > isotropyModel_
Exchange model.
autoPtr< DispersionModel< KinematicCloud< CloudType > > > dispersionModel_
Dispersion model.
void scaleSources()
Apply scaling to (transient) cloud sources.
autoPtr< PatchInteractionModel< KinematicCloud< CloudType > > > patchInteractionModel_
Patch interaction model.
List< DynamicList< parcelType * > > & cellOccupancy()
Return the cell occupancy information for each.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
void updateCellOccupancy()
Update (i.e. build) the cellOccupancy if it has.
const dictionary subModelProperties_
Sub-models dictionary.
const PackingModel< KinematicCloud< CloudType > > & packingModel() const
Return const access to the packing model.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
virtual autoPtr< Cloud< parcelType > > clone(const word &name)
Construct and return clone based on (this) with new name.
const DispersionModel< KinematicCloud< CloudType > > & dispersion() const
Return const-access to the dispersion model.
void evolve()
Evolve the cloud.
void postEvolve(const typename parcelType::trackingData &td)
Post-evolve.
volScalarField::Internal & UCoeff()
Return coefficient for carrier phase U equation.
functionType functions_
Optional cloud function objects.
const StochasticCollisionModel< KinematicCloud< CloudType > > & stochasticCollision() const
Return const-access to the stochastic collision model.
void evolveCloud(TrackCloudType &cloud, typename parcelType::trackingData &td)
Evolve the cloud.
autoPtr< List< DynamicList< parcelType * > > > cellOccupancyPtr_
Cell occupancy information for each parcel, (demand driven).
const DampingModel< KinematicCloud< CloudType > > & dampingModel() const
Return const access to the damping model.
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.
Random rndGen_
Random number generator - used by some injection routines.
void resetSourceTerms()
Reset the cloud source terms.
volScalarField::Internal & rhokTrans()
Return reference to mass for kinematic source.
autoPtr< PackingModel< KinematicCloud< CloudType > > > packingModel_
Packing model.
void buildCellOccupancy()
Build the cellOccupancy.
void updateMesh()
Update mesh.
bool log
Flag to write log into Info.
IOdictionary outputProperties_
Dictionary of output properties.
autoPtr< StochasticCollisionModel< KinematicCloud< CloudType > > > stochasticCollisionModel_
Stochastic collision model.
void solve(TrackCloudType &cloud, typename parcelType::trackingData &td)
Solve the cloud - calls all evolution functions.
virtual void writeObjects(objectRegistry &obr) const
Write particle fields as objects into the obr registry.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
Base class for packing models.
Templated patch interaction model class.
Inter-processor communications stream.
Definition Pstream.H:59
Templated stochastic collision model class.
Templated wall surface film model class.
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
virtual void readObjects(const objectRegistry &obr)
Read particle fields from objects in the obr registry.
Definition cloud.C:85
static const word prefix
The prefix to local: lagrangian.
Definition cloud.H:79
virtual void writeObjects(objectRegistry &obr) const
Write particle fields as objects into the obr registry.
Definition cloud.C:91
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static autoPtr< integrationScheme > New(const word &phiName, const dictionary &dict)
Select an integration scheme.
A class for handling keywords in dictionaries.
Definition keyType.H:69
Virtual abstract base class for templated KinematicCloud.
kinematicCloud()=default
Null constructor.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Registry of regIOobjects.
const Time & time() const noexcept
Return time registry.
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
UEqn relax()
regionModels::surfaceFilmModel & surfaceFilm
rDeltaTY field()
auto limits
Definition setRDeltaT.H:186
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
const cellShapeList & cells
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
word timeName
Definition getTimeIndex.H:3
#define Log_
Report write to Foam::Info if the class log switch is true.
Different types of constants.
Namespace for handling debugging switches.
Definition debug.C:45
ILList< DLListBase, T > IDLList
Definition IDLList.H:39
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
Tensor< scalar > tensor
Definition symmTensor.H:57
dimensionedScalar log(const dimensionedScalar &ds)
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
errorManip< error > abort(error &err)
Definition errorManip.H:139
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
labelList f(nPoints)
volScalarField & alpha
CEqn solve()
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299