tabulatedHeatTransfer.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
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace fv
37{
40}
41}
42
43
44// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45
47Foam::fv::tabulatedHeatTransfer::hTable()
48{
49 if (!hTable_)
50 {
51 hTable_.reset(new interpolation2DTable<scalar>(coeffs_));
52 }
53
54 return *hTable_;
55}
56
57
58const Foam::volScalarField& Foam::fv::tabulatedHeatTransfer::AoV()
59{
60 if (!AoV_)
61 {
62 AoV_.reset
63 (
65 (
66 IOobject
67 (
68 "AoV",
69 startTimeName_,
70 mesh_,
73 ),
74 mesh_
75 )
76 );
77 }
78
79 return *AoV_;
80}
81
82
83// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84
86(
87 const word& name,
88 const word& modelType,
89 const dictionary& dict,
90 const fvMesh& mesh
91)
92:
94 UName_(coeffs_.getOrDefault<word>("U", "U")),
95 UNbrName_(coeffs_.getOrDefault<word>("UNbr", "U")),
96 hTable_(),
97 AoV_(),
98 startTimeName_(mesh.time().timeName())
99{}
100
101
102// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103
105{
106 const auto& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName());
107
108 const auto& UNbr = nbrMesh.lookupObject<volVectorField>(UNbrName_);
109
110 const scalarField UMagNbr(mag(UNbr));
111
112 const scalarField UMagNbrMapped(interpolate(UMagNbr));
113
114 const auto& U = mesh_.lookupObject<volVectorField>(UName_);
115
116 scalarField& htcc = htc_.primitiveFieldRef();
117
118 forAll(htcc, i)
119 {
120 htcc[i] = hTable()(mag(U[i]), UMagNbrMapped[i]);
121 }
122
123 htcc = htcc*AoV();
124}
125
126
128{
130 {
131 return true;
132 }
133
134 return false;
135}
136
137
138// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
virtual bool read()
Re-read model coefficients if they have changed.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Intermediate class for handling inter-region heat exchanges.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:127
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:145
Applies a tabulated heat transfer model for inter-region heat exchanges.
virtual void calculateHtc()
Calculate the heat transfer coefficient.
2D table interpolation. The data must be in ascending order in both dimensions x and y.
const Type & lookupObject(const word &name, const bool recursive=false) const
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
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
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
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333