39template<
class CloudType>
41Foam::BrownianMotionForce<CloudType>::kModel()
const
59 <<
"Turbulence model not found in mesh database" <<
nl
60 <<
"Database objects include: " << obr.
sortedToc()
69template<
class CloudType>
79 lambda_(this->coeffs().getScalar(
"lambda")),
81 turbulence_(this->
coeffs().getBool(
"turbulence")),
83 useSpherical_(this->
coeffs().getOrDefault(
"spherical", true))
87template<
class CloudType>
97 turbulence_(bmf.turbulence_),
99 useSpherical_(bmf.useSpherical_)
105template<
class CloudType>
118template<
class CloudType>
150template<
class CloudType>
154 const typename CloudType::parcelType::trackingData&
td,
163 const scalar dp =
p.d();
164 const scalar Tc =
td.Tc();
166 const scalar
alpha = 2.0*lambda_/dp;
175 const label celli =
p.cell();
177 const scalar kc =
k[celli];
188 Random& rnd = this->owner().rndGen();
193 const scalar theta = rnd.
sample01<scalar>()*twoPi;
194 const scalar u = 2*rnd.
sample01<scalar>() - 1;
196 const scalar a =
sqrt(1 -
sqr(u));
compressible::turbulenceModel & turb
Calculates particle Brownian motion force.
BrownianMotionForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
virtual void cacheFields(const bool store)
Cache fields.
virtual ~BrownianMotionForce()
Destructor.
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.
Foam::List< Key > sortedToc(const Compare &comp) const
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Abstract base class for particle forces.
const fvMesh & mesh() const noexcept
Return the mesh database.
const CloudType & owner() const noexcept
Return const access to the cloud owner.
ParticleForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict, const word &forceType, const bool readCoeffs)
Construct from mesh.
const dictionary & coeffs() const noexcept
Return the force coefficients dictionary.
Type GaussNormal()
Return a sample whose components are normally distributed with zero mean and unity variance N(0,...
Type sample01()
Return a sample whose components lie in the range [0,1].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Helper container for force Su and Sp terms.
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Mesh data needed to do the Finite Volume discretisation.
Registry of regIOobjects.
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
A class for managing temporary objects.
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
bool movable() const noexcept
True if this is a non-null managed pointer with a unique ref-count.
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Fundamental dimensioned constants.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
constexpr const char *const group
Group name for atomic constants.
constexpr scalar pi(M_PI)
const dimensionedScalar k
Boltzmann constant.
Different types of constants.
dimensionedScalar pow5(const dimensionedScalar &ds)
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar sin(const dimensionedScalar &ds)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
static constexpr const zero Zero
Global zero (0).
scalarField Re(const UList< complex > &cmplx)
Extract real component.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void deleteDemandDrivenData(DataPtr &dataPtr)
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).