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 | Copyright (C) 2014-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 "virtualMassModel.H"
27 #include "phasePair.H"
28 #include "surfaceInterpolate.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  defineTypeNameAndDebug(virtualMassModel, 0);
35  defineRunTimeSelectionTable(virtualMassModel, dictionary);
36 }
37 
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& dict,
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 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
75 {
76  return Cvm()*pair_.continuous().rho();
77 }
78 
79 
81 {
82  return pair_.dispersed()*Ki();
83 }
84 
85 
87 {
88  return
90 }
91 
92 
93 bool Foam::virtualMassModel::writeData(Ostream& os) const
94 {
95  return os.good();
96 }
97 
98 
99 // ************************************************************************* //
virtual tmp< volScalarField > Ki() const
Return the phase-intensive virtual mass coefficient Ki.
virtual ~virtualMassModel()
Destructor.
dictionary dict
bool writeData(Ostream &os) const
Pure virtual writaData function.
virtualMassModel(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from a dictionary and a phase pair.
virtual const phaseModel & continuous() const
Continuous phase.
static const dimensionSet dimK
Coefficient dimensions.
virtual tmp< volScalarField > K() const
Return the virtual mass coefficient K.
Dimension set for the base types.
Definition: dimensionSet.H:118
dynamicFvMesh & mesh
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
word timeName
Definition: getTimeIndex.H:3
phaseModel & phase1
virtual const phaseModel & dispersed() const
Dispersed phase.
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.
const dimensionSet dimDensity
virtual tmp< surfaceScalarField > Kf() const
Return the virtual mass coefficient Kf.
virtual tmp< volScalarField > Cvm() const
Virtual mass coefficient.
A class for managing temporary objects.
Definition: PtrList.H:54
const phasePair & pair_
Phase pair.
Namespace for OpenFOAM.