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-2019 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 
40  (
41  diameterModel,
42  linearTsub,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const dictionary& diameterProperties,
54  const phaseModel& phase
55 )
56 :
57  diameterModel(diameterProperties, phase),
58  liquidPhaseName_(diameterProperties.lookup("liquidPhase")),
59  d2_("d2", dimLength, diameterProperties.lookupOrDefault("d2", 0.0015)),
60  Tsub2_
61  (
62  "Tsub2",
64  diameterProperties.lookupOrDefault("Tsub2", 0)
65  ),
66  d1_
67  (
68  "d1",
69  dimLength,
70  diameterProperties.lookupOrDefault("d1", 0.00015)
71  ),
72  Tsub1_
73  (
74  "Tsub1",
76  diameterProperties.lookupOrDefault("Tsub1", 13.5)
77  ),
78  d_
79  (
80  IOobject
81  (
82  IOobject::groupName("d", phase.name()),
83  phase_.time().timeName(),
84  phase_.mesh(),
85  IOobject::NO_READ,
86  IOobject::AUTO_WRITE
87  ),
88  phase_.mesh(),
89  d1_
90  )
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
95 
97 {}
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
103 {
104  return d_;
105 }
106 
107 
109 {
110  // Lookup the fluid model
111  const phaseSystem& fluid =
112  refCast<const phaseSystem>
113  (
114  phase_.mesh().lookupObject<phaseSystem>("phaseProperties")
115  );
116 
117  const phaseModel& liquid(fluid.phases()[liquidPhaseName_]);
118 
119  if (phase_.mesh().foundObject<saturationModel>("saturationModel"))
120  {
121  const saturationModel& satModel =
122  phase_.mesh().lookupObject<saturationModel>("saturationModel");
123 
124  const volScalarField Tsub
125  (
126  liquid.thermo().T() - satModel.Tsat(liquid.thermo().p())
127  );
128 
129  d_ = max
130  (
131  d1_,
132  min
133  (
134  d2_,
135  (d1_*(Tsub - Tsub2_) + d2_*(Tsub - Tsub1_))/(Tsub2_ - Tsub1_)
136  )
137  );
138  }
139 }
140 
141 
142 bool Foam::diameterModels::linearTsub::read(const dictionary& phaseProperties)
143 {
144  diameterModel::read(phaseProperties);
145  diameterProperties_.lookup("liquidPhase") >> liquidPhaseName_;
146  diameterProperties_.lookup("d2") >> d2_;
147  diameterProperties_.lookup("Tsub2") >> Tsub2_;
148  diameterProperties_.lookup("d1") >> d1_;
149  diameterProperties_.lookup("Tsub1") >> Tsub1_;
150 
151  return true;
152 }
153 
154 
155 // ************************************************************************* //
virtual tmp< volScalarField > d() const
Return the diameter field.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual bool read(const dictionary &phaseProperties)=0
Read phaseProperties dictionary.
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
word timeName
Definition: getTimeIndex.H:3
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual ~linearTsub()
Destructor.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
linearTsub(const dictionary &diameterProperties, const phaseModel &phase)
Construct from components.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
virtual void correct()
Correct the diameter field.