intersection with triangle plane for miss
[OpenFOAM-1.5.x.git] / applications / solvers / combustion / PDRFoam / PDRFoamAutoRefine.C
blobf4db015fe750fd31bf33bdca2faa517bd88733f4
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25 Application
26     PDRFoam
28 Description
29     Compressible premixed/partially-premixed combustion solver with turbulence 
30     modelling.
32     Combusting RANS code using the b-Xi two-equation model.
33     Xi may be obtained by either the solution of the Xi transport
34     equation or from an algebraic exression.  Both approaches are
35     based on Gulder's flame speed correlation which has been shown
36     to be appropriate by comparison with the results from the
37     spectral model.
39     Strain effects are encorporated directly into the Xi equation
40     but not in the algebraic approximation.  Further work need to be
41     done on this issue, particularly regarding the enhanced removal rate
42     caused by flame compression.  Analysis using results of the spectral
43     model will be required.
45     For cases involving very lean Propane flames or other flames which are
46     very strain-sensitive, a transport equation for the laminar flame
47     speed is present.  This equation is derived using heuristic arguments
48     involving the strain time scale and the strain-rate at extinction.
49     the transport velocity is the same as that for the Xi equation.
51     For large flames e.g. explosions additional modelling for the flame
52     wrinkling due to surface instabilities may be applied.
54     PDR (porosity/distributed resistance) modelling is included to handle
55     regions containing blockages which cannot be resolved by the mesh.
57 \*---------------------------------------------------------------------------*/
59 #include "fvCFD.H"
60 #include "dynamicFvMesh.H"
61 #include "hhuCombustionThermo.H"
62 #include "RASModel.H"
63 #include "laminarFlameSpeed.H"
64 #include "XiModel.H"
65 #include "PDRDragModel.H"
66 #include "ignition.H"
67 #include "Switch.H"
68 #include "bound.H"
69 #include "dynamicRefineFvMesh.H"
71 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 int main(int argc, char *argv[])
75 #   include "setRootCase.H"
77 #   include "createTime.H"
78 #   include "createDynamicFvMesh.H"
79 #   include "readCombustionProperties.H"
80 #   include "readEnvironmentalProperties.H"
81 #   include "createFields.H"
82 #   include "readPISOControls.H"
83 #   include "initContinuityErrs.H"
84 #   include "readTimeControls.H"
85 #   include "setInitialDeltaT.H"
87 scalar StCoNum = 0.0;
89 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91     Info<< "\nStarting time loop\n" << endl;
93     while (runTime.run())
94     {
95 #       include "readTimeControls.H"
96 #       include "readPISOControls.H"
97 #       include "CourantNo.H"
99 #       include "setDeltaT.H"
101         runTime++;
103         Info<< "\n\nTime = " << runTime.timeName() << endl;
105         // Indicators for refinement. Note: before runTime++
106         // only for postprocessing reasons.
107         tmp<volScalarField> tmagGradP = mag(fvc::grad(p));
108         volScalarField normalisedGradP
109         (
110             "normalisedGradP",
111             tmagGradP()/max(tmagGradP())
112         );
113         normalisedGradP.writeOpt() = IOobject::AUTO_WRITE;
114         tmagGradP.clear();
116         bool meshChanged = false;
117         {
118             // Make the fluxes absolute
119             fvc::makeAbsolute(phi, rho, U);
121             // Test : disable refinement for some cells
122             PackedList<1>& protectedCell =
123                 refCast<dynamicRefineFvMesh>(mesh).protectedCell();
124             if (protectedCell.size() == 0)
125             {
126                 protectedCell.setSize(mesh.nCells());
127                 protectedCell = 0;
128             }
130             forAll(betav, cellI)
131             {
132                 if (betav[cellI] < 0.99)
133                 {
134                     protectedCell[cellI] = 1;
135                 }
136             }
138             //volScalarField pIndicator("pIndicator", 
139             //    p*(fvc::laplacian(p))
140             //  / (
141             //        magSqr(fvc::grad(p))
142             //      + dimensionedScalar
143             //        (
144             //            "smallish",
145             //            sqr(p.dimensions()/dimLength),
146             //            1E-6
147             //        )
148             //    ));
149             //pIndicator.writeOpt() = IOobject::AUTO_WRITE;
151             // Flux estimate for introduced faces.
152             volVectorField rhoU("rhoU", rho*U);
154             // Do any mesh changes
155             meshChanged = mesh.update();
157 //        if (mesh.moving() || meshChanged)
158 //        {
159 //#           include "correctPhi.H"
160 //        }
162             // Make the fluxes relative to the mesh motion
163             fvc::makeRelative(phi, rho, U);
164         }
167 #       include "rhoEqn.H"
168 #       include "UEqn.H"
170         // --- PISO loop
171         for (int corr=1; corr<=nCorr; corr++)
172         {
173 #           include "bEqn.H"
174 #           include "ftEqn.H"
175 #           include "huEqn.H"
176 #           include "hEqn.H"
178             if (!ign.ignited())
179             {
180                 hu == h;
181             }
183 #           include "pEqn.H"
184         }
186         turbulence->correct();
188         runTime.write();
190         Info<< "\nExecutionTime = "
191              << runTime.elapsedCpuTime()
192              << " s\n" << endl;
193     }
195     Info<< "\n end\n";
197     return(0);
201 // ************************************************************************* //