virtualMassModel.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) 2014-2020 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"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(virtualMassModel, 0);
37  defineRunTimeSelectionTable(virtualMassModel, dictionary);
38 }
39 
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const dictionary& dict,
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 
68 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69 
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
77 {
78  return Cvm()*pair_.continuous().rho();
79 }
80 
81 
83 {
84  return pair_.dispersed()*Ki();
85 }
86 
87 
89 {
90  return
92 }
93 
94 
95 bool Foam::virtualMassModel::writeData(Ostream& os) const
96 {
97  return os.good();
98 }
99 
100 
101 // ************************************************************************* //
virtual ~virtualMassModel()
Destructor.
bool writeData(Ostream &os) const
Pure virtual writaData function.
dictionary dict
virtual tmp< surfaceScalarField > Kf() const
Return the virtual mass coefficient Kf.
virtualMassModel(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from a dictionary and a phase pair.
virtual tmp< volScalarField > rho() const =0
Return the density field.
virtual tmp< volScalarField > Ki() const
Return the phase-intensive virtual mass coefficient Ki.
static const dimensionSet dimK
Coefficient dimensions.
virtual tmp< volScalarField > Cvm() const
Virtual mass coefficient.
Dimension set for the base types.
Definition: dimensionSet.H:120
const phasePair & pair_
Phase pair.
dynamicFvMesh & mesh
virtual const phaseModel & continuous() const
Continuous phase.
word timeName
Definition: getTimeIndex.H:3
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
#define defineBlendedInterfacialModelTypeNameAndDebug(ModelType, DebugSwitch)
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > K() const
Return the virtual mass coefficient K.
Namespace for OpenFOAM.