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  for (phaseModel& phase : fluid.phases())
35  {
36  MRF.makeAbsolute(phase.phi().oldTime());
37  MRF.makeAbsolute(phase.phi());
38 
39  HbyAs.set(phasei, new volVectorField(phase.U()));
40  phiHbyAs.set(phasei, new surfaceScalarField(1.0*phase.phi()));
41 
42  ++phasei;
43  }
44 
46  (
47  IOobject
48  (
49  "phiHbyA",
50  runTime.timeName(),
51  mesh,
52  IOobject::NO_READ,
53  IOobject::AUTO_WRITE,
54  IOobject::REGISTER
55  ),
56  mesh,
58  );
59 
60  volScalarField rho("rho", fluid.rho());
62 
63  phasei = 0;
64  for (phaseModel& phase : fluid.phases())
65  {
66  const volScalarField& alpha = phase;
67 
68  alphafs.set(phasei, fvc::interpolate(alpha).ptr());
69  alphafs[phasei].rename("hmm" + alpha.name());
70 
71  volScalarField dragCoeffi
72  (
73  IOobject
74  (
75  "dragCoeffi",
76  runTime.timeName(),
77  mesh
78  ),
79  fluid.dragCoeff(phase, dragCoeffs())/phase.rho(),
81  );
82  dragCoeffi.correctBoundaryConditions();
83 
84  rAUs.set(phasei, (1.0/(UEqns[phasei].A() + dragCoeffi)).ptr());
85  rAlphaAUfs.set
86  (
87  phasei,
88  (
90  /fvc::interpolate(UEqns[phasei].A() + dragCoeffi)
91  ).ptr()
92  );
93 
95 
96  phiHbyAs[phasei] =
97  (
99  + MRF.zeroFilter
100  (
101  rAlphaAUfs[phasei]*fvc::ddtCorr(phase.U(), phase.phi())
102  )
103  );
104  MRF.makeRelative(phiHbyAs[phasei]);
105  MRF.makeRelative(phase.phi().oldTime());
106  MRF.makeRelative(phase.phi());
107 
108  phiHbyAs[phasei] +=
110  *(
111  fluid.surfaceTension(phase)*mesh.magSf()
112  + (phase.rho() - fvc::interpolate(rho))*(g & mesh.Sf())
113  - ghSnGradRho
114  )/phase.rho();
115 
116  auto dmIter = fluid.dragModels().cbegin();
117  auto dcIter = dragCoeffs().cbegin();
118 
119  for
120  (
121  ;
122  dmIter.good() && dcIter.good();
123  ++dmIter, ++dcIter
124  )
125  {
126  const phaseModel *phase2Ptr = nullptr;
127 
128  if (&phase == &dmIter()->phase1())
129  {
130  phase2Ptr = &dmIter()->phase2();
131  }
132  else if (&phase == &dmIter()->phase2())
133  {
134  phase2Ptr = &dmIter()->phase1();
135  }
136  else
137  {
138  continue;
139  }
140 
141  phiHbyAs[phasei] +=
142  fvc::interpolate((*dcIter())/phase.rho())
143  /fvc::interpolate(UEqns[phasei].A() + dragCoeffi)
144  *phase2Ptr->phi();
145 
146  HbyAs[phasei] +=
147  (1.0/phase.rho())*rAUs[phasei]*(*dcIter())
148  *phase2Ptr->U();
149 
150  // Alternative flux-pressure consistent drag
151  // but not momentum conservative
152  //
153  // HbyAs[phasei] += fvc::reconstruct
154  // (
155  // fvc::interpolate
156  // (
157  // (1.0/phase.rho())*rAUs[phasei]*(*dcIter())
158  // )*phase2Ptr->phi()
159  // );
160  }
161 
163 
164  ++phasei;
165  }
166 
168  (
169  IOobject
170  (
171  "rAUf",
172  runTime.timeName(),
173  mesh
174  ),
175  mesh,
176  dimensionedScalar(dimensionSet(-1, 3, 1, 0, 0), Zero)
177  );
178 
179  phasei = 0;
180  for (const phaseModel& phase : fluid.phases())
181  {
182  rAUf += mag(alphafs[phasei]*rAlphaAUfs[phasei])/phase.rho();
183 
184  ++phasei;
185  }
186 
187  // Update the fixedFluxPressure BCs to ensure flux consistency
188  {
189  surfaceScalarField::Boundary phib(phi.boundaryField());
190  phib = 0;
191  phasei = 0;
192  for (const phaseModel& phase : fluid.phases())
193  {
194  phib +=
195  alphafs[phasei].boundaryField()
196  *(mesh.Sf().boundaryField() & phase.U().boundaryField());
197 
198  ++phasei;
199  }
200 
202  (
203  p_rgh.boundaryFieldRef(),
204  (
205  phiHbyA.boundaryField() - MRF.relative(phib)
206  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
207  );
208  }
209 
210  while (pimple.correctNonOrthogonal())
211  {
212  fvScalarMatrix pEqnIncomp
213  (
216  );
217 
218  pEqnIncomp.setReference(pRefCell, pRefValue);
219 
220  solve
221  (
222  // (
223  // (alpha1/rho1)*pEqnComp1()
224  // + (alpha2/rho2)*pEqnComp2()
225  // ) +
226  pEqnIncomp,
227  p_rgh.select(pimple.finalInnerIter())
228  );
229 
230  if (pimple.finalNonOrthogonalIter())
231  {
232  surfaceScalarField mSfGradp("mSfGradp", pEqnIncomp.flux()/rAUf);
233 
234  phasei = 0;
235  phi = dimensionedScalar("phi", phi.dimensions(), Zero);
236 
237  for (phaseModel& phase : fluid.phases())
238  {
239  phase.phi() =
241  + rAlphaAUfs[phasei]*mSfGradp/phase.rho();
242 
243  phi +=
246  *mSfGradp/phase.rho();
247 
248  ++phasei;
249  }
250 
251  // dgdt =
252 
253  // (
254  // pos0(alpha2)*(pEqnComp2 & p)/rho2
255  // - pos0(alpha1)*(pEqnComp1 & p)/rho1
256  // );
257 
258  p_rgh.relax();
259 
260  p = p_rgh + rho*gh;
261 
262  mSfGradp = pEqnIncomp.flux()/rAUf;
263 
265 
266  phasei = 0;
267  for (phaseModel& phase : fluid.phases())
268  {
269  const volScalarField& alpha = phase;
270 
271  phase.U() =
272  HbyAs[phasei]
274  (
276  *(
277  (phase.rho() - fvc::interpolate(rho))
278  *(g & mesh.Sf())
279  - ghSnGradRho
280  + mSfGradp
281  )
282  )/phase.rho();
283 
284  //phase.U() = fvc::reconstruct(phase.phi());
285  phase.U().correctBoundaryConditions();
286 
287  U += alpha*phase.U();
288 
289  ++phasei;
290  }
291  }
292  }
293 
294  //p = max(p, pMin);
295 
296  #include "continuityErrs.H"
297 
298  // rho1 = rho10 + psi1*p_rgh;
299  // rho2 = rho20 + psi2*p_rgh;
300 
301  // Dp1Dt = fvc::DDt(phi1, p_rgh);
302  // Dp2Dt = fvc::DDt(phi2, p_rgh);
303 }
twoPhaseSystem & fluid
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
const word zeroGradientType
A zeroGradient patch field type.
rho
Definition: pEqn.H:1
p
Definition: pEqn.H:50
label phasei
Definition: pEqn.H:27
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
phiHbyA
Definition: pEqn.H:20
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:165
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
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:76
phaseModel & phase2
CEqn solve()
const surfaceScalarField & ghf
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
dynamicFvMesh & mesh
IOMRFZoneList & MRF
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
p_rgh
Definition: pEqn.H:139
const uniformDimensionedVectorField & g
pimpleControl & pimple
PtrList< volVectorField > HbyAs(fluid.phases().size())
const scalar pRefValue
const word & name() const noexcept
Return const reference to name.
const label pRefCell
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
PtrList< volScalarField > rAUs
Definition: pEqn.H:4
surfaceScalarField ghSnGradRho(ghf *fvc::snGrad(rho) *mesh.magSf())
PtrList< fvVectorMatrix > UEqns(fluid.phases().size())
PtrList< surfaceScalarField > alphafs(phases.size())
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
mesh interpolate(rAU)
setSnGrad< fixedFluxPressureFvPatchScalarField >(p_rgh.boundaryFieldRef(),(phiHbyA.boundaryField() - MRF.relative(phib))/(mesh.magSf().boundaryField() *rAUf.boundaryField()))
const volScalarField & gh
PtrList< surfaceScalarField > phiHbyAs(fluid.phases().size())
void correctBoundaryConditions()
Correct boundary field.
phi
Definition: pEqn.H:18
phib
Definition: pEqn.H:190
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
PtrList< surfaceScalarField > rAlphaAUfs(fluid.phases().size())
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity