34template<
class CloudType>
41 this->subModelProperties(),
48template<
class CloudType>
49template<
class TrackCloudType>
52 TrackCloudType&
cloud,
53 typename parcelType::trackingData&
td,
57 td.part() = parcelType::trackingData::tpVelocityHalfStep;
60 td.part() = parcelType::trackingData::tpLinearTrack;
66 this->updateCellOccupancy();
68 this->collision().collide();
70 td.part() = parcelType::trackingData::tpVelocityHalfStep;
76template<
class CloudType>
87template<
class CloudType>
88Foam::CollidingCloud<CloudType>::CollidingCloud
99 constProps_(this->particleProperties()),
100 collisionModel_(nullptr)
108 parcelType::readFields(*
this);
119 <<
"Collision modelling not currently available "
126template<
class CloudType>
127Foam::CollidingCloud<CloudType>::CollidingCloud
138template<
class CloudType>
139Foam::CollidingCloud<CloudType>::CollidingCloud
143 const CollidingCloud<CloudType>& c
153template<
class CloudType>
160template<
class CloudType>
173template<
class CloudType>
177 cloudCopyPtr_.clear();
181template<
class CloudType>
186 typename parcelType::trackingData
td(*
this);
193template<
class CloudType>
194template<
class TrackCloudType>
197 TrackCloudType& cloud,
198 typename parcelType::trackingData&
td
207 label nSubCycles = collision().nSubCycles();
211 Log_<<
" " << nSubCycles <<
" move-collide subCycles" <<
endl;
213 subCycleTime moveCollideSubCycle
215 const_cast<Time&
>(this->db().time()),
219 while(!(++moveCollideSubCycle).
end())
224 moveCollideSubCycle.endSubCycle();
233template<
class CloudType>
238 scalar rotationalKineticEnergy = rotationalKineticEnergyOfSystem();
239 reduce(rotationalKineticEnergy, sumOp<scalar>());
241 Log_<<
" Rotational kinetic energy = "
242 << rotationalKineticEnergy <<
nl;
const word cloudName(propsDict.get< word >("cloud"))
const uniformDimensionedVectorField & g
void deleteLostParticles()
Remove lost particles from cloud and delete.
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
Adds coolisions to kinematic clouds.
void setModels()
Set cloud sub-models.
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
void storeState()
Store the current cloud state.
autoPtr< CollisionModel< CollidingCloud< CloudType > > > collisionModel_
Collision model.
scalar rotationalKineticEnergyOfSystem() const
Total rotational kinetic energy in the system.
void cloudReset(CollidingCloud< CloudType > &c)
Reset state of cloud.
void moveCollide(TrackCloudType &cloud, typename parcelType::trackingData &td, const scalar deltaT)
Move-collide particles.
virtual autoPtr< Cloud< parcelType > > clone(const word &name)
Construct and return clone based on (this) with new name.
void evolve()
Evolve the cloud.
parcelType::constantProperties constProps_
Thermo parcel constant properties.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
virtual ~CollidingCloud()
Destructor.
const CollisionModel< CollidingCloud< CloudType > > & collision() const
Return const access to the collision model.
Templated collision model class.
const word & cloudName() const
const IOdictionary & particleProperties() const
const fvMesh & mesh() const
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Place holder for 'none' option.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A cloud is a registry collection of lagrangian particles.
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const noexcept
Return time registry.
Selector class for relaxation factors, solver type and solution.
A class for managing sub-cycling times.
void endSubCycle()
End the sub-cycling and reset the time-state.
A class for handling words, derived from Foam::string.
#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.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
bool isType(const U &obj)
Check if typeid of the object and Type are identical.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).