linearTsubDiameter.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) 2018-2022 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 "linearTsubDiameter.H"
27 #include "phaseSystem.H"
28 #include "saturationModel.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37  defineTypeNameAndDebug(linearTsub, 0);
38  addToRunTimeSelectionTable(diameterModel, linearTsub, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const dictionary& diameterProperties,
48  const phaseModel& phase
49 )
50 :
51  spherical(diameterProperties, phase),
52  liquidPhaseName_(diameterProperties.lookup("liquidPhase")),
53  d2_("d2", dimLength, diameterProperties.lookupOrDefault("d2", 0.0015)),
54  Tsub2_
55  (
56  "Tsub2",
58  diameterProperties.lookupOrDefault("Tsub2", 0)
59  ),
60  d1_("d1", dimLength, diameterProperties.lookupOrDefault("d1", 0.00015)),
61  Tsub1_
62  (
63  "Tsub1",
65  diameterProperties.lookupOrDefault("Tsub1", 13.5)
66  ),
67  d_
68  (
69  IOobject
70  (
71  IOobject::groupName("d", phase.name()),
72  phase.time().timeName(),
73  phase.mesh()
74  ),
75  phase.mesh(),
76  d1_
77  )
78 {
79  Info<< " d2: " << d2_.value() << endl
80  << " Tsub2: " << Tsub2_.value() << endl
81  << " d1: " << d1_.value() << endl
82  << " Tsub1: " << Tsub1_.value() << endl;
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
87 
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
95 {
96  return d_;
97 }
98 
99 
101 {
102  const phaseSystem& fluid = phase().fluid();
103  const phaseModel& liquid = fluid.phases()[liquidPhaseName_];
104  const phaseInterface interface(phase(), liquid);
105 
106  if (fluid.foundInterfacialModel<saturationModel>(interface))
107  {
108  const saturationModel& satModel =
109  fluid.lookupInterfacialModel<saturationModel>(interface);
110 
111  const volScalarField Tsub
112  (
113  satModel.Tsat(liquid.thermo().p()) - liquid.thermo().T()
114  );
115 
116  d_ = max
117  (
118  d1_,
119  min
120  (
121  d2_,
122  (d1_*(Tsub - Tsub2_) + d2_*(Tsub - Tsub1_))/(Tsub2_ - Tsub1_)
123  )
124  );
125  }
126 }
127 
128 
129 bool Foam::diameterModels::linearTsub::read(const dictionary& phaseProperties)
130 {
131  spherical::read(phaseProperties);
132 
133  diameterProperties().lookup("liquidPhase") >> liquidPhaseName_;
134 
135  d2_ =
137  (
138  dimLength,
139  diameterProperties().lookupOrDefault<scalar>("d2", 0.0015)
140  );
141  Tsub2_ =
143  (
145  diameterProperties().lookupOrDefault<scalar>("Tsub2", 0)
146  );
147  d1_ =
149  (
150  dimLength,
151  diameterProperties().lookupOrDefault("d1", 0.00015)
152  );
153  Tsub1_ =
155  (
157  diameterProperties().lookupOrDefault<scalar>("Tsub1", 13.5)
158  );
159 
160  return true;
161 }
162 
163 
164 // ************************************************************************* //
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual tmp< volScalarField > d() const
Get the diameter field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:35
stressControl lookup("compactNormalStress") >> compactNormalStress
const phaseModel & phase() const
Return the phase.
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
word timeName
Definition: getTimeIndex.H:3
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual ~linearTsub()
Destructor.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const phaseSystem & fluid() const
Return the system to which this phase belongs.
messageStream Info
linearTsub(const dictionary &diameterProperties, const phaseModel &phase)
Construct from dictionary and phase.
A class for managing temporary objects.
Definition: PtrList.H:53
const dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
const dimensionSet dimTemperature
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
virtual void correct()
Correct the model.