-
Notifications
You must be signed in to change notification settings - Fork 5
/
Copy pathcreateFields.H
166 lines (145 loc) · 3.69 KB
/
createFields.H
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
Info<< "Reading field p_rgh\n" << endl;
volScalarField p_rgh
(
IOobject
(
"p_rgh",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
Info<< "Reading field U\n" << endl;
volVectorField U
(
IOobject
(
"U",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
#include "createPhi.H"
Info<< "Constructing twoPhaseMixtureThermo\n" << endl;
twoPhaseMixtureThermo mixture(mesh);
volScalarField& alpha1(mixture.alpha1());
volScalarField& alpha2(mixture.alpha2());
Info<< "Reading thermophysical properties\n" << endl;
volScalarField& p = mixture.p();
volScalarField& T = mixture.T();
Info<< "Reading field rho_foam \n" << endl;
volScalarField rho_foam
(
IOobject
(
"rho_foam",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
Info<< "Reading field rho_gas \n" << endl;
volScalarField rho_gas
(
IOobject
(
"rho_gas",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
volScalarField rho
(
IOobject
(
"rho",
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
alpha1*rho_gas + (1-alpha1)*rho_foam
);
volScalarField Psi1
(
IOobject
(
"Psi1",
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
mesh,
dimensionedScalar("Psi1", dimensionSet(0,-2,2,0,0,0,0), 1e-5)
);
volScalarField Psi2
(
IOobject
(
"Psi2",
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
mesh,
dimensionedScalar("Psi2", dimensionSet(0,-2,2,0,0,0,0), 1e-2)
);
dimensionedScalar pMin(mixture.lookup("pMin"));
Info<< "Calculating field g.h\n" << endl;
volScalarField gh("gh", g & mesh.C());
surfaceScalarField ghf("ghf", g & mesh.Cf());
// Mass flux
// Initialisation does not matter because rhoPhi is reset after the
// alpha1 solution before it is used in the U equation.
surfaceScalarField rhoPhi
(
IOobject
(
"rhoPhi",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
fvc::interpolate(rho)*phi
);
volScalarField dgdt
(
pos(alpha2)*fvc::div(phi)/max(alpha2, scalar(0.0001))
);
// Construct interface from alpha1 distribution
interfaceProperties interface(alpha1, U, mixture);
#include "MomFields.H"
#include "KineticsFields.H"
#include "rheologyFields.H"
#include "thermalFields.H"
#include "simulationMode.H"
// Construct compressible turbulence model
autoPtr<compressible::turbulenceModel> turbulence
(
compressible::turbulenceModel::New(rho, U, rhoPhi, mixture)
);
Info<< "Creating field kinetic energy K\n" << endl;
volScalarField K("K", 0.5*magSqr(U));
dimensionedScalar densityDimension
(
"densityDimension", dimensionSet(1,-3,0,0,0,0,0), 1.0
);
dimensionedScalar initialFoamMass =
(
fvc::domainIntegrate(rho_foam*alpha2)
+ densityDimension*L0*rhoPoly*fvc::domainIntegrate(alpha2)
);
dimensionedScalar cumulativeContinuityError = 0.0;