pEqn.H
Go to the documentation of this file.
1 {
2  volScalarField rAU("rAU", 1.0/UEqn.A());
6  (
7  "phiHbyA",
10  );
11 
13  (
14  (
15  mixture.surfaceTensionForce()
17  )*rAUf*mesh.magSf()
18  );
19 
21 
22  // Update the pressure BCs to ensure flux consistency
24 
25  PtrList<fvScalarMatrix> p_rghEqnComps(mixture.phases().size());
26 
27  label phasei = 0;
28  forAllConstIters(mixture.phases(), phase)
29  {
30  const rhoThermo& thermo = phase().thermo();
31  const volScalarField& rho = thermo.rho()();
32 
33  p_rghEqnComps.set
34  (
35  phasei,
36  (
39  ).ptr()
40  );
41 
42  ++phasei;
43  }
44 
45  // Cache p_rgh prior to solve for density update
47 
48  while (pimple.correctNonOrthogonal())
49  {
50  fvScalarMatrix p_rghEqnIncomp
51  (
54  );
55 
56  tmp<fvScalarMatrix> p_rghEqnComp;
57 
58  phasei = 0;
59  forAllConstIters(mixture.phases(), phase)
60  {
61  tmp<fvScalarMatrix> hmm
62  (
63  (max(phase(), scalar(0))/phase().thermo().rho())
65  );
66 
67  if (phasei == 0)
68  {
69  p_rghEqnComp = hmm;
70  }
71  else
72  {
73  p_rghEqnComp.ref() += hmm;
74  }
75 
76  ++phasei;
77  }
78 
79  solve
80  (
81  p_rghEqnComp
82  + p_rghEqnIncomp,
83  mesh.solver(p_rgh.select(pimple.finalInnerIter()))
84  );
85 
86  if (pimple.finalNonOrthogonalIter())
87  {
88  phasei = 0;
89  for (phaseModel& phase : mixture.phases())
90  {
91  phase.dgdt() =
92  pos0(phase)
93  *(p_rghEqnComps[phasei] & p_rgh)/phase.thermo().rho();
94 
95  ++phasei;
96  }
97 
98  phi = phiHbyA + p_rghEqnIncomp.flux();
99 
100  U = HbyA
101  + rAU*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
102  U.correctBoundaryConditions();
103  }
104  }
105 
106  p = max(p_rgh + mixture.rho()*gh, pMin);
107 
108  // Update densities from change in p_rgh
109  mixture.correctRho(p_rgh - p_rgh_0);
110  rho = mixture.rho();
111 
112  // Correct p_rgh for consistency with p and the updated densities
113  p_rgh = p - rho*gh;
114  p_rgh.correctBoundaryConditions();
115 
116  K = 0.5*magSqr(U);
117 
118  Info<< "max(U) " << max(mag(U)).value() << endl;
119  Info<< "min(p_rgh) " << min(p_rgh).value() << endl;
120 }
Foam::fvc::reconstruct
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcReconstruct.C:56
Foam::constrainHbyA
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:35
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
ghf
const surfaceScalarField & ghf
Definition: setRegionFluidFields.H:18
pMin
const dimensionedScalar & pMin
Definition: setRegionFluidFields.H:58
phiHbyA
phiHbyA
Definition: pEqn.H:20
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
Sp
zeroField Sp
Definition: alphaSuSp.H:2
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
phasei
label phasei
Definition: pEqn.H:27
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::correction
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
gh
const volScalarField & gh
Definition: setRegionFluidFields.H:17
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
solve
CEqn solve()
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
HbyA
HbyA
Definition: pEqn.H:4
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
ddtCorr
ddtCorr
Definition: readControls.H:9
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
p
p
Definition: pEqn.H:50
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
rAU
volScalarField rAU(1.0/UEqn.A())
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
p_rghEqnComps
PtrList< fvScalarMatrix > p_rghEqnComps(mixture.phases().size())
U
U
Definition: pEqn.H:72
phig
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
p_rgh_0
volScalarField p_rgh_0(p_rgh)
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::fac::ddt
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:47
phi
phi
Definition: pEqn.H:18
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
p_rgh
p_rgh
Definition: pEqn.H:139
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
rho
rho
Definition: pEqn.H:1
mixture
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
interpolate
mesh interpolate(rAU)
constrainPressure
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
rAUf
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))