dragModel.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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "dragModel.H"
27 #include "phasePair.H"
28 #include "swarmCorrection.H"
29 #include "surfaceInterpolate.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(dragModel, 0);
36  defineRunTimeSelectionTable(dragModel, dictionary);
37 }
38 
39 const Foam::dimensionSet Foam::dragModel::dimK(1, -3, -1, 0, 0);
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const phasePair& pair,
47  const bool registerObject
48 )
49 :
50  regIOobject
51  (
52  IOobject
53  (
54  IOobject::groupName(typeName, pair.name()),
55  pair.phase1().mesh().time().timeName(),
56  pair.phase1().mesh(),
57  IOobject::NO_READ,
58  IOobject::NO_WRITE,
59  registerObject
60  )
61  ),
62  pair_(pair)
63 {}
64 
65 
67 (
68  const dictionary& dict,
69  const phasePair& pair,
70  const bool registerObject
71 )
72 :
73  regIOobject
74  (
75  IOobject
76  (
77  IOobject::groupName(typeName, pair.name()),
78  pair.phase1().mesh().time().timeName(),
79  pair.phase1().mesh(),
80  IOobject::NO_READ,
81  IOobject::NO_WRITE,
82  registerObject
83  )
84  ),
85  pair_(pair),
86  swarmCorrection_
87  (
88  swarmCorrection::New
89  (
90  dict.subDict("swarmCorrection"),
91  pair
92  )
93  )
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
98 
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
106 {
107  return
108  0.75
109  *CdRe()
110  *swarmCorrection_->Cs()
111  *pair_.continuous().rho()
112  *pair_.continuous().nu()
113  /sqr(pair_.dispersed().d());
114 }
115 
116 
118 {
119  return max(pair_.dispersed(), pair_.dispersed().residualAlpha())*Ki();
120 }
121 
122 
124 {
125  return
126  max
127  (
130  )*fvc::interpolate(Ki());
131 }
132 
133 
134 bool Foam::dragModel::writeData(Ostream& os) const
135 {
136  return os.good();
137 }
138 
139 
140 // ************************************************************************* //
dictionary dict
bool writeData(Ostream &os) const
Dummy write for regIOobject.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static const dimensionSet dimK
Coefficient dimensions.
Definition: dragModel.H:93
virtual const phaseModel & continuous() const
Continuous phase.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
const phasePair & pair_
Phase pair.
Definition: dragModel.H:62
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
virtual ~dragModel()
Destructor.
Dimension set for the base types.
Definition: dimensionSet.H:118
dynamicFvMesh & mesh
virtual tmp< volScalarField > K() const
Return the drag coefficient K.
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
word timeName
Definition: getTimeIndex.H:3
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
phaseModel & phase1
virtual const phaseModel & dispersed() const
Dispersed phase.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
tmp< volScalarField > d() const
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
dragModel(const dictionary &interfaceDict, const phaseModel &phase1, const phaseModel &phase2)
A class for managing temporary objects.
Definition: PtrList.H:54
autoPtr< swarmCorrection > swarmCorrection_
Swarm correction.
Definition: dragModel.H:65
virtual tmp< volScalarField > CdRe() const =0
Drag coefficient.
virtual tmp< volScalarField > Ki() const
Return the phase-intensive drag coefficient Ki.
virtual tmp< surfaceScalarField > Kf() const
Return the drag coefficient Kf.
Namespace for OpenFOAM.