variableHeatTransfer.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace fv
38{
41}
42}
43
44
45// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46
48(
49 const word& name,
50 const word& modelType,
51 const dictionary& dict,
52 const fvMesh& mesh
53)
54:
56 UNbrName_(coeffs_.getOrDefault<word>("UNbr", "U")),
57 a_(0),
58 b_(0),
59 c_(0),
60 ds_(0),
61 Pr_(0),
62 AoV_()
63{
64 if (master_)
65 {
66 coeffs_.readEntry("a", a_);
67 coeffs_.readEntry("b", b_);
68 coeffs_.readEntry("c", c_);
69 coeffs_.readEntry("ds", ds_);
70 coeffs_.readEntry("Pr", Pr_);
71
72 AoV_.reset
73 (
75 (
77 (
78 "AoV",
80 mesh_,
83 ),
84 mesh_
85 )
86 );
87 }
88}
89
90
91// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92
94{
95 const auto& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName());
96
97 const auto& nbrTurb =
99 (
101 );
102
103 const auto& nbrThermo =
104 nbrMesh.lookupObject<fluidThermo>(basicThermo::dictName);
105
106 const auto& UNbr = nbrMesh.lookupObject<volVectorField>(UNbrName_);
107
108 tmp<volScalarField> ReNbr(mag(UNbr)*ds_*nbrThermo.rho()/nbrTurb.mut());
109
110 tmp<volScalarField> NuNbr(a_*pow(ReNbr, b_)*pow(Pr_, c_));
111
112 scalarField htcNbr(NuNbr*nbrTurb.kappaEff()/ds_);
113
114 tmp<scalarField> htcNbrMapped(interpolate(htcNbr));
115
116 htc_.primitiveFieldRef() = htcNbrMapped*AoV_();
117}
118
119
121{
123 {
124 coeffs_.readIfPresent("UNbr", UNbrName_);
125
126 coeffs_.readIfPresent("a", a_);
127 coeffs_.readIfPresent("b", b_);
128 coeffs_.readIfPresent("c", c_);
129 coeffs_.readIfPresent("ds", ds_);
130 coeffs_.readIfPresent("Pr", Pr_);
131
132 return true;
133 }
134
135 return false;
136}
137
138
139// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual bool read()
Re-read model coefficients if they have changed.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
static const word dictName
Definition: basicThermo.H:256
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Intermediate class for handling inter-region heat exchanges.
bool master_
Master or slave region.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:127
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:145
Applies a variable heat transfer model depending on local values for inter-region heat exchanges.
virtual void calculateHtc()
Calculate the heat transfer coefficient.
const Type & lookupObject(const word &name, const bool recursive=false) const
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
Namespace for OpenFOAM.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
labelList fv(nPoints)
dictionary dict