Loading...
Searching...
No Matches
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-2015 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "isothermalDiameter.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace diameterModels
36{
38
40 (
41 diameterModel,
44 );
45}
46}
47
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
52(
53 const dictionary& diameterProperties,
54 const phaseModel& phase
55)
56:
57 diameterModel(diameterProperties, phase),
60{}
61
62
63// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
64
66{}
67
68
69// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70
72{
73 const volScalarField& p = phase_.U().db().lookupObject<volScalarField>
74 (
75 "p"
76 );
77
78 return d0_*cbrt(p0_/p);
79}
80
81
82bool Foam::diameterModels::isothermal::read(const dictionary& phaseProperties)
83{
84 diameterModel::read(phaseProperties);
85
86 diameterProperties_.readEntry("d0", d0_);
87 diameterProperties_.readEntry("p0", p0_);
88
89 return true;
90}
91
92
93// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
virtual bool read(const dictionary &phaseProperties)=0
Read phaseProperties dictionary.
dictionary diameterProperties_
Isothermal dispersed-phase particle diameter model.
isothermal(const dictionary &diameterProperties, const phaseModel &phase)
Construct from components.
virtual tmp< volScalarField > d() const
Return the diameter field.
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
virtual ~isothermal()=default
Destructor.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Helper class to manage multi-specie phase properties.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
Namespace for OpenFOAM.
const dimensionSet dimPressure
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar cbrt(const dimensionedScalar &ds)