Go to the documentation of this file.
35 template<
class CloudType>
42 this->subModelProperties(),
50 this->subModelProperties(),
58 this->subModelProperties(),
67 template<
class CloudType>
88 <<
"MPPIC modelling not available for steady state calculations"
102 template<
class CloudType>
110 packingModel_(
c.packingModel_->clone()),
111 dampingModel_(
c.dampingModel_->clone()),
112 isotropyModel_(
c.isotropyModel_->clone())
116 template<
class CloudType>
133 template<
class CloudType>
140 template<
class CloudType>
147 clone(this->
name() +
"Copy").ptr()
153 template<
class CloudType>
156 this->cloudReset(cloudCopyPtr_());
157 cloudCopyPtr_.clear();
161 template<
class CloudType>
166 typename parcelType::template
167 TrackingData<MPPICCloud<CloudType> > td(*
this);
174 template<
class CloudType>
175 template<
class TrackData>
182 td.part() = TrackData::tpLinearTrack;
183 CloudType::move(td, this->db().time().deltaTValue());
190 this->
forces().setCalcNonCoupled(
false);
191 this->
forces().setCalcCoupled(
false);
197 if (dampingModel_->active())
200 td.updateAverages(*
this);
203 dampingModel_->cacheFields(
true);
206 td.part() = TrackData::tpDampingNoTrack;
207 CloudType::move(td, this->db().time().deltaTValue());
210 td.part() = TrackData::tpCorrectTrack;
211 CloudType::move(td, this->db().time().deltaTValue());
214 dampingModel_->cacheFields(
false);
221 if (packingModel_->active())
224 td.updateAverages(*
this);
225 packingModel_->cacheFields(
true);
226 td.part() = TrackData::tpPackingNoTrack;
227 CloudType::move(td, this->db().time().deltaTValue());
228 td.part() = TrackData::tpCorrectTrack;
229 CloudType::move(td, this->db().time().deltaTValue());
230 packingModel_->cacheFields(
false);
237 if (isotropyModel_->active())
240 td.updateAverages(*
this);
243 isotropyModel_->calculate();
251 this->updateCellOccupancy();
254 this->
forces().setCalcNonCoupled(
true);
259 template<
class CloudType>
269 Info<<
" Min cell volume fraction = " << alphaMin <<
endl;
296 Info<<
" Min dense number of parcels = " << nMin <<
endl;
Selector class for relaxation factors, solver type and solution.
void setModels()
Set cloud sub-models.
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
const dimensionedScalar mu
Atomic mass unit.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
void storeState()
Store the current cloud state.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
const dimensionedVector & g
Ostream & endl(Ostream &os)
Add newline and flush stream.
Base class for collisional damping models.
This function object reads fields from the time directories and adds them to the mesh database for fu...
void evolve()
Evolve the cloud.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
void motion(TrackData &td)
Particle motion.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const List< DynamicList< molecule * > > & cellOccupancy
MPPICCloud(const MPPICCloud &)
Disallow default bitwise copy construct.
DSMCCloud< dsmcParcel > CloudType
Base class for packing models.
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
errorManipArg< error, int > exit(error &err, const int errNo=1)
This function object calculates the forces and moments by integrating the pressure and skin-friction ...
conserve internalField()+
void restoreState()
Reset the current cloud to the previously stored state.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
void readFields(const boolList &haveMesh, const fvMesh &mesh, const autoPtr< fvMeshSubset > &subsetterPtr, IOobjectList &allObjects, PtrList< GeoField > &fields)
const word cloudName(propsDict.lookup("cloudName"))
const dimensionedScalar c
Speed of light in a vacuum.
Base class for collisional return-to-isotropy models.
Type gMin(const FieldField< Field, Type > &f)
Generic GeometricField class.
Adds MPPIC modelling to kinematic clouds.
virtual ~MPPICCloud()
Destructor.
Type gMax(const FieldField< Field, Type > &f)
word name(const complex &)
Return a string representation of a complex.