Loading...
Searching...
No Matches
thermalBaffle.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "thermalBaffle.H"
30#include "fvm.H"
31#include "fvcDiv.H"
34#include "fvMatrices.H"
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39namespace Foam
40{
41namespace regionModels
42{
45
46// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47
48defineTypeNameAndDebug(thermalBaffle, 0);
49
50addToRunTimeSelectionTable(thermalBaffleModel, thermalBaffle, mesh);
51addToRunTimeSelectionTable(thermalBaffleModel, thermalBaffle, dictionary);
52
53// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54
56{
57 this->solution().readEntry("nNonOrthCorr", nNonOrthCorr_);
58 return regionModel1D::read();
59}
60
61
63{
64 this->solution().readEntry("nNonOrthCorr", nNonOrthCorr_);
66}
67
68
70{
72
74
75 auto tQ = volScalarField::New
76 (
77 "tQ",
79 regionMesh(),
81 );
82 auto& Q = tQ.ref();
83
84 volScalarField rho("rho", thermo_->rho());
85 volScalarField alpha("alpha", thermo_->alpha());
86
87
88 // If region is one-dimension variable thickness
90 {
91 // Scale K and rhoCp and fill Q in the internal baffle region.
92 const label patchi = intCoupledPatchIDs_[0];
93 const polyPatch& ppCoupled = rbm[patchi];
94
95 forAll(ppCoupled, localFacei)
96 {
97 const labelList& cells = boundaryFaceCells_[localFacei];
98 forAll(cells, i)
99 {
100 const label cellId = cells[i];
101
102 Q[cellId] =
103 qs_.boundaryField()[patchi][localFacei]
104 /thickness_[localFacei];
105
106 rho[cellId] *= delta_.value()/thickness_[localFacei];
107
108 alpha[cellId] *= delta_.value()/thickness_[localFacei];
109 }
110 }
111 }
112 else
113 {
114 Q = Q_;
115 }
116
117 fvScalarMatrix hEqn
118 (
119 fvm::ddt(rho, h_)
121 ==
122 Q
123 );
124
125 if (moveMesh_)
126 {
127 surfaceScalarField phiMesh
128 (
130 );
131
132 hEqn -= fvc::div(phiMesh);
133 }
134
135 hEqn.relax();
136 hEqn.solve();
137
138 thermo_->correct();
139
140 Info<< "T min/max = " << min(thermo_->T()) << ", "
141 << max(thermo_->T()) << endl;
142}
143
144
145// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
146
147thermalBaffle::thermalBaffle
148(
149 const word& modelType,
150 const fvMesh& mesh,
151 const dictionary& dict
152)
153:
154 thermalBaffleModel(modelType, mesh, dict),
155 nNonOrthCorr_(solution().get<label>("nNonOrthCorr")),
156 thermo_(solidThermo::New(regionMesh(), dict)),
157 h_(thermo_->he()),
158 qs_
159 (
161 (
162 "qs",
163 regionMesh().time().timeName(),
164 regionMesh().thisDb(),
165 IOobject::READ_IF_PRESENT,
166 IOobject::NO_WRITE
167 ),
168 regionMesh(),
170 ),
171 Q_
172 (
174 (
175 "Q",
176 regionMesh().time().timeName(),
177 regionMesh().thisDb(),
178 IOobject::READ_IF_PRESENT,
179 IOobject::AUTO_WRITE
180 ),
181 regionMesh(),
183 ),
184 radiation_
185 (
187 (
188 dict.subDict("radiation"),
189 thermo_->T()
190 )
192{
193 init();
194 thermo_->correct();
195}
196
197
198thermalBaffle::thermalBaffle
199(
200 const word& modelType,
201 const fvMesh& mesh
202)
203:
204 thermalBaffleModel(modelType, mesh),
205 nNonOrthCorr_(solution().get<label>("nNonOrthCorr")),
206 thermo_(solidThermo::New(regionMesh())),
207 h_(thermo_->he()),
208 qs_
209 (
211 (
212 "qs",
213 regionMesh().time().timeName(),
214 regionMesh().thisDb(),
215 IOobject::READ_IF_PRESENT,
216 IOobject::NO_WRITE
217 ),
218 regionMesh(),
220 ),
221 Q_
222 (
224 (
225 "Q",
226 regionMesh().time().timeName(),
227 regionMesh().thisDb(),
228 IOobject::READ_IF_PRESENT,
229 IOobject::NO_WRITE
230 ),
231 regionMesh(),
233 ),
234 radiation_
235 (
237 (
238 thermo_->T()
239 )
240 )
241{
242 init();
243 thermo_->correct();
244}
245
246
247// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
248
250{}
251
252
253// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
254
255void thermalBaffle::init()
256{
258 {
259 label patchi = intCoupledPatchIDs_[0];
260 const label qsb = qs_.boundaryField()[patchi].size();
261
262 if (qsb!= thickness_.size())
263 {
265 << "the boundary field of qs is "
266 << qsb << " and " << nl
267 << "the field 'thickness' is " << thickness_.size() << nl
268 << exit(FatalError);
269 }
270 }
272
273
275{}
276
277
279{
280 for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
281 {
282 solveEnergy();
283 }
284}
285
288{
289 return thermo_->Cp();
290}
291
294{
295 return radiation_->absorptionEmission().a();
296}
297
300{
301 return thermo_->rho();
302}
303
306{
307 return thermo_->kappa();
308}
309
311const volScalarField& thermalBaffle::T() const
312{
313 return thermo_->T();
314}
315
318{
319 return *thermo_;
320}
321
322
324{
325 const labelList& coupledPatches = intCoupledPatchIDs();
326
327 forAll(coupledPatches, i)
328 {
329 const label patchi = coupledPatches[i];
330 const fvPatchScalarField& ph = h_.boundaryField()[patchi];
331 const word patchName = regionMesh().boundary()[patchi].name();
332 Info<< indent << "Q : " << patchName << indent <<
333 gSum
334 (
335 mag(regionMesh().Sf().boundaryField()[patchi])
336 * ph.snGrad()
337 * thermo_->alpha().boundaryField()[patchi]
338 ) << endl;
339 }
340}
341
342
343// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
344
345} // End namespace thermalBaffleModels
346} // End namespace regionModels
347} // End namespace Foam
348
349// ************************************************************************* //
volScalarField & he
Definition YEEqn.H:52
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
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,...
const Type & value() const noexcept
Return const reference to value.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition fvMatrix.C:1094
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Top level model for radiation modelling.
Switch moveMesh_
Flag to allow mesh movement.
labelListList boundaryFaceCells_
Global cell IDs.
virtual bool read()
Read control parameters from dictionary.
const dictionary & solution() const
Return the solution dictionary.
const Time & time() const noexcept
Return the reference to the time database.
const fvMesh & regionMesh() const
Return the region mesh database.
const labelList & intCoupledPatchIDs() const noexcept
List of patch IDs internally coupled with the primary region.
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
static autoPtr< thermalBaffleModel > New(const fvMesh &mesh)
Return a reference to the selected model.
virtual const solidThermo & thermo() const
Return const reference to the solidThermo.
virtual const volScalarField & rho() const
Return density [Kg/m3].
volScalarField qs_
Surface energy source / [J/m2/s].
volScalarField & h_
Enthalpy/internal energy.
virtual const volScalarField & kappa() const
Return thermal conductivity [W/m/K].
virtual const volScalarField & T() const
Return temperature [K].
virtual const tmp< volScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
volScalarField Q_
Volumetric energy source / [J/m3/s].
tmp< scalarField > he(const scalarField &p, const scalarField &T, const label patchi) const
Return sensible enthalpy/internal energy.
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
virtual const volScalarField & kappaRad() const
Return solid absortivity [1/m].
autoPtr< solidThermo > thermo_
Solid thermo.
label nNonOrthCorr_
Number of non orthogonal correctors.
virtual bool read()
Read control parameters IO dictionary.
virtual void evolveRegion()
Evolve the thermal baffle.
Fundamental solid thermodynamic properties.
Definition solidThermo.H:51
Selector class for relaxation factors, solver type and solution.
Definition solution.H:95
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the divergence of the given field.
const cellShapeList & cells
label cellId
word timeName
Definition getTimeIndex.H:3
#define DebugInFunction
Report an information message using Foam::Info.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
Namespace for radiation modelling.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gSum(const FieldField< Field, Type > &f)
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimEnergy
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299