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-2023 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 heatTransferCoefficientModels
37 {
41 }
42 }
43 }
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::fv::heatTransferCoefficientModels::function1::readCoeffs
48 (
49  const dictionary& dict
50 )
51 {
52  UName_ = dict.lookupOrDefault<word>("U", "U");
53 
54  htcFunc_.reset(Function1<scalar>::New("htcFunc", dict).ptr());
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
61 (
62  const dictionary& dict,
63  const fvMesh& mesh
64 )
65 :
66  heatTransferCoefficientModel(typeName, dict, mesh),
67  UName_(word::null),
68  htcFunc_(nullptr)
69 {
70  readCoeffs(dict);
71 }
72 
73 
75 (
76  const dictionary& dict,
77  const interRegionModel& model
78 )
79 :
80  function1(dict, model.mesh())
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
85 
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
94 {
95  const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
96 
97  tmp<volScalarField> tHtc =
99  (
100  typedName("htc"),
101  mesh_,
103  zeroGradientFvPatchScalarField::typeName
104  );
105 
106  tHtc->primitiveFieldRef() = htcFunc_->value(mag(U));
107  tHtc->correctBoundaryConditions();
108 
109  return tHtc;
110 }
111 
112 
114 {}
115 
116 
118 (
119  const dictionary& dict
120 )
121 {
123  {
124  readCoeffs(dict);
125  return true;
126  }
127  else
128  {
129  return false;
130  }
131 }
132 
133 
134 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Base class for heat transfer coefficient modelling used in heat transfer fvModels.
virtual bool read(const dictionary &dict)
Read dictionary.
Function1 heat transfer model. The 1D function returns the heat transfer coefficient as a function of...
Definition: function1.H:68
virtual void correct()
Correct the heat transfer coefficient.
Definition: function1.C:113
function1(const dictionary &dict, const fvMesh &mesh)
Construct from dictionary and mesh.
Definition: function1.C:61
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: function1.C:118
virtual tmp< volScalarField > htc() const
Get the heat transfer coefficient.
Definition: function1.C:93
Base class for inter-region exchange.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
addToRunTimeSelectionTable(heatTransferCoefficientModel, constant, mesh)
Namespace for OpenFOAM.
const dimensionSet dimPower
const dimensionSet dimTemperature
dimensioned< scalar > mag(const dimensioned< Type > &)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:149
const dimensionSet dimArea
labelList fv(nPoints)
dictionary dict