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 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
49 {
50  this->solution().readEntry("nNonOrthCorr", nNonOrthCorr_);
51 
52  return true;
53 }
54 
55 
57 {
58  if (debug)
59  {
61  }
62 
63  const areaScalarField rhoCph(Cp()*rho()*h_);
64 
66  (
67  fam::ddt(rhoCph, T_)
69  ==
70  qs_
71  //+ q(T_) handled by faOption contactHeatFlux
72  + faOptions()(h_, rhoCph, T_)
73  );
74 
75  TEqn.relax();
76 
78 
79  TEqn.solve();
80 
81  faOptions().correct(T_);
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86 
88 (
89  const word& modelType,
90  const fvPatch& patch,
91  const dictionary& dict
92 )
93 :
94  thermalShellModel(modelType, patch, dict),
95  nNonOrthCorr_(1),
96  thermo_(dict.subDict("thermo")),
97  qs_
98  (
99  IOobject
100  (
101  "qs_" + regionName_,
102  primaryMesh().time().timeName(),
103  primaryMesh(),
106  ),
107  regionMesh(),
109  ),
110  h_
111  (
112  IOobject
113  (
114  "h_" + regionName_,
115  primaryMesh().time().timeName(),
116  primaryMesh(),
119  ),
120  regionMesh()
121  )
122 {
123  init();
124 }
125 
126 
127 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
128 
129 void thermalShell::init()
130 {}
131 
132 
134 {}
135 
136 
138 {
139  nNonOrthCorr_ = solution().get<label>("nNonOrthCorr");
140 
141  for (int nonOrth = 0; nonOrth <= nNonOrthCorr_; ++nonOrth)
142  {
143  solveEnergy();
144  }
145 
146  Info<< "T min/max = " << min(T_) << ", " << max(T_) << endl;
147 }
148 
149 
151 {
152  return tmp<areaScalarField>
153  (
154  new areaScalarField
155  (
156  IOobject
157  (
158  "Cps",
159  primaryMesh().time().timeName(),
160  primaryMesh(),
163  false
164  ),
165  regionMesh(),
167  zeroGradientFaPatchScalarField::typeName
168  )
169  );
170 }
171 
172 
174 {
175  return tmp<areaScalarField>
176  (
177  new areaScalarField
178  (
179  IOobject
180  (
181  "rhos",
182  primaryMesh().time().timeName(),
183  primaryMesh(),
186  false
187  ),
188  regionMesh(),
190  zeroGradientFaPatchScalarField::typeName
191  )
192  );
193 }
194 
195 
197 {
198  return tmp<areaScalarField>
199  (
200  new areaScalarField
201  (
202  IOobject
203  (
204  "kappas",
205  primaryMesh().time().timeName(),
206  primaryMesh(),
209  false
210  ),
211  regionMesh(),
213  (
215  thermo_.kappa()
216  ),
217  zeroGradientFaPatchScalarField::typeName
218  )
219  );
220 }
221 
222 
224 {}
225 
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 } // End namespace regionModels
230 } // End namespace Foam
231 
232 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::regionModels::regionFaModel::solution
const dictionary & solution() const
Return the solution dictionary.
Definition: regionFaModelI.H:106
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::regionModels::thermalShell::info
virtual void info()
Provide some feedback.
Definition: thermalShell.C:223
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:325
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::faMatrix
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition: faMatrix.H:59
Foam::regionModels::thermalShell::kappa
const tmp< areaScalarField > kappa() const
Return thermal conductivity [W/m/K].
Definition: thermalShell.C:196
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::regionModels::defineTypeNameAndDebug
defineTypeNameAndDebug(KirchhoffShell, 0)
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
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:140
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:129
Foam::regionModels::thermalShell::Cp
const tmp< areaScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
Definition: thermalShell.C:150
thermalShell.H
Foam::regionModels::thermalShell::read
virtual bool read(const dictionary &dict)
Read control parameters from dictionary.
Definition: thermalShell.C:48
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
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:81
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::regionModels::thermalShellModel::faOptions
Foam::fa::options & faOptions()
Return faOptions.
Definition: thermalShellModel.C:110
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:559
Foam::regionModels::thermalShell::rho
const tmp< areaScalarField > rho() const
Return density [Kg/m3].
Definition: thermalShell.C:173
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:528
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
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:224
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:314
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::regionModels::thermalShell::thermo_
solidProperties thermo_
Solid properties.
Definition: thermalShell.H:134
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:121
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:133
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:52
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:137
Foam::regionModels::thermalShell::thermalShell
thermalShell(const word &modelType, const fvPatch &patch, const dictionary &dict)
Construct from components and dict.
Definition: thermalShell.C:88
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:143
thermalShellModel
Intermediate class for thermal-shell finite-area models.
Foam::regionModels::thermalShellModel
Definition: thermalShellModel.H:108
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:301
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:55
Foam::regionModels::thermalShell::nNonOrthCorr_
label nNonOrthCorr_
Number of non orthogonal correctors.
Definition: thermalShell.H:128
Foam::GeometricField< scalar, faPatchField, areaMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::fa::optionList::correct
void correct(GeometricField< Type, faPatchField, areaMesh > &field)
Apply correction to field.
Definition: faOptionListTemplates.C:257
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:56
Foam::IOobject::MUST_READ
Definition: IOobject.H:120