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 | Copyright (C) 2011-2016 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 "tabulatedHeatTransfer.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace fv
34 {
35  defineTypeNameAndDebug(tabulatedHeatTransfer, 0);
37  (
38  option,
39  tabulatedHeatTransfer,
40  dictionary
41  );
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
49 Foam::fv::tabulatedHeatTransfer::hTable()
50 {
51  if (!hTable_.valid())
52  {
53  hTable_.reset(new interpolation2DTable<scalar>(coeffs_));
54  }
55 
56  return hTable_();
57 }
58 
59 
60 const Foam::volScalarField& Foam::fv::tabulatedHeatTransfer::AoV()
61 {
62  if (!AoV_.valid())
63  {
64  AoV_.reset
65  (
66  new volScalarField
67  (
68  IOobject
69  (
70  "AoV",
71  startTimeName_,
72  mesh_,
75  ),
76  mesh_
77  )
78  );
79  }
80 
81  return AoV_();
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86 
88 (
89  const word& name,
90  const word& modelType,
91  const dictionary& dict,
92  const fvMesh& mesh
93 )
94 :
95  interRegionHeatTransferModel(name, modelType, dict, mesh),
96  UName_(coeffs_.lookupOrDefault<word>("U", "U")),
97  UNbrName_(coeffs_.lookupOrDefault<word>("UNbr", "U")),
98  hTable_(),
99  AoV_(),
100  startTimeName_(mesh.time().timeName())
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
105 
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
113 {
114  const fvMesh& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName());
115 
116  const volVectorField& UNbr =
117  nbrMesh.lookupObject<volVectorField>(UNbrName_);
118 
119  const scalarField UMagNbr(mag(UNbr));
120 
121  const scalarField UMagNbrMapped(interpolate(UMagNbr));
122 
123  const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
124 
125  scalarField& htcc = htc_.primitiveFieldRef();
126 
127  forAll(htcc, i)
128  {
129  htcc[i] = hTable()(mag(U[i]), UMagNbrMapped[i]);
130  }
131 
132  htcc = htcc*AoV();
133 }
134 
135 
137 {
139  {
140  return true;
141  }
142  else
143  {
144  return false;
145  }
146 }
147 
148 
149 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
Base class for inter region heat exchange. The derived classes must provide the heat transfer coeffis...
U
Definition: pEqn.H:83
virtual bool read(const dictionary &dict)
Read dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
2D table interpolation. The data must be in ascending order in both dimensions x and y...
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
tabulatedHeatTransfer(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:644
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
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.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual ~tabulatedHeatTransfer()
Destructor.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual void calculateHtc()
Calculate the heat transfer coefficient.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual bool read(const dictionary &dict)
Read dictionary.
Namespace for OpenFOAM.