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

StCorr.H

Go to the documentation of this file.
00001     dimensionedScalar StCorr("StCorr", dimless, 1.0);
00002 
00003     if (ign.igniting())
00004     {
00005         // Calculate volume of ignition kernel
00006         dimensionedScalar Vk("Vk", dimVolume, gSum(c*mesh.V().field()));
00007         dimensionedScalar Ak("Ak", dimArea, 0.0);
00008 
00009         if (Vk.value() > SMALL)
00010         {
00011             // Calculate kernel area from its volume
00012             // and the dimensionality of the case
00013 
00014             switch(mesh.nGeometricD())
00015             {
00016                 case 3:
00017                 {
00018                     // Assume it is part-spherical
00019                     scalar sphereFraction
00020                     (
00021                         readScalar
00022                         (
00023                             combustionProperties.lookup
00024                             (
00025                                 "ignitionSphereFraction"
00026                             )
00027                         )
00028                     );
00029 
00030                     Ak = sphereFraction*4.0*mathematicalConstant::pi
00031                        *pow
00032                         (
00033                             3.0*Vk
00034                            /(sphereFraction*4.0*mathematicalConstant::pi),
00035                             2.0/3.0
00036                         );
00037                 }
00038                 break;
00039             
00040                 case 2:
00041                 {
00042                     // Assume it is part-circular
00043                     dimensionedScalar thickness
00044                     (
00045                         combustionProperties.lookup("ignitionThickness")
00046                     );
00047 
00048                     scalar circleFraction
00049                     (
00050                         readScalar
00051                         (
00052                             combustionProperties.lookup
00053                             (
00054                                 "ignitionCircleFraction"
00055                             )
00056                         )
00057                     );
00058 
00059                     Ak = circleFraction*mathematicalConstant::pi*thickness
00060                        *sqrt
00061                         (
00062                             4.0*Vk
00063                            /(circleFraction*thickness*mathematicalConstant::pi)
00064                         );
00065                 }
00066                 break;
00067 
00068                 case 1:
00069                     // Assume it is plane or two planes
00070                     Ak = dimensionedScalar
00071                     (
00072                         combustionProperties.lookup("ignitionKernelArea")
00073                     );
00074                 break;
00075             }
00076 
00077             // Calculate kernel area from b field consistent with the
00078             // discretisation of the b equation.
00079             volScalarField mgb =
00080                 fvc::div(nf, b, "div(phiSt,bprog)") - b*fvc::div(nf) + dMgb;
00081             dimensionedScalar AkEst = gSum(mgb*mesh.V().field());
00082 
00083             StCorr.value() = max(min((Ak/AkEst).value(), 10.0), 1.0);
00084 
00085             Info<< "StCorr = " << StCorr.value() << endl;
00086         }
00087     }
00088 
00089 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines