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

createFields.H

Go to the documentation of this file.
00001     Info<< "Reading thermophysical properties\n" << endl;
00002 
00003     autoPtr<hhuCombustionThermo> pThermo
00004     (
00005         hhuCombustionThermo::New(mesh)
00006     );
00007     hhuCombustionThermo& thermo = pThermo();
00008     basicMultiComponentMixture& composition = thermo.composition();
00009 
00010     volScalarField rho
00011     (
00012         IOobject
00013         (
00014             "rho",
00015             runTime.timeName(),
00016             mesh,
00017             IOobject::NO_READ,
00018             IOobject::AUTO_WRITE
00019         ),
00020         thermo.rho()
00021     );
00022 
00023     volScalarField& p = thermo.p();
00024     const volScalarField& psi = thermo.psi();
00025     volScalarField& h = thermo.h();
00026     volScalarField& hu = thermo.hu();
00027 
00028     volScalarField& b = composition.Y("bprog");
00029     Info<< "min(bprog) = " << min(b).value() << endl;
00030 
00031     const volScalarField& T = thermo.T();
00032 
00033 
00034     Info<< "\nReading field U\n" << endl;
00035     volVectorField U
00036     (
00037         IOobject
00038         (
00039             "U",
00040             runTime.timeName(),
00041             mesh,
00042             IOobject::MUST_READ,
00043             IOobject::AUTO_WRITE
00044         ),
00045         mesh
00046     );
00047 
00048 #   include <finiteVolume/compressibleCreatePhi.H>
00049 
00050 
00051     Info<< "Creating turbulence model\n" << endl;
00052     autoPtr<compressible::turbulenceModel> turbulence
00053     (
00054         compressible::turbulenceModel::New
00055         (
00056             rho,
00057             U,
00058             phi,
00059             thermo
00060         )
00061     );
00062 
00063     Info<< "Creating field DpDt\n" << endl;
00064     volScalarField DpDt =
00065         fvc::DDt(surfaceScalarField("phiU", phi/fvc::interpolate(rho)), p);
00066 
00067 
00068     Info<< "Creating field Xi\n" << endl;
00069     volScalarField Xi
00070     (
00071         IOobject
00072         (
00073             "Xi",
00074             runTime.timeName(),
00075             mesh,
00076             IOobject::MUST_READ,
00077             IOobject::AUTO_WRITE
00078         ),
00079         mesh
00080     );
00081 
00082 
00083     Info<< "Creating the unstrained laminar flame speed\n" << endl;
00084     autoPtr<laminarFlameSpeed> unstrainedLaminarFlameSpeed
00085     (
00086         laminarFlameSpeed::New(thermo)
00087     );
00088 
00089 
00090     Info<< "Reading strained laminar flame speed field Su\n" << endl;
00091     volScalarField Su
00092     (
00093         IOobject
00094         (
00095             "Su",
00096             runTime.timeName(),
00097             mesh,
00098             IOobject::MUST_READ,
00099             IOobject::AUTO_WRITE
00100         ),
00101         mesh
00102     );
00103 
00104     dimensionedScalar SuMin = 0.01*Su.average();
00105     dimensionedScalar SuMax = 4*Su.average();
00106 
00107     Info<< "Calculating turbulent flame speed field St\n" << endl;
00108     volScalarField St
00109     (
00110         IOobject
00111         (
00112             "St",
00113             runTime.timeName(),
00114             mesh,
00115             IOobject::NO_READ,
00116             IOobject::AUTO_WRITE
00117         ),
00118         Xi*Su
00119     );
00120 
00121 
00122     multivariateSurfaceInterpolationScheme<scalar>::fieldTable fields;
00123 
00124     if (composition.contains("ft"))
00125     {
00126         fields.add(composition.Y("ft"));
00127     }
00128 
00129     fields.add(b);
00130     fields.add(h);
00131     fields.add(hu);
00132 
00133 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines