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<< "Reading field U\n" << endl;
00030     volVectorField U
00031     (
00032         IOobject
00033         (
00034             "U",
00035             runTime.timeName(),
00036             mesh,
00037             IOobject::MUST_READ,
00038             IOobject::AUTO_WRITE
00039         ),
00040         mesh
00041     );
00042 
00043     #include <finiteVolume/createPhi.H>
00044 
00045     Info<< "Creating phaseChangeTwoPhaseMixture\n" << endl;
00046     autoPtr<phaseChangeTwoPhaseMixture> twoPhaseProperties =
00047         phaseChangeTwoPhaseMixture::New(U, phi);
00048 
00049     const dimensionedScalar& rho1 = twoPhaseProperties->rho1();
00050     const dimensionedScalar& rho2 = twoPhaseProperties->rho2();
00051     const dimensionedScalar& pSat = twoPhaseProperties->pSat();
00052 
00053     // Need to store rho for ddt(rho, U)
00054     volScalarField rho
00055     (
00056         IOobject
00057         (
00058             "rho",
00059             runTime.timeName(),
00060             mesh,
00061             IOobject::READ_IF_PRESENT
00062         ),
00063         alpha1*rho1 + (scalar(1) - alpha1)*rho2,
00064         alpha1.boundaryField().types()
00065     );
00066     rho.oldTime();
00067 
00068     // Construct interface from alpha1 distribution
00069     interfaceProperties interface(alpha1, U, twoPhaseProperties());
00070 
00071     // Construct incompressible turbulence model
00072     autoPtr<incompressible::turbulenceModel> turbulence
00073     (
00074         incompressible::turbulenceModel::New(U, phi, twoPhaseProperties())
00075     );
00076 
00077 
00078     Info<< "Calculating field g.h\n" << endl;
00079     volScalarField gh("gh", g & mesh.C());
00080     surfaceScalarField ghf("ghf", g & mesh.Cf());
00081 
00082     volScalarField p
00083     (
00084         IOobject
00085         (
00086             "p",
00087             runTime.timeName(),
00088             mesh,
00089             IOobject::NO_READ,
00090             IOobject::AUTO_WRITE
00091         ),
00092         p_rgh + rho*gh
00093     );
00094 
00095     label pRefCell = 0;
00096     scalar pRefValue = 0.0;
00097     setRefCell
00098     (
00099         p,
00100         p_rgh,
00101         mesh.solutionDict().subDict("PISO"),
00102         pRefCell,
00103         pRefValue
00104     );
00105 
00106     if (p_rgh.needReference())
00107     {
00108         p += dimensionedScalar
00109         (
00110             "p",
00111             p.dimensions(),
00112             pRefValue - getRefCellValue(p, pRefCell)
00113         );
00114         p_rgh = p - rho*gh;
00115     }
00116 
00117 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines