OpenFOAM logo
The Open Source CFD Toolbox
  Source Guide OpenCFD Solutions Contact OpenFOAM

PDRFoamAutoRefine.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software; you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by the
00013     Free Software Foundation; either version 2 of the License, or (at your
00014     option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM; if not, write to the Free Software Foundation,
00023     Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
00024 
00025 Application
00026     PDRFoam
00027 
00028 Description
00029     Compressible premixed/partially-premixed combustion solver with turbulence
00030     modelling.
00031 
00032     Combusting RANS code using the b-Xi two-equation model.
00033     Xi may be obtained by either the solution of the Xi transport
00034     equation or from an algebraic exression.  Both approaches are
00035     based on Gulder's flame speed correlation which has been shown
00036     to be appropriate by comparison with the results from the
00037     spectral model.
00038 
00039     Strain effects are encorporated directly into the Xi equation
00040     but not in the algebraic approximation.  Further work need to be
00041     done on this issue, particularly regarding the enhanced removal rate
00042     caused by flame compression.  Analysis using results of the spectral
00043     model will be required.
00044 
00045     For cases involving very lean Propane flames or other flames which are
00046     very strain-sensitive, a transport equation for the laminar flame
00047     speed is present.  This equation is derived using heuristic arguments
00048     involving the strain time scale and the strain-rate at extinction.
00049     the transport velocity is the same as that for the Xi equation.
00050 
00051     For large flames e.g. explosions additional modelling for the flame
00052     wrinkling due to surface instabilities may be applied.
00053 
00054     PDR (porosity/distributed resistance) modelling is included to handle
00055     regions containing blockages which cannot be resolved by the mesh.
00056 
00057 \*---------------------------------------------------------------------------*/
00058 
00059 #include "fvCFD.H"
00060 #include "dynamicFvMesh.H"
00061 #include "hhuCombustionThermo.H"
00062 #include "RASModel.H"
00063 #include "laminarFlameSpeed.H"
00064 #include "XiModel.H"
00065 #include "PDRDragModel.H"
00066 #include "ignition.H"
00067 #include "Switch.H"
00068 #include "bound.H"
00069 #include "dynamicRefineFvMesh.H"
00070 
00071 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00072 
00073 int main(int argc, char *argv[])
00074 {
00075     #include "setRootCase.H"
00076 
00077     #include "createTime.H"
00078     #include "createDynamicFvMesh.H"
00079     #include "readCombustionProperties.H"
00080     #include "readGravitationalAcceleration.H"
00081     #include "createFields.H"
00082     #include "readPISOControls.H"
00083     #include "initContinuityErrs.H"
00084     #include "readTimeControls.H"
00085     #include "setInitialDeltaT.H"
00086 
00087     scalar StCoNum = 0.0;
00088 
00089     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00090 
00091     Info<< "\nStarting time loop\n" << endl;
00092 
00093     while (runTime.run())
00094     {
00095         #include "readTimeControls.H"
00096         #include "readPISOControls.H"
00097         #include "CourantNo.H"
00098 
00099         #include "setDeltaT.H"
00100 
00101         runTime++;
00102 
00103         Info<< "\n\nTime = " << runTime.timeName() << endl;
00104 
00105         // Indicators for refinement. Note: before runTime++
00106         // only for postprocessing reasons.
00107         tmp<volScalarField> tmagGradP = mag(fvc::grad(p));
00108         volScalarField normalisedGradP
00109         (
00110             "normalisedGradP",
00111             tmagGradP()/max(tmagGradP())
00112         );
00113         normalisedGradP.writeOpt() = IOobject::AUTO_WRITE;
00114         tmagGradP.clear();
00115 
00116         bool meshChanged = false;
00117         {
00118             // Make the fluxes absolute
00119             fvc::makeAbsolute(phi, rho, U);
00120 
00121             // Test : disable refinement for some cells
00122             PackedBoolList& protectedCell =
00123                 refCast<dynamicRefineFvMesh>(mesh).protectedCell();
00124 
00125             if (protectedCell.empty())
00126             {
00127                 protectedCell.setSize(mesh.nCells());
00128                 protectedCell = 0;
00129             }
00130 
00131             forAll(betav, cellI)
00132             {
00133                 if (betav[cellI] < 0.99)
00134                 {
00135                     protectedCell[cellI] = 1;
00136                 }
00137             }
00138 
00139             //volScalarField pIndicator("pIndicator",
00140             //    p*(fvc::laplacian(p))
00141             //  / (
00142             //        magSqr(fvc::grad(p))
00143             //      + dimensionedScalar
00144             //        (
00145             //            "smallish",
00146             //            sqr(p.dimensions()/dimLength),
00147             //            1E-6
00148             //        )
00149             //    ));
00150             //pIndicator.writeOpt() = IOobject::AUTO_WRITE;
00151 
00152             // Flux estimate for introduced faces.
00153             volVectorField rhoU("rhoU", rho*U);
00154 
00155             // Do any mesh changes
00156             meshChanged = mesh.update();
00157 
00158 //        if (mesh.moving() || meshChanged)
00159 //        {
00160 //            #include "correctPhi.H"
00161 //        }
00162 
00163             // Make the fluxes relative to the mesh motion
00164             fvc::makeRelative(phi, rho, U);
00165         }
00166 
00167 
00168         #include "rhoEqn.H"
00169         #include "UEqn.H"
00170 
00171         // --- PISO loop
00172         for (int corr=1; corr<=nCorr; corr++)
00173         {
00174             #include "bEqn.H"
00175             #include "ftEqn.H"
00176             #include "huEqn.H"
00177             #include "hEqn.H"
00178 
00179             if (!ign.ignited())
00180             {
00181                 hu == h;
00182             }
00183 
00184             #include "pEqn.H"
00185         }
00186 
00187         turbulence->correct();
00188 
00189         runTime.write();
00190 
00191         Info<< "\nExecutionTime = "
00192              << runTime.elapsedCpuTime()
00193              << " s\n" << endl;
00194     }
00195 
00196     Info<< "\n end\n";
00197 
00198     return 0;
00199 }
00200 
00201 
00202 // ************************************************************************* //
Copyright © 2000-2009 OpenCFD Ltd