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-2020 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  IOobject io
64  (
65  "tqr",
67  primaryMesh()
68  );
69 
70  auto taqr =
72  (
73  io,
74  regionMesh(),
76  );
77 
78  if (!qrName_.empty() && qrName_ != "none")
79  {
80  vsm().mapToSurface<scalar>
81  (
83  taqr.ref().primitiveFieldRef()
84  );
85  }
86 
87  return taqr;
88 }
89 
90 
91 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
92 
94 {
95  if (debug)
96  {
98  }
99 
100  const areaScalarField rhoCph(Cp()*rho()*h_);
101 
103  (
104  fam::ddt(rhoCph, T_)
105  - fam::laplacian(kappa()*h_, T_)
106  ==
107  qs_
108  + qr()
109  + faOptions()(h_, rhoCph, T_)
110  );
111 
112  TEqn.relax();
113 
115 
116  TEqn.solve();
117 
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123 
125 (
126  const word& modelType,
127  const fvMesh& mesh,
128  const dictionary& dict
129 )
130 :
131  thermalShellModel(modelType, mesh, dict),
132  nNonOrthCorr_(1),
133  thermo_(dict.subDict("thermo")),
134  qs_
135  (
136  IOobject
137  (
138  "qs_" + regionName_,
139  primaryMesh().time().timeName(),
140  primaryMesh(),
141  IOobject::READ_IF_PRESENT,
142  IOobject::NO_WRITE
143  ),
144  regionMesh(),
146  ),
147  h_
148  (
149  IOobject
150  (
151  "h_" + regionName_,
152  primaryMesh().time().timeName(),
153  primaryMesh(),
154  IOobject::MUST_READ,
155  IOobject::AUTO_WRITE
156  ),
157  regionMesh()
158  ),
159  qrName_(dict.getOrDefault<word>("qr", "none")),
160  thickness_(dict.getOrDefault<scalar>("thickness", 0))
161 {
162  init(dict);
163 }
164 
166 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
167 
169 {}
170 
171 
173 {
174  nNonOrthCorr_ = solution().get<label>("nNonOrthCorr");
175 
176  for (int nonOrth = 0; nonOrth <= nNonOrthCorr_; ++nonOrth)
177  {
179  }
180 
181  Info<< T_.name() << " min/max = " << gMinMax(T_) << endl;
182 }
183 
184 
186 {
187  return tmp<areaScalarField>
188  (
189  new areaScalarField
190  (
191  IOobject
192  (
193  "Cps",
194  primaryMesh().time().timeName(),
195  primaryMesh(),
198  false
199  ),
200  regionMesh(),
202  zeroGradientFaPatchScalarField::typeName
203  )
204  );
205 }
206 
207 
209 {
210  return tmp<areaScalarField>
211  (
212  new areaScalarField
213  (
214  IOobject
215  (
216  "rhos",
217  primaryMesh().time().timeName(),
218  primaryMesh(),
221  false
222  ),
223  regionMesh(),
225  zeroGradientFaPatchScalarField::typeName
226  )
227  );
228 }
229 
230 
232 {
233  return tmp<areaScalarField>
234  (
235  new areaScalarField
236  (
237  IOobject
238  (
239  "kappas",
240  primaryMesh().time().timeName(),
241  primaryMesh(),
244  false
245  ),
246  regionMesh(),
248  (
250  thermo_.kappa()
251  ),
252  zeroGradientFaPatchScalarField::typeName
253  )
254  );
255 }
256 
257 
258 void thermalShell::info()
259 {}
260 
261 
262 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263 
264 } // End namespace regionModels
265 } // End namespace Foam
266 
267 // ************************************************************************* //
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:120
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
const volSurfaceMapping & vsm() const
Return mapping between surface and volume fields.
areaScalarField T_
Shell temperature.
virtual void info()
Provide some feedback.
Definition: thermalShell.C:251
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
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:201
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.
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:178
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:161
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
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
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:165
const faMesh & regionMesh() const
Return the region mesh database.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:203
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:1100
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:79
void solveEnergy()
Solve energy equation.
Definition: thermalShell.C:86
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:118
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
const tmp< areaScalarField > kappa() const
Return thermal conductivity [W/m/K].
Definition: thermalShell.C:224
areaScalarField h_
Thickness.
Definition: thermalShell.H:154
defineTypeNameAndDebug(KirchhoffShell, 0)
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:157