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 -------------------------------------------------------------------------------
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 "variableHeatTransfer.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace fv
38 {
39  defineTypeNameAndDebug(variableHeatTransfer, 0);
40  addToRunTimeSelectionTable(option, variableHeatTransfer, dictionary);
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  (
74  new volScalarField
75  (
76  IOobject
77  (
78  "AoV",
79  mesh_.time().timeName(),
80  mesh_,
81  IOobject::MUST_READ,
82  IOobject::AUTO_WRITE
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 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::fv::interRegionHeatTransferModel::nbrRegionName
const word & nbrRegionName() const
Return const access to the neighbour region name.
Definition: interRegionHeatTransferModelI.H:32
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::fv::interRegionHeatTransferModel
Intermediate class for handling inter-region heat exchanges.
Definition: interRegionHeatTransferModel.H:142
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::fluidThermo
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:52
variableHeatTransfer.H
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
Foam::Field< scalar >
Foam::fv::variableHeatTransfer::calculateHtc
virtual void calculateHtc()
Calculate the heat transfer coefficient.
Definition: variableHeatTransfer.C:93
Foam::fv::interRegionHeatTransferModel::interpolate
tmp< Field< Type > > interpolate(const interRegionHeatTransferModel &nbrModel, const Field< Type > &field) const
Interpolate field with nbrModel specified.
Foam::fv::variableHeatTransfer::variableHeatTransfer
variableHeatTransfer(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.
Definition: variableHeatTransfer.C:48
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::fv::interRegionHeatTransferModel::htc_
volScalarField htc_
Heat transfer coefficient [W/m2/k] times area/volume [1/m].
Definition: interRegionHeatTransferModel.H:166
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
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:123
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:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
fv
labelList fv(nPoints)
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:47
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::fv::interRegionHeatTransferModel::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: interRegionHeatTransferModel.C:264
Foam::fv::variableHeatTransfer::read
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: variableHeatTransfer.C:120
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
Foam::GeometricField< scalar, fvPatchField, volMesh >
turbulentFluidThermoModel.H
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:60