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;
29  (
30  PtrDictionary<phaseModel>,
31  mixture.phases(),
32  phase
33  )
34  {
35  const rhoThermo& thermo = phase().thermo();
36  const volScalarField& rho = thermo.rho()();
37 
38  p_rghEqnComps.set
39  (
40  phasei,
41  (
42  fvc::ddt(rho) + thermo.psi()*correction(fvm::ddt(p_rgh))
43  + fvc::div(phi, rho) - fvc::Sp(fvc::div(phi), rho)
44  ).ptr()
45  );
46 
47  phasei++;
48  }
49 
50  // Cache p_rgh prior to solve for density update
52 
53  while (pimple.correctNonOrthogonal())
54  {
55  fvScalarMatrix p_rghEqnIncomp
56  (
59  );
60 
61  tmp<fvScalarMatrix> p_rghEqnComp;
62 
63  phasei = 0;
65  (
66  PtrDictionary<phaseModel>,
67  mixture.phases(),
68  phase
69  )
70  {
71  tmp<fvScalarMatrix> p_rghEqnCompi
72  (
73  (max(phase(), scalar(0))/phase().thermo().rho())
74  *p_rghEqnComps[phasei]
75  );
76 
77  if (phasei == 0)
78  {
79  p_rghEqnComp = p_rghEqnCompi;
80  }
81  else
82  {
83  p_rghEqnComp.ref() += p_rghEqnCompi;
84  }
85 
86  phasei++;
87  }
88 
89  solve(p_rghEqnComp + p_rghEqnIncomp);
90 
91  if (pimple.finalNonOrthogonalIter())
92  {
93  p = max(p_rgh + mixture.rho()*gh, pMin);
94  p_rgh = p - rho*gh;
95 
96  phasei = 0;
98  (
99  PtrDictionary<phaseModel>,
100  mixture.phases(),
101  phase
102  )
103  {
104  phase().dgdt() =
105  pos0(phase())
106  *(p_rghEqnComps[phasei] & p_rgh)/phase().thermo().rho();
107 
108  phasei++;
109  }
110 
111  phi = phiHbyA + p_rghEqnIncomp.flux();
112 
113  U = HbyA
114  + rAU*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
115  U.correctBoundaryConditions();
116  }
117  }
118 
119  // Update densities from change in p_rgh
120  mixture.correctRho(p_rgh - p_rgh_0);
121  mixture.correct();
122 
123  // Correct p_rgh for consistency with p and the updated densities
124  p_rgh = p - rho*gh;
125  p_rgh.correctBoundaryConditions();
126 
127  K = 0.5*magSqr(U);
128 
129  Info<< "max(U) " << max(mag(U)).value() << endl;
130  Info<< "min(p_rgh) " << min(p_rgh).value() << endl;
131 }
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
rAU
Definition: pEqn.H:1
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
fluidReactionThermo & thermo
Definition: createFields.H:28
pimpleNoLoopControl & pimple
label phasei
Definition: pEqn.H:27
p_rgh
Definition: pEqn.H:159
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
dimensionedScalar pMin("pMin", dimPressure, mixture)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
CGAL::Exact_predicates_exact_constructions_kernel K
phi
Definition: pEqn.H:99
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
rho
Definition: pEqn.H:1
volScalarField p_rgh_0(p_rgh)
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedScalar pos0(const dimensionedScalar &ds)
const surfaceScalarField & ghf
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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.
surfaceScalarField phig(-rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
phiHbyA
Definition: pEqn.H:32
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
rhoEqn solve()
U
Definition: pEqn.H:72
volVectorField & HbyA
Definition: pEqn.H:13
const volScalarField & gh
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
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
PtrList< fvScalarMatrix > p_rghEqnComps(mixture.phases().size())
p
Definition: pEqn.H:101
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF)
fvVectorMatrix & UEqn
Definition: UEqn.H:13