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-2015 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 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(dragModel, 0);
39  defineRunTimeSelectionTable(dragModel, dictionary);
40 }
41 
42 const Foam::dimensionSet Foam::dragModel::dimK(1, -3, -1, 0, 0);
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const phasePair& pair,
50  const bool registerObject
51 )
52 :
53  regIOobject
54  (
55  IOobject
56  (
57  IOobject::groupName(typeName, pair.name()),
58  pair.phase1().mesh().time().timeName(),
59  pair.phase1().mesh(),
60  IOobject::NO_READ,
61  IOobject::NO_WRITE,
62  registerObject
63  )
64  ),
65  pair_(pair)
66 {}
67 
68 
70 (
71  const dictionary& dict,
72  const phasePair& pair,
73  const bool registerObject
74 )
75 :
76  regIOobject
77  (
78  IOobject
79  (
80  IOobject::groupName(typeName, pair.name()),
81  pair.phase1().mesh().time().timeName(),
82  pair.phase1().mesh(),
83  IOobject::NO_READ,
84  IOobject::NO_WRITE,
85  registerObject
86  )
87  ),
88  pair_(pair),
89  swarmCorrection_
90  (
91  swarmCorrection::New
92  (
93  dict.subDict("swarmCorrection"),
94  pair
95  )
96  )
97 {}
98 
99 
100 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
101 
104 (
105  const dictionary& dict,
106  const phasePair& pair
107 )
108 {
109  const word modelType(dict.get<word>("type"));
110 
111  Info<< "Selecting dragModel for "
112  << pair << ": " << modelType << endl;
113 
114  auto* ctorPtr = dictionaryConstructorTable(modelType);
115 
116  if (!ctorPtr)
117  {
119  (
120  dict,
121  "dragModel",
122  modelType,
123  *dictionaryConstructorTablePtr_
124  ) << exit(FatalIOError);
125  }
126 
127  return ctorPtr(dict, pair, true);
128 }
129 
130 
131 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
132 
134 {}
135 
136 
137 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138 
140 {
141  return
142  0.75
143  *CdRe()
144  *swarmCorrection_->Cs()
145  *pair_.continuous().rho()
146  *pair_.continuous().nu()
147  /sqr(pair_.dispersed().d());
148 }
149 
150 
152 {
153  return max(pair_.dispersed(), pair_.dispersed().residualAlpha())*Ki();
154 }
155 
156 
158 {
159  return
160  max
161  (
162  fvc::interpolate(pair_.dispersed()),
163  pair_.dispersed().residualAlpha()
164  )*fvc::interpolate(Ki());
165 }
166 
167 
168 bool Foam::dragModel::writeData(Ostream& os) const
169 {
170  return os.good();
171 }
172 
173 
174 // ************************************************************************* //
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
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:487
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.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
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 INVALID.
Definition: exprTraits.C:52
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
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:274
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
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:615
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 ...