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"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace regionModels
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 defineTypeNameAndDebug(thermalShell, 0);
43 
44 addToRunTimeSelectionTable(thermalShellModel, thermalShell, dictionary);
45 
46 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 bool thermalShell::init(const dictionary& dict)
49 {
50  if (thickness_ > 0)
51  {
52  h_ = dimensionedScalar("thickness", dimLength, thickness_);
53  }
54 
55  this->solution().readEntry("nNonOrthCorr", nNonOrthCorr_);
56 
57  return true;
58 }
59 
60 
61 tmp<areaScalarField> thermalShell::qr()
62 {
63  auto taqr =
65  (
66  IOobject
67  (
68  "tqr",
69  regionMesh().time().timeName(),
70  regionMesh().thisDb()
71  ),
72  regionMesh(),
74  );
75 
76  if (!qrName_.empty() && qrName_ != "none")
77  {
78  vsm().mapToSurface<scalar>
79  (
81  taqr.ref().primitiveFieldRef()
82  );
83  }
84 
85  return taqr;
86 }
87 
88 
89 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
90 
92 {
93  if (debug)
94  {
96  }
97 
98  const areaScalarField rhoCph(Cp()*rho()*h_);
99 
101  (
102  fam::ddt(rhoCph, T_)
103  - fam::laplacian(kappa()*h_, T_)
104  ==
105  qs_
106  + qr()
107  + faOptions()(h_, rhoCph, T_)
108  );
109 
110  TEqn.relax();
111 
113 
114  TEqn.solve();
115 
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
123 (
124  const word& modelType,
125  const fvMesh& mesh,
126  const dictionary& dict
127 )
128 :
129  thermalShellModel(modelType, mesh, dict),
130  nNonOrthCorr_(1),
131  thermo_(dict.subDict("thermo")),
132  qs_
133  (
134  IOobject
135  (
136  "qs_" + regionName_,
137  regionMesh().time().timeName(),
138  regionMesh().thisDb(),
139  IOobject::READ_IF_PRESENT,
140  IOobject::NO_WRITE
141  ),
142  regionMesh(),
144  ),
145  h_
146  (
147  IOobject
148  (
149  "h_" + regionName_,
150  regionMesh().time().timeName(),
151  regionMesh().thisDb(),
152  IOobject::MUST_READ,
153  IOobject::AUTO_WRITE
154  ),
155  regionMesh()
156  ),
157  qrName_(dict.getOrDefault<word>("qr", "none")),
158  thickness_(dict.getOrDefault<scalar>("thickness", 0))
159 {
160  init(dict);
161 }
162 
164 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
165 
167 {}
168 
169 
171 {
172  nNonOrthCorr_ = solution().get<label>("nNonOrthCorr");
173 
174  for (int nonOrth = 0; nonOrth <= nNonOrthCorr_; ++nonOrth)
175  {
177  }
178 
179  Info<< T_.name() << " min/max = " << gMinMax(T_) << endl;
180 }
181 
182 
184 {
185  return tmp<areaScalarField>
186  (
187  new areaScalarField
188  (
189  IOobject
190  (
191  "Cps",
192  regionMesh().time().timeName(),
193  regionMesh().thisDb(),
197  ),
198  regionMesh(),
201  )
202  );
203 }
204 
205 
207 {
208  return tmp<areaScalarField>
209  (
210  new areaScalarField
211  (
212  IOobject
213  (
214  "rhos",
215  regionMesh().time().timeName(),
216  regionMesh().thisDb(),
220  ),
221  regionMesh(),
224  )
225  );
226 }
227 
228 
230 {
231  return tmp<areaScalarField>
232  (
233  new areaScalarField
234  (
235  IOobject
236  (
237  "kappas",
238  regionMesh().time().timeName(),
239  regionMesh().thisDb(),
243  ),
244  regionMesh(),
246  (
248  thermo_.kappa()
249  ),
251  )
252  );
253 }
254 
255 
256 void thermalShell::info()
257 {}
258 
259 
260 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
261 
262 } // End namespace regionModels
263 } // End namespace Foam
264 
265 // ************************************************************************* //
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:133
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:249
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:159
const tmp< areaScalarField > rho() const
Return density [Kg/m3].
Definition: thermalShell.C:199
Ignore writing from objectRegistry::writeObject()
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:149
void correct(GeometricField< Type, faPatchField, areaMesh > &field)
Apply correction to field.
solidProperties thermo_
Solid properties.
Definition: thermalShell.H:141
const tmp< areaScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
Definition: thermalShell.C:176
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:159
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 thickness.
Definition: thermalShell.H:164
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
const fvMesh & primaryMesh() const noexcept
Return the reference to the primary mesh database.
virtual void evolveRegion()
Evolve the thermal baffle.
Definition: thermalShell.C:163
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)
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.
int debug
Static debugging option.
scalar kappa() const
Thermal conductivity [W/(m.K)].
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1101
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:84
Nothing to be read.
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:116
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
const tmp< areaScalarField > kappa() const
Return thermal conductivity [W/m/K].
Definition: thermalShell.C:222
areaScalarField h_
Thickness.
Definition: thermalShell.H:154
defineTypeNameAndDebug(KirchhoffShell, 0)
Do not request registration (bool: false)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:80
#define InfoInFunction
Report an information message using Foam::Info.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127