variableHeatTransfer.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-2020 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 "variableHeatTransfer.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
36  defineTypeNameAndDebug(variableHeatTransfer, 0);
38  (
39  option,
40  variableHeatTransfer,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const word& name,
52  const word& modelType,
53  const dictionary& dict,
54  const fvMesh& mesh
55 )
56 :
57  interRegionHeatTransferModel(name, modelType, dict, mesh),
58  UNbrName_(coeffs_.lookupOrDefault<word>("UNbr", "U")),
59  a_(0),
60  b_(0),
61  c_(0),
62  ds_(0),
63  Pr_(0),
64  AoV_()
65 {
66  if (master_)
67  {
68  a_ = coeffs_.lookup<scalar>("a");
69  b_ = coeffs_.lookup<scalar>("b");
70  c_ = coeffs_.lookup<scalar>("c");
71  ds_ = coeffs_.lookup<scalar>("ds");
72  Pr_ = coeffs_.lookup<scalar>("Pr");
73  AoV_.reset
74  (
75  new volScalarField
76  (
77  IOobject
78  (
79  "AoV",
80  mesh_.time().timeName(),
81  mesh_,
84  ),
85  mesh_
86  )
87  );
88  }
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
93 
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 {
102  const fvMesh& nbrMesh =
103  mesh_.time().lookupObject<fvMesh>(nbrRegionName());
104 
105  const thermophysicalTransportModel& nbrTtm =
107  (
108  thermophysicalTransportModel::typeName
109  );
110 
111  const compressibleMomentumTransportModel& nbrTurb =
112  nbrTtm.momentumTransport();
113 
114  const fluidThermo& nbrThermo = nbrTtm.thermo();
115 
116  const volVectorField& UNbr =
117  nbrMesh.lookupObject<volVectorField>(UNbrName_);
118 
119  const volScalarField ReNbr(mag(UNbr)*ds_*nbrThermo.rho()/nbrTurb.mut());
120 
121  const volScalarField NuNbr(a_*pow(ReNbr, b_)*pow(Pr_, c_));
122 
123  const scalarField htcNbr(NuNbr*nbrTtm.kappaEff()/ds_);
124 
125  const scalarField htcNbrMapped(interpolate(htcNbr));
126 
127  htc_.primitiveFieldRef() = htcNbrMapped*AoV_;
128 }
129 
130 
132 {
134  {
135  coeffs_.readIfPresent("UNbr", UNbrName_);
136 
137  coeffs_.readIfPresent("a", a_);
138  coeffs_.readIfPresent("b", b_);
139  coeffs_.readIfPresent("c", c_);
140  coeffs_.readIfPresent("ds", ds_);
141  coeffs_.readIfPresent("Pr", Pr_);
142 
143  return true;
144  }
145  else
146  {
147  return false;
148  }
149 }
150 
151 
152 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
Base class for inter region heat exchange. The derived classes must provide the heat transfer coeffic...
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
virtual ~variableHeatTransfer()
Destructor.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:48
virtual void calculateHtc()
Calculate the heat transfer coefficient.
variableHeatTransfer(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.
Abstract base class for thermophysical transport models (RAS, LES and laminar).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual bool read(const dictionary &dict)
Read dictionary.
virtual bool read(const dictionary &dict)
Read dictionary.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Abstract base class for turbulence models (RAS, LES and laminar).
Namespace for OpenFOAM.