All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-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 "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 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 
46 {
47  return d_;
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const dictionary& diameterProperties,
56  const phaseModel& phase
57 )
58 :
59  spherical(diameterProperties, phase),
60  liquidPhaseName_(diameterProperties.lookup("liquidPhase")),
61  d2_("d2", dimLength, diameterProperties.lookupOrDefault("d2", 0.0015)),
62  Tsub2_
63  (
64  "Tsub2",
66  diameterProperties.lookupOrDefault("Tsub2", 0)
67  ),
68  d1_
69  (
70  "d1",
71  dimLength,
72  diameterProperties.lookupOrDefault("d1", 0.00015)
73  ),
74  Tsub1_
75  (
76  "Tsub1",
78  diameterProperties.lookupOrDefault("Tsub1", 13.5)
79  ),
80  d_(dRef())
81 {
82  d_ = d1_;
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
87 
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
95 {
96  // Lookup the fluid model
97  const phaseSystem& fluid =
98  refCast<const phaseSystem>
99  (
100  phase().mesh().lookupObject<phaseSystem>("phaseProperties")
101  );
102 
103  const phaseModel& liquid(fluid.phases()[liquidPhaseName_]);
104 
105  if (phase().mesh().foundObject<saturationModel>("saturationModel"))
106  {
107  const saturationModel& satModel =
108  phase().mesh().lookupObject<saturationModel>("saturationModel");
109 
110  const volScalarField Tsub
111  (
112  liquid.thermo().T() - satModel.Tsat(liquid.thermo().p())
113  );
114 
115  d_ = max
116  (
117  d1_,
118  min
119  (
120  d2_,
121  (d1_*(Tsub - Tsub2_) + d2_*(Tsub - Tsub1_))/(Tsub2_ - Tsub1_)
122  )
123  );
124  }
125 }
126 
127 
128 bool Foam::diameterModels::linearTsub::read(const dictionary& phaseProperties)
129 {
130  spherical::read(phaseProperties);
131 
132  diameterProperties().lookup("liquidPhase") >> liquidPhaseName_;
133  diameterProperties().lookup("d2") >> d2_;
134  diameterProperties().lookup("Tsub2") >> Tsub2_;
135  diameterProperties().lookup("d1") >> d1_;
136  diameterProperties().lookup("Tsub1") >> Tsub1_;
137 
138  return true;
139 }
140 
141 
142 // ************************************************************************* //
virtual tmp< volScalarField > calcD() const
Get the diameter field.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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:57
stressControl lookup("compactNormalStress") >> compactNormalStress
const phaseModel & phase() const
Return the phase.
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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
const dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812
virtual void correct()
Correct the diameter field.