reactingOneDim.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-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2020 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 "reactingOneDim.H"
31 #include "fvm.H"
32 #include "fvcDiv.H"
33 #include "fvcVolumeIntegrate.H"
34 #include "fvcLaplacian.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace regionModels
42 {
43 namespace pyrolysisModels
44 {
45 
46 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47 
48 defineTypeNameAndDebug(reactingOneDim, 0);
49 
50 addToRunTimeSelectionTable(pyrolysisModel, reactingOneDim, mesh);
51 addToRunTimeSelectionTable(pyrolysisModel, reactingOneDim, dictionary);
52 
53 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54 
55 void reactingOneDim::readReactingOneDimControls()
56 {
57  const dictionary& solution = this->solution().subDict("SIMPLE");
58  solution.readEntry("nNonOrthCorr", nNonOrthCorr_);
59  time().controlDict().readEntry("maxDi", maxDiff_);
60  coeffs().readEntry("minimumDelta", minimumDelta_);
61  gasHSource_ = coeffs().getOrDefault("gasHSource", false);
62  coeffs().readEntry("qrHSource", qrHSource_);
64  coeffs().getOrDefault("useChemistrySolvers", true);
65 }
66 
67 
69 {
71  {
72  readReactingOneDimControls();
73  return true;
74  }
75 
76  return false;
77 }
78 
79 
81 {
83  {
84  readReactingOneDimControls();
85  return true;
86  }
87 
88  return false;
89 }
90 
91 
93 {
94  // Update local qr from coupled qr field
95  qr_ == dimensionedScalar(qr_.dimensions(), Zero);
96 
97  // Retrieve field from coupled region using mapped boundary conditions
99 
100  volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
101 
103  {
104  const label patchi = intCoupledPatchIDs_[i];
105 
106  // qr is positive going in the solid
107  // If the surface is emitting the radiative flux is set to zero
108  qrBf[patchi] = max(qrBf[patchi], scalar(0));
109  }
110 
111  const vectorField& cellC = regionMesh().cellCentres();
112 
114 
115  // Propagate qr through 1-D regions
116  label localPyrolysisFacei = 0;
118  {
119  const label patchi = intCoupledPatchIDs_[i];
120 
121  const scalarField& qrp = qr_.boundaryField()[patchi];
122  const vectorField& Cf = regionMesh().Cf().boundaryField()[patchi];
123 
124  forAll(qrp, facei)
125  {
126  const scalar qr0 = qrp[facei];
127  point Cf0 = Cf[facei];
128  const labelList& cells = boundaryFaceCells_[localPyrolysisFacei++];
129  scalar kappaInt = 0.0;
130  forAll(cells, k)
131  {
132  const label celli = cells[k];
133  const point& Cf1 = cellC[celli];
134  const scalar delta = mag(Cf1 - Cf0);
135  kappaInt += kappa()[celli]*delta;
136  qr_[celli] = qr0*exp(-kappaInt);
137  Cf0 = Cf1;
138  }
139  }
140  }
141 }
142 
143 
145 {
146  phiHsGas_ == dimensionedScalar(phiHsGas_.dimensions(), Zero);
147  phiGas_ == dimensionedScalar(phiGas_.dimensions(), Zero);
148 
149  const speciesTable& gasTable = solidChemistry_->gasTable();
150 
151  forAll(gasTable, gasI)
152  {
153  tmp<volScalarField> tHsiGas =
154  solidChemistry_->gasHs(solidThermo_->p(), solidThermo_->T(), gasI);
155 
156  const volScalarField& HsiGas = tHsiGas();
157 
158  const volScalarField::Internal& RRiGas = solidChemistry_->RRg(gasI);
159 
160  surfaceScalarField::Boundary& phiGasBf = phiGas_.boundaryFieldRef();
161 
162  label totalFaceId = 0;
164  {
165  const label patchi = intCoupledPatchIDs_[i];
166 
167  scalarField& phiGasp = phiGasBf[patchi];
168  const scalarField& cellVol = regionMesh().V();
169 
170  forAll(phiGasp, facei)
171  {
172  const labelList& cells = boundaryFaceCells_[totalFaceId++];
173  scalar massInt = 0.0;
175  {
176  const label celli = cells[k];
177  massInt += RRiGas[celli]*cellVol[celli];
178  phiHsGas_[celli] += massInt*HsiGas[celli];
179  }
180 
181  phiGasp[facei] += massInt;
182 
183  if (debug)
184  {
185  Info<< " Gas : " << gasTable[gasI]
186  << " on patch : " << patchi
187  << " mass produced at face(local) : "
188  << facei
189  << " is : " << massInt
190  << " [kg/s] " << endl;
191  }
192  }
193  }
194  }
195 }
196 
197 
199 {
200  if (qrHSource_)
201  {
202  updateqr();
203  }
204 
205  //Note: Commented out as the sensible gas energy is included in energy eq.
206  //updatePhiGas();
207 }
208 
209 
211 {
212  Info<< "Initial/final volumes = " << gSum(deltaV) << endl;
213 
214  // Move the mesh
215  const labelList moveMap = moveMesh(deltaV, minimumDelta_);
216 
217  // Flag any cells that have not moved as non-reacting
218  forAll(moveMap, i)
219  {
220  if (moveMap[i] == 1)
221  {
222  solidChemistry_->setCellReacting(i, false);
223  }
224  }
225 }
226 
227 
229 {
231 
232  if (!moveMesh_)
233  {
234  fvScalarMatrix rhoEqn
235  (
236  fvm::ddt(rho_) == -solidChemistry_->RRg()
237  );
238 
239  rhoEqn.solve();
240  }
241  else
242  {
243  const scalarField deltaV
244  (
246  );
247 
248  updateMesh(deltaV);
249  }
250 }
251 
252 
254 {
256 
257  volScalarField Yt(0.0*Ys_[0]);
258 
259  for (label i=0; i<Ys_.size()-1; i++)
260  {
261  volScalarField& Yi = Ys_[i];
262 
263  fvScalarMatrix YiEqn
264  (
265  fvm::ddt(rho_, Yi) == solidChemistry_->RRs(i)
266  );
267 
268  if (regionMesh().moving())
269  {
270  surfaceScalarField phiYiRhoMesh
271  (
273  );
274 
275  YiEqn -= fvc::div(phiYiRhoMesh);
276 
277  }
278 
279  YiEqn.solve(regionMesh().solver("Yi"));
280  Yi.max(0.0);
281  Yt += Yi;
282  }
283 
284  Ys_[Ys_.size() - 1] = 1.0 - Yt;
285 
286 }
287 
288 
290 {
292 
294 
295  fvScalarMatrix hEqn
296  (
297  fvm::ddt(rho_, h_)
300  - fvc::laplacian(kappa(), T())
301  ==
303  + solidChemistry_->RRsHs()
304  );
305 
306 /*
307  NOTE: gas Hs is included in hEqn
308 
309  if (gasHSource_)
310  {
311  const surfaceScalarField phiGas(fvc::interpolate(phiHsGas_));
312  hEqn += fvc::div(phiGas);
313  }
314 */
315 
316  if (qrHSource_)
317  {
319  hEqn += fvc::div(phiqr);
320  }
321 
322 /*
323  NOTE: The moving mesh option is only correct for reaction such as
324  Solid -> Gas, thus the ddt term is compensated exactly by chemistrySh and
325  the mesh flux is not necessary.
326 
327  if (regionMesh().moving())
328  {
329  surfaceScalarField phihMesh
330  (
331  fvc::interpolate(rho_*h_)*regionMesh().phi()
332  );
333 
334  hEqn -= fvc::div(phihMesh);
335  }
336 */
337  hEqn.relax();
338  hEqn.solve();
339 }
340 
341 
343 {
344  /*
345  totalGasMassFlux_ = 0;
346  forAll(intCoupledPatchIDs_, i)
347  {
348  const label patchi = intCoupledPatchIDs_[i];
349  totalGasMassFlux_ += gSum(phiGas_.boundaryField()[patchi]);
350  }
351  */
352 
353  if (infoOutput_)
354  {
356 
357  addedGasMass_ +=
359  lostSolidMass_ +=
361  }
362 }
363 
364 
365 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
366 
367 reactingOneDim::reactingOneDim
368 (
369  const word& modelType,
370  const fvMesh& mesh,
371  const word& regionType
372 )
373 :
374  pyrolysisModel(modelType, mesh, regionType),
375  solidThermo_(solidReactionThermo::New(regionMesh())),
376  solidChemistry_(basicSolidChemistryModel::New(solidThermo_())),
377  radiation_(radiation::radiationModel::New(solidThermo_->T())),
378  rho_
379  (
380  IOobject
381  (
382  "rho",
383  regionMesh().time().timeName(),
384  regionMesh(),
387  ),
388  solidThermo_->rho()
389  ),
390  Ys_(solidThermo_->composition().Y()),
391  h_(solidThermo_->he()),
392  nNonOrthCorr_(-1),
393  maxDiff_(10),
394  minimumDelta_(1e-4),
395 
396  phiGas_
397  (
398  IOobject
399  (
400  "phiGas",
401  time().timeName(),
402  regionMesh(),
405  ),
406  regionMesh(),
408  ),
409 
410  phiHsGas_
411  (
412  IOobject
413  (
414  "phiHsGas",
415  time().timeName(),
416  regionMesh(),
419  ),
420  regionMesh(),
422  ),
423 
424  chemistryQdot_
425  (
426  IOobject
427  (
428  "chemistryQdot",
429  time().timeName(),
430  regionMesh(),
433  ),
434  regionMesh(),
436  ),
437 
438  qr_
439  (
440  IOobject
441  (
442  "qr",
443  time().timeName(),
444  regionMesh(),
447  ),
448  regionMesh()
449  ),
450 
451  lostSolidMass_(dimensionedScalar(dimMass, Zero)),
452  addedGasMass_(dimensionedScalar(dimMass, Zero)),
453  totalGasMassFlux_(0.0),
454  totalHeatRR_(dimensionedScalar(dimEnergy/dimTime, Zero)),
455  gasHSource_(false),
456  qrHSource_(false),
457  useChemistrySolvers_(true)
458 {
459  if (active_)
460  {
461  read();
462  }
463 }
464 
465 
466 reactingOneDim::reactingOneDim
467 (
468  const word& modelType,
469  const fvMesh& mesh,
470  const dictionary& dict,
471  const word& regionType
472 )
473 :
474  pyrolysisModel(modelType, mesh, dict, regionType),
475  solidThermo_(solidReactionThermo::New(regionMesh())),
476  solidChemistry_(basicSolidChemistryModel::New(solidThermo_())),
477  radiation_(radiation::radiationModel::New(solidThermo_->T())),
478  rho_
479  (
480  IOobject
481  (
482  "rho",
483  regionMesh().time().timeName(),
484  regionMesh(),
487  ),
488  solidThermo_->rho()
489  ),
490  Ys_(solidThermo_->composition().Y()),
491  h_(solidThermo_->he()),
492  nNonOrthCorr_(-1),
493  maxDiff_(10),
494  minimumDelta_(1e-4),
495 
496  phiGas_
497  (
498  IOobject
499  (
500  "phiGas",
501  time().timeName(),
502  regionMesh(),
505  ),
506  regionMesh(),
508  ),
509 
510  phiHsGas_
511  (
512  IOobject
513  (
514  "phiHsGas",
515  time().timeName(),
516  regionMesh(),
519  ),
520  regionMesh(),
522  ),
523 
524  chemistryQdot_
525  (
526  IOobject
527  (
528  "chemistryQdot",
529  time().timeName(),
530  regionMesh(),
533  ),
534  regionMesh(),
536  ),
537 
538  qr_
539  (
540  IOobject
541  (
542  "qr",
543  time().timeName(),
544  regionMesh(),
547  ),
548  regionMesh()
549  ),
550 
551  lostSolidMass_(dimensionedScalar(dimMass, Zero)),
552  addedGasMass_(dimensionedScalar(dimMass, Zero)),
553  totalGasMassFlux_(0.0),
554  totalHeatRR_(dimensionedScalar(dimEnergy/dimTime, Zero)),
555  gasHSource_(false),
556  qrHSource_(false),
557  useChemistrySolvers_(true)
558 {
559  if (active_)
560  {
561  read(dict);
562  }
563 }
564 
565 
566 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
567 
569 {}
570 
571 
572 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
573 
574 scalar reactingOneDim::addMassSources(const label patchi, const label facei)
575 {
576  label index = 0;
578  {
579  if (primaryPatchIDs_[i] == patchi)
580  {
581  index = i;
582  break;
583  }
584  }
585 
586  const label localPatchId = intCoupledPatchIDs_[index];
587 
588  const scalar massAdded = phiGas_.boundaryField()[localPatchId][facei];
589 
590  if (debug)
591  {
592  Info<< "\nPyrolysis region: " << type() << "added mass : "
593  << massAdded << endl;
594  }
595 
596  return massAdded;
597 }
598 
599 
601 {
602  scalar DiNum = -GREAT;
603 
604  if (regionMesh().nInternalFaces() > 0)
605  {
606  surfaceScalarField KrhoCpbyDelta
607  (
611  );
612 
613  DiNum = max(KrhoCpbyDelta.primitiveField())*time().deltaTValue();
614  }
615 
616  return returnReduce(DiNum, maxOp<scalar>());
617 }
618 
619 
621 {
622  return maxDiff_;
623 }
624 
625 
627 {
628  return rho_;
629 }
630 
631 
633 {
634  return solidThermo_->T();
635 }
636 
637 
639 {
640  return solidThermo_->Cp();
641 }
642 
643 
645 {
646  return radiation_->absorptionEmission().a();
647 }
648 
649 
651 {
652  return solidThermo_->kappa();
653 }
654 
655 
657 {
658  return phiGas_;
659 }
660 
661 
663 {
665 }
666 
667 
669 {
670  Info<< "\nEvolving pyrolysis in region: " << regionMesh().name() << endl;
671 
673  {
674  solidChemistry_->solve(time().deltaTValue());
675  }
676  else
677  {
678  solidChemistry_->calculate();
679  }
680 
681  solveContinuity();
682 
683  chemistryQdot_ = solidChemistry_->Qdot()();
684 
685  updateFields();
686 
688 
689  for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
690  {
691  solveEnergy();
692  }
693 
695 
696  solidThermo_->correct();
697 
698  Info<< "pyrolysis min/max(T) = "
699  << gMin(solidThermo_->T().primitiveField())
700  << ", "
701  << gMax(solidThermo_->T().primitiveField())
702  << endl;
703 }
704 
705 
707 {
708  Info<< "\nPyrolysis in region: " << regionMesh().name() << endl;
709 
710  Info<< indent << "Total gas mass produced [kg] = "
711  << addedGasMass_.value() << nl
712  << indent << "Total solid mass lost [kg] = "
713  << lostSolidMass_.value() << nl
714  //<< indent << "Total pyrolysis gases [kg/s] = "
715  //<< totalGasMassFlux_ << nl
716  << indent << "Total heat release rate [J/s] = "
717  << totalHeatRR_.value() << nl;
718 }
719 
720 
721 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
722 
723 } // End namespace Foam
724 } // End namespace regionModels
725 } // End namespace pyrolysisModels
726 
727 // ************************************************************************* //
Foam::regionModels::pyrolysisModels::pyrolysisModel
Base class for pyrolysis models.
Definition: pyrolysisModel.H:62
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::regionModels::pyrolysisModels::reactingOneDim::solidThermo_
autoPtr< solidReactionThermo > solidThermo_
Reference to solid thermo.
Definition: reactingOneDim.H:83
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::maxOp
Definition: ops.H:223
Foam::regionModels::pyrolysisModels::reactingOneDim::~reactingOneDim
virtual ~reactingOneDim()
Destructor.
Definition: reactingOneDim.C:568
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::regionModels::pyrolysisModels::reactingOneDim::kappa
virtual tmp< volScalarField > kappa() const
Return the region thermal conductivity [W/m/k].
Definition: reactingOneDim.C:650
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::regionModels::pyrolysisModels::reactingOneDim::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve region.
Definition: reactingOneDim.C:662
Foam::regionModels::pyrolysisModels::reactingOneDim::phiGas
virtual const surfaceScalarField & phiGas() const
Return the total gas mass flux to primary region [kg/m2/s].
Definition: reactingOneDim.C:656
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::regionModels::pyrolysisModels::reactingOneDim::solveSpeciesMass
void solveSpeciesMass()
Solve solid species mass conservation.
Definition: reactingOneDim.C:253
Foam::regionModels::pyrolysisModels::reactingOneDim::maxDiff_
scalar maxDiff_
Maximum diffusivity.
Definition: reactingOneDim.H:110
Foam::TimeState::deltaT
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::regionModels::pyrolysisModels::reactingOneDim::useChemistrySolvers_
bool useChemistrySolvers_
Use chemistry solvers (ode or sequential)
Definition: reactingOneDim.H:162
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::polyMesh::moving
bool moving() const
Is mesh moving.
Definition: polyMesh.H:522
Foam::fvc::domainIntegrate
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcVolumeIntegrate.C:88
Foam::dictionary::dictionary
dictionary()
Default construct, a top-level empty dictionary.
Definition: dictionary.C:81
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::regionModels::pyrolysisModels::reactingOneDim::updatePhiGas
void updatePhiGas()
Update enthalpy flux for pyrolysis gases.
Definition: reactingOneDim.C:144
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Definition: GeometricFieldI.H:53
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::regionModels::pyrolysisModels::reactingOneDim::h_
volScalarField & h_
Definition: reactingOneDim.H:101
fvcDiv.H
Calculate the divergence of the given field.
Foam::regionModels::regionModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:40
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::regionModels::pyrolysisModels::reactingOneDim::solveContinuity
void solveContinuity()
Solve continuity equation.
Definition: reactingOneDim.C:228
Foam::regionModels::pyrolysisModels::reactingOneDim::solidRegionDiffNo
virtual scalar solidRegionDiffNo() const
Mean diffusion number of the solid region.
Definition: reactingOneDim.C:600
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::regionModels::pyrolysisModels::reactingOneDim::addedGasMass_
dimensionedScalar addedGasMass_
Cumulative mass generation of the gas phase [kg].
Definition: reactingOneDim.H:144
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::regionModels::pyrolysisModels::reactingOneDim::rho_
volScalarField rho_
Density [kg/m3].
Definition: reactingOneDim.H:95
Foam::regionModels::pyrolysisModels::reactingOneDim::gasHSource_
bool gasHSource_
Add gas enthalpy source term.
Definition: reactingOneDim.H:156
Foam::surfaceInterpolation::deltaCoeffs
virtual const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
Definition: surfaceInterpolation.C:113
Foam::regionModels::pyrolysisModels::reactingOneDim::updateFields
void updateFields()
Update submodels.
Definition: reactingOneDim.C:198
Foam::regionModels::pyrolysisModels::reactingOneDim::chemistryQdot_
volScalarField chemistryQdot_
Heat release rate [J/s/m3].
Definition: reactingOneDim.H:125
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::regionModels::regionModel1D::boundaryFaceCells_
labelListList boundaryFaceCells_
Global cell IDs.
Definition: regionModel1D.H:83
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:47
Foam::regionModels::pyrolysisModels::defineTypeNameAndDebug
defineTypeNameAndDebug(noPyrolysis, 0)
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::regionModels::pyrolysisModels::reactingOneDim::calculateMassTransfer
void calculateMassTransfer()
Mass check.
Definition: reactingOneDim.C:342
Foam::Field< vector >
Foam::regionModels::regionModel1D::moveMesh
tmp< labelField > moveMesh(const scalarList &deltaV, const scalar minDelta=0.0)
Move mesh points according to change in cell volumes.
Definition: regionModel1D.C:195
Foam::regionModels::pyrolysisModels::reactingOneDim::maxDiff
virtual scalar maxDiff() const
Return max diffusivity allowed in the solid.
Definition: reactingOneDim.C:620
Foam::hashedWordList
A wordList with hashed named lookup, which can be faster in some situations than using the normal lis...
Definition: hashedWordList.H:54
Yt
volScalarField Yt(0.0 *Y[0])
Foam::radiation::radiationModel::New
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
Definition: radiationModelNew.C:36
Foam::solver
Base class for solution control classes.
Definition: solver.H:51
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:365
Foam::regionModels::pyrolysisModels::reactingOneDim::nNonOrthCorr_
label nNonOrthCorr_
Number of non-orthogonal correctors.
Definition: reactingOneDim.H:107
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::regionModels::pyrolysisModels::reactingOneDim::T
virtual const volScalarField & T() const
Return const temperature [K].
Definition: reactingOneDim.C:632
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
DiNum
scalar DiNum
Definition: solidRegionDiffusionNo.H:1
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:528
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
Foam::regionModels::pyrolysisModels::reactingOneDim::Ys_
PtrList< volScalarField > & Ys_
List of solid components.
Definition: reactingOneDim.H:98
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::regionModels::pyrolysisModels::reactingOneDim::read
bool read()
Read control parameters from dictionary.
Definition: reactingOneDim.C:68
Foam::regionModels::pyrolysisModels::reactingOneDim::rho
const volScalarField & rho() const
Fields.
Definition: reactingOneDim.C:626
Foam::regionModels::regionModel::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:500
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
Foam::regionModels::pyrolysisModels::pyrolysisModel::read
virtual bool read()
Read control parameters.
Definition: pyrolysisModel.C:52
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::regionModels::pyrolysisModels::reactingOneDim::minimumDelta_
scalar minimumDelta_
Minimum delta for combustion.
Definition: reactingOneDim.H:113
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::pyrolysisModels::reactingOneDim::kappaRad
virtual tmp< volScalarField > kappaRad() const
Return the region absorptivity [1/m].
Definition: reactingOneDim.C:644
Foam::regionModels::regionModel::solution
const dictionary & solution() const
Return the solution dictionary.
Definition: regionModelI.H:99
Foam::regionModels::pyrolysisModels::reactingOneDim::totalHeatRR_
dimensionedScalar totalHeatRR_
Total heat release rate [J/s].
Definition: reactingOneDim.H:150
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::regionModels::regionModel1D::nMagSf
const surfaceScalarField & nMagSf() const
Return the face area magnitudes / [m2].
Definition: regionModel1DI.H:55
fvm.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::Time::controlDict
const dictionary & controlDict() const
Return read access to the controlDict dictionary.
Definition: Time.H:364
Foam::regionModels::regionModel::infoOutput_
Switch infoOutput_
Active information output.
Definition: regionModel.H:96
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::regionModels::pyrolysisModels::reactingOneDim::addMassSources
virtual scalar addMassSources(const label patchi, const label facei)
External hook to add mass to the primary region.
Definition: reactingOneDim.C:574
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
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::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:320
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::regionModels::pyrolysisModels::reactingOneDim::evolveRegion
virtual void evolveRegion()
Evolve the pyrolysis equations.
Definition: reactingOneDim.C:668
Foam::basicSolidChemistryModel::New
static autoPtr< basicSolidChemistryModel > New(solidReactionThermo &thermo)
Selector.
Definition: basicSolidChemistryModelNew.C:34
Foam::regionModels::pyrolysisModels::reactingOneDim::phiGas_
surfaceScalarField phiGas_
Total gas mass flux to the primary region [kg/m2/s].
Definition: reactingOneDim.H:119
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::regionModels::regionModel::time_
const Time & time_
Reference to the time database.
Definition: regionModel.H:90
fvcLaplacian.H
Calculate the laplacian of the given field.
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::regionModels::pyrolysisModels::reactingOneDim::phiHsGas_
volScalarField phiHsGas_
Sensible enthalpy gas flux [J/m2/s].
Definition: reactingOneDim.H:122
absorptionEmissionModel.H
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::Vector< scalar >
Foam::regionModels::pyrolysisModels::reactingOneDim::updateMesh
void updateMesh(const scalarField &mass0)
Update/move mesh based on change in mass.
Definition: reactingOneDim.C:210
Foam::TimeState::deltaTValue
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:42
Foam::List< label >
Foam::regionModels::pyrolysisModels::reactingOneDim::qrHSource_
bool qrHSource_
Add in depth radiation source term.
Definition: reactingOneDim.H:159
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::regionModels::pyrolysisModels::reactingOneDim::lostSolidMass_
dimensionedScalar lostSolidMass_
Cumulative lost mass of the condensed phase [kg].
Definition: reactingOneDim.H:141
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::regionModels::pyrolysisModels::reactingOneDim::info
virtual void info()
Provide some feedback.
Definition: reactingOneDim.C:706
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
fvcVolumeIntegrate.H
Volume integrate volField creating a volField.
Foam::regionModels::pyrolysisModels::reactingOneDim::radiation_
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
Definition: reactingOneDim.H:89
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:301
Foam::regionModels::regionModel1D::moveMesh_
Switch moveMesh_
Flag to allow mesh movement.
Definition: regionModel1D.H:98
Foam::regionModels::pyrolysisModels::reactingOneDim::qr_
volScalarField qr_
Coupled region radiative heat flux [W/m2].
Definition: reactingOneDim.H:135
Foam::regionModels::regionModel::intCoupledPatchIDs_
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:114
Foam::GeometricField::max
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
Definition: GeometricField.C:1143
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:309
Foam::regionModels::regionModel::coeffs
const dictionary & coeffs() const
Return the model coefficients dictionary.
Definition: regionModelI.H:92
Foam::regionModels::pyrolysisModels::reactingOneDim::solveEnergy
void solveEnergy()
Solve energy.
Definition: reactingOneDim.C:289
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::regionModels::pyrolysisModels::reactingOneDim::solidChemistry_
autoPtr< basicSolidChemistryModel > solidChemistry_
Reference to the solid chemistry model.
Definition: reactingOneDim.H:86
Foam::solidReactionThermo::New
static autoPtr< solidReactionThermo > New(const fvMesh &, const word &phaseName=word::null)
Standard selection based on fvMesh.
Definition: solidReactionThermo.C:67
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
reactingOneDim.H
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::fvMesh::Cf
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
Definition: fvMeshGeometry.C:352
Foam::regionModels::pyrolysisModels::reactingOneDim::updateqr
void updateqr()
Update radiative flux in pyrolysis region.
Definition: reactingOneDim.C:92
Foam::regionModels::pyrolysisModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(pyrolysisModel, noPyrolysis, mesh)
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::regionModels::pyrolysisModels::reactingOneDim::Cp
virtual const tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
Definition: reactingOneDim.C:638
Foam::fvMesh::name
const word & name() const
Return reference to name.
Definition: fvMesh.H:295
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:179
Foam::regionModels::regionModel::primaryPatchIDs_
labelList primaryPatchIDs_
List of patch IDs on the primary region coupled to this region.
Definition: regionModel.H:111
Foam::IOobject::MUST_READ
Definition: IOobject.H:120