function1.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 "function1.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
36 namespace heatTransferModels
37 {
38  defineTypeNameAndDebug(function1, 0);
39  addToRunTimeSelectionTable(heatTransferModel, function1, mesh);
40  addToRunTimeSelectionTable(heatTransferModel, function1, model);
41 }
42 }
43 }
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::fv::heatTransferModels::function1::readCoeffs()
48 {
49  UName_ = coeffs().lookupOrDefault<word>("U", "U");
50 
51  htcFunc_.reset(Function1<scalar>::New("htcFunc", coeffs()).ptr());
52 }
53 
54 
55 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56 
58 (
59  const dictionary& dict,
60  const fvMesh& mesh
61 )
62 :
63  heatTransferModel(typeName, dict, mesh),
64  UName_(word::null),
65  htcFunc_(nullptr)
66 {
67  readCoeffs();
68 }
69 
70 
72 (
73  const dictionary& dict,
74  const interRegionModel& model
75 )
76 :
77  function1(dict, model.mesh())
78 {}
79 
80 
81 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
82 
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
91 {
92  const volVectorField& U = mesh().lookupObject<volVectorField>(UName_);
93 
94  tmp<volScalarField> tHtc =
96  (
97  type() + ":htc",
98  mesh(),
100  zeroGradientFvPatchScalarField::typeName
101  );
102 
103  tHtc->primitiveFieldRef() = htcFunc_->value(mag(U));
104  tHtc->correctBoundaryConditions();
105 
106  return tHtc;
107 }
108 
109 
111 {}
112 
113 
115 {
116  if (heatTransferModel::read(dict))
117  {
118  readCoeffs();
119  return true;
120  }
121  else
122  {
123  return false;
124  }
125 }
126 
127 
128 // ************************************************************************* //
const dimensionSet dimArea
const dimensionSet dimPower
function1(const dictionary &dict, const fvMesh &mesh)
Construct from dictionary and mesh.
Definition: function1.C:58
U
Definition: pEqn.H:72
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.
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.
virtual ~function1()
Destructor.
Definition: function1.C:83
Base class for inter-region exchange.
labelList fv(nPoints)
static const word null
An empty word.
Definition: word.H:77
Function1 heat transfer model. The 1D function returns the heat transfer coefficient as a function of...
Definition: function1.H:67
defineTypeNameAndDebug(constant, 0)
virtual tmp< volScalarField > htc() const
Get the heat transfer coefficient.
Definition: function1.C:90
virtual void correct()
Correct the heat transfer coefficient.
Definition: function1.C:110
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: function1.C:114
virtual bool read(const dictionary &dict)
Read dictionary.
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)
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionSet dimTemperature
Namespace for OpenFOAM.
static autoPtr< Function1< scalar > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32