dragModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 #include "BlendedInterfacialModel.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(dragModel, 0);
38  defineRunTimeSelectionTable(dragModel, dictionary);
39 }
40 
41 const Foam::dimensionSet Foam::dragModel::dimK(1, -3, -1, 0, 0);
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const phasePair& pair,
49  const bool registerObject
50 )
51 :
52  regIOobject
53  (
54  IOobject
55  (
56  IOobject::groupName(typeName, pair.name()),
57  pair.phase1().mesh().time().timeName(),
58  pair.phase1().mesh(),
59  IOobject::NO_READ,
60  IOobject::NO_WRITE,
61  registerObject
62  )
63  ),
64  pair_(pair)
65 {}
66 
67 
69 (
70  const dictionary& dict,
71  const phasePair& pair,
72  const bool registerObject
73 )
74 :
75  regIOobject
76  (
77  IOobject
78  (
79  IOobject::groupName(typeName, pair.name()),
80  pair.phase1().mesh().time().timeName(),
81  pair.phase1().mesh(),
82  IOobject::NO_READ,
83  IOobject::NO_WRITE,
84  registerObject
85  )
86  ),
87  pair_(pair),
88  swarmCorrection_
89  (
90  swarmCorrection::New
91  (
92  dict.subDict("swarmCorrection"),
93  pair
94  )
95  )
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
100 
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 {
109  return
110  0.75
111  *CdRe()
112  *swarmCorrection_->Cs()
113  *pair_.continuous().rho()
114  *pair_.continuous().nu()
115  /sqr(pair_.dispersed().d());
116 }
117 
118 
120 {
121  return max(pair_.dispersed(), pair_.dispersed().residualAlpha())*Ki();
122 }
123 
124 
126 {
127  return
128  max
129  (
132  )*fvc::interpolate(Ki());
133 }
134 
135 
136 bool Foam::dragModel::writeData(Ostream& os) const
137 {
138  return os.good();
139 }
140 
141 
142 // ************************************************************************* //
dictionary dict
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define defineBlendedInterfacialModelTypeNameAndDebug(ModelType, DebugSwitch)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static const dimensionSet dimK
Coefficient dimensions.
Definition: dragModel.H:93
tmp< volScalarField > d() const
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.
virtual tmp< surfaceScalarField > Kf() const
Return the drag coefficient Kf.
virtual tmp< volScalarField > K() const
Return the drag coefficient K.
Dimension set for the base types.
Definition: dimensionSet.H:120
dynamicFvMesh & mesh
virtual const phaseModel & continuous() const
Continuous phase.
bool writeData(Ostream &os) const
Dummy write for regIOobject.
word timeName
Definition: getTimeIndex.H:3
virtual tmp< volScalarField > Ki() const
Return the phase-intensive drag coefficient Ki.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
phaseModel & phase1
virtual const phaseModel & dispersed() const
Dispersed phase.
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
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)
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
A class for managing temporary objects.
Definition: PtrList.H:53
autoPtr< swarmCorrection > swarmCorrection_
Swarm correction.
Definition: dragModel.H:65
virtual tmp< volScalarField > CdRe() const =0
Drag coefficient.
Namespace for OpenFOAM.