39 template<
class CloudType>
42 dispersionModel_.reset
51 patchInteractionModel_.reset
60 stochasticCollisionModel_.reset
83 solution_.integrationSchemes()
89 template<
class CloudType>
90 template<
class TrackCloudType>
93 TrackCloudType&
cloud,
94 typename parcelType::trackingData& td
97 this->changeTimeStep();
99 if (solution_.steadyState())
106 if (solution_.coupled())
108 cloud.resetSourceTerms();
111 if (solution_.transient())
113 label preInjectionSize = this->size();
115 this->surfaceFilm().inject(
cloud);
119 if (preInjectionSize != this->size())
121 updateCellOccupancy();
122 preInjectionSize = this->size();
125 injectors_.inject(
cloud, td);
131 stochasticCollision().update(td);
135 injectors_.injectSteadyState(
cloud, td);
137 CloudType::move(
cloud, td);
140 if (solution_.coupled() && solution_.transient())
142 cloud.scaleSources();
145 if (solution_.coupled() && solution_.steadyState())
154 if (solution_.steadyState())
156 cloud.restoreState();
161 template<
class CloudType>
164 if (cellOccupancyPtr_.empty())
166 cellOccupancyPtr_.reset
171 else if (cellOccupancyPtr_().size() != this->mesh().nCells())
176 cellOccupancyPtr_().setSize(this->mesh().nCells());
193 template<
class CloudType>
199 if (cellOccupancyPtr_.valid())
201 buildCellOccupancy();
206 template<
class CloudType>
213 this->writePositions();
216 this->dispersion().cacheFields(
false);
218 forces_.cacheFields(
false);
220 functions_.postEvolve();
222 solution_.nextIter();
224 if (this->db().time().writeTime())
226 outputProperties_.writeObject
229 IOstream::currentVersion,
230 this->db().time().writeCompression(),
237 template<
class CloudType>
240 CloudType::cloudReset(
c);
242 forces_.transfer(
c.forces_);
244 functions_.transfer(
c.functions_);
246 injectors_.transfer(
c.injectors_);
248 dispersionModel_.reset(
c.dispersionModel_.ptr());
249 patchInteractionModel_.reset(
c.patchInteractionModel_.ptr());
250 stochasticCollisionModel_.reset(
c.stochasticCollisionModel_.ptr());
251 filmModel_.reset(
c.filmModel_.ptr());
253 UIntegrator_.reset(
c.UIntegrator_.ptr());
259 template<
class CloudType>
271 cloudCopyPtr_(nullptr),
277 this->mesh().time().constant(),
288 this->mesh().time().
name(),
295 solution_(this->mesh(), particleProperties_.subDict(
"solution")),
296 constProps_(particleProperties_),
299 particleProperties_.subOrEmptyDict(
"subModels", true)
301 cpuLoad_(particleProperties_.lookupOrDefault(
"cpuLoad", false)),
303 stdNormal_(rndGen_.generator()),
305 cellLengthScale_(
mag(
cbrt(this->mesh().V()))),
315 subModelProperties_.subOrEmptyDict
325 particleProperties_.subOrEmptyDict(
"cloudFunctions")
329 subModelProperties_.subOrEmptyDict(
"injectionModels"),
332 dispersionModel_(nullptr),
333 patchInteractionModel_(nullptr),
334 stochasticCollisionModel_(nullptr),
336 UIntegrator_(nullptr),
343 this->
name() +
":UTrans",
344 this->db().time().
name(),
359 this->
name() +
":UCoeff",
360 this->db().time().
name(),
385 template<
class CloudType>
400 template<
class CloudType>
408 cloudCopyPtr_(nullptr),
409 particleProperties_(
c.particleProperties_),
410 outputProperties_(
c.outputProperties_),
411 solution_(
c.solution_),
412 constProps_(
c.constProps_),
413 subModelProperties_(
c.subModelProperties_),
414 cpuLoad_(
c.cpuLoad_),
416 stdNormal_(
c.stdNormal_),
417 cellOccupancyPtr_(nullptr),
418 cellLengthScale_(
c.cellLengthScale_),
423 pAmbient_(
c.pAmbient_),
425 functions_(
c.functions_),
426 injectors_(
c.injectors_),
427 dispersionModel_(
c.dispersionModel_->
clone()),
428 patchInteractionModel_(
c.patchInteractionModel_->
clone()),
429 stochasticCollisionModel_(
c.stochasticCollisionModel_->
clone()),
430 filmModel_(
c.filmModel_->
clone()),
431 UIntegrator_(
c.UIntegrator_->
clone()),
438 this->
name() +
":UTrans",
439 this->db().time().
name(),
455 this->db().time().
name(),
467 template<
class CloudType>
476 cloudCopyPtr_(nullptr),
482 mesh.time().constant(),
493 name +
"OutputProperties",
494 this->mesh().time().
name(),
505 cpuLoad_(
c.cpuLoad_),
507 stdNormal_(rndGen_.generator()),
508 cellOccupancyPtr_(nullptr),
509 cellLengthScale_(
c.cellLengthScale_),
514 pAmbient_(
c.pAmbient_),
515 forces_(*this, mesh),
518 dispersionModel_(nullptr),
519 patchInteractionModel_(nullptr),
520 stochasticCollisionModel_(nullptr),
522 UIntegrator_(nullptr),
530 template<
class CloudType>
537 template<
class CloudType>
543 parcel.rho() = constProps_.rho0();
547 template<
class CloudType>
551 const label injectori
554 if (parcel.typeId() == -1)
556 parcel.typeId() = constProps_.parcelTypeId();
561 template<
class CloudType>
574 template<
class CloudType>
577 cloudReset(cloudCopyPtr_());
578 cloudCopyPtr_.clear();
582 template<
class CloudType>
585 UTransRef().primitiveFieldRef() =
Zero;
586 UCoeffRef().primitiveFieldRef() = 0.0;
590 template<
class CloudType>
599 const scalar coeff = solution_.relaxCoeff(
name);
600 field = field0 + coeff*(field - field0);
604 template<
class CloudType>
612 const scalar coeff = solution_.relaxCoeff(
name);
617 template<
class CloudType>
628 template<
class CloudType>
631 this->scale(UTrans_(),
"U");
632 this->scale(UCoeff_(),
"U");
636 template<
class CloudType>
641 label nGeometricD = this->mesh().nGeometricD();
643 Info<<
nl <<
"Solving " << nGeometricD <<
"-D cloud " << this->
name()
646 this->dispersion().cacheFields(
true);
647 forces_.cacheFields(
true);
648 updateCellOccupancy();
650 pAmbient_ = constProps_.dict().template
651 lookupOrDefault<scalar>(
"pAmbient", pAmbient_);
653 functions_.preEvolve();
657 template<
class CloudType>
660 if (solution_.canEvolve())
662 typename parcelType::trackingData td(*
this);
669 template<
class CloudType>
670 template<
class TrackCloudType>
673 TrackCloudType&
cloud,
674 typename parcelType::trackingData& td
679 updateCellOccupancy();
683 template<
class CloudType>
692 p.patchData(this->mesh(), nw, Up);
693 Up /= this->mesh().time().deltaTValue();
699 if (!this->mesh().moving() && isA<wallPolyPatch>(pp))
708 if (U_.boundaryField()[
patchi].fixesValue())
710 const vector Uw1 = U_.boundaryField()[
patchi][patchFacei];
712 U_.oldTime().boundaryField()[
patchi][patchFacei];
714 const vector Uw = Uw0 +
p.stepFraction()*(Uw1 - Uw0);
718 Up = (nnw & Up) + Uw - (nnw & Uw);
724 template<
class CloudType>
729 updateCellOccupancy();
731 injectors_.topoChange();
733 cellLengthScale_ =
mag(
cbrt(this->mesh().V()));
737 template<
class CloudType>
742 updateCellOccupancy();
744 injectors_.topoChange();
746 cellLengthScale_ =
mag(
cbrt(this->mesh().V()));
750 template<
class CloudType>
755 updateCellOccupancy();
757 injectors_.topoChange();
759 cellLengthScale_ =
mag(
cbrt(this->mesh().V()));
763 template<
class CloudType>
766 vector linearMomentum = linearMomentumOfSystem();
769 scalar linearKineticEnergy = linearKineticEnergyOfSystem();
773 <<
" Current number of parcels = "
775 <<
" Current mass in system = "
777 <<
" Linear momentum = "
778 << linearMomentum <<
nl
779 <<
" |Linear momentum| = "
780 <<
mag(linearMomentum) <<
nl
781 <<
" Linear kinetic energy = "
782 << linearKineticEnergy <<
nl;
784 injectors_.info(
Info);
785 this->surfaceFilm().info(
Info);
786 this->patchInteraction().info(
Info);
#define forAll(list, i)
Loop across all elements in list.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
const List< DynamicList< molecule * > > & cellOccupancy
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
void deleteLostParticles()
Remove lost particles from cloud and delete.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td)
Move the particles.
Templated base class for dsmc cloud.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
InfoProxy< IOobject > info() const
Return info proxy.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Templated base class for momentum cloud.
void postEvolve()
Post-evolve.
cloudSolution solution_
Solution properties.
autoPtr< volVectorField::Internal > UTrans_
Momentum.
void setModels()
Set cloud sub-models.
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
virtual ~MomentumCloud()
Destructor.
void scale(DimensionedField< Type, volMesh > &field, const word &name) const
Scale field.
void storeState()
Store the current cloud state.
void patchData(const parcelType &p, const polyPatch &pp, vector &normal, vector &Up) const
Calculate the patch normal and velocity to interact with,.
void relax(DimensionedField< Type, volMesh > &field, const DimensionedField< Type, volMesh > &field0, const word &name) const
Relax field.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
void scaleSources()
Apply scaling to (transient) cloud sources.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
MomentumCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g, const bool readFields=true)
Construct given carrier fields.
void updateCellOccupancy()
Update (i.e. build) the cellOccupancy if it has.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
void evolve()
Evolve the cloud.
void cloudReset(MomentumCloud< CloudType > &c)
Reset state of cloud.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
void checkParcelProperties(parcelType &parcel, const label injectori)
Check parcel properties.
void preEvolve()
Pre-evolve.
void resetSourceTerms()
Reset the cloud source terms.
void buildCellOccupancy()
Build the cellOccupancy.
autoPtr< volScalarField::Internal > UCoeff_
Coefficient for carrier phase U equation.
void setParcelThermoProperties(parcelType &parcel)
Set parcel thermo properties.
void relaxSources(const MomentumCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
void solve(TrackCloudType &cloud, typename parcelType::trackingData &td)
Solve the cloud - calls all evolution functions.
Templated patch interaction model class.
Templated stochastic collision model class.
Templated wall surface film model class.
const Switch resetSourcesOnStartup() const
Return const access to the reset sources flag.
A cloud is a collection of lagrangian particles.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Base-class for fluid thermodynamic properties.
Mesh data needed to do the Finite Volume discretisation.
label index() const
Return the index of this patch in the boundaryMesh.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
A patch is a list of labels that address the faces in the global face list.
label whichFace(const label l) const
Return label of face in patch from global face label.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar mu
Atomic mass unit.
const dimensionedScalar c
Speed of light in a vacuum.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
word name(const bool)
Return a word representation of a bool.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the specified type.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionSet dimMass
const dimensionSet dimVelocity
dimensionedScalar cbrt(const dimensionedScalar &ds)
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
const word cloudName(propsDict.lookup("cloudName"))