Go to the documentation of this file.
35 template<
class CloudType>
38 const scalar a = 0.147;
40 scalar
h =
log(1.0 -
y*
y)/a;
54 template<
class CloudType>
75 <<
"Turbulence model not found in mesh database" <<
nl
76 <<
"Database objects include: " << obr.
sortedToc()
86 template<
class CloudType>
103 template<
class CloudType>
120 template<
class CloudType>
127 template<
class CloudType>
142 kPtr_ = tk.operator->();
158 template<
class CloudType>
170 const scalar dp =
p.d();
171 const scalar Tc =
p.Tc();
173 const scalar eta = rndGen_.sample01<scalar>();
174 const scalar
alpha = 2.0*lambda_/dp;
182 const label cellI =
p.cell();
184 const scalar kc =
k[cellI];
190 const scalar rhoRatio =
p.rho()/
p.rhoc();
196 const scalar sqrt2 =
sqrt(2.0);
197 for (
label i = 0; i < 3; i++)
199 const scalar
x = rndGen_.sample01<scalar>();
200 const scalar eta = sqrt2*erfInv(2*
x - 1.0);
201 value.
Su()[i] = mass*
f*eta;
static word groupName(Name name, const word &group)
A class for handling words, derived from string.
BrownianMotionForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
A class for managing temporary objects.
Template functions to aid in the implementation of demand driven data.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual ~BrownianMotionForce()
Destructor.
scalarField Re(const UList< complex > &cf)
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual void cacheFields(const bool store)
Cache fields.
const char *const group
Group name for atomic constants.
scalar erfInv(const scalar y) const
Inverse erf for Gaussian distribution.
const Type & value() const
Return const reference to value.
dimensionedScalar exp(const dimensionedScalar &ds)
Random & rndGen()
Return refernce to the random object.
bool isTmp() const
Return true if this is really a temporary object.
Registry of regIOobjects.
void deleteDemandDrivenData(DataPtr &dataPtr)
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 dimensionedScalar h
Planck constant.
const fvMesh & mesh() const
Return refernce to the mesh.
Helper container for force Su and Sp terms.
Abstract base class for particle forces.
Templated base class for dsmc cloud.
Abstract base class for turbulence models (RAS, LES and laminar).
A list of keyword definitions, which are a keyword followed by any number of values (e....
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Mesh data needed to do the Finite Volume discretisation.
T * ptr() const
Return tmp pointer for reuse.
errorManip< error > abort(error &err)
List< Key > sortedToc() const
Return the table of contents as a sorted list.
Calculates particle Brownian motion force.
virtual forceSuSp calcCoupled(const typename CloudType::parcelType &p, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the non-coupled force.
bool foundObject(const word &name) const
Is the named Type found?
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const scalar lambda_
Molecular free path length [m].
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
ParcelType parcelType
Type of parcel the cloud was instantiated for.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Generic GeometricField class.
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
bool turbulence_
Turbulence flag.
cachedRandom & rndGen_
Reference to the cloud random number generator.
tmp< volScalarField > kModel() const
Return the k field from the turbulence model.
stressControl lookup("compactNormalStress") >> compactNormalStress