thermalShell.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2019-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "thermalShell.H"
30 #include "fvPatchFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace regionModels
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(thermalShell, 0);
42 
43 addToRunTimeSelectionTable(thermalShellModel, thermalShell, dictionary);
44 
45 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 bool thermalShell::init(const dictionary& dict)
48 {
49  if (thickness_ > 0)
50  {
51  h_ = dimensionedScalar("thickness", dimLength, thickness_);
52  }
53 
54  this->solution().readEntry("nNonOrthCorr", nNonOrthCorr_);
55 
56  return true;
57 }
58 
59 
60 tmp<areaScalarField> thermalShell::qr()
61 {
62  auto taqr =
64  (
65  IOobject
66  (
67  "tqr",
68  regionMesh().time().timeName(),
69  regionMesh().thisDb()
70  ),
71  regionMesh(),
73  );
74 
75  if (!qrName_.empty() && qrName_ != "none")
76  {
77  vsm().mapToSurface<scalar>
78  (
80  taqr.ref().primitiveFieldRef()
81  );
82  }
83 
84  return taqr;
85 }
86 
87 
88 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
89 
91 {
93 
94  const areaScalarField rhoCph(Cp()*rho()*h_);
95 
97  (
98  fam::ddt(rhoCph, T_)
100  ==
101  qs_
102  + qr()
103  + faOptions()(h_, rhoCph, T_)
104  );
105 
106  TEqn.relax();
107 
109 
110  TEqn.solve();
111 
113 }
114 
115 
116 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
117 
119 (
120  const word& modelType,
121  const fvMesh& mesh,
122  const dictionary& dict
123 )
124 :
125  thermalShellModel(modelType, mesh, dict),
126  nNonOrthCorr_(1),
127  thermo_(dict.subDict("thermo")),
128  qs_
129  (
130  IOobject
131  (
132  "qs_" + regionName_,
133  regionMesh().time().timeName(),
134  regionMesh().thisDb(),
135  IOobject::READ_IF_PRESENT,
136  IOobject::NO_WRITE
137  ),
138  regionMesh(),
140  ),
141  h_
142  (
143  IOobject
144  (
145  "h_" + regionName_,
146  regionMesh().time().timeName(),
147  regionMesh().thisDb(),
148  IOobject::MUST_READ,
149  IOobject::AUTO_WRITE
150  ),
151  regionMesh()
152  ),
153  qrName_(dict.getOrDefault<word>("qr", "none")),
154  thickness_(dict.getOrDefault<scalar>("thickness", 0))
155 {
156  init(dict);
157 }
158 
160 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
161 
163 {}
164 
165 
167 {
168  nNonOrthCorr_ = solution().get<label>("nNonOrthCorr");
169 
170  for (int nonOrth = 0; nonOrth <= nNonOrthCorr_; ++nonOrth)
171  {
173  }
174 
175  Info<< T_.name() << " min/max = " << gMinMax(T_) << endl;
176 }
177 
178 
180 {
181  return areaScalarField::New
182  (
183  "Cps",
188  );
189 }
190 
191 
193 {
194  return areaScalarField::New
195  (
196  "rhos",
201  );
202 }
203 
204 
206 {
207  return areaScalarField::New
208  (
209  "kappas",
211  regionMesh(),
213  (
216  ),
218  );
219 }
220 
221 
222 void thermalShell::info()
223 {}
224 
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 } // End namespace regionModels
229 } // End namespace Foam
230 
231 // ************************************************************************* //
void mapToSurface(const GeometricBoundaryField< Type, fvPatchField, volMesh > &, Field< Type > &result) const
Map volume boundary fields as area field.
label nNonOrthCorr_
Number of non orthogonal correctors.
Definition: thermalShell.H:145
dictionary dict
const dictionary & solution() const
Return the solution dictionary.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
scalar Cp() const
Specific heat capacity [J/(kg.K)].
scalar rho() const
Density [kg/m3].
void constrain(faMatrix< Type > &eqn)
Apply constraints to equation.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
const volSurfaceMapping & vsm() const
Return mapping between surface and volume fields.
areaScalarField T_
Shell temperature.
virtual void info()
Provide some feedback.
Definition: thermalShell.C:215
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const word qrName_
Name of the primary region radiative flux.
Definition: thermalShell.H:171
const tmp< areaScalarField > rho() const
Return density [kg/m3].
Definition: thermalShell.C:185
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: faPatchField.H:191
Macros for easy insertion into run-time selection tables.
areaScalarField qs_
External surface energy source [J/m2/s].
Definition: thermalShell.H:161
void correct(GeometricField< Type, faPatchField, areaMesh > &field)
Apply correction to field.
solidProperties thermo_
Solid properties.
Definition: thermalShell.H:153
const tmp< areaScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
Definition: thermalShell.C:172
tmp< faMatrix< Type > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famLaplacian.C:42
Foam::fa::options & faOptions() noexcept
Return faOptions.
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
Definition: thermalShell.C:155
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
scalar thickness_
Uniform film thickness [m].
Definition: thermalShell.H:176
A class for handling words, derived from Foam::string.
Definition: word.H:63
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
#define DebugInFunction
Report an information message using Foam::Info.
const fvMesh & primaryMesh() const noexcept
Return the reference to the primary mesh database.
virtual void evolveRegion()
Evolve the thermal baffle.
Definition: thermalShell.C:159
const faMesh & regionMesh() const
Return the region mesh database.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static tmp< GeometricField< scalar, faPatchField, areaMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=faPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
const dimensionSet dimPower
tmp< faMatrix< Type > > ddt(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famDdt.C:41
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
scalar kappa() const
Thermal conductivity [W/(m.K)].
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1095
const dimensionSet dimEnergy
const dimensionSet dimDensity
addToRunTimeSelectionTable(vibrationShellModel, KirchhoffShell, dictionary)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void solveEnergy()
Solve energy equation.
Definition: thermalShell.C:83
Template specialisation for scalar faMatrix.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const Time & time() const noexcept
Return the reference to the time database.
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
thermalShell(const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components and dict.
Definition: thermalShell.C:112
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
const tmp< areaScalarField > kappa() const
Return thermal conductivity [W/m/K].
Definition: thermalShell.C:198
areaScalarField h_
Film thickness [m].
Definition: thermalShell.H:166
defineTypeNameAndDebug(KirchhoffShell, 0)
Do not request registration (bool: false)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127