thermalBaffle.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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  (
79  (
80  "tQ",
81  regionMesh(),
83  )
84  );
85 
86  volScalarField& Q = tQ.ref();
87 
88  volScalarField rho("rho", thermo_->rho());
89  volScalarField alpha("alpha", thermo_->alpha());
90 
91 
92  // If region is one-dimension variable thickness
93  if (oneD_ && !constantThickness_)
94  {
95  // Scale K and rhoCp and fill Q in the internal baffle region.
96  const label patchi = intCoupledPatchIDs_[0];
97  const polyPatch& ppCoupled = rbm[patchi];
98 
99  forAll(ppCoupled, localFacei)
100  {
101  const labelList& cells = boundaryFaceCells_[localFacei];
102  forAll(cells, i)
103  {
104  const label cellId = cells[i];
105 
106  Q[cellId] =
107  Qs_.boundaryField()[patchi][localFacei]
108  /thickness_[localFacei];
109 
110  rho[cellId] *= delta_.value()/thickness_[localFacei];
111 
112  alpha[cellId] *= delta_.value()/thickness_[localFacei];
113  }
114  }
115  }
116  else
117  {
118  Q = Q_;
119  }
120 
121  fvScalarMatrix hEqn
122  (
123  fvm::ddt(rho, h_)
125  ==
126  Q
127  );
128 
129  if (moveMesh_)
130  {
131  surfaceScalarField phiMesh
132  (
134  );
135 
136  hEqn -= fvc::div(phiMesh);
137  }
138 
139  hEqn.relax();
140  hEqn.solve();
141 
142  thermo_->correct();
143 
144  Info<< "T min/max = " << min(thermo_->T()) << ", "
145  << max(thermo_->T()) << endl;
146 }
147 
148 
149 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
150 
152 (
153  const word& modelType,
154  const fvMesh& mesh,
155  const dictionary& dict
156 )
157 :
158  thermalBaffleModel(modelType, mesh, dict),
159  nNonOrthCorr_(readLabel(solution().lookup("nNonOrthCorr"))),
161  h_(thermo_->he()),
162  Qs_
163  (
164  IOobject
165  (
166  "Qs",
167  regionMesh().time().timeName(),
168  regionMesh(),
171  ),
172  regionMesh(),
174  ),
175  Q_
176  (
177  IOobject
178  (
179  "Q",
180  regionMesh().time().timeName(),
181  regionMesh(),
184  ),
185  regionMesh(),
187  ),
188  radiation_
189  (
191  (
192  dict.subDict("radiation"),
193  thermo_->T()
194  )
195  )
196 {
197  init();
198  thermo_->correct();
199 }
200 
201 
203 (
204  const word& modelType,
205  const fvMesh& mesh
206 )
207 :
208  thermalBaffleModel(modelType, mesh),
209  nNonOrthCorr_(readLabel(solution().lookup("nNonOrthCorr"))),
211  h_(thermo_->he()),
212  Qs_
213  (
214  IOobject
215  (
216  "Qs",
217  regionMesh().time().timeName(),
218  regionMesh(),
221  ),
222  regionMesh(),
224  ),
225  Q_
226  (
227  IOobject
228  (
229  "Q",
230  regionMesh().time().timeName(),
231  regionMesh(),
234  ),
235  regionMesh(),
237  ),
238  radiation_
239  (
241  (
242  thermo_->T()
243  )
244  )
245 {
246  init();
247  thermo_->correct();
248 }
249 
250 
251 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
252 
254 {}
255 
256 
257 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
258 
259 void thermalBaffle::init()
260 {
261  if (oneD_ && !constantThickness_)
262  {
264  const label Qsb = Qs_.boundaryField()[patchi].size();
265 
266  if (Qsb!= thickness_.size())
267  {
269  << "the boundary field of Qs is "
270  << Qsb << " and " << nl
271  << "the field 'thickness' is " << thickness_.size() << nl
272  << exit(FatalError);
273  }
274  }
275 }
276 
277 
279 {}
280 
281 
283 {
284  for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
285  {
286  solveEnergy();
287  }
288 }
289 
290 
292 {
293  return thermo_->Cp();
294 }
295 
296 
298 {
299  return radiation_->absorptionEmission().a();
300 }
301 
302 
304 {
305  return thermo_->rho();
306 }
307 
308 
310 {
311  return thermo_->kappa();
312 }
313 
314 
316 {
317  return thermo_->T();
318 }
319 
320 
322 {
323  return thermo_;
324 }
325 
326 
328 {
329  const labelList& coupledPatches = intCoupledPatchIDs();
330 
331  forAll(coupledPatches, i)
332  {
333  const label patchi = coupledPatches[i];
334  const fvPatchScalarField& ph = h_.boundaryField()[patchi];
335  const word patchName = regionMesh().boundary()[patchi].name();
336  Info<< indent << "Q : " << patchName << indent <<
337  gSum
338  (
339  mag(regionMesh().Sf().boundaryField()[patchi])
340  * ph.snGrad()
341  * thermo_->alpha().boundaryField()[patchi]
342  ) << endl;
343  }
344 }
345 
346 
347 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
348 
349 } // end namespace thermalBaffleModels
350 } // end namespace regionModels
351 } // end namespace Foam
352 
353 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
autoPtr< solidThermo > thermo_
Solid thermo.
Definition: thermalBaffle.H:79
thermalBaffle(const word &modelType, const fvMesh &mesh)
Construct from components.
volScalarField Qs_
Surface energy source / [J/m2/s].
Definition: thermalBaffle.H:88
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
thermalBaffleModel(const fvMesh &mesh)
Construct null from mesh.
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
autoPtr< radiationModel > radiation_
Pointer to radiation model.
Definition: thermalBaffle.H:97
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:226
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:186
surfaceScalarField & phi
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:158
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:73
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:256
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
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:699
volScalarField & h_
Enthalpy/internal energy.
Definition: thermalBaffle.H:82
Switch moveMesh_
Flag to allow mesh movement.
Definition: regionModel1D.H:96
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
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:91
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:97
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:58
static const char nl
Definition: Ostream.H:265
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:48
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:521
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
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 > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
mesh Sf()
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:111
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.
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/m^3].
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:583
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:540
#define InfoInFunction
Report an information message using Foam::Info.