Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / applications / solvers / multiphase / interMixingFoam / createFields.H
blob6fdd341758add6b15133b85e7962e6683cc95ff0
1     Info<< "Reading field p\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     );
30     Info<< "Reading field alpha2\n" << endl;
31     volScalarField alpha2
32     (
33         IOobject
34         (
35             "alpha2",
36             runTime.timeName(),
37             mesh,
38             IOobject::MUST_READ,
39             IOobject::AUTO_WRITE
40         ),
41         mesh
42     );
45     Info<< "Reading field alpha3\n" << endl;
46     volScalarField alpha3
47     (
48         IOobject
49         (
50             "alpha3",
51             runTime.timeName(),
52             mesh,
53             IOobject::MUST_READ,
54             IOobject::AUTO_WRITE
55         ),
56         mesh
57     );
59     alpha3 == 1.0 - alpha1 - alpha2;
62     Info<< "Reading field U\n" << endl;
63     volVectorField U
64     (
65         IOobject
66         (
67             "U",
68             runTime.timeName(),
69             mesh,
70             IOobject::MUST_READ,
71             IOobject::AUTO_WRITE
72         ),
73         mesh
74     );
76 #   include "createPhi.H"
78     threePhaseMixture threePhaseProperties(U, phi);
80     const dimensionedScalar& rho1 = threePhaseProperties.rho1();
81     const dimensionedScalar& rho2 = threePhaseProperties.rho2();
82     const dimensionedScalar& rho3 = threePhaseProperties.rho3();
84     dimensionedScalar D23(threePhaseProperties.lookup("D23"));
86     // Need to store rho for ddt(rho, U)
87     volScalarField rho
88     (
89         IOobject
90         (
91             "rho",
92             runTime.timeName(),
93             mesh,
94             IOobject::READ_IF_PRESENT
95         ),
96         alpha1*rho1 + alpha2*rho2 + alpha3*rho3,
97         alpha1.boundaryField().types()
98     );
99     rho.oldTime();
102     // Mass flux
103     // Initialisation does not matter because rhoPhi is reset after the
104     // alpha solution before it is used in the U equation.
105     surfaceScalarField rhoPhi
106     (
107         IOobject
108         (
109             "rho*phi",
110             runTime.timeName(),
111             mesh,
112             IOobject::NO_READ,
113             IOobject::NO_WRITE
114         ),
115         rho1*phi
116     );
118     Info<< "Calculating field g.h\n" << endl;
119     volScalarField gh("gh", g & mesh.C());
120     surfaceScalarField ghf("gh", g & mesh.Cf());
123     volScalarField p
124     (
125         IOobject
126         (
127             "p",
128             runTime.timeName(),
129             mesh,
130             IOobject::NO_READ,
131             IOobject::AUTO_WRITE
132         ),
133         pd + rho*gh
134     );
136     label pdRefCell = 0;
137     scalar pdRefValue = 0.0;
138     setRefCell(pd, mesh.solutionDict().subDict("PISO"), pdRefCell, pdRefValue);
140     scalar pRefValue = 0.0;
142     if (pd.needReference())
143     {
144         pRefValue = readScalar
145         (
146             mesh.solutionDict().subDict("PISO").lookup("pRefValue")
147         );
149         p += dimensionedScalar
150         (
151             "p",
152             p.dimensions(),
153             pRefValue - getRefCellValue(p, pdRefCell)
154         );
155     }
157     // Construct interface from alpha distribution
158     threePhaseInterfaceProperties interface(threePhaseProperties);
161     // Construct incompressible turbulence model
162     autoPtr<incompressible::turbulenceModel> turbulence
163     (
164         incompressible::turbulenceModel::New(U, phi, threePhaseProperties)
165     );