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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2018 OpenFOAM Foundation
9  Copyright (C) 2019-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "dragModel.H"
30 #include "phasePair.H"
31 #include "swarmCorrection.H"
32 #include "surfaceInterpolate.H"
33 #include "BlendedInterfacialModel.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(dragModel, 0);
41  defineRunTimeSelectionTable(dragModel, dictionary);
42 }
43 
44 const Foam::dimensionSet Foam::dragModel::dimK(1, -3, -1, 0, 0);
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const phasePair& pair,
52  const bool registerObject
53 )
54 :
56  (
57  IOobject
58  (
59  IOobject::groupName(typeName, pair.name()),
60  pair.phase1().mesh().time().timeName(),
61  pair.phase1().mesh(),
62  IOobject::NO_READ,
63  IOobject::NO_WRITE,
64  registerObject
65  )
66  ),
67  pair_(pair)
68 {}
69 
70 
72 (
73  const dictionary& dict,
74  const phasePair& pair,
75  const bool registerObject
76 )
77 :
79  (
80  IOobject
81  (
82  IOobject::groupName(typeName, pair.name()),
83  pair.phase1().mesh().time().timeName(),
84  pair.phase1().mesh(),
85  IOobject::NO_READ,
86  IOobject::NO_WRITE,
87  registerObject
88  )
89  ),
90  pair_(pair),
91  swarmCorrection_
92  (
94  (
95  dict.subDict("swarmCorrection"),
96  pair
97  )
98  )
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
103 
106 (
107  const dictionary& dict,
108  const phasePair& pair
109 )
110 {
111  const word modelType(dict.get<word>("type"));
112 
113  Info<< "Selecting dragModel for "
114  << pair << ": " << modelType << endl;
115 
116  auto* ctorPtr = dictionaryConstructorTable(modelType);
117 
118  if (!ctorPtr)
119  {
121  (
122  dict,
123  "dragModel",
124  modelType,
125  *dictionaryConstructorTablePtr_
126  ) << abort(FatalIOError);
127  }
129  return ctorPtr(dict, pair, true);
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
136 {}
137 
138 
139 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140 
142 {
143  return
144  0.75
145  *CdRe()
146  *swarmCorrection_->Cs()
147  *pair_.continuous().rho()
148  *pair_.continuous().nu()
149  /sqr(pair_.dispersed().d());
150 }
151 
154 {
155  return max(pair_.dispersed(), pair_.dispersed().residualAlpha())*Ki();
156 }
157 
158 
160 {
161  return
162  max
163  (
164  fvc::interpolate(pair_.dispersed()),
165  pair_.dispersed().residualAlpha()
166  )*fvc::interpolate(Ki());
167 }
168 
169 
171 {
172  return os.good();
173 }
174 
175 
176 // ************************************************************************* //
dictionary dict
defineBlendedInterfacialModelTypeNameAndDebug(massTransferModel, 0)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
static const dimensionSet dimK
Coefficient dimensions.
Definition: dragModel.H:100
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
word timeName
Definition: getTimeIndex.H:3
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
bool writeData(Ostream &os) const
Dummy write for regIOobject.
Definition: dragModel.C:163
virtual tmp< volScalarField > K() const
Return the drag coefficient K.
Definition: dragModel.C:146
dragModel(const phasePair &pair, const bool registerObject)
Definition: dragModel.C:43
errorManip< error > abort(error &err)
Definition: errorManip.H:139
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
phaseModel & phase1
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
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.
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:281
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:66
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual ~dragModel()
Destructor.
Definition: dragModel.C:128
virtual tmp< surfaceScalarField > Kf() const
Return the drag coefficient Kf.
Definition: dragModel.C:152
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static autoPtr< dragModel > New(const dictionary &dict, const phasePair &pair)
Definition: dragModel.C:99
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:635
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:49
virtual tmp< volScalarField > Ki() const
Return the phase-intensive drag coefficient Ki.
Definition: dragModel.C:134
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...