FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

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-2010 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
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your 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, see <http://www.gnu.org/licenses/>.
00023     Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
00024 
00025 Application
00026     PDRFoam
00027 
00028 Description
00029     Solver for compressible premixed/partially-premixed combustion with
00030     turbulence 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 incorporated 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 <finiteVolume/fvCFD.H>
00060 #include <dynamicFvMesh/dynamicFvMesh.H>
00061 #include <reactionThermophysicalModels/hhuCombustionThermo.H>
00062 #include <compressibleRASModels/RASModel.H>
00063 #include <laminarFlameSpeedModels/laminarFlameSpeed.H>
00064 #include "XiModels/XiModel/XiModel.H"
00065 #include "PDRModels/dragModels/PDRDragModel/PDRDragModel.H"
00066 #include <engine/ignition.H>
00067 #include <OpenFOAM/Switch.H>
00068 #include <finiteVolume/bound.H>
00069 #include <dynamicFvMesh/dynamicRefineFvMesh.H>
00070 
00071 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00072 
00073 int main(int argc, char *argv[])
00074 {
00075     #include <OpenFOAM/setRootCase.H>
00076 
00077     #include <OpenFOAM/createTime.H>
00078     #include <dynamicFvMesh/createDynamicFvMesh.H>
00079     #include "readCombustionProperties.H"
00080     #include <finiteVolume/readGravitationalAcceleration.H>
00081     #include "createFields.H"
00082     #include <finiteVolume/initContinuityErrs.H>
00083     #include <finiteVolume/readTimeControls.H>
00084     #include <finiteVolume/CourantNo.H>
00085     #include <finiteVolume/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 <finiteVolume/readTimeControls.H>
00096         #include <finiteVolume/readPISOControls.H>
00097         #include <finiteVolume/CourantNo.H>
00098 
00099         #include "setDeltaT.H"
00100 
00101         // Indicators for refinement. Note: before runTime++
00102         // only for postprocessing reasons.
00103         tmp<volScalarField> tmagGradP = mag(fvc::grad(p));
00104         volScalarField normalisedGradP
00105         (
00106             "normalisedGradP",
00107             tmagGradP()/max(tmagGradP())
00108         );
00109         normalisedGradP.writeOpt() = IOobject::AUTO_WRITE;
00110         tmagGradP.clear();
00111 
00112         runTime++;
00113 
00114         Info<< "\n\nTime = " << runTime.timeName() << endl;
00115 
00116 
00117         bool meshChanged = false;
00118         {
00119             // Make the fluxes absolute
00120             fvc::makeAbsolute(phi, rho, U);
00121 
00122             // Test : disable refinement for some cells
00123             PackedBoolList& protectedCell =
00124                 refCast<dynamicRefineFvMesh>(mesh).protectedCell();
00125 
00126             if (protectedCell.empty())
00127             {
00128                 protectedCell.setSize(mesh.nCells());
00129                 protectedCell = 0;
00130             }
00131 
00132             forAll(betav, cellI)
00133             {
00134                 if (betav[cellI] < 0.99)
00135                 {
00136                     protectedCell[cellI] = 1;
00137                 }
00138             }
00139 
00140             //volScalarField pIndicator("pIndicator",
00141             //    p*(fvc::laplacian(p))
00142             //  / (
00143             //        magSqr(fvc::grad(p))
00144             //      + dimensionedScalar
00145             //        (
00146             //            "smallish",
00147             //            sqr(p.dimensions()/dimLength),
00148             //            1E-6
00149             //        )
00150             //    ));
00151             //pIndicator.writeOpt() = IOobject::AUTO_WRITE;
00152 
00153             // Flux estimate for introduced faces.
00154             volVectorField rhoU("rhoU", rho*U);
00155 
00156             // Do any mesh changes
00157             meshChanged = mesh.update();
00158 
00159 //        if (mesh.moving() || meshChanged)
00160 //        {
00161 //            #include "correctPhi.H"
00162 //        }
00163 
00164             // Make the fluxes relative to the mesh motion
00165             fvc::makeRelative(phi, rho, U);
00166         }
00167 
00168 
00169         #include "rhoEqn.H"
00170         #include "UEqn.H"
00171 
00172         // --- PISO loop
00173         for (int corr=1; corr<=nCorr; corr++)
00174         {
00175             #include "bEqn.H"
00176             #include "ftEqn.H"
00177             #include "huEqn.H"
00178             #include "hEqn.H"
00179 
00180             if (!ign.ignited())
00181             {
00182                 hu == h;
00183             }
00184 
00185             #include "pEqn.H"
00186         }
00187 
00188         turbulence->correct();
00189 
00190         runTime.write();
00191 
00192         Info<< "\nExecutionTime = "
00193              << runTime.elapsedCpuTime()
00194              << " s\n" << endl;
00195     }
00196 
00197     Info<< "\n end\n";
00198 
00199     return 0;
00200 }
00201 
00202 
00203 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines