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-2021 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 heatTransferModels
37 {
41 }
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::fv::heatTransferModels::constant::readCoeffs()
49 {
51  (
52  "htc",
53  mesh().time().constant(),
54  mesh(),
57  );
58 
59  if (coeffs().found("htc"))
60  {
61  htc_ =
63  (
64  "htc",
66  coeffs()
67  );
68  htcPtr_.clear();
69  }
70  else if (htcIO.headerOk())
71  {
73  htcPtr_.set(new volScalarField(htcIO, mesh()));
74  }
75  else
76  {
77  FatalIOErrorInFunction(coeffs())
78  << "Heat transfer coefficient (htc) not found. A uniform htc "
79  << "value should be specified, or a non-uniform field should "
80  << "exist at " << htcIO.objectPath()
81  << exit(FatalIOError);
82  }
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
89 (
90  const dictionary& dict,
91  const fvMesh& mesh
92 )
93 :
94  heatTransferModel(typeName, dict, mesh),
95  htc_("htc", dimPower/dimTemperature/dimArea, NaN),
96  htcPtr_(nullptr)
97 {
98  readCoeffs();
99 }
100 
101 
103 (
104  const dictionary& dict,
105  const interRegionModel& model
106 )
107 :
108  constant(dict, model.mesh())
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
113 
115 {}
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
122 {
123  if (!htcPtr_.valid())
124  {
125  return volScalarField::New(type() + ":htc", mesh(), htc_);
126  }
127  else
128  {
129  return htcPtr_();
130  }
131 }
132 
133 
135 {}
136 
137 
139 {
140  if (heatTransferModel::read(dict))
141  {
142  readCoeffs();
143  return true;
144  }
145  else
146  {
147  return false;
148  }
149 }
150 
151 
152 // ************************************************************************* //
const dimensionSet dimArea
const dimensionSet dimPower
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
fileName objectPath() const
Return the object path for this Type.
Definition: IOobject.H:564
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Base class for heat transfer coefficient modelling used in heat transfer fvModels. Area per unit volume [1/m] (AoV) must be provided as a value in the coefficients dictionary or as a field in constant.
Templated form of IOobject providing type information for file reading and header type checking...
Definition: IOobject.H:537
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
fvMesh & mesh
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
virtual void correct()
Correct the heat transfer coefficient.
Definition: constant.C:134
Base class for inter-region exchange.
labelList fv(nPoints)
constant(const dictionary &dict, const fvMesh &mesh)
Construct from dictionary and mesh.
Definition: constant.C:89
defineTypeNameAndDebug(constant, 0)
Constant heat transfer model. The heat transfer coefficient [W/m^2/K] (htc) must be provided as a val...
Definition: constant.H:63
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: constant.C:138
virtual bool read(const dictionary &dict)
Read dictionary.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
addToRunTimeSelectionTable(heatTransferModel, constant, mesh)
virtual ~constant()
Destructor.
Definition: constant.C:114
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > htc() const
Get the heat transfer coefficient.
Definition: constant.C:121
bool found
const dimensionSet dimTemperature
Namespace for OpenFOAM.
IOerror FatalIOError