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 field p_rgh\n" << endl;
00002     volScalarField p_rgh
00003     (
00004         IOobject
00005         (
00006             "p_rgh",
00007             runTime.timeName(),
00008             mesh,
00009             IOobject::MUST_READ,
00010             IOobject::AUTO_WRITE
00011         ),
00012         mesh
00013     );
00014 
00015     Info<< "Reading field alpha1\n" << endl;
00016     volScalarField alpha1
00017     (
00018         IOobject
00019         (
00020             "alpha1",
00021             runTime.timeName(),
00022             mesh,
00023             IOobject::MUST_READ,
00024             IOobject::AUTO_WRITE
00025         ),
00026         mesh
00027     );
00028 
00029     Info<< "Calculating field alpha1\n" << endl;
00030     volScalarField alpha2("alpha2", scalar(1) - alpha1);
00031 
00032     Info<< "Reading field U\n" << endl;
00033     volVectorField U
00034     (
00035         IOobject
00036         (
00037             "U",
00038             runTime.timeName(),
00039             mesh,
00040             IOobject::MUST_READ,
00041             IOobject::AUTO_WRITE
00042         ),
00043         mesh
00044     );
00045 
00046     #include <finiteVolume/createPhi.H>
00047 
00048 
00049     Info<< "Reading transportProperties\n" << endl;
00050     twoPhaseMixture twoPhaseProperties(U, phi);
00051 
00052     dimensionedScalar rho10
00053     (
00054         twoPhaseProperties.subDict
00055         (
00056             twoPhaseProperties.phase1Name()
00057         ).lookup("rho0")
00058     );
00059 
00060     dimensionedScalar rho20
00061     (
00062         twoPhaseProperties.subDict
00063         (
00064             twoPhaseProperties.phase2Name()
00065         ).lookup("rho0")
00066     );
00067 
00068     dimensionedScalar psi1
00069     (
00070         twoPhaseProperties.subDict
00071         (
00072             twoPhaseProperties.phase1Name()
00073         ).lookup("psi")
00074     );
00075 
00076     dimensionedScalar psi2
00077     (
00078         twoPhaseProperties.subDict
00079         (
00080             twoPhaseProperties.phase2Name()
00081         ).lookup("psi")
00082     );
00083 
00084     dimensionedScalar pMin(twoPhaseProperties.lookup("pMin"));
00085 
00086     Info<< "Calculating field g.h\n" << endl;
00087     volScalarField gh("gh", g & mesh.C());
00088     surfaceScalarField ghf("ghf", g & mesh.Cf());
00089 
00090     volScalarField p
00091     (
00092         IOobject
00093         (
00094             "p",
00095             runTime.timeName(),
00096             mesh,
00097             IOobject::NO_READ,
00098             IOobject::AUTO_WRITE
00099         ),
00100         max
00101         (
00102             (p_rgh + gh*(alpha1*rho10 + alpha2*rho20))
00103            /(1.0 - gh*(alpha1*psi1 + alpha2*psi2)),
00104             pMin
00105         )
00106     );
00107 
00108     volScalarField rho1 = rho10 + psi1*p;
00109     volScalarField rho2 = rho20 + psi2*p;
00110 
00111     volScalarField rho
00112     (
00113         IOobject
00114         (
00115             "rho",
00116             runTime.timeName(),
00117             mesh,
00118             IOobject::READ_IF_PRESENT,
00119             IOobject::AUTO_WRITE
00120         ),
00121         alpha1*rho1 + alpha2*rho2
00122     );
00123 
00124 
00125     // Mass flux
00126     // Initialisation does not matter because rhoPhi is reset after the
00127     // alpha1 solution before it is used in the U equation.
00128     surfaceScalarField rhoPhi
00129     (
00130         IOobject
00131         (
00132             "rho*phi",
00133             runTime.timeName(),
00134             mesh,
00135             IOobject::NO_READ,
00136             IOobject::NO_WRITE
00137         ),
00138         fvc::interpolate(rho)*phi
00139     );
00140 
00141     volScalarField dgdt =
00142         pos(alpha2)*fvc::div(phi)/max(alpha2, scalar(0.0001));
00143 
00144     // Construct interface from alpha1 distribution
00145     interfaceProperties interface(alpha1, U, twoPhaseProperties);
00146 
00147     // Construct incompressible turbulence model
00148     autoPtr<incompressible::turbulenceModel> turbulence
00149     (
00150         incompressible::turbulenceModel::New(U, phi, twoPhaseProperties)
00151     );
00152 
00153 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines