heatTransferAv.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) 2021-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 "heatTransferAv.H"
27 #include "heatTransfer.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 void Foam::fv::heatTransferAv::readCoeffs(const dictionary& dict)
32 {
33  typeIOobject<volScalarField> AvIO
34  (
35  "Av",
36  mesh_.time().constant(),
37  mesh_,
40  );
41  typeIOobject<volScalarField> AoVIO
42  (
43  "AoV",
44  mesh_.time().constant(),
45  mesh_,
48  );
49 
50  if (dict.found("Av"))
51  {
53  AvPtr_.clear();
54  }
55  else if (dict.found("AoV"))
56  {
57  Av_ = dimensionedScalar("AoV", dimless/dimLength, dict);
58  AvPtr_.clear();
59  }
60  else if (AvIO.headerOk())
61  {
62  Av_ = dimensionedScalar("Av", dimless/dimLength, NaN);
63  AvPtr_.set(new volScalarField(AvIO, mesh_));
64  }
65  else if (AoVIO.headerOk())
66  {
67  Av_ = dimensionedScalar("AoV", dimless/dimLength, NaN);
68  AvPtr_.set(new volScalarField(AoVIO, mesh_));
69  }
70  else
71  {
73  << "Area per unit volume (Av) not found. A uniform Av value "
74  << "should be specified, or a non-uniform field should exist at "
75  << AvIO.objectPath() << exit(FatalIOError);
76  }
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
83 (
84  const dictionary& dict,
85  const fvMesh& mesh
86 )
87 :
88  mesh_(mesh),
89  Av_("Av", dimless/dimLength, NaN),
90  AvPtr_(nullptr)
91 {
92  readCoeffs(dict);
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97 
99 {}
100 
101 
102 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
103 
105 {
106  if (!AvPtr_.valid())
107  {
108  return volScalarField::New(typedName<heatTransfer>("Av"), mesh_, Av_);
109  }
110  else
111  {
112  return AvPtr_();
113  }
114 }
115 
116 
118 {
119  readCoeffs(dict);
120 
121  return true;
122 }
123 
124 
125 // ************************************************************************* //
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
virtual ~heatTransferAv()
Destructor.
virtual bool read(const dictionary &dict)
Read dictionary.
heatTransferAv(const dictionary &dict, const fvMesh &mesh)
Construct from dictionary and mesh.
tmp< volScalarField > Av() const
Get the area per unit volume.
A class for managing temporary objects.
Definition: tmp.H:55
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const dimensionSet dimless
const dimensionSet dimLength
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
IOerror FatalIOError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict