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-2019 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34namespace Foam
35{
36namespace diameterModels
37{
39
41 (
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 (
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 * * * * * * * * * * * * * //
84
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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Abstract base-class for dispersed-phase particle diameter models.
virtual bool read(const dictionary &phaseProperties)=0
Read phaseProperties dictionary.
const phaseModel & phase() const
Return the phase.
diameterModel(const dictionary &dict, const phaseModel &phase)
const phaseModel & phase_
dictionary diameterProperties_
const dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
Isothermal dispersed-phase particle diameter model.
virtual void correct()
Correct the diameter field.
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.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
Helper class to manage multi-specie phase properties.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
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
dynamicFvMesh & mesh
word timeName
Definition getTimeIndex.H:3
Namespace for OpenFOAM.
const dimensionSet dimPressure
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127