variable.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-2024 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 "variable.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
37 namespace heatTransferCoefficientModels
38 {
42 }
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::fv::heatTransferCoefficientModels::variable::readCoeffs
50 (
51  const dictionary& dict
52 )
53 {
54  UName_ = dict.lookupOrDefault<word>("U", "U");
55 
56  a_ = dict.lookup<scalar>("a");
57  b_ = dict.lookup<scalar>("b");
58  c_ = dict.lookup<scalar>("c");
59  L_ = dimensionedScalar("L", dimLength, dict);
60  Pr_ = dimensionedScalar("Pr", dimless, dict);
61 }
62 
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
67 (
68  const dictionary& dict,
69  const fvMesh& mesh
70 )
71 :
72  heatTransferCoefficientModel(typeName, dict, mesh),
73  UName_(word::null),
74  a_(NaN),
75  b_(NaN),
76  c_(NaN),
77  L_("L", dimLength, NaN),
78  Pr_("Pr", dimless, NaN),
79  htc_
80  (
81  IOobject
82  (
83  typedName("htc"),
84  mesh.time().name(),
85  mesh,
86  IOobject::NO_READ,
87  IOobject::NO_WRITE
88  ),
89  mesh,
91  zeroGradientFvPatchScalarField::typeName
92  )
93 {
94  readCoeffs(dict);
95 }
96 
97 
99 (
100  const dictionary& dict,
101  const interRegionModel& model
102 )
103 :
104  variable(dict, model.mesh())
105 {
106  readCoeffs(dict);
107 }
108 
109 
110 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
111 
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
121  mesh_.lookupType<fluidThermophysicalTransportModel>();
122 
124  ttm.momentumTransport();
125 
126  const volVectorField& U =
127  mesh_.lookupObject<volVectorField>(UName_);
128 
129  const volScalarField Re(mag(U)*L_/mtm.nuEff());
130  const volScalarField Nu(a_*pow(Re, b_)*pow(Pr_, c_));
131 
132  htc_ = Nu*ttm.kappaEff()/L_;
133  htc_.correctBoundaryConditions();
134 }
135 
136 
138 (
139  const dictionary& dict
140 )
141 {
143  {
144  readCoeffs(dict);
145  return true;
146  }
147  else
148  {
149  return false;
150  }
151 }
152 
153 
154 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Base class for single-phase compressible turbulence models.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for fluid thermophysical transport models RAS, LES and laminar.
const compressibleMomentumTransportModel & momentumTransport() const
Access function to momentum transport model.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Base class for heat transfer coefficient modelling used in heat transfer fvModels.
virtual bool read(const dictionary &dict)
Read dictionary.
Variable heat transfer model depending on local values. The Nusselt number is calculated as:
Definition: variable.H:92
virtual void correct()
Correct the heat transfer coefficient.
Definition: variable.C:118
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: variable.C:138
variable(const dictionary &dict, const fvMesh &mesh)
Construct from dictionary and mesh.
Definition: variable.C:67
Base class for inter-region exchange.
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent conductivity.
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
addToRunTimeSelectionTable(heatTransferCoefficientModel, constant, mesh)
Namespace for OpenFOAM.
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimPower
const dimensionSet dimless
const dimensionSet dimLength
const dimensionSet dimTemperature
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > mag(const dimensioned< Type > &)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
const dimensionSet dimArea
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
labelList fv(nPoints)
dictionary dict