Go to the documentation of this file.
34 template<
class CloudType>
41 this->subModelProperties(),
48 template<
class CloudType>
49 template<
class TrackCloudType>
52 TrackCloudType&
cloud,
53 typename parcelType::trackingData& td,
57 td.part() = parcelType::trackingData::tpVelocityHalfStep;
58 CloudType::move(
cloud, td, deltaT);
60 td.part() = parcelType::trackingData::tpLinearTrack;
61 CloudType::move(
cloud, td, deltaT);
66 this->updateCellOccupancy();
68 this->collision().collide();
70 td.part() = parcelType::trackingData::tpVelocityHalfStep;
71 CloudType::move(
cloud, td, deltaT);
76 template<
class CloudType>
79 CloudType::cloudReset(
c);
81 collisionModel_.reset(
c.collisionModel_.ptr());
87 template<
class CloudType>
99 constProps_(this->particleProperties()),
100 collisionModel_(nullptr)
109 this->deleteLostParticles();
119 <<
"Collision modelling not currently available "
126 template<
class CloudType>
129 CollidingCloud<CloudType>&
c,
134 collisionModel_(
c.collisionModel_->clone())
138 template<
class CloudType>
147 collisionModel_(nullptr)
153 template<
class CloudType>
160 template<
class CloudType>
167 clone(this->
name() +
"Copy").ptr()
173 template<
class CloudType>
176 cloudReset(cloudCopyPtr_());
177 cloudCopyPtr_.clear();
181 template<
class CloudType>
186 typename parcelType::trackingData td(*
this);
188 this->
solve(*
this, td);
193 template<
class CloudType>
194 template<
class TrackCloudType>
197 TrackCloudType& cloud,
198 typename parcelType::trackingData& td
207 label nSubCycles = collision().nSubCycles();
211 Info<<
" " << nSubCycles <<
" move-collide subCycles" <<
endl;
213 subCycleTime moveCollideSubCycle
215 const_cast<Time&
>(this->db().time()),
219 while(!(++moveCollideSubCycle).
end())
221 moveCollide(cloud, td, this->db().time().deltaTValue());
224 moveCollideSubCycle.endSubCycle();
228 moveCollide(
cloud, td, this->db().time().deltaTValue());
233 template<
class CloudType>
238 scalar rotationalKineticEnergy = rotationalKineticEnergyOfSystem();
239 reduce(rotationalKineticEnergy, sumOp<scalar>());
241 Info<<
" Rotational kinetic energy = "
242 << rotationalKineticEnergy <<
nl;
Adds coolisions to kinematic clouds.
void cloudReset(CollidingCloud< CloudType > &c)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Selector class for relaxation factors, solver type and solution.
const word cloudName(propsDict.get< word >("cloud"))
A class for handling words, derived from Foam::string.
const dimensionedScalar mu
void moveCollide(TrackCloudType &cloud, typename parcelType::trackingData &td, const scalar deltaT)
Ostream & endl(Ostream &os)
A class for managing sub-cycling times.
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
DSMCCloud< dsmcParcel > CloudType
Templated base class for dsmc cloud.
reduce(hasMovingMesh, orOp< bool >())
Generic dimensioned Type class.
constexpr auto end(C &c) -> decltype(c.end())
Mesh data needed to do the Finite Volume discretisation.
const uniformDimensionedVectorField & g
errorManipArg< error, int > exit(error &err, const int errNo=1)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
bool isType(const Type &t)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
A cloud is a registry collection of lagrangian particles.
#define FatalErrorInFunction
virtual ~CollidingCloud()
Place holder for 'none' option.
const dimensionedScalar c
word name(const expressions::valueTypeCode typeCode)
Generic GeometricField class.
Templated collision model class.