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 | Copyright (C) 2011-2014 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 "reactingOneDim.H"
29 #include "surfaceInterpolate.H"
30 #include "fvm.H"
31 #include "fvcDiv.H"
32 #include "fvcVolumeIntegrate.H"
33 #include "fvMatrices.H"
35 #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.lookup("nNonOrthCorr") >> nNonOrthCorr_;
59  time().controlDict().lookup("maxDi") >> maxDiff_;
60 
61  coeffs().lookup("minimumDelta") >> minimumDelta_;
62 
63  coeffs().lookup("gasHSource") >> gasHSource_;
64  coeffs().lookup("QrHSource") >> QrHSource_;
66  coeffs().lookupOrDefault<bool>("useChemistrySolvers", true);
67 
68 }
69 
70 
72 {
74  {
75  readReactingOneDimControls();
76  return true;
77  }
78  else
79  {
80  return false;
81  }
82 }
83 
84 
86 {
87  if (pyrolysisModel::read(dict))
88  {
89  readReactingOneDimControls();
90  return true;
91  }
92  else
93  {
94  return false;
95  }
96 }
97 
98 
100 {
101  // Update local Qr from coupled Qr field
102  Qr_ == dimensionedScalar("zero", Qr_.dimensions(), 0.0);
103 
104  // Retrieve field from coupled region using mapped boundary conditions
106 
108  {
109  const label patchI = intCoupledPatchIDs_[i];
110 
111  scalarField& Qrp = Qr_.boundaryField()[patchI];
112 
113  // Qr is positive going in the solid
114  // If the surface is emitting the radiative flux is set to zero
115  Qrp = max(Qrp, scalar(0.0));
116  }
117 
118  const vectorField& cellC = regionMesh().cellCentres();
119 
121 
122  // Propagate Qr through 1-D regions
123  label localPyrolysisFaceI = 0;
125  {
126  const label patchI = intCoupledPatchIDs_[i];
127 
128  const scalarField& Qrp = Qr_.boundaryField()[patchI];
129  const vectorField& Cf = regionMesh().Cf().boundaryField()[patchI];
130 
131  forAll(Qrp, faceI)
132  {
133  const scalar Qr0 = Qrp[faceI];
134  point Cf0 = Cf[faceI];
135  const labelList& cells = boundaryFaceCells_[localPyrolysisFaceI++];
136  scalar kappaInt = 0.0;
137  forAll(cells, k)
138  {
139  const label cellI = cells[k];
140  const point& Cf1 = cellC[cellI];
141  const scalar delta = mag(Cf1 - Cf0);
142  kappaInt += kappa()[cellI]*delta;
143  Qr_[cellI] = Qr0*exp(-kappaInt);
144  Cf0 = Cf1;
145  }
146  }
147  }
148 }
149 
150 
152 {
153  phiHsGas_ == dimensionedScalar("zero", phiHsGas_.dimensions(), 0.0);
154  phiGas_ == dimensionedScalar("zero", phiGas_.dimensions(), 0.0);
155 
156  const speciesTable& gasTable = solidChemistry_->gasTable();
157 
158  forAll(gasTable, gasI)
159  {
160  tmp<volScalarField> tHsiGas =
161  solidChemistry_->gasHs(solidThermo_.p(), solidThermo_.T(), gasI);
162 
163  const volScalarField& HsiGas = tHsiGas();
164 
165  const DimensionedField<scalar, volMesh>& RRiGas =
166  solidChemistry_->RRg(gasI);
167 
168  label totalFaceId = 0;
170  {
171  const label patchI = intCoupledPatchIDs_[i];
172 
173  scalarField& phiGasp = phiGas_.boundaryField()[patchI];
174  const scalarField& cellVol = regionMesh().V();
175 
176  forAll(phiGasp, faceI)
177  {
178  const labelList& cells = boundaryFaceCells_[totalFaceId++];
179  scalar massInt = 0.0;
180  forAllReverse(cells, k)
181  {
182  const label cellI = cells[k];
183  massInt += RRiGas[cellI]*cellVol[cellI];
184  phiHsGas_[cellI] += massInt*HsiGas[cellI];
185  }
186 
187  phiGasp[faceI] += massInt;
188 
189  if (debug)
190  {
191  Info<< " Gas : " << gasTable[gasI]
192  << " on patch : " << patchI
193  << " mass produced at face(local) : "
194  << faceI
195  << " is : " << massInt
196  << " [kg/s] " << endl;
197  }
198  }
199  }
200  tHsiGas().clear();
201  }
202 }
203 
204 
206 {
207  if (QrHSource_)
208  {
209  updateQr();
210  }
211 
212  updatePhiGas();
213 }
214 
215 
217 {
218  if (!moveMesh_)
219  {
220  return;
221  }
222 
223  const scalarField newV(mass0/rho_);
224 
225  Info<< "Initial/final volumes = " << gSum(regionMesh().V()) << ", "
226  << gSum(newV) << " [m3]" << endl;
227 
228  // move the mesh
229  const labelList moveMap = moveMesh(regionMesh().V() - newV, minimumDelta_);
230 
231  // flag any cells that have not moved as non-reacting
232  forAll(moveMap, i)
233  {
234  if (moveMap[i] == 0)
235  {
236  solidChemistry_->setCellReacting(i, false);
237  }
238  }
239 }
240 
241 
243 {
244  if (debug)
245  {
246  Info<< "reactingOneDim::solveContinuity()" << endl;
247  }
248 
249  const scalarField mass0 = rho_*regionMesh().V();
250 
251  fvScalarMatrix rhoEqn
252  (
253  fvm::ddt(rho_)
254  ==
255  - solidChemistry_->RRg()
256  );
257 
258  if (regionMesh().moving())
259  {
260  surfaceScalarField phiRhoMesh
261  (
263  );
264 
265  rhoEqn += fvc::div(phiRhoMesh);
266  }
267 
268  rhoEqn.solve();
269 
270  updateMesh(mass0);
271 }
272 
273 
275 {
276  if (debug)
277  {
278  Info<< "reactingOneDim::solveSpeciesMass()" << endl;
279  }
280 
281  volScalarField Yt(0.0*Ys_[0]);
282 
283  for (label i=0; i<Ys_.size()-1; i++)
284  {
285  volScalarField& Yi = Ys_[i];
286 
287  fvScalarMatrix YiEqn
288  (
289  fvm::ddt(rho_, Yi)
290  ==
291  solidChemistry_->RRs(i)
292  );
293 
294  if (regionMesh().moving())
295  {
296  surfaceScalarField phiYiRhoMesh
297  (
299  );
300 
301  YiEqn += fvc::div(phiYiRhoMesh);
302 
303  }
304 
305  YiEqn.solve(regionMesh().solver("Yi"));
306  Yi.max(0.0);
307  Yt += Yi;
308  }
309 
310  Ys_[Ys_.size() - 1] = 1.0 - Yt;
311 
312 }
313 
314 
316 {
317  if (debug)
318  {
319  Info<< "reactingOneDim::solveEnergy()" << endl;
320  }
321 
323 
324  fvScalarMatrix hEqn
325  (
326  fvm::ddt(rho_, h_)
329  - fvc::laplacian(kappa(), T())
330  ==
332  - fvm::Sp(solidChemistry_->RRg(), h_)
333  );
334 
335  if (gasHSource_)
336  {
338  hEqn += fvc::div(phiGas);
339  }
340 
341  if (QrHSource_)
342  {
344  hEqn += fvc::div(phiQr);
345  }
346 
347  if (regionMesh().moving())
348  {
349  surfaceScalarField phihMesh
350  (
352  );
353 
354  hEqn += fvc::div(phihMesh);
355  }
356 
357  hEqn.relax();
358  hEqn.solve();
359 }
360 
361 
363 {
364  totalGasMassFlux_ = 0;
366  {
367  const label patchI = intCoupledPatchIDs_[i];
369  }
370 
371  if (infoOutput_)
372  {
374 
375  addedGasMass_ +=
377  lostSolidMass_ +=
379  }
380 }
381 
382 
383 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
384 
385 reactingOneDim::reactingOneDim
386 (
387  const word& modelType,
388  const fvMesh& mesh,
389  const word& regionType
390 )
391 :
392  pyrolysisModel(modelType, mesh, regionType),
394  solidThermo_(solidChemistry_->solidThermo()),
396  rho_
397  (
398  IOobject
399  (
400  "rho",
401  regionMesh().time().timeName(),
402  regionMesh(),
405  ),
406  solidThermo_.rho()
407  ),
409  h_(solidThermo_.he()),
410  nNonOrthCorr_(-1),
411  maxDiff_(10),
412  minimumDelta_(1e-4),
413 
414  phiGas_
415  (
416  IOobject
417  (
418  "phiGas",
419  time().timeName(),
420  regionMesh(),
423  ),
424  regionMesh(),
425  dimensionedScalar("zero", dimMass/dimTime, 0.0)
426  ),
427 
428  phiHsGas_
429  (
430  IOobject
431  (
432  "phiHsGas",
433  time().timeName(),
434  regionMesh(),
437  ),
438  regionMesh(),
439  dimensionedScalar("zero", dimEnergy/dimTime, 0.0)
440  ),
441 
443  (
444  IOobject
445  (
446  "chemistrySh",
447  time().timeName(),
448  regionMesh(),
451  ),
452  regionMesh(),
454  ),
455 
456  Qr_
457  (
458  IOobject
459  (
460  "Qr",
461  time().timeName(),
462  regionMesh(),
465  ),
466  regionMesh()
467  //dimensionedScalar("zero", dimEnergy/dimArea/dimTime, 0.0),
468  //zeroGradientFvPatchVectorField::typeName
469  ),
470 
472  addedGasMass_(dimensionedScalar("zero", dimMass, 0.0)),
473  totalGasMassFlux_(0.0),
475  gasHSource_(false),
476  QrHSource_(false),
478 {
479  if (active_)
480  {
481  read();
482  }
483 }
484 
485 
486 reactingOneDim::reactingOneDim
487 (
488  const word& modelType,
489  const fvMesh& mesh,
490  const dictionary& dict,
491  const word& regionType
492 )
493 :
494  pyrolysisModel(modelType, mesh, dict, regionType),
496  solidThermo_(solidChemistry_->solidThermo()),
498  rho_
499  (
500  IOobject
501  (
502  "rho",
503  regionMesh().time().timeName(),
504  regionMesh(),
507  ),
508  solidThermo_.rho()
509  ),
511  h_(solidThermo_.he()),
512  nNonOrthCorr_(-1),
513  maxDiff_(10),
514  minimumDelta_(1e-4),
515 
516  phiGas_
517  (
518  IOobject
519  (
520  "phiGas",
521  time().timeName(),
522  regionMesh(),
525  ),
526  regionMesh(),
527  dimensionedScalar("zero", dimMass/dimTime, 0.0)
528  ),
529 
530  phiHsGas_
531  (
532  IOobject
533  (
534  "phiHsGas",
535  time().timeName(),
536  regionMesh(),
539  ),
540  regionMesh(),
541  dimensionedScalar("zero", dimEnergy/dimTime, 0.0)
542  ),
543 
545  (
546  IOobject
547  (
548  "chemistrySh",
549  time().timeName(),
550  regionMesh(),
553  ),
554  regionMesh(),
556  ),
557 
558  Qr_
559  (
560  IOobject
561  (
562  "Qr",
563  time().timeName(),
564  regionMesh(),
567  ),
568  regionMesh()
569  ),
570 
572  addedGasMass_(dimensionedScalar("zero", dimMass, 0.0)),
573  totalGasMassFlux_(0.0),
575  gasHSource_(false),
576  QrHSource_(false),
578 {
579  if (active_)
580  {
581  read(dict);
582  }
583 }
584 
585 
586 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
587 
589 {}
590 
591 
592 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
593 
594 scalar reactingOneDim::addMassSources(const label patchI, const label faceI)
595 {
596  label index = 0;
598  {
599  if (primaryPatchIDs_[i] == patchI)
600  {
601  index = i;
602  break;
603  }
604  }
605 
606  const label localPatchId = intCoupledPatchIDs_[index];
607 
608  const scalar massAdded = phiGas_.boundaryField()[localPatchId][faceI];
609 
610  if (debug)
611  {
612  Info<< "\nPyrolysis region: " << type() << "added mass : "
613  << massAdded << endl;
614  }
615 
616  return massAdded;
617 }
618 
619 
621 {
622  scalar DiNum = -GREAT;
623 
624  if (regionMesh().nInternalFaces() > 0)
625  {
626  surfaceScalarField KrhoCpbyDelta
627  (
631  );
632 
633  DiNum = max(KrhoCpbyDelta.internalField())*time().deltaTValue();
634  }
635 
636  return DiNum;
637 }
638 
639 
641 {
642  return maxDiff_;
643 }
644 
645 
647 {
648  return rho_;
649 }
650 
651 
653 {
654  return solidThermo_.T();
655 }
656 
657 
659 {
660  return solidThermo_.Cp();
661 }
662 
663 
665 {
666  return radiation_->absorptionEmission().a();
667 }
668 
669 
671 {
672  return solidThermo_.kappa();
673 }
674 
675 
677 {
678  return phiGas_;
679 }
680 
681 
683 {
685 
686  // Initialise all cells as able to react
687  forAll(h_, cellI)
688  {
689  solidChemistry_->setCellReacting(cellI, true);
690  }
691 }
692 
693 
695 {
696  Info<< "\nEvolving pyrolysis in region: " << regionMesh().name() << endl;
697 
699  {
700  solidChemistry_->solve(time().deltaTValue());
701  }
702  else
703  {
704  solidChemistry_->calculate();
705  }
706 
707  solveContinuity();
708 
709  chemistrySh_ = solidChemistry_->Sh()();
710 
711  updateFields();
712 
714 
715  for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
716  {
717  solveEnergy();
718  }
719 
721 
723 
724  Info<< "pyrolysis min/max(T) = "
726  << ", "
728  << endl;
729 }
730 
731 
733 {
734  Info<< "\nPyrolysis in region: " << regionMesh().name() << endl;
735 
736  Info<< indent << "Total gas mass produced [kg] = "
737  << addedGasMass_.value() << nl
738  << indent << "Total solid mass lost [kg] = "
739  << lostSolidMass_.value() << nl
740  << indent << "Total pyrolysis gases [kg/s] = "
741  << totalGasMassFlux_ << nl
742  << indent << "Total heat release rate [J/s] = "
743  << totalHeatRR_.value() << nl;
744 }
745 
746 
747 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
748 
749 } // End namespace Foam
750 } // End namespace regionModels
751 } // End namespace pyrolysisModels
752 
753 // ************************************************************************* //
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
Definition: solidThermo.C:120
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
labelListList boundaryFaceCells_
Global cell IDs.
Definition: regionModel1D.H:84
dictionary()
Construct top-level dictionary null.
Definition: dictionary.C:113
virtual void info()
Provide some feedback.
virtual scalar solidRegionDiffNo() const
Mean diffusion number of the solid region.
tmp< labelField > moveMesh(const scalarList &deltaV, const scalar minDelta=0.0)
Move mesh points according to change in cell volumes.
PtrList< volScalarField > & Ys_
List of solid components.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void updateQr()
Update radiative flux in pyrolysis region.
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:478
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
label nNonOrthCorr_
Number of non-orthogonal correctors.
virtual tmp< volScalarField > kappa() const
Return the region thermal conductivity [W/m/k].
virtual const tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:174
scalar totalGasMassFlux_
Total mass gas flux at the pyrolysing walls [kg/s].
solidReactionThermo & solidThermo_
Reference to solid thermo.
static autoPtr< basicSolidChemistryModel > New(const fvMesh &mesh, const word &phaseName=word::null)
Selector.
virtual const surfaceScalarField & phiGas() const
Return the total gas mass flux to primary region [kg/m2/s].
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
bool gasHSource_
Add gas enthalpy source term.
Switch active_
Active flag.
Definition: regionModel.H:93
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
virtual void correct()=0
Update properties.
labelList primaryPatchIDs_
List of patch IDs on the primary region coupled to this region.
Definition: regionModel.H:114
virtual tmp< volScalarField > a(const label bandI=0) const
Absorption coefficient (net)
addToRunTimeSelectionTable(pyrolysisModel, noPyrolysis, mesh)
virtual tmp< volScalarField > kappa() const =0
Thermal diffusivity for temperature of mixture [J/m/s/K].
volScalarField Yt(0.0 *Y[0])
virtual scalar maxDiff() const
Return max diffusivity allowed in the solid.
A class for handling words, derived from string.
Definition: word.H:59
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
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
const surfaceScalarField & nMagSf() const
Return the face area magnitudes / [m2].
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:466
scalar DiNum
InternalField & internalField()
Return internal field.
Surface Interpolation.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimEnergy
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61
const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
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:68
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
void solveContinuity()
Solve continuity equation.
const vectorField & cellCentres() const
virtual tmp< volScalarField > kappaRad() const
Return the region absorptivity [1/m].
messageStream Info
dimensionedScalar lostSolidMass_
Cumulative lost mass of the condensed phase [kg].
virtual basicSpecieMixture & composition()=0
Return the composition of the multi-component mixture.
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
void updatePhiGas()
Update enthalpy flux for pyrolysis gases.
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:564
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensionedScalar addedGasMass_
Cumulative mass generation of the gas phase [kg].
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:638
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:117
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
Namespace for OpenFOAM.
Calculate the laplacian of the given field.
surfaceScalarField phiGas_
Total gas mass flux to the primary region [kg/m2/s].
autoPtr< basicSolidChemistryModel > solidChemistry_
Reference to the solid chemistry model.
Type gSum(const FieldField< Field, Type > &f)
bool QrHSource_
Add in depth radiation source term.
dimensionedScalar totalHeatRR_
Total heat release rate [J/s].
static const char nl
Definition: Ostream.H:260
const double e
Elementary charge.
Definition: doubleFloat.H:78
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
scalar deltaTValue() const
Return time step value.
Definition: TimeState.H:100
volScalarField Qr_
Coupled region radiative heat flux [W/m2].
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Definition: UList.H:421
label k
Boltzmann constant.
virtual void preEvolveRegion()
Pre-evolve region.
virtual void evolveRegion()
Evolve the pyrolysis equations.
scalar delta
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
const volScalarField & rho() const
Fields.
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
const absorptionEmissionModel & absorptionEmission() const
Access to absorptionEmission model.
bool useChemistrySolvers_
Use chemistry solvers (ode or sequential)
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
const dictionary & solution() const
Return the solution dictionary.
Definition: regionModelI.H:109
virtual const volScalarField & T() const
Return const temperature [K].
surfaceScalarField & phi
const cellShapeList & cells
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
const dictionary & coeffs() const
Return the model coefficients dictionary.
Definition: regionModelI.H:102
scalar minimumDelta_
Minimum delta for combustion.
#define forAllReverse(list, i)
Definition: UList.H:424
Calculate the divergence of the given field.
volScalarField phiHsGas_
Sensible enthalpy gas flux [J/m2/s].
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
void solveSpeciesMass()
Solve solid species mass conservation.
void correctBoundaryConditions()
Correct boundary field.
Switch moveMesh_
Flag to allow mesh movement.
Definition: regionModel1D.H:99
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
volScalarField chemistrySh_
Heat release [J/s/m3].
virtual scalar addMassSources(const label patchI, const label faceI)
External hook to add mass to the primary region.
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
virtual const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: basicThermo.C:490
const Time & time_
Reference to the time database.
Definition: regionModel.H:90
void updateMesh(const scalarField &mass0)
Update/move mesh based on change in mass.
bool read()
Read control parameters from dictionary.
A wordList with hashed indices for faster lookup by name.
void max(const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
dimensionedScalar deltaT() const
Return time step.
Definition: TimeState.C:79
virtual bool read()
Read control parameters.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
bool moving() const
Is mesh moving.
Definition: polyMesh.H:493
Volume integrate volField creating a volField.
const dictionary & controlDict() const
Definition: Time.H:286
Switch infoOutput_
Active information output.
Definition: regionModel.H:96
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:118
const Type & value() const
Return const reference to value.
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
word timeName
Definition: getTimeIndex.H:3
A special matrix type and solver, designed for finite volume solutions of scalar equations.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].