pEqn.H
Go to the documentation of this file.
1 {
2  // rho1 = rho10 + psi1*p_rgh;
3  // rho2 = rho20 + psi2*p_rgh;
4 
5  // tmp<fvScalarMatrix> pEqnComp1;
6  // tmp<fvScalarMatrix> pEqnComp2;
7 
8  // // if (transonic)
9  // //{
10  // //}
11  // // else
12  // {
13  // surfaceScalarField phid1("phid1", fvc::interpolate(psi1)*phi1);
14  // surfaceScalarField phid2("phid2", fvc::interpolate(psi2)*phi2);
15 
16  // pEqnComp1 =
17  // fvc::ddt(rho1) + psi1*correction(fvm::ddt(p_rgh))
18  // + fvc::div(phid1, p_rgh)
19  // - fvc::Sp(fvc::div(phid1), p_rgh);
20 
21  // pEqnComp2 =
22  // fvc::ddt(rho2) + psi2*correction(fvm::ddt(p_rgh))
23  // + fvc::div(phid2, p_rgh)
24  // - fvc::Sp(fvc::div(phid2), p_rgh);
25  // }
26 
27  PtrList<surfaceScalarField> alphafs(fluid.phases().size());
28  PtrList<volVectorField> HbyAs(fluid.phases().size());
29  PtrList<surfaceScalarField> phiHbyAs(fluid.phases().size());
30  PtrList<volScalarField> rAUs(fluid.phases().size());
31  PtrList<surfaceScalarField> rAlphaAUfs(fluid.phases().size());
32 
33  phasei = 0;
34  forAllIter(PtrDictionary<phaseModel>, fluid.phases(), iter)
35  {
36  phaseModel& phase = iter();
37 
38  MRF.makeAbsolute(phase.phi().oldTime());
39  MRF.makeAbsolute(phase.phi());
40 
41  HbyAs.set(phasei, new volVectorField(phase.U()));
42  phiHbyAs.set(phasei, new surfaceScalarField(1.0*phase.phi()));
43 
44  phasei++;
45  }
46 
48  (
49  IOobject
50  (
51  "phiHbyA",
52  runTime.timeName(),
53  mesh,
54  IOobject::NO_READ,
55  IOobject::AUTO_WRITE
56  ),
57  mesh,
58  dimensionedScalar("phiHbyA", dimArea*dimVelocity, 0)
59  );
60 
61  volScalarField rho("rho", fluid.rho());
63 
64  phasei = 0;
65  forAllIter(PtrDictionary<phaseModel>, fluid.phases(), iter)
66  {
67  phaseModel& phase = iter();
68  const volScalarField& alpha = phase;
69 
70  alphafs.set(phasei, fvc::interpolate(alpha).ptr());
71  alphafs[phasei].rename("hmm" + alpha.name());
72 
73  volScalarField dragCoeffi
74  (
75  IOobject
76  (
77  "dragCoeffi",
78  runTime.timeName(),
79  mesh
80  ),
81  fluid.dragCoeff(phase, dragCoeffs())/phase.rho(),
82  zeroGradientFvPatchScalarField::typeName
83  );
84  dragCoeffi.correctBoundaryConditions();
85 
86  rAUs.set(phasei, (1.0/(UEqns[phasei].A() + dragCoeffi)).ptr());
87  rAlphaAUfs.set
88  (
89  phasei,
90  (
91  alphafs[phasei]
92  /fvc::interpolate(UEqns[phasei].A() + dragCoeffi)
93  ).ptr()
94  );
95 
97 
98  phiHbyAs[phasei] =
99  (
100  fvc::flux(HbyAs[phasei])
101  + MRF.zeroFilter
102  (
103  rAlphaAUfs[phasei]*fvc::ddtCorr(phase.U(), phase.phi())
104  )
105  );
106  MRF.makeRelative(phiHbyAs[phasei]);
107  MRF.makeRelative(phase.phi().oldTime());
108  MRF.makeRelative(phase.phi());
109 
110  phiHbyAs[phasei] +=
112  *(
113  fluid.surfaceTension(phase)*mesh.magSf()
114  + (phase.rho() - fvc::interpolate(rho))*(g & mesh.Sf())
115  - ghSnGradRho
116  )/phase.rho();
117 
118  multiphaseSystem::dragModelTable::const_iterator dmIter =
119  fluid.dragModels().begin();
120  multiphaseSystem::dragCoeffFields::const_iterator dcIter =
121  dragCoeffs().begin();
122  for
123  (
124  ;
125  dmIter != fluid.dragModels().end() && dcIter != dragCoeffs().end();
126  ++dmIter, ++dcIter
127  )
128  {
129  const phaseModel *phase2Ptr = nullptr;
130 
131  if (&phase == &dmIter()->phase1())
132  {
133  phase2Ptr = &dmIter()->phase2();
134  }
135  else if (&phase == &dmIter()->phase2())
136  {
137  phase2Ptr = &dmIter()->phase1();
138  }
139  else
140  {
141  continue;
142  }
143 
144  phiHbyAs[phasei] +=
145  fvc::interpolate((*dcIter())/phase.rho())
146  /fvc::interpolate(UEqns[phasei].A() + dragCoeffi)
147  *phase2Ptr->phi();
148 
149  HbyAs[phasei] +=
150  (1.0/phase.rho())*rAUs[phasei]*(*dcIter())
151  *phase2Ptr->U();
152 
153  // Alternative flux-pressure consistent drag
154  // but not momentum conservative
155  //
156  // HbyAs[phasei] += fvc::reconstruct
157  // (
158  // fvc::interpolate
159  // (
160  // (1.0/phase.rho())*rAUs[phasei]*(*dcIter())
161  // )*phase2Ptr->phi()
162  // );
163  }
164 
166 
167  phasei++;
168  }
169 
171  (
172  IOobject
173  (
174  "rAUf",
175  runTime.timeName(),
176  mesh
177  ),
178  mesh,
179  dimensionedScalar("rAUf", dimensionSet(-1, 3, 1, 0, 0), 0)
180  );
181 
182  phasei = 0;
183  forAllIter(PtrDictionary<phaseModel>, fluid.phases(), iter)
184  {
185  phaseModel& phase = iter();
186  rAUf += mag(alphafs[phasei]*rAlphaAUfs[phasei])/phase.rho();
187 
188  phasei++;
189  }
190 
191  // Update the fixedFluxPressure BCs to ensure flux consistency
192  {
193  surfaceScalarField::Boundary phib(phi.boundaryField());
194  phib = 0;
195  phasei = 0;
196  forAllIter(PtrDictionary<phaseModel>, fluid.phases(), iter)
197  {
198  phaseModel& phase = iter();
199 
200  phib +=
201  alphafs[phasei].boundaryField()
202  *(mesh.Sf().boundaryField() & phase.U().boundaryField());
203 
204  phasei++;
205  }
206 
208  (
209  p_rgh.boundaryFieldRef(),
210  (
211  phiHbyA.boundaryField() - MRF.relative(phib)
212  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
213  );
214  }
215 
216  while (pimple.correctNonOrthogonal())
217  {
218  fvScalarMatrix pEqnIncomp
219  (
221  - fvm::laplacian(rAUf, p_rgh)
222  );
223 
224  pEqnIncomp.setReference(pRefCell, pRefValue);
225 
226  solve
227  (
228  // (
229  // (alpha1/rho1)*pEqnComp1()
230  // + (alpha2/rho2)*pEqnComp2()
231  // ) +
232  pEqnIncomp,
233  mesh.solver(p_rgh.select(pimple.finalInnerIter()))
234  );
235 
236  if (pimple.finalNonOrthogonalIter())
237  {
238  surfaceScalarField mSfGradp("mSfGradp", pEqnIncomp.flux()/rAUf);
239 
240  phasei = 0;
241  phi = dimensionedScalar("phi", phi.dimensions(), 0);
242  forAllIter(PtrDictionary<phaseModel>, fluid.phases(), iter)
243  {
244  phaseModel& phase = iter();
245 
246  phase.phi() =
248  + rAlphaAUfs[phasei]*mSfGradp/phase.rho();
249  phi +=
251  + mag(alphafs[phasei]*rAlphaAUfs[phasei])
252  *mSfGradp/phase.rho();
253 
254  phasei++;
255  }
256 
257  // dgdt =
258 
259  // (
260  // pos0(alpha2)*(pEqnComp2 & p)/rho2
261  // - pos0(alpha1)*(pEqnComp1 & p)/rho1
262  // );
263 
264  p_rgh.relax();
265 
266  p = p_rgh + rho*gh;
267 
268  mSfGradp = pEqnIncomp.flux()/rAUf;
269 
270  U = dimensionedVector("U", dimVelocity, Zero);
271 
272  phasei = 0;
273  forAllIter(PtrDictionary<phaseModel>, fluid.phases(), iter)
274  {
275  phaseModel& phase = iter();
276  const volScalarField& alpha = phase;
277 
278  phase.U() =
279  HbyAs[phasei]
281  (
282  rAlphaAUfs[phasei]
283  *(
284  (phase.rho() - fvc::interpolate(rho))
285  *(g & mesh.Sf())
286  - ghSnGradRho
287  + mSfGradp
288  )
289  )/phase.rho();
290 
291  // phase.U() = fvc::reconstruct(phase.phi());
292  phase.U().correctBoundaryConditions();
293 
294  U += alpha*phase.U();
295 
296  phasei++;
297  }
298  }
299  }
300 
301  // p = max(p, pMin);
302 
303  #include "continuityErrs.H"
304 
305  // rho1 = rho10 + psi1*p_rgh;
306  // rho2 = rho20 + psi2*p_rgh;
307 
308  // Dp1Dt = fvc::DDt(phi1, p_rgh);
309  // Dp2Dt = fvc::DDt(phi2, p_rgh);
310 }
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
rho
Definition: pEqn.H:1
pimpleNoLoopControl & pimple
multiphaseSystem & fluid
Definition: createFields.H:11
p
Definition: pEqn.H:50
label phasei
Definition: pEqn.H:27
PtrList< volScalarField > rAUs
Definition: pEqn.H:4
IOMRFZoneList & MRF
phiHbyA
Definition: pEqn.H:20
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:170
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
autoPtr< multiphaseSystem::dragCoeffFields > dragCoeffs(fluid.dragCoeffs())
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
rhoEqn solve()
static const zero Zero
Definition: zero.H:97
forAllIter(PtrDictionary< phaseModel >, fluid.phases(), iter)
Definition: pEqn.H:34
p_rgh
Definition: pEqn.H:152
const surfaceScalarField & ghf
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
PtrList< volVectorField > HbyAs(fluid.phases().size())
PtrList< surfaceScalarField > alphafs(phases.size())
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
U
Definition: pEqn.H:72
surfaceScalarField ghSnGradRho(ghf *fvc::snGrad(rho) *mesh.magSf())
PtrList< fvVectorMatrix > UEqns(fluid.phases().size())
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label pRefCell
Definition: createFields.H:106
setSnGrad< fixedFluxPressureFvPatchScalarField >(p_rgh.boundaryFieldRef(),(phiHbyA.boundaryField() - MRF.relative(phib))/(mesh.magSf().boundaryField() *rAUf.boundaryField()))
PtrList< surfaceScalarField > phiHbyAs(fluid.phases().size())
void correctBoundaryConditions()
Correct boundary field.
const volScalarField & gh
dimensioned< scalar > mag(const dimensioned< Type > &)
phi
Definition: pEqn.H:18
phib
Definition: pEqn.H:194
phaseModel & phase2
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const dimensionedVector & g
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
PtrList< surfaceScalarField > rAlphaAUfs(fluid.phases().size())
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
scalar pRefValue
Definition: createFields.H:107