Loading...
Searching...
No Matches
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-2015 OpenFOAM Foundation
9 Copyright (C) 2020-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
39}
40
42
43
44// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45
47(
48 const dictionary& dict,
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
69// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
70
73(
74 const dictionary& dict,
75 const phasePair& pair
76)
77{
78 const word modelType(dict.get<word>("type"));
79
80 Info<< "Selecting virtualMassModel for "
81 << pair << ": " << modelType << endl;
82
83 auto* ctorPtr = dictionaryConstructorTable(modelType);
84
85 if (!ctorPtr)
86 {
88 (
89 dict,
90 "virtualMassModel",
91 modelType,
92 *dictionaryConstructorTablePtr_
93 ) << exit(FatalIOError);
94 }
95
96 return ctorPtr(dict, pair, true);
97}
98
99
100// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101
102Foam::tmp<Foam::volScalarField> Foam::virtualMassModel::Ki() const
103{
104 return Cvm()*pair_.continuous().rho();
105}
106
107
108Foam::tmp<Foam::volScalarField> Foam::virtualMassModel::K() const
109{
110 return pair_.dispersed()*Ki();
111}
112
113
114Foam::tmp<Foam::surfaceScalarField> Foam::virtualMassModel::Kf() const
115{
116 return
117 fvc::interpolate(pair_.dispersed())*fvc::interpolate(Ki());
118}
119
120
122{
123 return os.good();
124}
125
126
127// ************************************************************************* //
phaseModel & phase1
bool good() const noexcept
True if next operation might succeed.
Definition IOstream.H:281
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
virtual tmp< surfaceScalarField > Kf() const
Return the virtual mass coefficient Kf.
static const dimensionSet dimK
Coefficient dimensions.
bool writeData(Ostream &os) const
Pure virtual writeData function.
static autoPtr< virtualMassModel > New(const dictionary &dict, const phasePair &pair)
virtual tmp< volScalarField > Ki() const
Return the phase-intensive virtual mass coefficient Ki.
virtual tmp< volScalarField > K() const
Return the virtual mass coefficient K.
virtualMassModel(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from a dictionary and a phase pair.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
OBJstream os(runTime.globalPath()/outputName)
auto & name
word timeName
Definition getTimeIndex.H:3
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.
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict