PDRFoamAutoRefine.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-2021 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 Application
25  PDRFoam
26 
27 Description
28  Solver for compressible premixed/partially-premixed combustion with
29  turbulence modelling.
30 
31  Combusting RANS code using the b-Xi two-equation model.
32  Xi may be obtained by either the solution of the Xi transport
33  equation or from an algebraic expression. Both approaches are
34  based on Gulder's flame speed correlation which has been shown
35  to be appropriate by comparison with the results from the
36  spectral model.
37 
38  Strain effects are incorporated directly into the Xi equation
39  but not in the algebraic approximation. Further work need to be
40  done on this issue, particularly regarding the enhanced removal rate
41  caused by flame compression. Analysis using results of the spectral
42  model will be required.
43 
44  For cases involving very lean Propane flames or other flames which are
45  very strain-sensitive, a transport equation for the laminar flame
46  speed is present. This equation is derived using heuristic arguments
47  involving the strain time scale and the strain-rate at extinction.
48  the transport velocity is the same as that for the Xi equation.
49 
50  For large flames e.g. explosions additional modelling for the flame
51  wrinkling due to surface instabilities may be applied.
52 
53  PDR (porosity/distributed resistance) modelling is included to handle
54  regions containing blockages which cannot be resolved by the mesh.
55 
56 \*---------------------------------------------------------------------------*/
57 
58 #include "fvCFD.H"
59 #include "psiuReactionThermo.H"
62 #include "laminarFlameSpeed.H"
63 #include "XiModel.H"
64 #include "PDRDragModel.H"
65 #include "ignition.H"
66 #include "Switch.H"
67 #include "bound.H"
68 #include "pimpleControl.H"
69 #include "fvModels.H"
70 #include "fvConstraints.H"
71 
72 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
74 int main(int argc, char *argv[])
75 {
76  #include "postProcess.H"
77 
78  #include "setRootCaseLists.H"
79  #include "createTime.H"
80  #include "createMesh.H"
81  #include "createControl.H"
82  #include "readCombustionProperties.H"
83  #include "readGravitationalAcceleration.H"
84  #include "createFields.H"
85  #include "createFieldRefs.H"
86  #include "initContinuityErrs.H"
87  #include "createTimeControls.H"
88  #include "compressibleCourantNo.H"
89  #include "setInitialDeltaT.H"
90 
91  turbulence->validate();
92  scalar StCoNum = 0.0;
93 
94  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
95 
96  Info<< "\nStarting time loop\n" << endl;
97 
98  bool hasChanged = false;
99 
100  while (pimple.run(runTime))
101  {
102  #include "readTimeControls.H"
103  #include "compressibleCourantNo.H"
104  #include "setDeltaT.H"
105 
106  // Indicators for refinement. Note: before runTime++
107  // only for postprocessing reasons.
108  tmp<volScalarField> tmagGradP = mag(fvc::grad(p));
109  volScalarField normalisedGradP
110  (
111  "normalisedGradP",
112  tmagGradP()/max(tmagGradP())
113  );
114  normalisedGradP.writeOpt() = IOobject::AUTO_WRITE;
115  tmagGradP.clear();
116 
117  runTime++;
118 
119  Info<< "\n\nTime = " << runTime.timeName() << endl;
120 
121  {
122  // Make the fluxes absolute
124 
125  // Test : disable refinement for some cells
126  // PackedBoolList& protectedCell =
127  // refCast<dynamicRefineFvMesh>(mesh).protectedCell();
128 
129  // if (protectedCell.empty())
130  // {
131  // protectedCell.setSize(mesh.nCells());
132  // protectedCell = 0;
133  // }
134 
135  // forAll(betav, celli)
136  // {
137  // if (betav[celli] < 0.99)
138  // {
139  // protectedCell[celli] = 1;
140  // }
141  // }
142 
143  // Flux estimate for introduced faces.
144  volVectorField rhoU("rhoU", rho*U);
145 
146  // Do any mesh changes
147  bool meshChanged = mesh.update();
148 
149 
150  if (meshChanged)
151  {
152  hasChanged = true;
153  }
154 
155  if (runTime.write() && hasChanged)
156  {
157  betav.write();
158  Lobs.write();
159  CT.write();
160  drag->writeFields();
161  flameWrinkling->writeFields();
162  hasChanged = false;
163  }
164 
165  // Make the fluxes relative to the mesh motion
167  }
168 
169 
170  #include "rhoEqn.H"
171 
172  // --- Pressure-velocity PIMPLE corrector loop
173  while (pimple.loop())
174  {
175  fvModels.correct();
176 
177  #include "UEqn.H"
178 
179  // --- Pressure corrector loop
180  while (pimple.correct())
181  {
182  #include "bEqn.H"
183  #include "ftEqn.H"
184  #include "EauEqn.H"
185  #include "EaEqn.H"
186 
187  if (!ign.ignited())
188  {
189  thermo.heu() == thermo.he();
190  }
191 
192  #include "pEqn.H"
193  }
194 
195  if (pimple.turbCorr())
196  {
197  turbulence->correct();
198  thermophysicalTransport->correct();
199  }
200  }
201 
202  runTime.write();
203 
204  Info<< "\nExecutionTime = "
205  << runTime.elapsedCpuTime()
206  << " s\n" << endl;
207  }
208 
209  Info<< "\n end\n";
210 
211  return 0;
212 }
213 
214 
215 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
fluidReactionThermo & thermo
Definition: createFields.H:28
pimpleNoLoopControl & pimple
U
Definition: pEqn.H:72
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:353
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Info<< "Reading strained laminar flame speed field Su\"<< endl;volScalarField Su(IOobject("Su", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field betav\"<< endl;volScalarField betav(IOobject("betav", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Lobs\"<< endl;volScalarField Lobs(IOobject("Lobs", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field CT\"<< endl;volSymmTensorField CT(IOobject("CT", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Nv\"<< endl;volScalarField Nv(IOobject("Nv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field nsv\"<< endl;volSymmTensorField nsv(IOobject("nsv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);IOdictionary PDRProperties(IOobject("PDRProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));autoPtr< PDRDragModel > drag
Definition: createFields.H:192
Read the control parameters used by setDeltaT.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
autoPtr< XiModel > flameWrinkling
Create the flame-wrinkling model.
Definition: createFields.H:202
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
phi
Definition: correctPhi.H:3
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Bound the given scalar field if it has gone unbounded.
fluidReactionThermophysicalTransportModel & thermophysicalTransport
Foam::fvModels & fvModels
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:128
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:89
Read the control parameters used by setDeltaT.