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 -------------------------------------------------------------------------------
11 License
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 
29 #include "tabulatedHeatTransfer.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
38  defineTypeNameAndDebug(tabulatedHeatTransfer, 0);
40  (
41  option,
42  tabulatedHeatTransfer,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
52 Foam::fv::tabulatedHeatTransfer::hTable()
53 {
54  if (!hTable_.valid())
55  {
56  hTable_.reset(new interpolation2DTable<scalar>(coeffs_));
57  }
58 
59  return *hTable_;
60 }
61 
62 
63 const Foam::volScalarField& Foam::fv::tabulatedHeatTransfer::AoV()
64 {
65  if (!AoV_.valid())
66  {
67  AoV_.reset
68  (
69  new volScalarField
70  (
71  IOobject
72  (
73  "AoV",
74  startTimeName_,
75  mesh_,
78  ),
79  mesh_
80  )
81  );
82  }
83 
84  return *AoV_;
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89 
91 (
92  const word& name,
93  const word& modelType,
94  const dictionary& dict,
95  const fvMesh& mesh
96 )
97 :
99  UName_(coeffs_.getOrDefault<word>("U", "U")),
100  UNbrName_(coeffs_.getOrDefault<word>("UNbr", "U")),
101  hTable_(),
102  AoV_(),
103  startTimeName_(mesh.time().timeName())
104 {}
105 
106 
107 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
108 
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 {
117  const fvMesh& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName());
118 
119  const volVectorField& UNbr =
120  nbrMesh.lookupObject<volVectorField>(UNbrName_);
121 
122  const scalarField UMagNbr(mag(UNbr));
123 
124  const scalarField UMagNbrMapped(interpolate(UMagNbr));
125 
126  const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
127 
128  scalarField& htcc = htc_.primitiveFieldRef();
129 
130  forAll(htcc, i)
131  {
132  htcc[i] = hTable()(mag(U[i]), UMagNbrMapped[i]);
133  }
134 
135  htcc = htcc*AoV();
136 }
137 
138 
140 {
142  {
143  return true;
144  }
145 
146  return false;
147 }
148 
149 
150 // ************************************************************************* //
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::interRegionHeatTransferModel
Base class for inter region heat exchange. The derived classes must provide the heat transfer coeffis...
Definition: interRegionHeatTransferModel.H:59
tabulatedHeatTransfer.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::fv::tabulatedHeatTransfer::read
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: tabulatedHeatTransfer.C:139
Foam::interpolation2DTable
2D table interpolation. The data must be in ascending order in both dimensions x and y.
Definition: interpolation2DTable.H:53
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::Field< scalar >
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::fv::option::coeffs_
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:88
Foam::fv::tabulatedHeatTransfer::tabulatedHeatTransfer
tabulatedHeatTransfer(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.
Definition: tabulatedHeatTransfer.C:91
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::fv::tabulatedHeatTransfer::calculateHtc
virtual void calculateHtc()
Calculate the heat transfer coefficient.
Definition: tabulatedHeatTransfer.C:115
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
fv
labelList fv(nPoints)
U
U
Definition: pEqn.H:72
Foam::fv::tabulatedHeatTransfer::~tabulatedHeatTransfer
virtual ~tabulatedHeatTransfer()
Destructor.
Definition: tabulatedHeatTransfer.C:109
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::fv::interRegionHeatTransferModel::read
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: interRegionHeatTransferModel.C:266
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:248
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::MUST_READ
Definition: IOobject.H:120