49template<
class CloudType>
126template<
class CloudType>
127template<
class TrackCloudType>
130 TrackCloudType&
cloud,
131 typename parcelType::trackingData&
td
136 log = solution_.log();
138 if (solution_.steadyState())
146 if (solution_.coupled())
155 evolveCloud(cloud,
td);
157 if (solution_.coupled())
159 cloud.scaleSources();
165 cloud.postEvolve(
td);
167 if (solution_.steadyState())
169 cloud.restoreState();
174template<
class CloudType>
177 if (!cellOccupancyPtr_)
179 cellOccupancyPtr_.reset
181 new List<DynamicList<parcelType*>>(mesh_.nCells())
184 else if (cellOccupancyPtr_().size() != mesh_.nCells())
189 cellOccupancyPtr_().setSize(mesh_.nCells());
192 List<DynamicList<parcelType*>>&
cellOccupancy = cellOccupancyPtr_();
199 for (parcelType&
p : *
this)
206template<
class CloudType>
212 if (cellOccupancyPtr_)
219template<
class CloudType>
220template<
class TrackCloudType>
223 TrackCloudType& cloud,
224 typename parcelType::trackingData&
td
227 if (solution_.coupled())
229 cloud.resetSourceTerms();
232 if (solution_.transient())
234 label preInjectionSize = this->size();
240 if (preInjectionSize != this->size())
242 updateCellOccupancy();
243 preInjectionSize = this->size();
246 injectors_.inject(cloud,
td);
250 cloud.motion(cloud,
td);
252 stochasticCollision().update(
td, solution_.trackTime());
258 injectors_.injectSteadyState(cloud,
td, solution_.trackTime());
260 td.part() = parcelType::trackingData::tpLinearTrack;
266template<
class CloudType>
269 const typename parcelType::trackingData&
td
276 this->writePositions();
279 this->dispersion().cacheFields(
false);
281 this->patchInteraction().postEvolve();
283 forces_.cacheFields(
false);
285 functions_.postEvolve(
td);
287 solution_.nextIter();
289 if (this->db().time().writeTime())
291 outputProperties_.writeObject
295 IOstreamOption::ASCII,
296 this->db().time().writeCompression()
302 if (this->dampingModel().active())
304 this->dampingModel().cacheFields(
false);
306 if (this->packingModel().active())
313template<
class CloudType>
316 CloudType::cloudReset(c);
320 forces_.transfer(c.forces_);
322 functions_.transfer(c.functions_);
324 injectors_.transfer(c.injectors_);
326 dispersionModel_.reset(c.dispersionModel_.ptr());
327 patchInteractionModel_.reset(c.patchInteractionModel_.ptr());
328 stochasticCollisionModel_.reset(c.stochasticCollisionModel_.ptr());
329 surfaceFilmModel_.reset(c.surfaceFilmModel_.ptr());
331 packingModel_.reset(c.packingModel_.ptr());
332 dampingModel_.reset(c.dampingModel_.ptr());
333 isotropyModel_.reset(c.isotropyModel_.ptr());
341template<
class CloudType>
342Foam::KinematicCloud<CloudType>::KinematicCloud
354 cloudCopyPtr_(nullptr),
492 parcelType::readFields(*this);
493 this->deleteLostParticles();
497 if (solution_.resetSourcesOnStartup())
504template<
class CloudType>
505Foam::KinematicCloud<CloudType>::KinematicCloud
513 cloudCopyPtr_(nullptr),
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_),
527 pAmbient_(c.pAmbient_),
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()),
536 packingModel_(c.packingModel_->clone()),
537 dampingModel_(c.dampingModel_->clone()),
538 isotropyModel_(c.isotropyModel_->clone()),
540 UIntegrator_(c.UIntegrator_->clone()),
593template<
class CloudType>
594Foam::KinematicCloud<CloudType>::KinematicCloud
603 cloudCopyPtr_(nullptr),
621 name +
"OutputProperties",
632 subModelProperties_(),
634 cellOccupancyPtr_(nullptr),
635 cellLengthScale_(c.cellLengthScale_),
640 pAmbient_(c.pAmbient_),
641 forces_(*this,
mesh),
644 dispersionModel_(nullptr),
645 patchInteractionModel_(nullptr),
646 stochasticCollisionModel_(nullptr),
647 surfaceFilmModel_(nullptr),
649 packingModel_(nullptr),
650 dampingModel_(nullptr),
651 isotropyModel_(nullptr),
653 UIntegrator_(nullptr),
663template<
class CloudType>
667 const scalar lagrangianDt
671 if (constProps_.rho0() != -1)
678template<
class CloudType>
682 const scalar lagrangianDt,
683 const bool fullyDescribed
686 const scalar carrierDt = mesh_.time().deltaTValue();
687 parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
689 if (parcel.typeId() == -1)
691 parcel.typeId() = constProps_.parcelTypeId();
694 if (parcel.rho() == -1)
697 <<
"The kinematic cloud needs rho0 in the constantProperties "
698 <<
" dictionary. " <<
nl
704template<
class CloudType>
709 static_cast<KinematicCloud<CloudType>*
>
717template<
class CloudType>
721 cloudCopyPtr_.clear();
725template<
class CloudType>
728 rhokTrans().field() =
Zero;
734template<
class CloudType>
748template<
class CloudType>
761template<
class CloudType>
773template<
class CloudType>
782template<
class CloudType>
785 const typename parcelType::trackingData&
td
790 label nGeometricD = mesh_.nGeometricD();
792 Log_<<
"\nSolving" << nGeometricD <<
"-D cloud " << this->
name() << endl;
794 this->dispersion().cacheFields(
true);
795 forces_.cacheFields(
true);
797 pAmbient_ = constProps_.dict().template
798 getOrDefault<scalar>(
"pAmbient", pAmbient_);
800 if (this->dampingModel().active() || this->packingModel().active())
802 const_cast<typename parcelType::trackingData&
>(
td).updateAverages(*
this);
805 if (this->dampingModel().active())
807 this->dampingModel().cacheFields(
true);
809 if (this->packingModel().active())
811 this->packingModel().cacheFields(
true);
820template<
class CloudType>
823 if (solution_.canEvolve())
825 typename parcelType::trackingData
td(*
this);
831template<
class CloudType>
832template<
class TrackCloudType>
835 TrackCloudType& cloud,
836 typename parcelType::trackingData&
td
839 td.part() = parcelType::trackingData::tpLinearTrack;
848 updateCellOccupancy();
852template<
class CloudType>
869 const label patchi =
pp.index();
870 const label patchFacei =
pp.whichFace(
p.face());
876 if (
U_.boundaryField()[patchi].fixesValue())
878 const vector Uw1(
U_.boundaryField()[patchi][patchFacei]);
880 U_.oldTime().boundaryField()[patchi][patchFacei];
882 const scalar
f =
p.currentTimeFraction();
884 const vector Uw(Uw0 +
f*(Uw1 - Uw0));
888 Up = (nnw & Up) + Uw - (nnw & Uw);
894template<
class CloudType>
903template<
class CloudType>
912template<
class CloudType>
918 const scalar linearKineticEnergy =
923 const scalar particlePerParcel =
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;
948 if (this->
db().
time().writeTime())
955 Log_<<
" Min cell volume fraction = " <<
limits.min() <<
nl
956 <<
" Max cell volume fraction = " <<
limits.max() <<
endl;
982 Log_<<
" Min dense number of parcels = " << nMin <<
endl;
987template<
class CloudType>
994template<
class CloudType>
997 parcelType::writeObjects(*
this, obr);
const List< DynamicList< molecule * > > & cellOccupancy
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
void writePositions() const
Write positions to <cloudName>_positions.obj file.
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
label size() const noexcept
The number of elements in list.
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.
@ 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,...
const objectRegistry & db() const noexcept
Return the local objectRegistry.
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
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.
autoPtr< volVectorField::Internal > UTrans_
void setModels()
Set cloud sub-models.
scalar massInSystem() const
scalar totalParticlePerParcel() const
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 tmp< volScalarField > theta() const
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.
const tmp< volScalarField > alpha() const
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.
scalarField cellLengthScale_
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.
vector linearMomentumOfSystem() const
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.
const volVectorField & U_
List< DynamicList< parcelType * > > & cellOccupancy()
Return the cell occupancy information for each.
const volScalarField & mu_
const dimensionedVector & g_
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.
parcelType::constantProperties constProps_
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.
const volScalarField & rho_
IOdictionary particleProperties_
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.
const dimensionedVector & g() const
const fvMesh & mesh() const
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.
autoPtr< volScalarField::Internal > UCoeff_
void updateMesh()
Update mesh.
bool log
Flag to write log into Info.
IOdictionary outputProperties_
Dictionary of output properties.
scalar linearKineticEnergyOfSystem() const
const volScalarField & mu() const
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...
Base class for packing models.
Templated patch interaction model class.
Inter-processor communications stream.
Templated stochastic collision model class.
Templated wall surface film model class.
A cloud is a registry collection of lagrangian particles.
virtual void readObjects(const objectRegistry &obr)
Read particle fields from objects in the obr registry.
static const word prefix
The prefix to local: lagrangian.
virtual void writeObjects(objectRegistry &obr) const
Write particle fields as objects into the obr registry.
Mesh data needed to do the Finite Volume discretisation.
static autoPtr< integrationScheme > New(const word &phiName, const dictionary &dict)
Select an integration scheme.
A class for handling keywords in dictionaries.
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.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
regionModels::surfaceFilmModel & surfaceFilm
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
#define Log_
Report write to Foam::Info if the class log switch is true.
Different types of constants.
Namespace for handling debugging switches.
ILList< DLListBase, T > IDLList
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
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.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
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)
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
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.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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).
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
#define forAll(list, i)
Loop across all elements in list.