heatTransferCoeff.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) 2017-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
28 #include "heatTransferCoeff.H"
29 #include "heatTransferCoeffModel.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace functionObjects
37 {
38  defineTypeNameAndDebug(heatTransferCoeff, 0);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
45 
47 {
49 
50  htcModelPtr_->calc(htc, htcModelPtr_->q());
51 
52  htc *= L_/kappa_;
53 
54  return true;
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
61 (
62  const word& name,
63  const Time& runTime,
64  const dictionary& dict
65 )
66 :
68  L_(1),
69  kappa_(1),
70  htcModelPtr_(heatTransferCoeffModel::New(dict, mesh_, fieldName_))
71 {
72  read(dict);
73 
74  setResultName(typeName, "htc:" + htcModelPtr_->type());
75 
76  auto* heatTransferCoeffPtr =
77  new volScalarField
78  (
79  IOobject
80  (
82  mesh_.time().timeName(),
83  mesh_,
87  ),
88  mesh_,
90  );
91 
92  mesh_.objectRegistry::store(heatTransferCoeffPtr);
93 }
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
99 {
100  if (!fieldExpression::read(dict) || !htcModelPtr_->read(dict))
101  {
102  return false;
103  }
104 
105  L_ = dict.getCheckOrDefault<scalar>("L", 1, scalarMinMax::ge(0));
106  kappa_ =
107  dict.getCheckOrDefault<scalar>("kappa", 1, scalarMinMax::ge(SMALL));
108 
109  return true;
110 }
111 
112 
113 // ************************************************************************* //
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
word resultName_
Name of result field.
virtual bool calc()
Calculate the heat transfer coefficient field and return true if successful.
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Lookup and return non-const reference to the object of the given Type. Fatal if not found or the wron...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
engineTime & runTime
virtual bool read(const dictionary &dict)
Read the fieldExpression data.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Computes the heat transfer coefficient [W/(m^2 K)] as a volScalarField for a given set of patches...
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Abstract base-class for Time/database function objects.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
A base class for heat transfer coefficient models.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
static MinMax< scalar > ge(const scalar &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:24
heatTransferCoeff(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
const dimensionSet dimPower
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc...
Nothing to be read.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
virtual bool read(const dictionary &dict)
Read the top-level dictionary.
const fvMesh & mesh_
Reference to the fvMesh.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127