All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
buoyancyEnergy.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) 2015-2021 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 "buoyancyEnergy.H"
27 #include "fvMatrices.H"
28 #include "basicThermo.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
37  defineTypeNameAndDebug(buoyancyEnergy, 0);
38 
40  (
41  fvModel,
42  buoyancyEnergy,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::fv::buoyancyEnergy::readCoeffs()
52 {
53  phaseName_ = coeffs().lookupOrDefault<word>("phase", word::null);
54 
55  UName_ =
56  coeffs().lookupOrDefault<word>
57  (
58  "U",
59  IOobject::groupName("U", phaseName_)
60  );
61 }
62 
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
67 (
68  const word& name,
69  const word& modelType,
70  const dictionary& dict,
71  const fvMesh& mesh
72 )
73 :
74  fvModel(name, modelType, dict, mesh),
75  phaseName_(word::null),
76  UName_(word::null)
77 {
78  readCoeffs();
79 }
80 
81 
82 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83 
85 {
86  const basicThermo& thermo =
87  mesh().lookupObject<basicThermo>
88  (
90  );
91 
92  return wordList(1, thermo.he().name());
93 }
94 
95 
97 (
98  const volScalarField& rho,
99  fvMatrix<scalar>& eqn,
100  const word& fieldName
101 ) const
102 {
104  mesh().lookupObject<uniformDimensionedVectorField>("g");
105 
106  const volVectorField& U = mesh().lookupObject<volVectorField>(UName_);
107 
108  eqn += rho*(U&g);
109 }
110 
111 
113 (
114  const volScalarField& alpha,
115  const volScalarField& rho,
116  fvMatrix<scalar>& eqn,
117  const word& fieldName
118 ) const
119 {
121  mesh().lookupObject<uniformDimensionedVectorField>("g");
122 
123  const volVectorField& U = mesh().lookupObject<volVectorField>(UName_);
124 
125  eqn += alpha*rho*(U&g);
126 }
127 
128 
130 {
131  if (fvModel::read(dict))
132  {
133  readCoeffs();
134  return true;
135  }
136  else
137  {
138  return false;
139  }
140 }
141 
142 
143 // ************************************************************************* //
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
const word & name() const
Return name.
Definition: IOobject.H:303
fluidReactionThermo & thermo
Definition: createFields.H:28
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:77
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:158
Finite volume model abstract base class.
Definition: fvModel.H:55
Macros for easy insertion into run-time selection tables.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
static word groupName(Name name, const word &group)
static const word null
An empty word.
Definition: word.H:77
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
buoyancyEnergy(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
static const word dictName
Name of the thermophysical properties dictionary.
Definition: basicThermo.H:129
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit contribution to compressible energy equation.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A special matrix type and solver, designed for finite volume solutions of scalar equations.
virtual bool read(const dictionary &dict)
Read source dictionary.
const dimensionedVector & g
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Namespace for OpenFOAM.