fixedInterfacialAreaDiameter.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) 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 
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  diameterModel(diameterProperties, phase),
50  AvbyAlpha_(inv(dimLength), -1),
51  AvbyAlphaFieldPtr_()
52 {
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58 
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
64 
66 const
67 {
68  if (AvbyAlphaFieldPtr_.valid())
69  {
70  const volScalarField& AvbyAlpha = AvbyAlphaFieldPtr_;
71  return 6/AvbyAlpha;
72  }
73 
74  return volScalarField::New
75  (
76  IOobject::groupName("d", phase().name()),
77  phase().mesh(),
78  6/AvbyAlpha_
79  );
80 }
81 
82 
84 const
85 {
86  if (AvbyAlphaFieldPtr_.valid())
87  {
88  const volScalarField& AvbyAlpha = AvbyAlphaFieldPtr_;
89  return phase()*AvbyAlpha;
90  }
91 
92  return phase()*AvbyAlpha_;
93 }
94 
95 
98 {
100 
101  AvbyAlpha_ = dimensionedScalar
102  (
103  inv(dimLength),
104  diameterProperties().lookupOrDefault<scalar>("AvbyAlpha", -1)
105  );
106 
107  if (AvbyAlpha_.value() < 0 && !AvbyAlphaFieldPtr_.valid())
108  {
109  Info<< "fixedInterfacialArea: Uniform AvbyAlpha not provided. "
110  << "Looking up " << IOobject::groupName("AvbyAlpha", phase().name())
111  << endl;
112 
113  AvbyAlphaFieldPtr_ =
114  new volScalarField
115  (
116  IOobject
117  (
119  (
120  "AvbyAlpha",
121  phase().name()
122  ),
123  phase().mesh().time().name(),
124  phase().mesh(),
127  ),
128  phase().mesh()
129  );
130 
131  AvbyAlphaFieldPtr_->max
132  (
133  phaseProperties.lookupOrDefault<scalar>
134  (
135  "residualAvbyAlpha",
136  rootSmall
137  )
138  );
139  }
140 
141  return true;
142 }
143 
144 
145 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:52
const dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
Definition: diameterModel.C:62
fixedInterfacialArea dispersed-phase diameter model. The interfacial are is set by providing phase su...
fixedInterfacialArea(const dictionary &diameterProperties, const phaseModel &phase)
Construct from dictionary and phase.
virtual bool read(const dictionary &diameterProperties)
Read diameterProperties dictionary.
virtual tmp< volScalarField > d() const
Get the diameter field.
virtual tmp< volScalarField > Av() const
Get the surface area per unit volume field.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
const dimensionSet dimLength
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47