Go to the documentation of this file.
34 template<
class CloudType>
50 template<
class CloudType>
68 Info<<
"Constructing particle forces" <<
endl;
75 const word& modelName = dEntry.keyword();
115 template<
class CloudType>
130 template<
class CloudType>
137 template<
class CloudType>
142 this->operator[](i).cacheFields(store);
147 template<
class CloudType>
151 const typename CloudType::parcelType::trackingData& td,
164 value += this->operator[](i).calcCoupled(
p, td, dt, mass,
Re, muc);
172 template<
class CloudType>
176 const typename CloudType::parcelType::trackingData& td,
190 this->operator[](i).calcNonCoupled(
p, td, dt, mass,
Re, muc);
198 template<
class CloudType>
202 const typename CloudType::parcelType::trackingData& td,
206 scalar massEff = mass;
209 massEff += this->operator[](i).massAdd(
p, td, mass);
A keyword and a list of tokens is an 'entry'.
void set(List< bool > &bools, const labelRange &range)
Set the specified range 'on' in a boolList.
virtual void cacheFields(const bool store)
Cache fields.
virtual forceSuSp calcNonCoupled(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the non-coupled force.
A class for handling words, derived from Foam::string.
static constexpr const zero Zero
Global zero (0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual ~ParticleForceList()
Destructor.
#define forAll(list, i)
Loop across all elements in list.
patchWriters resize(patchIds.size())
messageStream Info
Information stream (stdout output on master, null elsewhere)
Helper container for force Su and Sp terms.
Abstract base class for particle forces.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Templated base class for dsmc cloud.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
ParticleForceList(CloudType &owner, const fvMesh &mesh)
Null constructor.
Mesh data needed to do the Finite Volume discretisation.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
virtual forceSuSp calcCoupled(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the coupled force.
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
virtual scalar massEff(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar mass) const
Return the effective mass.
scalarField Re(const UList< complex > &cf)
Extract real component.
ParcelType parcelType
Type of parcel the cloud was instantiated for.