constantDiameter.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 "constantDiameter.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  d_("d", dimLength, diameterProperties)
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
55 
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
63 {
64  return volScalarField::New
65  (
66  IOobject::groupName("d", phase().name()),
67  phase().mesh(),
68  d_
69  );
70 }
71 
72 
74 {
76 
77  diameterProperties().lookup("d") >> d_;
78 
79  return true;
80 }
81 
82 
83 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
static word groupName(Name name, const word &group)
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
Constant dispersed-phase particle diameter model.
constant(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 ~constant()
Destructor.
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: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.
const dimensionSet dimLength
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47