constant.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 "constant.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
36 namespace heatTransferCoefficientModels
37 {
41 }
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::fv::heatTransferCoefficientModels::constant::readCoeffs
49 (
50  const dictionary& dict
51 )
52 {
53  typeIOobject<volScalarField> htcIO
54  (
55  "htc",
56  mesh_.time().constant(),
57  mesh_,
60  );
61 
62  if (dict.found("htc"))
63  {
64  htc_ =
66  (
67  "htc",
69  dict
70  );
71  htcPtr_.clear();
72  }
73  else if (htcIO.headerOk())
74  {
76  htcPtr_.set(new volScalarField(htcIO, mesh_));
77  }
78  else
79  {
81  << "Heat transfer coefficient (htc) not found. A uniform htc "
82  << "value should be specified, or a non-uniform field should "
83  << "exist at " << htcIO.objectPath()
84  << exit(FatalIOError);
85  }
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
92 (
93  const dictionary& dict,
94  const fvMesh& mesh
95 )
96 :
97  heatTransferCoefficientModel(typeName, dict, mesh),
98  htc_("htc", dimPower/dimTemperature/dimArea, NaN),
99  htcPtr_(nullptr)
100 {
101  readCoeffs(dict);
102 }
103 
104 
106 (
107  const dictionary& dict,
108  const interRegionModel& model
109 )
110 :
111  constant(dict, model.mesh())
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
116 
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
125 {
126  if (!htcPtr_.valid())
127  {
128  return volScalarField::New(typedName("htc"), mesh_, htc_);
129  }
130  else
131  {
132  return htcPtr_();
133  }
134 }
135 
136 
138 {}
139 
140 
142 (
143  const dictionary& dict
144 )
145 {
147  {
148  readCoeffs(dict);
149  return true;
150  }
151  else
152  {
153  return false;
154  }
155 }
156 
157 
158 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
Base class for heat transfer coefficient modelling used in heat transfer fvModels.
const fvMesh & mesh_
Reference to the mesh.
virtual bool read(const dictionary &dict)
Read dictionary.
Constant heat transfer model. The heat transfer coefficient [W/m^2/K] (htc) must be provided as a val...
Definition: constant.H:64
virtual void correct()
Correct the heat transfer coefficient.
Definition: constant.C:137
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: constant.C:142
constant(const dictionary &dict, const fvMesh &mesh)
Construct from dictionary and mesh.
Definition: constant.C:92
virtual tmp< volScalarField > htc() const
Get the heat transfer coefficient.
Definition: constant.C:124
Base class for inter-region exchange.
A class for managing temporary objects.
Definition: tmp.H:55
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
addToRunTimeSelectionTable(heatTransferCoefficientModel, constant, mesh)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const dimensionSet dimPower
const dimensionSet dimTemperature
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
IOerror FatalIOError
const dimensionSet dimArea
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList fv(nPoints)
dictionary dict