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 -------------------------------------------------------------------------------
11 License
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 
39 namespace Foam
40 {
41 namespace regionModels
42 {
43 namespace thermalBaffleModels
44 {
45 
46 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47 
49 
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 {
72 
73  const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
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
89  if (oneD_ && !constantThickness_)
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 
147 thermalBaffle::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  (
160  IOobject
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  (
173  IOobject
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  )
191  )
192 {
193  init();
194  thermo_->correct();
195 }
196 
197 
198 thermalBaffle::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  (
210  IOobject
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  (
223  IOobject
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 
255 void thermalBaffle::init()
256 {
257  if (oneD_ && !constantThickness_)
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  }
271 }
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 
299 const volScalarField& thermalBaffle::rho() const
300 {
301  return thermo_->rho();
302 }
303 
306 {
307  return thermo_->kappa();
308 }
309 
311 const volScalarField& thermalBaffle::T() const
312 {
313  return thermo_->T();
314 }
315 
317 const solidThermo& thermalBaffle::thermo() const
318 {
319  return *thermo_;
320 }
321 
322 
323 void thermalBaffle::info()
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 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
volScalarField & he
Definition: YEEqn.H:52
autoPtr< solidThermo > thermo_
Solid thermo.
Definition: thermalBaffle.H:95
const Type & value() const noexcept
Return const reference to value.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:493
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:213
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
label nNonOrthCorr_
Number of non orthogonal correctors.
Definition: thermalBaffle.H:87
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
virtual const solidThermo & thermo() const
Return const reference to the solidThermo.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const labelList & intCoupledPatchIDs() const noexcept
List of patch IDs internally coupled with the primary region.
Definition: regionModelI.H:97
volScalarField & h_
Enthalpy/internal energy.
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...
Switch moveMesh_
Flag to allow mesh movement.
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
virtual void evolveRegion()
Evolve the thermal baffle.
dynamicFvMesh & mesh
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].
addToRunTimeSelectionTable(thermalBaffleModel, noThermo, mesh)
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
#define DebugInFunction
Report an information message using Foam::Info.
labelListList boundaryFaceCells_
Global cell IDs.
Definition: regionModel1D.H:91
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const dictionary & solution() const
Return the solution dictionary.
Definition: regionModelI.H:55
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:26
Top level model for radiation modelling.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
virtual const tmp< volScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
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())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
Calculate the divergence of the given field.
virtual bool read()
Read control parameters IO dictionary.
Definition: thermalBaffle.C:48
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:48
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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.
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:41
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
Information stream (stdout output on master, null elsewhere)
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:92
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:135
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual const volScalarField & kappa() const
Return thermal conductivity [W/m/K].
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
virtual void info()
Provide some feedback.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
virtual const volScalarField & rho() const
Return density [Kg/m3].
Do not request registration (bool: false)
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
volScalarField qs_
Surface energy source / [J/m2/s].
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127