virtualMassModel.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) 2014-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 "virtualMassModel.H"
30 #include "phasePair.H"
31 #include "surfaceInterpolate.H"
32 #include "BlendedInterfacialModel.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(virtualMassModel, 0);
40  defineRunTimeSelectionTable(virtualMassModel, dictionary);
41 }
42 
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const dictionary& dict,
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 
71 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
72 
75 (
76  const dictionary& dict,
77  const phasePair& pair
78 )
79 {
80  const word modelType(dict.get<word>("type"));
81 
82  Info<< "Selecting virtualMassModel for "
83  << pair << ": " << modelType << endl;
84 
85  auto* ctorPtr = dictionaryConstructorTable(modelType);
86 
87  if (!ctorPtr)
88  {
90  (
91  dict,
92  "virtualMassModel",
93  modelType,
94  *dictionaryConstructorTablePtr_
95  ) << abort(FatalIOError);
96  }
97 
98  return ctorPtr(dict, pair, true);
99 }
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 {
106  return Cvm()*pair_.continuous().rho();
107 }
108 
111 {
112  return pair_.dispersed()*Ki();
113 }
114 
115 
117 {
118  return
119  fvc::interpolate(pair_.dispersed())*fvc::interpolate(Ki());
120 }
121 
122 
124 {
125  return os.good();
126 }
127 
128 
129 // ************************************************************************* //
bool writeData(Ostream &os) const
Pure virtual writeData function.
dictionary dict
defineBlendedInterfacialModelTypeNameAndDebug(massTransferModel, 0)
virtualMassModel(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from a dictionary and a phase pair.
virtual tmp< volScalarField > K() const
Return the virtual mass coefficient K.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static const dimensionSet dimK
Coefficient dimensions.
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
virtual tmp< surfaceScalarField > Kf() const
Return the virtual mass coefficient Kf.
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.
const dimensionSet dimDensity
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)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
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
static autoPtr< virtualMassModel > New(const dictionary &dict, const phasePair &pair)
virtual tmp< volScalarField > Ki() const
Return the phase-intensive virtual mass coefficient Ki.
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...