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 transportProperties\n" << endl;
00002 
00003     IOdictionary transportProperties
00004     (
00005         IOobject
00006         (
00007             "transportProperties",
00008             runTime.constant(),
00009             mesh,
00010             IOobject::MUST_READ,
00011             IOobject::NO_WRITE
00012         )
00013     );
00014 
00015     autoPtr<phaseModel> phasea = phaseModel::New
00016     (
00017         mesh,
00018         transportProperties,
00019         "a"
00020     );
00021 
00022     autoPtr<phaseModel> phaseb = phaseModel::New
00023     (
00024         mesh,
00025         transportProperties,
00026         "b"
00027     );
00028 
00029     volVectorField& Ua = phasea->U();
00030     surfaceScalarField& phia = phasea->phi();
00031     const dimensionedScalar& rhoa = phasea->rho();
00032     const dimensionedScalar& nua = phasea->nu();
00033 
00034     volVectorField& Ub = phaseb->U();
00035     surfaceScalarField& phib = phaseb->phi();
00036     const dimensionedScalar& rhob = phaseb->rho();
00037     const dimensionedScalar& nub = phaseb->nu();
00038 
00039     Info<< "Reading field alpha\n" << endl;
00040     volScalarField alpha
00041     (
00042         IOobject
00043         (
00044             "alpha",
00045             runTime.timeName(),
00046             mesh,
00047             IOobject::MUST_READ,
00048             IOobject::AUTO_WRITE
00049         ),
00050         mesh
00051     );
00052 
00053     volScalarField beta
00054     (
00055         IOobject
00056         (
00057             "beta",
00058             runTime.timeName(),
00059             mesh,
00060             IOobject::NO_READ,
00061             IOobject::NO_WRITE
00062         ),
00063         scalar(1) - alpha
00064         //,alpha.boundaryField().types()
00065     );
00066 
00067     Info<< "Reading field p\n" << endl;
00068     volScalarField p
00069     (
00070         IOobject
00071         (
00072             "p",
00073             runTime.timeName(),
00074             mesh,
00075             IOobject::MUST_READ,
00076             IOobject::AUTO_WRITE
00077         ),
00078         mesh
00079     );
00080 
00081     volVectorField U
00082     (
00083         IOobject
00084         (
00085             "U",
00086             runTime.timeName(),
00087             mesh,
00088             IOobject::NO_READ,
00089             IOobject::AUTO_WRITE
00090         ),
00091         alpha*Ua + beta*Ub
00092     );
00093 
00094     dimensionedScalar Cvm
00095     (
00096         transportProperties.lookup("Cvm")
00097     );
00098 
00099     dimensionedScalar Cl
00100     (
00101         transportProperties.lookup("Cl")
00102     );
00103 
00104     dimensionedScalar Ct
00105     (
00106         transportProperties.lookup("Ct")
00107     );
00108 
00109     surfaceScalarField phi
00110     (
00111         IOobject
00112         (
00113             "phi",
00114             runTime.timeName(),
00115             mesh
00116         ),
00117         fvc::interpolate(alpha)*phia + fvc::interpolate(beta)*phib
00118     );
00119 
00120     volScalarField rho
00121     (
00122         IOobject
00123         (
00124             "rho",
00125             runTime.timeName(),
00126             mesh
00127         ),
00128         alpha*rhoa + beta*rhob
00129     );
00130 
00131     #include "../bubbleFoam/createRASTurbulence.H"
00132 
00133     Info<< "Calculating field DDtUa and DDtUb\n" << endl;
00134 
00135     volVectorField DDtUa =
00136         fvc::ddt(Ua)
00137       + fvc::div(phia, Ua)
00138       - fvc::div(phia)*Ua;
00139 
00140     volVectorField DDtUb =
00141         fvc::ddt(Ub)
00142       + fvc::div(phib, Ub)
00143       - fvc::div(phib)*Ub;
00144 
00145 
00146     Info<< "Calculating field g.h\n" << endl;
00147     volScalarField gh("gh", g & mesh.C());
00148 
00149     IOdictionary interfacialProperties
00150     (
00151         IOobject
00152         (
00153             "interfacialProperties",
00154             runTime.constant(),
00155             mesh,
00156             IOobject::MUST_READ,
00157             IOobject::NO_WRITE
00158         )
00159     );
00160 
00161     autoPtr<dragModel> draga = dragModel::New
00162     (
00163         interfacialProperties,
00164         alpha,
00165         phasea,
00166         phaseb
00167     );
00168 
00169     autoPtr<dragModel> dragb = dragModel::New
00170     (
00171         interfacialProperties,
00172         beta,
00173         phaseb,
00174         phasea
00175     );
00176 
00177     word dragPhase("blended");
00178     if (interfacialProperties.found("dragPhase"))
00179     {
00180         dragPhase = word(interfacialProperties.lookup("dragPhase"));
00181 
00182         bool validDrag =
00183             dragPhase == "a" || dragPhase == "b" || dragPhase == "blended";
00184 
00185         if (!validDrag)
00186         {
00187             FatalErrorIn(args.executable())
00188                 << "invalid dragPhase " << dragPhase
00189                 << exit(FatalError);
00190         }
00191     }
00192 
00193     Info << "dragPhase is " << dragPhase << endl;
00194     kineticTheoryModel kineticTheory
00195     (
00196         phasea,
00197         Ub,
00198         alpha,
00199         draga
00200     );
00201 
00202     surfaceScalarField rUaAf
00203     (
00204         IOobject
00205         (
00206             "rUaAf",
00207             runTime.timeName(),
00208             mesh,
00209             IOobject::NO_READ,
00210             IOobject::NO_WRITE
00211         ),
00212         mesh,
00213         dimensionedScalar("zero", dimensionSet(0, 0, 1, 0, 0), 0.0)
00214     );
00215 
00216     surfaceScalarField ppMagf
00217     (
00218         IOobject
00219         (
00220             "ppMagf",
00221             runTime.timeName(),
00222             mesh,
00223             IOobject::NO_READ,
00224             IOobject::NO_WRITE
00225         ),
00226         mesh,
00227         dimensionedScalar("zero", dimensionSet(0, 2, -1, 0, 0), 0.0)
00228     );
00229 
00230 
00231     label pRefCell = 0;
00232     scalar pRefValue = 0.0;
00233     setRefCell(p, mesh.solutionDict().subDict("PISO"), pRefCell, pRefValue);
00234 
00235 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines