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-2025 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"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& diameterProperties,
47  const phaseModel& phase
48 )
49 :
50  spherical(diameterProperties, phase),
51  liquidPhaseName_(diameterProperties.lookup("liquidPhase")),
52  d2_("d2", dimLength, diameterProperties.lookupOrDefault("d2", 0.0015)),
53  Tsub2_
54  (
55  "Tsub2",
57  diameterProperties.lookupOrDefault("Tsub2", 0)
58  ),
59  d1_("d1", dimLength, diameterProperties.lookupOrDefault("d1", 0.00015)),
60  Tsub1_
61  (
62  "Tsub1",
64  diameterProperties.lookupOrDefault("Tsub1", 13.5)
65  ),
66  saturationModelPtr_
67  (
69  (
70  "saturationTemperature",
71  diameterProperties
72  ).ptr()
73  ),
74  d_
75  (
76  IOobject
77  (
78  IOobject::groupName("d", phase.name()),
79  phase.time().name(),
80  phase.mesh()
81  ),
82  phase.mesh(),
83  d1_
84  )
85 {
86  Info<< " d2: " << d2_.value() << endl
87  << " Tsub2: " << Tsub2_.value() << endl
88  << " d1: " << d1_.value() << endl
89  << " Tsub1: " << Tsub1_.value() << endl;
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
94 
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
102 {
103  return d_;
104 }
105 
106 
108 {
109  const phaseModel& liquid = phase().fluid().phases()[liquidPhaseName_];
110 
111  const volScalarField Tsub
112  (
113  saturationModelPtr_->Tsat(liquid.fluidThermo().p())
114  - liquid.thermo().T()
115  );
116 
117  d_ = max
118  (
119  d1_,
120  min
121  (
122  d2_,
123  (d1_*(Tsub - Tsub2_) + d2_*(Tsub - Tsub1_))/(Tsub2_ - Tsub1_)
124  )
125  );
126 }
127 
128 
130 {
132 
133  diameterProperties().lookup("liquidPhase") >> liquidPhaseName_;
134 
135  d2_.readIfPresent(diameterProperties());
136  Tsub2_.readIfPresent(diameterProperties());
137  d1_.readIfPresent(diameterProperties());
138  Tsub1_.readIfPresent(diameterProperties());
139 
140  saturationModelPtr_.reset
141  (
143  (
144  "saturationTemperature",
145  diameterProperties()
146  ).ptr()
147  );
148 
149  return true;
150 }
151 
152 
153 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:52
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
Definition: diameterModel.C:62
linearTsub(const dictionary &diameterProperties, const phaseModel &phase)
Construct from dictionary and phase.
virtual void correct()
Correct the model.
virtual tmp< volScalarField > d() const
Get the diameter field.
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
Base class for models which represent spherical diameter models, providing a common implementation of...
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:53
Helper class to manage multi-specie phase properties.
Model to describe the dependence of saturation temperature on pressure.
static autoPtr< saturationTemperatureModel > New(const dictionary &dict)
Select with dictionary.
A class for managing temporary objects.
Definition: tmp.H:55
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
addToRunTimeSelectionTable(diameterModel, constant, dictionary)
defineTypeNameAndDebug(constant, 0)
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
const dimensionSet dimLength
const dimensionSet dimTemperature
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.