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-2018 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 {
35  defineTypeNameAndDebug(constant, 0);
36 
38  (
39  diameterModel,
40  constant,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const dictionary& diameterProperties,
52  const phaseModel& phase
53 )
54 :
55  diameterModel(diameterProperties, phase),
56  d_("d", dimLength, diameterProperties_)
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 {
70  return volScalarField::New
71  (
72  "d",
73  phase_.U().mesh(),
74  d_
75  );
76 }
77 
78 
79 bool Foam::diameterModels::constant::read(const dictionary& phaseProperties)
80 {
81  diameterModel::read(phaseProperties);
82 
83  diameterProperties_.lookup("d") >> d_;
84 
85  return true;
86 }
87 
88 
89 // ************************************************************************* //
tmp< volScalarField > d() const
Return the phase mean diameter field.
virtual bool read(const dictionary &phaseProperties)=0
Read phaseProperties dictionary.
virtual ~constant()
Destructor.
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
const volVectorField & U() const
Definition: phaseModel.H:170
virtual bool read(const dictionary &diameterProperties)
Read diameterProperties dictionary.
constant(const dictionary &dict, const phaseModel &phase)
Construct from components.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
dictionary diameterProperties_
Definition: diameterModel.H:58
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
const phaseModel & phase_
Definition: diameterModel.H:58
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583