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-------------------------------------------------------------------------------
10License
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
35namespace Foam
36{
37namespace regionModels
38{
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
43
45
46// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47
49{
50 if (thickness_ > 0)
51 {
53 }
54
55 this->solution().readEntry("nNonOrthCorr", nNonOrthCorr_);
56
57 return true;
58}
59
60
61tmp<areaScalarField> thermalShell::qr()
62{
64 (
65 "tqr",
68 );
69
70 auto taqr =
72 (
73 io,
74 regionMesh(),
76 );
77
78 if (qrName_ != "none")
79 {
80 auto& aqr = taqr.ref();
81
83
84 const volScalarField::Boundary& vqr = qr.boundaryField();
85
86 aqr.primitiveFieldRef() = vsm().mapToSurface<scalar>(vqr);
87 }
88
89 return taqr;
90}
91
92
93// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
94
96{
97 if (debug)
98 {
100 }
101
102 const areaScalarField rhoCph(Cp()*rho()*h_);
103
105 (
106 fam::ddt(rhoCph, T_)
108 ==
109 qs_
110 + qr()
111 + faOptions()(h_, rhoCph, T_)
112 );
113
114 TEqn.relax();
115
117
118 TEqn.solve();
119
121}
122
123
124// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
125
127(
128 const word& modelType,
129 const fvPatch& patch,
130 const dictionary& dict
131)
132:
133 thermalShellModel(modelType, patch, dict),
134 nNonOrthCorr_(1),
135 thermo_(dict.subDict("thermo")),
136 qs_
137 (
139 (
140 "qs_" + regionName_,
141 primaryMesh().time().timeName(),
142 primaryMesh(),
143 IOobject::READ_IF_PRESENT,
144 IOobject::NO_WRITE
145 ),
146 regionMesh(),
148 ),
149 h_
150 (
152 (
153 "h_" + regionName_,
154 primaryMesh().time().timeName(),
155 primaryMesh(),
156 IOobject::MUST_READ,
157 IOobject::AUTO_WRITE
158 ),
159 regionMesh()
160 ),
161 qrName_(dict.getOrDefault<word>("qr", "none")),
162 thickness_(dict.getOrDefault<scalar>("thickness", 0))
163{
164 init(dict);
165}
166
167
168// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
169
171{}
172
173
175{
176 nNonOrthCorr_ = solution().get<label>("nNonOrthCorr");
177
178 for (int nonOrth = 0; nonOrth <= nNonOrthCorr_; ++nonOrth)
179 {
180 solveEnergy();
181 }
182
183 Info<< "T min/max = " << min(T_) << ", " << max(T_) << endl;
184}
185
186
188{
190 (
192 (
194 (
195 "Cps",
197 primaryMesh(),
200 false
201 ),
202 regionMesh(),
204 zeroGradientFaPatchScalarField::typeName
205 )
206 );
207}
208
209
211{
213 (
215 (
217 (
218 "rhos",
220 primaryMesh(),
223 false
224 ),
225 regionMesh(),
227 zeroGradientFaPatchScalarField::typeName
228 )
229 );
230}
231
232
234{
236 (
238 (
240 (
241 "kappas",
243 primaryMesh(),
246 false
247 ),
248 regionMesh(),
250 (
252 thermo_.kappa()
253 ),
254 zeroGradientFaPatchScalarField::typeName
255 )
256 );
257}
258
259
261{}
262
263
264// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265
266} // End namespace regionModels
267} // End namespace Foam
268
269// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition: faMatrix.H:76
void correct(GeometricField< Type, faPatchField, areaMesh > &field)
Apply correction to field.
void constrain(faMatrix< Type > &eqn)
Apply constraints to equation.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1092
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
const Type & lookupObject(const word &name, const bool recursive=false) const
const Time & time() const
Return the reference to the time database.
const dictionary & solution() const
Return the solution dictionary.
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
const volSurfaceMapping & vsm() const
Return mapping between surface and volume fields.
const faMesh & regionMesh() const
Return the region mesh database.
areaScalarField T_
Shell temperature.
Foam::fa::options & faOptions() noexcept
Return faOptions.
areaScalarField h_
Thickness.
Definition: thermalShell.H:146
const tmp< areaScalarField > rho() const
Return density [Kg/m3].
Definition: thermalShell.C:210
scalar thickness_
Uniform thickness.
Definition: thermalShell.H:152
void solveEnergy()
Solve energy equation.
Definition: thermalShell.C:95
areaScalarField qs_
External surface energy source / [J/m2/s].
Definition: thermalShell.H:143
const word qrName_
Name of the primary region radiative flux.
Definition: thermalShell.H:149
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
Definition: thermalShell.C:170
solidProperties thermo_
Solid properties.
Definition: thermalShell.H:137
const tmp< areaScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
Definition: thermalShell.C:187
virtual void info()
Provide some feedback.
Definition: thermalShell.C:260
label nNonOrthCorr_
Number of non orthogonal correctors.
Definition: thermalShell.H:131
const tmp< areaScalarField > kappa() const
Return thermal conductivity [W/m/K].
Definition: thermalShell.C:233
virtual void evolveRegion()
Evolve the thermal baffle.
Definition: thermalShell.C:174
scalar Cp() const
Specific heat capacity [J/(kg.K)].
scalar kappa() const
Thermal conductivity [W/(m.K)].
scalar rho() const
Density [kg/m3].
A class for managing temporary objects.
Definition: tmp.H:65
tmp< Field< Type > > mapToSurface(const GeometricBoundaryField< Type, fvPatchField, volMesh > &df) const
Map volume boundary field to surface.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
word timeName
Definition: getTimeIndex.H:3
#define InfoInFunction
Report an information message using Foam::Info.
tmp< faMatrix< Type > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famLaplacian.C:49
tmp< faMatrix< Type > > ddt(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famDdt.C:48
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimPower
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const dimensionSet dimDensity
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dictionary dict