Removed unneeded lib dependency from mdInitialise
[foam-extend-3.2.git] / applications / solvers / multiphase / interFoam / createFields.H
blob8012f798c402b2f9a923f0bb6efe5c44250af446
1     Info<< "Reading field pd\n" << endl;
2     volScalarField pd
3     (
4         IOobject
5         (
6             "pd",
7             runTime.timeName(),
8             mesh,
9             IOobject::MUST_READ,
10             IOobject::AUTO_WRITE
11         ),
12         mesh
13     );
15     Info<< "Reading field alpha1\n" << endl;
16     volScalarField alpha1
17     (
18         IOobject
19         (
20             "alpha1",
21             runTime.timeName(),
22             mesh,
23             IOobject::MUST_READ,
24             IOobject::AUTO_WRITE
25         ),
26         mesh
27     );
29     Info<< "Reading field U\n" << endl;
30     volVectorField U
31     (
32         IOobject
33         (
34             "U",
35             runTime.timeName(),
36             mesh,
37             IOobject::MUST_READ,
38             IOobject::AUTO_WRITE
39         ),
40         mesh
41     );
43 #   include "createPhi.H"
46     Info<< "Reading transportProperties\n" << endl;
47     twoPhaseMixture twoPhaseProperties(U, phi, "alpha1");
49     const dimensionedScalar& rho1 = twoPhaseProperties.rho1();
50     const dimensionedScalar& rho2 = twoPhaseProperties.rho2();
53     // Need to store rho for ddt(rho, U)
54     volScalarField rho
55     (
56         IOobject
57         (
58             "rho",
59             runTime.timeName(),
60             mesh,
61             IOobject::READ_IF_PRESENT
62         ),
63         alpha1*rho1 + (scalar(1) - alpha1)*rho2,
64         alpha1.boundaryField().types()
65     );
66     rho.oldTime();
69     // Mass flux
70     // Initialisation does not matter because rhoPhi is reset after the
71     // alpha1 solution before it is used in the U equation.
72     surfaceScalarField rhoPhi
73     (
74         IOobject
75         (
76             "rho*phi",
77             runTime.timeName(),
78             mesh,
79             IOobject::NO_READ,
80             IOobject::NO_WRITE
81         ),
82         rho1*phi
83     );
86     Info<< "Calculating field g.h\n" << endl;
87     volScalarField gh("gh", g & mesh.C());
88     surfaceScalarField ghf("gh", g & mesh.Cf());
90     volScalarField p
91     (
92         IOobject
93         (
94             "p",
95             runTime.timeName(),
96             mesh,
97             IOobject::NO_READ,
98             IOobject::AUTO_WRITE
99         ),
100         pd + rho*gh
101     );
104     label pdRefCell = 0;
105     scalar pdRefValue = 0.0;
106     setRefCell
107     (
108         pd,
109         mesh.solutionDict().subDict("PIMPLE"),
110         pdRefCell,
111         pdRefValue
112     );
114     scalar pRefValue = 0.0;
116     if (pd.needReference())
117     {
118         pRefValue = readScalar
119         (
120             mesh.solutionDict().subDict("PIMPLE").lookup("pRefValue")
121         );
123         p += dimensionedScalar
124         (
125             "p",
126             p.dimensions(),
127             pRefValue - getRefCellValue(p, pdRefCell)
128         );
129     }
131     // Construct interface from alpha1 distribution
132     interfaceProperties interface(alpha1, U, twoPhaseProperties);
134      // Construct incompressible turbulence model
135     autoPtr<incompressible::turbulenceModel> turbulence
136     (
137         incompressible::turbulenceModel::New(U, phi, twoPhaseProperties)
138     );