Skip to content
GitLab
Projects
Groups
Snippets
Help
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
O
OpenFOAM-plus
Project overview
Project overview
Details
Activity
Releases
Repository
Repository
Files
Commits
Branches
Tags
Contributors
Graph
Compare
Issues
339
Issues
339
List
Boards
Labels
Service Desk
Milestones
Merge Requests
0
Merge Requests
0
Operations
Operations
Incidents
Analytics
Analytics
Repository
Value Stream
Wiki
Wiki
Members
Members
Collapse sidebar
Close sidebar
Activity
Graph
Create a new issue
Commits
Issue Boards
Open sidebar
Development
OpenFOAM-plus
Commits
eda41a68
Commit
eda41a68
authored
Feb 01, 2019
by
sergio
Browse files
Options
Browse Files
Download
Email Patches
Plain Diff
ENH: Changing energy Eq in interCondensatingEvaporatingFoam
from e to T. T proved to be more generic solution.
parent
7dc44ae0
Changes
4
Hide whitespace changes
Inline
Side-by-side
Showing
4 changed files
with
33 additions
and
40 deletions
+33
-40
applications/solvers/multiphase/interCondensatingEvaporatingFoam/UEqn.H
...olvers/multiphase/interCondensatingEvaporatingFoam/UEqn.H
+0
-25
applications/solvers/multiphase/interCondensatingEvaporatingFoam/createFields.H
...ultiphase/interCondensatingEvaporatingFoam/createFields.H
+14
-3
applications/solvers/multiphase/interCondensatingEvaporatingFoam/eEqn.H
...olvers/multiphase/interCondensatingEvaporatingFoam/eEqn.H
+18
-9
applications/solvers/multiphase/interCondensatingEvaporatingFoam/interCondensatingEvaporatingFoam.C
...nsatingEvaporatingFoam/interCondensatingEvaporatingFoam.C
+1
-3
No files found.
applications/solvers/multiphase/interCondensatingEvaporatingFoam/UEqn.H
deleted
100644 → 0
View file @
7dc44ae0
fvVectorMatrix
UEqn
(
fvm
::
ddt
(
rho
,
U
)
+
fvm
::
div
(
rhoPhi
,
U
)
+
turbulence
->
divDevRhoReff
(
rho
,
U
)
);
UEqn
.
relax
();
if
(
pimple
.
momentumPredictor
())
{
solve
(
UEqn
==
fvc
::
reconstruct
(
(
interface
.
surfaceTensionForce
()
-
ghf
*
fvc
::
snGrad
(
rho
)
-
fvc
::
snGrad
(
p_rgh
)
)
*
mesh
.
magSf
()
)
);
}
applications/solvers/multiphase/interCondensatingEvaporatingFoam/createFields.H
View file @
eda41a68
...
...
@@ -40,9 +40,6 @@ autoPtr<temperaturePhaseChangeTwoPhaseMixture> mixture =
temperaturePhaseChangeTwoPhaseMixture
::
New
(
thermo
(),
mesh
);
// Correct e from T and alpha
//thermo->correct();
volScalarField
&
alpha1
(
thermo
->
alpha1
());
volScalarField
&
alpha2
(
thermo
->
alpha2
());
...
...
@@ -139,3 +136,17 @@ volScalarField pDivU
dimensionedScalar
(
p
.
dimensions
()
/
dimTime
,
Zero
)
);
// Need to store rho for ddt(rhoCp, U)
volScalarField
rhoCp
(
IOobject
(
"rhoCp"
,
runTime
.
timeName
(),
mesh
,
IOobject
::
NO_READ
,
IOobject
::
NO_WRITE
),
rho
*
thermo
->
Cp
()
);
rhoCp
.
oldTime
();
applications/solvers/multiphase/interCondensatingEvaporatingFoam/eEqn.H
View file @
eda41a68
{
tmp
<
volScalarField
>
tcp
(
thermo
->
Cp
());
const
volScalarField
&
cp
=
tcp
();
rhoCp
=
rho
*
cp
;
kappaEff
=
thermo
->
kappa
()
+
rho
*
cp
*
turbulence
->
nut
()
/
Prt
;
...
...
@@ -11,18 +12,26 @@
pDivU
=
(
p
*
fvc
::
div
(
rhoPhi
/
fvc
::
interpolate
(
rho
)));
}
fvScalarMatrix
eEqn
const
surfaceScalarField
rhoCpPhi
(
fvc
::
interpolate
(
cp
)
*
rhoPhi
);
Pair
<
tmp
<
volScalarField
>>
vDotAlphal
=
mixture
->
mDot
();
const
volScalarField
&
vDotcAlphal
=
vDotAlphal
[
0
]();
const
volScalarField
&
vDotvAlphal
=
vDotAlphal
[
1
]();
const
volScalarField
vDotvmcAlphal
(
vDotvAlphal
-
vDotcAlphal
);
fvScalarMatrix
TEqn
(
fvm
::
ddt
(
rho
,
e
)
+
fvm
::
div
(
rhoPhi
,
e
)
-
fvm
::
laplacian
(
kappaEff
/
cp
,
e
)
+
pDivU
fvm
::
ddt
(
rhoCp
,
T
)
+
fvm
::
div
(
rhoCpPhi
,
T
,
"div(rhoCpPhi,T)"
)
-
fvm
::
Sp
(
fvc
::
ddt
(
rhoCp
)
+
fvc
::
div
(
rhoCpPhi
),
T
)
-
fvm
::
laplacian
(
kappaEff
,
T
)
+
thermo
->
hc
()
*
vDotvmcAlphal
+
pDivU
);
eEqn
.
relax
();
eEqn
.
solve
();
thermo
->
correct
();
TEqn
.
relax
();
TEqn
.
solve
();
Info
<<
"min/max(T) = "
<<
min
(
T
).
value
()
<<
", "
<<
max
(
T
).
value
()
<<
endl
;
...
...
applications/solvers/multiphase/interCondensatingEvaporatingFoam/interCondensatingEvaporatingFoam.C
View file @
eda41a68
...
...
@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016 OpenCFD Ltd.
\\ / A nd | Copyright (C) 2016
-2019
OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
...
...
@@ -79,8 +79,6 @@ int main(int argc, char *argv[])
#include "setInitialDeltaT.H"
volScalarField
&
T
=
thermo
->
T
();
volScalarField
&
e
=
thermo
->
he
();
e
.
oldTime
();
turbulence
->
validate
();
...
...
Write
Preview
Markdown
is supported
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment