DAFoam
v4.0.2
Discrete Adjoint with OpenFOAM
dafoam
src
adjoint
DASolver
DARhoPimpleFoam
pEqnRhoPimple.H
Go to the documentation of this file.
1
rho
=
thermo
.rho();
2
DAUtility::boundVar(
allOptions
,
rho
, pimplePrintToScreen);
3
rho
.relax();
4
5
volScalarField
rAU
(1.0 /
UEqn
.A());
6
surfaceScalarField
rhorAUf
(
"rhorAUf"
, fvc::interpolate(
rho
*
rAU
));
7
//***************** NOTE *******************
8
// constrainHbyA has been used since OpenFOAM-v1606; however, it may degrade the accuracy of derivatives
9
// because constraining variables will create discontinuity. Here we have a option to use the old
10
// implementation in OpenFOAM-3.0+ and before (no constraint for HbyA)
11
autoPtr<volVectorField>
HbyAPtr
=
nullptr
;
12
label
useConstrainHbyA
= daOptionPtr_->getOption<label>(
"useConstrainHbyA"
);
13
if
(
useConstrainHbyA
)
14
{
15
HbyAPtr
.reset(
new
volVectorField(constrainHbyA(
rAU
*
UEqn
.H(),
U
,
p
)));
16
}
17
else
18
{
19
HbyAPtr
.reset(
new
volVectorField(
"HbyA"
,
U
));
20
HbyAPtr
() =
rAU
*
UEqn
.H();
21
}
22
volVectorField&
HbyA
=
HbyAPtr
();
23
24
if
(
pimple
.nCorrPISO() <= 1)
25
{
26
tUEqn
.clear();
27
}
28
29
surfaceScalarField
phiHbyA
(
"phiHbyA"
, fvc::interpolate(
rho
) * fvc::flux(
HbyA
));
30
31
// NOTE: we don't support transonic = true
32
33
while
(
pimple
.correctNonOrthogonal())
34
{
35
fvScalarMatrix pEqn(
36
fvm::ddt(
psi
,
p
)
37
+ fvc::div(
phiHbyA
)
38
- fvm::laplacian(
rhorAUf
,
p
));
39
40
// get the solver performance info such as initial
41
// and final residuals
42
SolverPerformance<scalar> solverP = pEqn.solve(
mesh
.solver(
p
.select(
pimple
.finalInnerIter())));
43
44
DAUtility::primalResidualControl(solverP, pimplePrintToScreen,
"p"
, daGlobalVarPtr_->primalMaxRes);
45
46
if
(
pimple
.finalNonOrthogonalIter())
47
{
48
phi
=
phiHbyA
+ pEqn.flux();
49
}
50
}
51
52
#include "
rhoEqnRhoPimple.H
"
53
54
if
(pimplePrintToScreen)
55
{
56
#include "
continuityErrsPython.H
"
57
}
58
59
// Explicitly relax pressure for momentum corrector
60
p
.relax();
61
62
// bound p
63
DAUtility::boundVar(
allOptions
,
p
, pimplePrintToScreen);
64
65
rho
=
thermo
.rho();
66
DAUtility::boundVar(
allOptions
,
rho
, pimplePrintToScreen);
67
rho
.relax();
68
69
U
=
HbyA
-
rAU
* fvc::grad(
p
);
70
// bound U
71
DAUtility::boundVar(
allOptions
,
U
, pimplePrintToScreen);
72
U
.correctBoundaryConditions();
73
74
K
= 0.5 * magSqr(
U
);
75
76
if
(
thermo
.dpdt())
77
{
78
dpdt
= fvc::ddt(
p
);
79
}
pimple
pimpleControlDF & pimple
Definition:
createRefsPimpleDyM.H:5
HbyAPtr
autoPtr< volVectorField > HbyAPtr
Definition:
pEqnRhoPimple.H:11
rhorAUf
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
UEqn
fvVectorMatrix & UEqn
Definition:
UEqnPimpleDyM.H:15
thermo
fluidThermo & thermo
Definition:
createRefsRhoPimple.H:6
continuityErrsPython.H
rho
rho
Definition:
pEqnRhoPimple.H:1
HbyA
volVectorField & HbyA
Definition:
pEqnRhoPimple.H:22
mesh
fvMesh & mesh
Definition:
createRefsHeatTransfer.H:4
dpdt
volScalarField & dpdt
Definition:
createRefsRhoPimple.H:11
U
U
Definition:
pEqnRhoPimple.H:69
p
volScalarField & p
Definition:
createRefsPimpleDyM.H:6
phiHbyA
surfaceScalarField phiHbyA("phiHbyA", fvc::interpolate(rho) *fvc::flux(HbyA))
phi
surfaceScalarField & phi
Definition:
createRefsPimpleDyM.H:8
tUEqn
tmp< fvVectorMatrix > tUEqn(fvm::ddt(U)+fvm::div(phi, U)+turbulencePtr_->divDevReff(U) - fvSource)
allOptions
const dictionary & allOptions
Definition:
createRefsPimpleDyM.H:14
K
K
Definition:
pEqnRhoPimple.H:74
rAU
volScalarField rAU(1.0/UEqn.A())
rhoEqnRhoPimple.H
useConstrainHbyA
label useConstrainHbyA
Definition:
pEqnRhoPimple.H:12
psi
const volScalarField & psi
Definition:
createRefsRhoPimple.H:14
Generated by
1.8.17