Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
applications
solvers
heatTransfer
buoyantBoussinesqSimpleFoam
pEqn.H
Go to the documentation of this file.
1
{
2
volScalarField
rUA
(
"rUA"
, 1.0/
UEqn
().
A
());
3
surfaceScalarField
rUAf
(
"(1|A(U))"
,
fvc::interpolate
(
rUA
));
4
5
U
=
rUA
*
UEqn
().H();
6
UEqn
.clear();
7
8
phi
=
fvc::interpolate
(
U
) &
mesh
.Sf();
9
adjustPhi
(
phi
,
U
,
p_rgh
);
10
11
surfaceScalarField
buoyancyPhi =
rUAf
*
ghf
*
fvc::snGrad
(
rhok
)*
mesh
.magSf();
12
phi
-= buoyancyPhi;
13
14
for
(
int
nonOrth=0; nonOrth<=
nNonOrthCorr
; nonOrth++)
15
{
16
fvScalarMatrix
p_rghEqn
17
(
18
fvm::laplacian
(
rUAf
,
p_rgh
) ==
fvc::div
(
phi
)
19
);
20
21
p_rghEqn.setReference(
pRefCell
,
getRefCellValue
(
p_rgh
,
pRefCell
));
22
23
// retain the residual from the first iteration
24
if
(nonOrth == 0)
25
{
26
eqnResidual
= p_rghEqn.solve().initialResidual();
27
maxResidual
=
max
(
eqnResidual
,
maxResidual
);
28
}
29
else
30
{
31
p_rghEqn.solve();
32
}
33
34
if
(nonOrth ==
nNonOrthCorr
)
35
{
36
// Calculate the conservative fluxes
37
phi
-= p_rghEqn.flux();
38
39
// Explicitly relax pressure for momentum corrector
40
p_rgh
.relax();
41
42
// Correct the momentum source with the pressure gradient flux
43
// calculated from the relaxed pressure
44
U
-=
rUA
*
fvc::reconstruct
((buoyancyPhi + p_rghEqn.flux())/
rUAf
);
45
U
.correctBoundaryConditions();
46
}
47
}
48
49
#include <
finiteVolume/continuityErrs.H
>
50
51
p
=
p_rgh
+
rhok
*
gh
;
52
53
if
(
p_rgh
.needReference())
54
{
55
p
+=
dimensionedScalar
56
(
57
"p"
,
58
p
.dimensions(),
59
pRefValue
-
getRefCellValue
(
p
,
pRefCell
)
60
);
61
p_rgh
=
p
-
rhok
*
gh
;
62
}
63
}