Loading...
Searching...
No Matches
tabulatedNTUHeatTransfer.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) 2015-2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "surfaceFields.H"
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace fv
37{
40}
41}
42
43const Foam::Enum
44<
46>
48({
49 { geometryModeType::gmCalculated, "calculated" },
50 { geometryModeType::gmUser, "user" },
51});
52
53
54// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55
57Foam::fv::tabulatedNTUHeatTransfer::ntuTable()
58{
59 if (!ntuTable_)
60 {
61 ntuTable_.reset(new interpolation2DTable<scalar>(coeffs_));
62 }
63
64 return *ntuTable_;
65}
66
67
68const Foam::basicThermo& Foam::fv::tabulatedNTUHeatTransfer::thermo
69(
70 const fvMesh& mesh
71) const
72{
73 if (!mesh.foundObject<basicThermo>(basicThermo::dictName))
74 {
76 << " on mesh " << mesh.name()
77 << " could not find " << basicThermo::dictName
78 << exit(FatalError);
79 }
80
81 return mesh.lookupObject<basicThermo>(basicThermo::dictName);
82}
83
84
85void Foam::fv::tabulatedNTUHeatTransfer::initialiseGeometry()
86{
87 if (Ain_ < 0)
88 {
89 geometryMode_ = geometryModelNames_.get("geometryMode", coeffs_);
90
91 Info<< "Region " << mesh_.name() << " " << type() << " " << name_
92 << " " << geometryModelNames_[geometryMode_] << " geometry:" << nl;
93
94 switch (geometryMode_)
95 {
96 case gmCalculated:
97 {
98 const auto& nbrMesh =
99 mesh_.time().lookupObject<fvMesh>(nbrRegionName());
100
101 word inletPatchName(coeffs_.get<word>("inletPatch"));
102 word inletPatchNbrName(coeffs_.get<word>("inletPatchNbr"));
103
104 Info<< " Inlet patch : " << inletPatchName << nl
105 << " Inlet patch neighbour : " << inletPatchNbrName
106 << nl;
107
108 label patchI = mesh_.boundary().findPatchID(inletPatchName);
109 label patchINbr =
110 nbrMesh.boundary().findPatchID(inletPatchNbrName);
111
112 scalar alpha(coeffs_.get<scalar>("inletBlockageRatio"));
113
115 {
117 << "Inlet patch blockage ratio must be between 0 and 1"
118 << ". Current value: " << alpha
119 << abort(FatalError);
120 }
121
122 scalar alphaNbr(coeffs_.get<scalar>("inletBlockageRatioNbr"));
123
124 if (alphaNbr < 0 || alphaNbr > 1)
125 {
127 << "Inlet patch neighbour blockage ratio must be "
128 << "between 0 and 1. Current value: " << alphaNbr
129 << abort(FatalError);
130 }
131
132 Info<< " Inlet blockage ratio : " << alpha << nl
133 << " Inlet blockage ratio neighbour : " << alphaNbr
134 << nl;
135
136 Ain_ =
137 (scalar(1) - alpha)
138 *gSum(mesh_.magSf().boundaryField()[patchI]);
139
140 AinNbr_ =
141 (scalar(1) - alphaNbr)
142 *gSum(nbrMesh.magSf().boundaryField()[patchINbr]);
143
144 scalar beta(coeffs_.get<scalar>("coreBlockageRatio"));
145
146 if (beta < 0 || beta > 1)
147 {
149 << "Core volume blockage ratio must be between 0 and 1"
150 << ". Current value: " << beta
151 << abort(FatalError);
152 }
153
154 Info<< " Core volume blockage ratio : " << beta << nl;
155
156 Vcore_ = (scalar(1) - beta)*meshInterp().V();
157
158 break;
159 }
160 case gmUser:
161 {
162 coeffs_.readEntry("Ain", Ain_);
163 coeffs_.readEntry("AinNbr", AinNbr_);
164
165 if (!coeffs_.readIfPresent("Vcore", Vcore_))
166 {
167 Vcore_ = meshInterp().V();
168 }
169
170 break;
171 }
172 default:
173 {
175 << "Unhandled enumeration " << geometryMode_
176 << abort(FatalError);
177 }
178 }
179
180 Info<< " Inlet area local : " << Ain_ << nl
181 << " Inlet area neighbour : " << AinNbr_ << nl
182 << " Core volume : " << Vcore_ << nl
184 }
185}
186
187
188// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
189
191(
192 const word& name,
193 const word& modelType,
194 const dictionary& dict,
195 const fvMesh& mesh
196)
197:
198 interRegionHeatTransferModel(name, modelType, dict, mesh),
199 UName_(coeffs_.getOrDefault<word>("U", "U")),
200 UNbrName_(coeffs_.getOrDefault<word>("UNbr", "U")),
201 rhoName_(coeffs_.getOrDefault<word>("rho", "rho")),
202 rhoNbrName_(coeffs_.getOrDefault<word>("rhoNbr", "rho")),
203 ntuTable_(),
204 geometryMode_(gmCalculated),
205 Ain_(-1),
206 AinNbr_(-1),
207 Vcore_(-1)
208{}
209
210
211// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
212
214{
215 initialiseGeometry();
216
217 const auto& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName());
218
219 const basicThermo& thermo = this->thermo(mesh_);
220 const basicThermo& thermoNbr = this->thermo(nbrMesh);
221 const volScalarField Cp(thermo.Cp());
222 const volScalarField CpNbr(thermoNbr.Cp());
223
224 // Calculate scaled mass flow for primary region
225 const auto& U = mesh_.lookupObject<volVectorField>(UName_);
226 const auto& rho = mesh_.lookupObject<volScalarField>(rhoName_);
227 const scalarField mDot(mag(U)*rho*Ain_);
228
229 // Calculate scaled mass flow for neighbour region
230 const auto& UNbr = nbrMesh.lookupObject<volVectorField>(UNbrName_);
231 const scalarField UMagNbr(mag(UNbr));
232 const scalarField UMagNbrMapped(interpolate(UMagNbr));
233 const auto& rhoNbr =
234 nbrMesh.lookupObject<volScalarField>(rhoNbrName_).internalField();
235 const scalarField rhoNbrMapped(interpolate(rhoNbr));
236 const scalarField mDotNbr(UMagNbrMapped*rhoNbrMapped*AinNbr_);
237
238 scalarField& htcc = htc_.primitiveFieldRef();
239 const interpolation2DTable<Foam::scalar>& ntuTable = this->ntuTable();
240
241 forAll(htcc, cellI)
242 {
243 scalar Cpc = Cp[cellI];
244 scalar CpcNbr = CpNbr[cellI];
245 scalar mDotc = mDot[cellI];
246 scalar mDotcNbr = mDotNbr[cellI];
247 scalar Cmin = min(Cpc*mDotc, CpcNbr*mDotcNbr);
248 scalar ntu = ntuTable(mDotc, mDotcNbr);
249
250 htcc[cellI] = Cmin*ntu/Vcore_;
251 }
252}
253
254
256{
258 {
259 coeffs_.readIfPresent("U", UName_);
260 coeffs_.readIfPresent("UNbr", UNbrName_);
261 coeffs_.readIfPresent("rho", rhoName_);
262 coeffs_.readIfPresent("rhoNbr", rhoNbrName_);
263
264 // Force geometry re-initialisation
265 Ain_ = -1;
266 initialiseGeometry();
267
268 return true;
269 }
270
271 return false;
272}
273
274
275// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
static const word dictName
The dictionary name ("thermophysicalProperties").
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
interRegionHeatTransferModel(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.
volScalarField htc_
Heat transfer coefficient [W/m2/k] times area/volume [1/m].
const word & nbrRegionName() const
Return const access to the neighbour region name.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition fvOptionIO.C:48
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
Applies a tabulated heat transfer model for inter-region heat exchanges. The heat flux is calculated ...
static const Enum< geometryModeType > geometryModelNames_
Names for geometryModeType.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void calculateHtc()
Calculate the heat transfer coefficient.
tabulatedNTUHeatTransfer(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
geometryModeType
Options for the geometry mode type.
2D table interpolation. The data must be in ascending order in both dimensions x and y.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
const volScalarField & Cp
Definition EEqn.H:7
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for finite-volume.
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dictionary dict
psiReactionThermo & thermo
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.