Loading...
Searching...
No Matches
PlessisMasliyahDragForce.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
29#include "volFields.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class CloudType>
35(
37 const fvMesh& mesh,
38 const dictionary& dict
39)
40:
42 alphac_
43 (
44 this->mesh().template lookupObject<volScalarField>
45 (
46 this->coeffs().getWord("alphac")
47 )
48 )
49{}
50
51
52template<class CloudType>
54(
56)
57:
59 alphac_
60 (
61 this->mesh().template lookupObject<volScalarField>
62 (
63 this->coeffs().getWord("alphac")
64 )
65 )
66{}
67
68
69// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70
71template<class CloudType>
73(
74 const typename CloudType::parcelType& p,
75 const typename CloudType::parcelType::trackingData& td,
76 const scalar dt,
77 const scalar mass,
78 const scalar Re,
79 const scalar muc
80) const
81{
82 const scalar alphac = alphac_[p.cell()];
83
84 const scalar cbrtAlphap = cbrt(1.0 - alphac);
85
86 const scalar A =
87 26.8*pow3(alphac)
88 /(
89 sqr(cbrtAlphap)
90 *(1.0 - cbrtAlphap)
91 *sqr(1.0 - sqr(cbrtAlphap))
92 + SMALL
93 );
94
95 // (P:Eq. 36)
96 const scalar B = sqr(alphac)/sqr(1.0 - sqr(cbrtAlphap));
97
98 return forceSuSp
99 (
100 Zero,
101 (mass/p.rho())
102 *(A*(1.0 - alphac)/alphac + B*Re)*muc/(alphac*sqr(p.d()))
103 );
104}
105
106
107// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
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.
Particle-drag model wherein drag forces (per unit carrier-fluid velocity) are dynamically computed ba...
PlessisMasliyahDragForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
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.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Helper container for force Su and Sp terms.
Definition forceSuSp.H:63
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
volScalarField & p
dynamicFvMesh & mesh
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
DSMCCloud< dsmcParcel > CloudType
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
scalarField Re(const UList< complex > &cmplx)
Extract real component.
dimensionedScalar cbrt(const dimensionedScalar &ds)
dictionary dict