isothermalDiameter.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) 2011-2023 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 "isothermalDiameter.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace diameterModels
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& diameterProperties,
46  const phaseModel& phase
47 )
48 :
49  spherical(diameterProperties, phase),
50  d0_("d0", dimLength, diameterProperties),
51  p0_("p0", dimPressure, diameterProperties),
52  d_
53  (
54  IOobject
55  (
56  IOobject::groupName("d", phase.name()),
57  phase.time().name(),
58  phase.mesh()
59  ),
60  phase.mesh(),
61  d0_
62  )
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
75 {
76  return d_;
77 }
78 
79 
81 {
82  const volScalarField& p = phase().db().lookupObject<volScalarField>("p");
83 
84  d_ = d0_*pow(p0_/p, 1.0/3.0);
85 }
86 
87 
89 {
91 
92  diameterProperties().lookup("d0") >> d0_;
93  diameterProperties().lookup("p0") >> p0_;
94 
95  return true;
96 }
97 
98 
99 // ************************************************************************* //
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
Isothermal dispersed-phase particle diameter model.
virtual void correct()
Correct the model.
isothermal(const dictionary &diameterProperties, const phaseModel &phase)
Construct from dictionary and phase.
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 keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Helper class to manage multi-specie phase properties.
A class for managing temporary objects.
Definition: tmp.H:55
addToRunTimeSelectionTable(diameterModel, constant, dictionary)
defineTypeNameAndDebug(constant, 0)
Namespace for OpenFOAM.
const dimensionSet dimPressure
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimLength
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
volScalarField & p