isothermalDiameter.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2019 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "isothermalDiameter.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace diameterModels
37 {
38  defineTypeNameAndDebug(isothermal, 0);
39 
41  (
42  diameterModel,
43  isothermal,
44  dictionary
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const dictionary& diameterProperties,
55  const phaseModel& phase
56 )
57 :
58  diameterModel(diameterProperties, phase),
59  d0_("d0", dimLength, diameterProperties_),
60  p0_("p0", dimPressure, diameterProperties_),
61  d_
62  (
63  IOobject
64  (
65  IOobject::groupName("d", phase.name()),
66  phase_.time().timeName(),
67  phase_.mesh(),
68  IOobject::NO_READ,
69  IOobject::AUTO_WRITE
70  ),
71  phase_.mesh(),
72  d0_
73  )
74 {}
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
80 {
81  return d_;
82 }
83 
84 
86 {
87  const volScalarField& p = phase_.db().lookupObject<volScalarField>("p");
88 
89  d_ = d0_*pow(p0_/p, 1.0/3.0);
90 }
91 
92 
94 {
96 
97  diameterProperties_.readEntry("d0", d0_);
98  diameterProperties_.readEntry("p0", p0_);
99 
100  return true;
101 }
102 
103 
104 // ************************************************************************* //
Foam::dimPressure
const dimensionSet dimPressure
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::diameterModels::isothermal::read
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
Definition: isothermalDiameter.C:93
Foam::phaseProperties
Helper class to manage multi-specie phase properties.
Definition: phaseProperties.H:60
Foam::diameterModels::defineTypeNameAndDebug
defineTypeNameAndDebug(constant, 0)
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::diameterModels::isothermal::correct
virtual void correct()
Correct the diameter field.
Definition: isothermalDiameter.C:85
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
isothermalDiameter.H
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::diameterModel::read
virtual bool read(const dictionary &phaseProperties)=0
Read phaseProperties dictionary.
Definition: diameterModel.C:96
Foam::diameterModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(diameterModel, constant, dictionary)
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::diameterModels::isothermal::isothermal
isothermal(const dictionary &dict, const phaseModel &phase)
Construct from components.
Definition: isothermalDiameter.C:52
Foam::diameterModels::isothermal::d
tmp< volScalarField > d() const
Return the phase mean diameter field.
Definition: isothermalDiameter.C:65