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 
43 
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_ != "none")
79  {
80  auto& aqr = taqr.ref();
81 
82  const auto qr = primaryMesh().lookupObject<volScalarField>(qrName_);
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_)
107  - fam::laplacian(kappa()*h_, T_)
108  ==
109  qs_
110  + qr()
111  + faOptions()(h_, rhoCph, T_)
112  );
113 
114  TEqn.relax();
115 
117 
118  TEqn.solve();
119 
120  faOptions().correct(T_);
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  (
138  IOobject
139  (
140  "qs_" + regionName_,
141  primaryMesh().time().timeName(),
142  primaryMesh(),
145  ),
146  regionMesh(),
148  ),
149  h_
150  (
151  IOobject
152  (
153  "h_" + regionName_,
154  primaryMesh().time().timeName(),
155  primaryMesh(),
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 {
189  return tmp<areaScalarField>
190  (
191  new areaScalarField
192  (
193  IOobject
194  (
195  "Cps",
196  primaryMesh().time().timeName(),
197  primaryMesh(),
200  false
201  ),
202  regionMesh(),
204  zeroGradientFaPatchScalarField::typeName
205  )
206  );
207 }
208 
209 
211 {
212  return tmp<areaScalarField>
213  (
214  new areaScalarField
215  (
216  IOobject
217  (
218  "rhos",
219  primaryMesh().time().timeName(),
220  primaryMesh(),
223  false
224  ),
225  regionMesh(),
227  zeroGradientFaPatchScalarField::typeName
228  )
229  );
230 }
231 
232 
234 {
235  return tmp<areaScalarField>
236  (
237  new areaScalarField
238  (
239  IOobject
240  (
241  "kappas",
242  primaryMesh().time().timeName(),
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 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::regionModels::regionFaModel::solution
const dictionary & solution() const
Return the solution dictionary.
Definition: regionFaModelI.H:132
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::regionModels::thermalShell::info
virtual void info()
Provide some feedback.
Definition: thermalShell.C:260
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::faMatrix
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition: faMatricesFwd.H:43
Foam::regionModels::thermalShell::qrName_
const word qrName_
Name of the primary region radiative flux.
Definition: thermalShell.H:149
Foam::regionModels::thermalShell::kappa
const tmp< areaScalarField > kappa() const
Return thermal conductivity [W/m/K].
Definition: thermalShell.C:233
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::regionModels::defineTypeNameAndDebug
defineTypeNameAndDebug(KirchhoffShell, 0)
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::regionModels::regionFaModel::primaryMesh
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionFaModelI.H:33
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::regionModels::thermalShell::qs_
areaScalarField qs_
External surface energy source / [J/m2/s].
Definition: thermalShell.H:143
Foam::fam::laplacian
tmp< faMatrix< Type > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famLaplacian.C:49
Foam::dimDensity
const dimensionSet dimDensity
Foam::regionModels::thermalShellModel::T_
areaScalarField T_
Shell temperature.
Definition: thermalShellModel.H:123
Foam::regionModels::thermalShell::Cp
const tmp< areaScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
Definition: thermalShell.C:187
thermalShell.H
Foam::volSurfaceMapping::mapToSurface
tmp< Field< Type > > mapToSurface(const typename GeometricField< Type, fvPatchField, volMesh >::Boundary &df) const
Map Boundary field to surface.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::regionModels::thermalShell::thickness_
scalar thickness_
Uniform thickness.
Definition: thermalShell.H:152
Foam::solidProperties::rho
scalar rho() const
Density [kg/m3].
Definition: solidPropertiesI.H:33
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::regionModels::regionFaModel::vsm
const volSurfaceMapping & vsm() const
Return mapping between surface and volume fields.
Definition: regionFaModel.C:110
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::fam::ddt
tmp< faMatrix< Type > > ddt(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famDdt.C:48
TEqn
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1183
Foam::regionModels::thermalShell::rho
const tmp< areaScalarField > rho() const
Return density [Kg/m3].
Definition: thermalShell.C:210
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
init
mesh init(true)
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::solidProperties::kappa
scalar kappa() const
Thermal conductivity [W/(m.K)].
Definition: solidPropertiesI.H:45
Foam::fa::optionList::constrain
void constrain(faMatrix< Type > &eqn)
Apply constraints to equation.
Definition: faOptionListTemplates.C:249
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::regionModels::thermalShell::thermo_
solidProperties thermo_
Solid properties.
Definition: thermalShell.H:137
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::regionModels::thermalShellModel::faOptions
Foam::fa::options & faOptions() noexcept
Return faOptions.
Definition: thermalShellModel.H:199
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
thermalShell
Thermal-shell finite-area model. It solves the energy equation in 2D. The coupling with the 3D region...
Foam::dimPower
const dimensionSet dimPower
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::regionModels::regionFaModel::regionMesh
const faMesh & regionMesh() const
Return the region mesh database.
Definition: regionFaModelI.H:63
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::regionModels::thermalShell::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
Definition: thermalShell.C:170
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::regionModels::regionFaModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionFaModelI.H:39
zeroGradientFaPatchFields.H
Foam::regionModels::thermalShell::evolveRegion
virtual void evolveRegion()
Evolve the thermal baffle.
Definition: thermalShell.C:174
Foam::regionModels::thermalShell::thermalShell
thermalShell(const word &modelType, const fvPatch &patch, const dictionary &dict)
Construct from components and dict.
Definition: thermalShell.C:127
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
fvPatchFields.H
Foam::regionModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(vibrationShellModel, KirchhoffShell, dictionary)
Foam::regionModels::thermalShell::h_
areaScalarField h_
Thickness.
Definition: thermalShell.H:146
thermalShellModel
Intermediate class for thermal-shell finite-area models.
Foam::regionModels::thermalShellModel
Definition: thermalShellModel.H:108
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:319
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::regionModels::thermalShell::nNonOrthCorr_
label nNonOrthCorr_
Number of non orthogonal correctors.
Definition: thermalShell.H:131
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::fa::optionList::correct
void correct(GeometricField< Type, faPatchField, areaMesh > &field)
Apply correction to field.
Definition: faOptionListTemplates.C:290
Foam::solidProperties::Cp
scalar Cp() const
Specific heat capacity [J/(kg.K)].
Definition: solidPropertiesI.H:39
Foam::regionModels::thermalShell::solveEnergy
void solveEnergy()
Solve energy equation.
Definition: thermalShell.C:95
Foam::IOobject::MUST_READ
Definition: IOobject.H:185