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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "thermalBaffle.H"
27 
28 #include "fvm.H"
29 #include "fvcDiv.H"
32 #include "fvMatrices.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace regionModels
40 {
41 namespace thermalBaffleModels
42 {
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 defineTypeNameAndDebug(thermalBaffle, 0);
47 
48 addToRunTimeSelectionTable(thermalBaffleModel, thermalBaffle, mesh);
49 addToRunTimeSelectionTable(thermalBaffleModel, thermalBaffle, dictionary);
50 
51 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
52 
54 {
55  this->solution().lookup("nNonOrthCorr") >> nNonOrthCorr_;
56  return regionModel1D::read();
57 }
58 
59 
61 {
62  this->solution().lookup("nNonOrthCorr") >> nNonOrthCorr_;
63  return regionModel1D::read(dict);
64 }
65 
66 
68 {
69  if (debug)
70  {
72  }
73 
74  const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
75 
77  (
78  new volScalarField
79  (
80  IOobject
81  (
82  "tQ",
83  regionMesh().time().timeName(),
84  regionMesh(),
87  false
88  ),
89  regionMesh(),
91  )
92  );
93 
94  volScalarField& Q = tQ.ref();
95 
96  volScalarField rho("rho", thermo_->rho());
97  volScalarField alpha("alpha", thermo_->alpha());
98 
99 
100  //If region is one-dimension variable thickness
101  if (oneD_ && !constantThickness_)
102  {
103  // Scale K and rhoCp and fill Q in the internal baffle region.
104  const label patchi = intCoupledPatchIDs_[0];
105  const polyPatch& ppCoupled = rbm[patchi];
106 
107  forAll(ppCoupled, localFacei)
108  {
109  const labelList& cells = boundaryFaceCells_[localFacei];
110  forAll(cells, i)
111  {
112  const label cellId = cells[i];
113 
114  Q[cellId] =
115  Qs_.boundaryField()[patchi][localFacei]
116  /thickness_[localFacei];
117 
118  rho[cellId] *= delta_.value()/thickness_[localFacei];
119 
120  alpha[cellId] *= delta_.value()/thickness_[localFacei];
121  }
122  }
123  }
124  else
125  {
126  Q = Q_;
127  }
128 
129  fvScalarMatrix hEqn
130  (
131  fvm::ddt(rho, h_)
133  ==
134  Q
135  );
136 
137  if (moveMesh_)
138  {
139  surfaceScalarField phiMesh
140  (
142  );
143 
144  hEqn -= fvc::div(phiMesh);
145  }
146 
147  hEqn.relax();
148  hEqn.solve();
149 
150  thermo_->correct();
151 
152  Info<< "T min/max = " << min(thermo_->T()) << ", "
153  << max(thermo_->T()) << endl;
154 }
155 
156 
157 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
158 
159 thermalBaffle::thermalBaffle
160 (
161  const word& modelType,
162  const fvMesh& mesh,
163  const dictionary& dict
164 )
165 :
166  thermalBaffleModel(modelType, mesh, dict),
167  nNonOrthCorr_(readLabel(solution().lookup("nNonOrthCorr"))),
169  h_(thermo_->he()),
170  Qs_
171  (
172  IOobject
173  (
174  "Qs",
175  regionMesh().time().timeName(),
176  regionMesh(),
179  ),
180  regionMesh(),
182  (
183  "zero",
185  Zero
186  )
187  ),
188  Q_
189  (
190  IOobject
191  (
192  "Q",
193  regionMesh().time().timeName(),
194  regionMesh(),
197  ),
198  regionMesh(),
200  (
201  "zero",
203  Zero
204  )
205  ),
206  radiation_
207  (
209  (
210  dict.subDict("radiation"),
211  thermo_->T()
212  )
213  )
214 {
215  init();
216  thermo_->correct();
217 }
218 
219 
220 thermalBaffle::thermalBaffle
221 (
222  const word& modelType,
223  const fvMesh& mesh
224 )
225 :
226  thermalBaffleModel(modelType, mesh),
227  nNonOrthCorr_(readLabel(solution().lookup("nNonOrthCorr"))),
229  h_(thermo_->he()),
230  Qs_
231  (
232  IOobject
233  (
234  "Qs",
235  regionMesh().time().timeName(),
236  regionMesh(),
239  ),
240  regionMesh(),
242  (
243  "zero",
245  Zero
246  )
247  ),
248  Q_
249  (
250  IOobject
251  (
252  "Q",
253  regionMesh().time().timeName(),
254  regionMesh(),
257  ),
258  regionMesh(),
260  (
261  "zero",
263  Zero
264  )
265  ),
266  radiation_
267  (
269  (
270  thermo_->T()
271  )
272  )
273 {
274  init();
275  thermo_->correct();
276 }
277 
278 
279 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
280 
282 {}
283 
284 
285 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
286 
287 void thermalBaffle::init()
288 {
289  if (oneD_ && !constantThickness_)
290  {
292  const label Qsb = Qs_.boundaryField()[patchi].size();
293 
294  if (Qsb!= thickness_.size())
295  {
297  << "the boundary field of Qs is "
298  << Qsb << " and " << nl
299  << "the field 'thickness' is " << thickness_.size() << nl
300  << exit(FatalError);
301  }
302  }
303 }
304 
305 
307 {}
308 
309 
311 {
312  for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
313  {
314  solveEnergy();
315  }
316 }
317 
318 
320 {
321  return thermo_->Cp();
322 }
323 
324 
326 {
327  return radiation_->absorptionEmission().a();
328 }
329 
330 
332 {
333  return thermo_->rho();
334 }
335 
336 
338 {
339  return thermo_->kappa();
340 }
341 
342 
344 {
345  return thermo_->T();
346 }
347 
348 
350 {
351  return thermo_;
352 }
353 
354 
356 {
357  const labelList& coupledPatches = intCoupledPatchIDs();
358 
359  forAll(coupledPatches, i)
360  {
361  const label patchi = coupledPatches[i];
362  const fvPatchScalarField& ph = h_.boundaryField()[patchi];
363  const word patchName = regionMesh().boundary()[patchi].name();
364  Info<< indent << "Q : " << patchName << indent <<
365  gSum
366  (
367  mag(regionMesh().Sf().boundaryField()[patchi])
368  * ph.snGrad()
369  * thermo_->alpha().boundaryField()[patchi]
370  ) << endl;
371  }
372 }
373 
374 
375 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
376 
377 } // end namespace thermalBaffleModels
378 } // end namespace regionModels
379 } // end namespace Foam
380 
381 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
surfaceScalarField & phi
autoPtr< solidThermo > thermo_
Solid thermo.
Definition: thermalBaffle.H:87
volScalarField Qs_
Surface energy source / [J/m2/s].
Definition: thermalBaffle.H:96
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:59
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:223
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:216
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
label nNonOrthCorr_
Number of non orthogonal correctors.
Definition: thermalBaffle.H:81
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual const solidThermo & thermo() const
Return const reference to the solidThermo.
static autoPtr< solidThermo > New(const fvMesh &, const word &phaseName=word::null)
Return a pointer to a new solidThermo created from.
Definition: solidThermo.C:92
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
Macros for easy insertion into run-time selection tables.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
volScalarField & h_
Enthalpy/internal energy.
Definition: thermalBaffle.H:90
Switch moveMesh_
Flag to allow mesh movement.
Definition: regionModel1D.H:96
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
virtual void evolveRegion()
Evolve the thermal baffle.
Type gSum(const FieldField< Field, Type > &f)
const cellShapeList & cells
virtual const volScalarField & kappaRad() const
Return solid absortivity [1/m].
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
volScalarField Q_
Volumetric energy source / [J/m3/s].
Definition: thermalBaffle.H:99
addToRunTimeSelectionTable(thermalBaffleModel, noThermo, mesh)
A class for handling words, derived from string.
Definition: word.H:59
labelListList boundaryFaceCells_
Global cell IDs.
Definition: regionModel1D.H:81
const dictionary & solution() const
Return the solution dictionary.
Definition: regionModelI.H:103
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
const Type & value() const
Return const reference to value.
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:91
label readLabel(Istream &is)
Definition: label.H:64
virtual const tmp< volScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
Foam::polyBoundaryMesh.
Calculate the divergence of the given field.
virtual bool read()
Read control parameters IO dictionary.
Definition: thermalBaffle.C:53
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:57
static const char nl
Definition: Ostream.H:262
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:54
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:523
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionSet dimEnergy
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.
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
virtual bool read()
Read control parameters from dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label cellId
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual const volScalarField & T() const
Return temperature [K].
A special matrix type and solver, designed for finite volume solutions of scalar equations.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:114
A class for managing temporary objects.
Definition: PtrList.H:53
virtual const volScalarField & kappa() const
Return thermal conductivity [W/m/K].
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:179
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual void info()
Provide some feedback.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual const volScalarField & rho() const
Return density [Kg/m3].
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545
#define InfoInFunction
Report an information message using Foam::Info.