Loading...
Searching...
No Matches
createFields.H
Go to the documentation of this file.
1#include "createRDeltaT.H"
2
3Info<< "Reading field p_rgh\n" << endl;
5(
6 IOobject
7 (
8 "p_rgh",
9 runTime.timeName(),
10 mesh,
11 IOobject::MUST_READ,
12 IOobject::AUTO_WRITE
13 ),
14 mesh
15);
16
17Info<< "Reading field U\n" << endl;
19(
20 IOobject
21 (
22 "U",
23 runTime.timeName(),
24 mesh,
25 IOobject::MUST_READ,
26 IOobject::AUTO_WRITE
27 ),
28 mesh
29);
30
31#include "createPhi.H"
32
33Info<< "Reading transportProperties\n" << endl;
34immiscibleIncompressibleTwoPhaseMixture mixture(U, phi);
35
38
39const dimensionedScalar& rho1 = mixture.rho1();
40const dimensionedScalar& rho2 = mixture.rho2();
41
42// Need to store rho for ddt(rho, U)
44(
45 IOobject
46 (
47 "rho",
48 runTime.timeName(),
49 mesh,
50 IOobject::READ_IF_PRESENT,
51 IOobject::AUTO_WRITE
52 ),
54);
55rho.oldTime();
56
57// Need to store mu as incompressibleTwoPhaseMixture does not store it
59(
60 IOobject
61 (
62 "mu",
63 runTime.timeName(),
64 mesh,
65 IOobject::READ_IF_PRESENT
66 ),
67 mixture.mu()
68);
69
70
71// Mass flux
73(
74 IOobject
75 (
76 "rhoPhi",
77 runTime.timeName(),
78 mesh,
79 IOobject::NO_READ,
80 IOobject::NO_WRITE
81 ),
82 fvc::interpolate(rho)*phi
83);
84
85#include "readGravitationalAcceleration.H"
86#include "readhRef.H"
87#include "gh.H"
88
90(
91 IOobject
92 (
93 "p",
94 runTime.timeName(),
95 mesh,
96 IOobject::NO_READ,
97 IOobject::AUTO_WRITE
98 ),
99 p_rgh + rho*gh
100);
101
102label pRefCell = 0;
103scalar pRefValue = 0.0;
105(
106 p,
107 p_rgh,
108 mesh.solutionDict().subDict("PIMPLE"),
109 pRefCell,
111);
112
113if (p_rgh.needReference())
114{
116 (
117 "p",
118 p.dimensions(),
120 );
121 p_rgh = p - rho*gh;
122}
123
124mesh.setFluxRequired(p_rgh.name());
125mesh.setFluxRequired(alpha1.name());
126
127// alphac must be constructed before the cloud
128// so that the drag-models can find it
129volScalarField alphac
130(
131 IOobject
132 (
133 "alphac",
134 runTime.timeName(),
135 mesh,
136 IOobject::READ_IF_PRESENT,
137 IOobject::AUTO_WRITE
138 ),
139 mesh,
140 dimensionedScalar(dimless, Zero),
141 fvPatchFieldBase::zeroGradientType()
142);
143alphac.oldTime();
144
145volScalarField alphacRho(alphac*rho);
146alphacRho.oldTime();
147
148Info<< "Constructing kinematicCloud " << endl;
149basicKinematicCloud kinematicCloud
150(
151 "kinematicCloud",
152 rho,
153 U,
154 mu,
155 g
156);
157
158// Particle fraction upper limit
159scalar alphacMin
160(
161 1.0
162 - (
163 kinematicCloud.particleProperties().subDict("constantProperties")
164 .get<scalar>("alphaMax")
165 )
166);
167
168// Update alphac from the particle locations
169alphac = max(1.0 - kinematicCloud.theta(), alphacMin);
170alphac.correctBoundaryConditions();
171
172surfaceScalarField alphacf("alphacf", fvc::interpolate(alphac));
173
174// Phase mass flux
175surfaceScalarField alphaRhoPhic("alphaRhoPhic", alphacf*rhoPhi);
176
177// Volumetric phase flux
178surfaceScalarField alphaPhic("alphaPhic", alphacf*phi);
179
180autoPtr
181<
182 PhaseCompressibleTurbulenceModel
183 <
184 immiscibleIncompressibleTwoPhaseMixture
185 >
187(
188 PhaseCompressibleTurbulenceModel
189 <
190 immiscibleIncompressibleTwoPhaseMixture
191 >::New
192 (
193 alphac,
194 rho,
195 U,
196 alphaRhoPhic,
197 rhoPhi,
198 mixture
199 )
200);
201
202#include "createMRF.H"
rhoPhi
Definition rhoEqn.H:10
const uniformDimensionedVectorField & g
volScalarField & p_rgh
const scalar pRefValue
const label pRefCell
const volScalarField & gh
const volScalarField & alpha1
volScalarField & rho2
const volScalarField & alpha2
volScalarField & rho1
Cloud class to introduce kinematic parcels.
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
engineTime & runTime
compressible::turbulenceModel & turbulence
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
scalar getRefCellValue(const volScalarField &field, const label refCelli)
Return the current value of field in the reference cell.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
setRefCell(p, pimple.dict(), pRefCell, pRefValue)
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture