function2.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 "function2.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
36 namespace heatTransferModels
37 {
38  defineTypeNameAndDebug(function2, 0);
39  addToRunTimeSelectionTable(heatTransferModel, function2, model);
40 }
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::fv::heatTransferModels::function2::readCoeffs()
48 {
49  UName_ = coeffs().lookupOrDefault<word>("U", "U");
50  UNbrName_ = coeffs().lookupOrDefault<word>("UNbr", "U");
51 
52  htcFunc_.reset(Function2<scalar>::New("htcFunc", coeffs()).ptr());
53 }
54 
55 
56 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57 
59 (
60  const dictionary& dict,
61  const interRegionModel& model
62 )
63 :
64  heatTransferModel(typeName, dict, model),
65  model_(model),
66  UName_(word::null),
67  UNbrName_(word::null),
68  htcFunc_(),
69  htc_
70  (
71  IOobject
72  (
73  type() + ":htc",
74  model.mesh().time().timeName(),
75  model.mesh(),
78  ),
79  model.mesh(),
81  zeroGradientFvPatchScalarField::typeName
82  )
83 {
84  readCoeffs();
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
89 
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
97 {
98  const fvMesh& mesh = model_.mesh();
99  const fvMesh& nbrMesh = model_.nbrMesh();
100 
101  const volVectorField& U = mesh.lookupObject<volVectorField>(UName_);
102 
103  const volVectorField& UNbr =
104  nbrMesh.lookupObject<volVectorField>(UNbrName_);
105  const scalarField UMagNbr(mag(UNbr));
106  const scalarField UMagNbrMapped(model_.interpolate(UMagNbr));
107 
108  htc_.primitiveFieldRef() = htcFunc_->value(mag(U()), UMagNbrMapped);
109  htc_.correctBoundaryConditions();
110 }
111 
112 
114 {
115  if (heatTransferModel::read(dict))
116  {
117  readCoeffs();
118  return true;
119  }
120  else
121  {
122  return false;
123  }
124 }
125 
126 
127 // ************************************************************************* //
virtual ~function2()
Destructor.
Definition: function2.C:90
const dimensionSet dimArea
const dimensionSet dimPower
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 autoPtr< Function2< scalar > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function2New.C:32
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
fvMesh & mesh
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: function2.C:113
Macros for easy insertion into run-time selection tables.
Base class for inter-region exchange.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
labelList fv(nPoints)
virtual void correct()
Correct the heat transfer coefficient.
Definition: function2.C:96
static const word null
An empty word.
Definition: word.H:77
defineTypeNameAndDebug(constant, 0)
function2(const dictionary &dict, const interRegionModel &model)
Construct from dictionary and model.
Definition: function2.C:59
virtual bool read(const dictionary &dict)
Read dictionary.
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)
dimensioned< scalar > mag(const dimensioned< Type > &)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
const dimensionSet dimTemperature
Namespace for OpenFOAM.