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 -------------------------------------------------------------------------------
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 "thermalBaffle.H"
29 
30 #include "fvm.H"
31 #include "fvcDiv.H"
34 #include "fvMatrices.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace regionModels
42 {
43 namespace thermalBaffleModels
44 {
45 
46 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47 
48 defineTypeNameAndDebug(thermalBaffle, 0);
49 
50 addToRunTimeSelectionTable(thermalBaffleModel, thermalBaffle, mesh);
51 addToRunTimeSelectionTable(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_);
65  return regionModel1D::read(dict);
66 }
67 
68 
70 {
71  if (debug)
72  {
74  }
75 
76  const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
77 
79  (
80  new volScalarField
81  (
82  IOobject
83  (
84  "tQ",
85  regionMesh().time().timeName(),
86  regionMesh(),
89  false
90  ),
91  regionMesh(),
93  )
94  );
95 
96  volScalarField& Q = tQ.ref();
97 
98  volScalarField rho("rho", thermo_->rho());
99  volScalarField alpha("alpha", thermo_->alpha());
100 
101 
102  // If region is one-dimension variable thickness
103  if (oneD_ && !constantThickness_)
104  {
105  // Scale K and rhoCp and fill Q in the internal baffle region.
106  const label patchi = intCoupledPatchIDs_[0];
107  const polyPatch& ppCoupled = rbm[patchi];
108 
109  forAll(ppCoupled, localFacei)
110  {
111  const labelList& cells = boundaryFaceCells_[localFacei];
112  forAll(cells, i)
113  {
114  const label cellId = cells[i];
115 
116  Q[cellId] =
117  qs_.boundaryField()[patchi][localFacei]
118  /thickness_[localFacei];
119 
120  rho[cellId] *= delta_.value()/thickness_[localFacei];
121 
122  alpha[cellId] *= delta_.value()/thickness_[localFacei];
123  }
124  }
125  }
126  else
127  {
128  Q = Q_;
129  }
130 
131  fvScalarMatrix hEqn
132  (
133  fvm::ddt(rho, h_)
135  ==
136  Q
137  );
138 
139  if (moveMesh_)
140  {
141  surfaceScalarField phiMesh
142  (
144  );
145 
146  hEqn -= fvc::div(phiMesh);
147  }
148 
149  hEqn.relax();
150  hEqn.solve();
151 
152  thermo_->correct();
153 
154  Info<< "T min/max = " << min(thermo_->T()) << ", "
155  << max(thermo_->T()) << endl;
156 }
157 
158 
159 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
160 
161 thermalBaffle::thermalBaffle
162 (
163  const word& modelType,
164  const fvMesh& mesh,
165  const dictionary& dict
166 )
167 :
168  thermalBaffleModel(modelType, mesh, dict),
169  nNonOrthCorr_(solution().get<label>("nNonOrthCorr")),
170  thermo_(solidThermo::New(regionMesh(), dict)),
171  h_(thermo_->he()),
172  qs_
173  (
174  IOobject
175  (
176  "qs",
177  regionMesh().time().timeName(),
178  regionMesh(),
181  ),
182  regionMesh(),
184  ),
185  Q_
186  (
187  IOobject
188  (
189  "Q",
190  regionMesh().time().timeName(),
191  regionMesh(),
194  ),
195  regionMesh(),
197  ),
198  radiation_
199  (
201  (
202  dict.subDict("radiation"),
203  thermo_->T()
204  )
205  )
206 {
207  init();
208  thermo_->correct();
209 }
210 
211 
212 thermalBaffle::thermalBaffle
213 (
214  const word& modelType,
215  const fvMesh& mesh
216 )
217 :
218  thermalBaffleModel(modelType, mesh),
219  nNonOrthCorr_(solution().get<label>("nNonOrthCorr")),
220  thermo_(solidThermo::New(regionMesh())),
221  h_(thermo_->he()),
222  qs_
223  (
224  IOobject
225  (
226  "qs",
227  regionMesh().time().timeName(),
228  regionMesh(),
231  ),
232  regionMesh(),
234  ),
235  Q_
236  (
237  IOobject
238  (
239  "Q",
240  regionMesh().time().timeName(),
241  regionMesh(),
244  ),
245  regionMesh(),
247  ),
248  radiation_
249  (
251  (
252  thermo_->T()
253  )
254  )
255 {
256  init();
257  thermo_->correct();
258 }
259 
260 
261 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
262 
264 {}
265 
266 
267 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
268 
269 void thermalBaffle::init()
270 {
271  if (oneD_ && !constantThickness_)
272  {
273  label patchi = intCoupledPatchIDs_[0];
274  const label qsb = qs_.boundaryField()[patchi].size();
275 
276  if (qsb!= thickness_.size())
277  {
279  << "the boundary field of qs is "
280  << qsb << " and " << nl
281  << "the field 'thickness' is " << thickness_.size() << nl
282  << exit(FatalError);
283  }
284  }
285 }
286 
287 
289 {}
290 
291 
293 {
294  for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
295  {
296  solveEnergy();
297  }
298 }
299 
300 
302 {
303  return thermo_->Cp();
304 }
305 
306 
308 {
309  return radiation_->absorptionEmission().a();
310 }
311 
312 
314 {
315  return thermo_->rho();
316 }
317 
318 
320 {
321  return thermo_->kappa();
322 }
323 
324 
326 {
327  return thermo_->T();
328 }
329 
330 
332 {
333  return *thermo_;
334 }
335 
336 
338 {
339  const labelList& coupledPatches = intCoupledPatchIDs();
340 
341  forAll(coupledPatches, i)
342  {
343  const label patchi = coupledPatches[i];
344  const fvPatchScalarField& ph = h_.boundaryField()[patchi];
345  const word patchName = regionMesh().boundary()[patchi].name();
346  Info<< indent << "Q : " << patchName << indent <<
347  gSum
348  (
349  mag(regionMesh().Sf().boundaryField()[patchi])
350  * ph.snGrad()
351  * thermo_->alpha().boundaryField()[patchi]
352  ) << endl;
353  }
354 }
355 
356 
357 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
358 
359 } // End namespace thermalBaffleModels
360 } // End namespace regionModels
361 } // End namespace Foam
362 
363 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< scalar >
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::oneD_
bool oneD_
Is it one dimension.
Definition: thermalBaffleModel.H:89
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::regionModels::regionModel1D::read
virtual bool read()
Read control parameters from dictionary.
Definition: regionModel1D.C:170
Foam::regionModels::thermalBaffleModels::thermalBaffle::T
virtual const volScalarField & T() const
Return temperature [K].
Definition: thermalBaffle.C:325
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:217
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::regionModels::thermalBaffleModels::thermalBaffle::kappa
virtual const volScalarField & kappa() const
Return thermal conductivity [W/m/K].
Definition: thermalBaffle.C:319
Foam::solidThermo::New
static autoPtr< solidThermo > New(const fvMesh &, const word &phaseName=word::null)
Return a pointer to a new solidThermo created from.
Definition: solidThermo.C:82
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::delta_
dimensionedScalar delta_
Baffle mesh thickness.
Definition: thermalBaffleModel.H:86
Foam::tmp< volScalarField >
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::constantThickness_
bool constantThickness_
Is thickness constant.
Definition: thermalBaffleModel.H:92
Foam::regionModels::thermalBaffleModels::thermalBaffle::qs_
volScalarField qs_
Surface energy source / [J/m2/s].
Definition: thermalBaffle.H:98
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::regionModels::thermalBaffleModels::thermalBaffle::h_
volScalarField & h_
Enthalpy/internal energy.
Definition: thermalBaffle.H:92
thermalBaffle.H
fvcDiv.H
Calculate the divergence of the given field.
Foam::regionModels::thermalBaffleModels::thermalBaffle::~thermalBaffle
virtual ~thermalBaffle()
Destructor.
Definition: thermalBaffle.C:263
Foam::regionModels::regionModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:39
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::regionModels::thermalBaffleModels::thermalBaffle::solveEnergy
void solveEnergy()
Solve energy equation.
Definition: thermalBaffle.C:69
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:404
Foam::regionModels::thermalBaffleModels::thermalBaffle::Q_
volScalarField Q_
Volumetric energy source / [J/m3/s].
Definition: thermalBaffle.H:101
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
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::regionModels::thermalBaffleModels::thermalBaffle::thermo
virtual const solidThermo & thermo() const
Return const reference to the solidThermo.
Definition: thermalBaffle.C:331
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::thickness_
scalarField thickness_
Baffle physical thickness.
Definition: thermalBaffleModel.H:83
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::regionModels::regionModel1D::boundaryFaceCells_
labelListList boundaryFaceCells_
Global cell IDs.
Definition: regionModel1D.H:83
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::solidThermo
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:52
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::regionModels::thermalBaffleModels::thermalBaffleModel
Definition: thermalBaffleModel.H:60
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::radiation::radiationModel::New
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
Definition: radiationModelNew.C:36
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:574
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:63
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
Foam::regionModels::thermalBaffleModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(thermalBaffleModel, noThermo, mesh)
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
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
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
Foam::regionModels::regionModel::solution
const dictionary & solution() const
Return the solution dictionary.
Definition: regionModelI.H:106
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
fvm.H
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvc::interpolate
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.
Foam::regionModels::thermalBaffleModels::thermalBaffle::evolveRegion
virtual void evolveRegion()
Evolve the thermal baffle.
Definition: thermalBaffle.C:292
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:307
Foam::regionModels::thermalBaffleModels::thermalBaffle::radiation_
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
Definition: thermalBaffle.H:107
cellId
label cellId
Definition: interrogateWallPatches.H:67
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::regionModels::thermalBaffleModels::thermalBaffle::kappaRad
virtual const volScalarField & kappaRad() const
Return solid absortivity [1/m].
Definition: thermalBaffle.C:307
Foam::regionModels::thermalBaffleModels::thermalBaffle::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
Definition: thermalBaffle.C:288
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:589
Foam::regionModels::thermalBaffleModels::defineTypeNameAndDebug
defineTypeNameAndDebug(noThermo, 0)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::regionModels::thermalBaffleModels::thermalBaffle::rho
virtual const volScalarField & rho() const
Return density [Kg/m3].
Definition: thermalBaffle.C:313
Foam::regionModels::thermalBaffleModels::thermalBaffle::info
virtual void info()
Provide some feedback.
Definition: thermalBaffle.C:337
Foam::regionModels::thermalBaffleModels::thermalBaffle::Cp
virtual const tmp< volScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
Definition: thermalBaffle.C:301
absorptionEmissionModel.H
Foam::List< label >
Foam::regionModels::thermalBaffleModels::thermalBaffle::read
virtual bool read()
Read control parameters IO dictionary.
Definition: thermalBaffle.C:55
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:298
Foam::regionModels::regionModel1D::moveMesh_
Switch moveMesh_
Flag to allow mesh movement.
Definition: regionModel1D.H:98
Foam::regionModels::regionModel::intCoupledPatchIDs_
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:116
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::regionModels::thermalBaffleModels::thermalBaffle::thermo_
autoPtr< solidThermo > thermo_
Solid thermo.
Definition: thermalBaffle.H:89
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
zeroGradientFvPatchFields.H
Foam::regionModels::regionModel::intCoupledPatchIDs
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:182
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::regionModels::thermalBaffleModels::thermalBaffle::nNonOrthCorr_
label nNonOrthCorr_
Number of non orthogonal correctors.
Definition: thermalBaffle.H:83