MoDeNa  1.0
Software framework facilitating sequential multi-scale modelling
pEqn.H
Go to the documentation of this file.
1 
7 {
9  Psi1 = rho_gas/p_rgh;
10  Psi2 = rho_foam/p_rgh;
11 
12  volScalarField rAU("rAU", 1.0/UEqn.A());
13  surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU));
14 
15  volVectorField HbyA("HbyA", U);
16  HbyA = rAU*UEqn.H();
17 
18  surfaceScalarField phiHbyA
19  (
20  "phiHbyA",
21  (fvc::interpolate(HbyA) & mesh.Sf())
22  + fvc::interpolate(rho*rAU)*fvc::ddtCorr(U, phi)
23  );
24 
25  surfaceScalarField phig
26  (
27  (
28  interface.surfaceTensionForce()
29  //fvc::interpolate(interface.sigmaK())*fvc::snGrad(alpha1)
30  - ghf*fvc::snGrad(rho)
31  )*rAUf*mesh.magSf()
32  );
33 
34  phiHbyA += phig;
35  // Update the pressure BCs to ensure flux consistency
36  constrainPressure(p_rgh, U, phiHbyA, rAUf);
37 
38  tmp<fvScalarMatrix> p_rghEqnComp1;
39  tmp<fvScalarMatrix> p_rghEqnComp2;
40 
41  if (pimple.transonic())
42  {
43  surfaceScalarField phid1("phid1", fvc::interpolate(Psi1)*phi);
44  surfaceScalarField phid2("phid2", fvc::interpolate(Psi2)*phi);
45 
46  p_rghEqnComp1 =
47  fvc::ddt(rho_gas) + fvc::div(phi, rho_gas)
48  - fvc::Sp(fvc::div(phi), rho_gas)
49  + correction
50  (
51  Psi1*fvm::ddt(p_rgh)
52  + fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
53  );
54  deleteDemandDrivenData(p_rghEqnComp1.ref().faceFluxCorrectionPtr());
55  p_rghEqnComp1.ref().relax();
56 
57  p_rghEqnComp2 =
58  fvc::ddt(rho_foam) + fvc::div(phi, rho_foam)
59  - fvc::Sp(fvc::div(phi), rho_foam)
60  + correction
61  (
62  Psi2*fvm::ddt(p_rgh)
63  + fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
64  );
65  deleteDemandDrivenData(p_rghEqnComp2.ref().faceFluxCorrectionPtr());
66  p_rghEqnComp2.ref().relax();
67  }
68 
69  else
70  {
71  p_rghEqnComp1 =
72  (
73  fvc::ddt(alpha1, rho_gas)
74  + fvc::div(fvc::interpolate(alpha1)*phi, rho_gas)
75  - fvc::Sp(fvc::ddt(alpha1)
76  + fvc::div(fvc::interpolate(alpha1)*phi), rho_gas)
77  )/rho_gas
78  + (alpha1*Psi1/rho_gas)*correction(fvm::ddt(p_rgh));
79 
80  p_rghEqnComp2 =
81  (
82  fvc::ddt(alpha2, rho_foam)
83  + fvc::div(fvc::interpolate(alpha2)*phi, rho_foam)
84  - fvc::Sp(fvc::ddt(alpha2)
85  + fvc::div(fvc::interpolate(alpha2)*phi), rho_foam)
86  )/rho_foam
87  + (alpha2*Psi2/rho_foam)*correction(fvm::ddt(p_rgh));
88  }
89 
90  // Cache p_rgh prior to solve for density update
91  volScalarField p_rgh_0(p_rgh);
92 
93  while (pimple.correctNonOrthogonal())
94  {
95  fvScalarMatrix p_rghEqnIncomp
96  (
97  fvc::div(phiHbyA)
98  - fvm::laplacian(rAUf, p_rgh)
99  );
100  solve
101  (
102  (
103  p_rghEqnComp1() + p_rghEqnComp2()
104  )
105  + p_rghEqnIncomp,
106  mesh.solver(p_rgh.select(pimple.finalInnerIter()))
107  );
108  if (pimple.finalNonOrthogonalIter())
109  {
110  p = max(p_rgh + (alpha1*rho_gas + (1.0 - alpha1)*rho_foam)*gh, pMin);
111  p_rgh = p - (alpha1*rho_gas + (1.0 - alpha1)*rho_foam)*gh;
112  dimensionedScalar dumdgdt ("dumdgdt",dimensionSet(1,-3,0,0,0,0,0), 1.0);
113  dgdt =
114  dumdgdt*
115  (
116  pos((1.0 - alpha1))*(p_rghEqnComp2 & p_rgh)/rho_foam
117  - pos(alpha1)*(p_rghEqnComp1 & p_rgh)/rho_gas
118  );
119 
120  phi = phiHbyA + p_rghEqnIncomp.flux();
121  U = HbyA
122  + rAU*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
123  U.correctBoundaryConditions();
124  }
125  }
126  p = max(p_rgh + (alpha1*rho_gas + (1.0 - alpha1)*rho_foam)*gh, pMin);
127 
128  // Update densities from change in p_rgh
129  rho_gas += Psi1*(p_rgh - p_rgh_0);
130  rho_foam += Psi2*(p_rgh - p_rgh_0);
131  rho = alpha1*rho_gas + (1 - alpha1)*rho_foam;
132  K = 0.5*magSqr(U);
133  Info<< "max(U) " << max(mag(U)).value() << endl;
134  Info<< "min(p_rgh) " << min(p_rgh).value() << endl;
135 }
136 ///@endcond