Skip to content
GitLab
Projects
Groups
Snippets
/
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
Menu
Open sidebar
Development
openfoam
Commits
370bccf1
Commit
370bccf1
authored
Nov 12, 2013
by
mattijs
Browse files
Merge branch 'master' of /home/dm4/OpenFOAM/OpenFOAM-dev
parents
13848fcc
75c192fb
Changes
9
Hide whitespace changes
Inline
Side-by-side
src/lagrangian/intermediate/clouds/Templates/KinematicCloud/KinematicCloud.C
View file @
370bccf1
...
...
@@ -321,6 +321,7 @@ Foam::KinematicCloud<CloudType>::KinematicCloud
:
-
1
),
cellOccupancyPtr_
(),
cellLengthScale_
(
cbrt
(
mesh_
.
V
())),
rho_
(
rho
),
U_
(
U
),
mu_
(
mu
),
...
...
@@ -421,6 +422,7 @@ Foam::KinematicCloud<CloudType>::KinematicCloud
subModelProperties_
(
c
.
subModelProperties_
),
rndGen_
(
c
.
rndGen_
,
true
),
cellOccupancyPtr_
(
NULL
),
cellLengthScale_
(
c
.
cellLengthScale_
),
rho_
(
c
.
rho_
),
U_
(
c
.
U_
),
mu_
(
c
.
mu_
),
...
...
@@ -511,6 +513,7 @@ Foam::KinematicCloud<CloudType>::KinematicCloud
subModelProperties_
(
dictionary
::
null
),
rndGen_
(
0
,
0
),
cellOccupancyPtr_
(
NULL
),
cellLengthScale_
(
c
.
cellLengthScale_
),
rho_
(
c
.
rho_
),
U_
(
c
.
U_
),
mu_
(
c
.
mu_
),
...
...
@@ -842,7 +845,9 @@ void Foam::KinematicCloud<CloudType>::patchData
template
<
class
CloudType
>
void
Foam
::
KinematicCloud
<
CloudType
>::
updateMesh
()
{
updateCellOccupancy
();
injectors_
.
updateMesh
();
cellLengthScale_
=
cbrt
(
mesh_
.
V
());
}
...
...
src/lagrangian/intermediate/clouds/Templates/KinematicCloud/KinematicCloud.H
View file @
370bccf1
...
...
@@ -165,6 +165,9 @@ protected:
//- Cell occupancy information for each parcel, (demand driven)
autoPtr
<
List
<
DynamicList
<
parcelType
*>
>
>
cellOccupancyPtr_
;
//- Cell length scale
scalarField
cellLengthScale_
;
// References to the carrier gas fields
...
...
@@ -368,6 +371,9 @@ public:
// if particles are removed or created.
inline
List
<
DynamicList
<
parcelType
*>
>&
cellOccupancy
();
//- Return the cell length scale
inline
const
scalarField
&
cellLengthScale
()
const
;
// References to the carrier gas fields
...
...
src/lagrangian/intermediate/clouds/Templates/KinematicCloud/KinematicCloudI.H
View file @
370bccf1
...
...
@@ -517,6 +517,14 @@ Foam::KinematicCloud<CloudType>::cellOccupancy()
}
template
<
class
CloudType
>
inline
const
Foam
::
scalarField
&
Foam
::
KinematicCloud
<
CloudType
>::
cellLengthScale
()
const
{
return
cellLengthScale_
;
}
template
<
class
CloudType
>
inline
Foam
::
DimensionedField
<
Foam
::
vector
,
Foam
::
volMesh
>&
Foam
::
KinematicCloud
<
CloudType
>::
UTrans
()
...
...
src/lagrangian/intermediate/parcels/Templates/KinematicParcel/KinematicParcel.C
View file @
370bccf1
...
...
@@ -265,7 +265,7 @@ bool Foam::KinematicParcel<ParcelType>::move
const
polyMesh
&
mesh
=
td
.
cloud
().
pMesh
();
const
polyBoundaryMesh
&
pbMesh
=
mesh
.
boundaryMesh
();
const
scalarField
&
V
=
mesh
.
cellVolumes
();
const
scalarField
&
cellLengthScale
=
td
.
cloud
().
cellLengthScale
();
const
scalar
maxCo
=
td
.
cloud
().
solution
().
maxCo
();
scalar
tEnd
=
(
1
.
0
-
p
.
stepFraction
())
*
trackTime
;
...
...
@@ -290,7 +290,7 @@ bool Foam::KinematicParcel<ParcelType>::move
if
(
p
.
active
()
&&
moving
&&
(
magU
>
ROOTVSMALL
))
{
const
scalar
d
=
dt
*
magU
;
const
scalar
dCorr
=
min
(
d
,
maxCo
*
c
brt
(
V
[
cellI
])
)
;
const
scalar
dCorr
=
min
(
d
,
maxCo
*
c
ellLengthScale
[
cellI
]);
dt
*=
dCorr
/
d
*
p
.
trackToFace
(
p
.
position
()
+
dCorr
*
U_
/
magU
,
td
);
...
...
src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InflationInjection/InflationInjection.C
View file @
370bccf1
...
...
@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-201
2
OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-201
3
OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
...
...
@@ -213,8 +213,8 @@ Foam::label Foam::InflationInjection<CloudType>::parcelsToInject
if
((
time0
>=
0
.
0
)
&&
(
time0
<
duration_
))
{
volumeAccumulator_
+=
fraction_
*
flowRateProfile_
.
integrate
(
time0
,
time1
);
volumeAccumulator_
+=
fraction_
*
flowRateProfile_
.
integrate
(
time0
,
time1
);
}
labelHashSet
cellCentresUsed
;
...
...
src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.C
View file @
370bccf1
...
...
@@ -586,27 +586,25 @@ void Foam::InjectionModel<CloudType>::inject(TrackData& td)
pPtr
->
rho
()
);
const
scalar
mParcel0
=
pPtr
->
nParticle
()
*
pPtr
->
mass
();
if
(
!
pPtr
->
move
(
td
,
dt
))
{
massAdded
+=
mParcel0
;
delete
pPtr
;
}
else
if
(
pPtr
->
nParticle
()
>=
1
.
0
)
{
if
(
pPtr
->
nParticle
()
>=
1
.
0
)
parcelsAdded
++
;
massAdded
+=
pPtr
->
nParticle
()
*
pPtr
->
mass
();
if
(
pPtr
->
move
(
td
,
dt
))
{
td
.
cloud
().
addParticle
(
pPtr
);
massAdded
+=
mParcel0
;
parcelsAdded
++
;
}
else
{
delayedVolume
+=
pPtr
->
nParticle
()
*
pPtr
->
volume
();
delete
pPtr
;
}
}
else
{
delayedVolume
+=
pPtr
->
nParticle
()
*
pPtr
->
volume
();
delete
pPtr
;
}
}
}
}
...
...
src/lagrangian/spray/parcels/Templates/SprayParcel/SprayParcel.C
View file @
370bccf1
...
...
@@ -93,13 +93,11 @@ void Foam::SprayParcel<ParcelType>::calc
td
.
cloud
().
constProps
().
setTMax
(
TMax
);
// store the parcel properties
const
scalarField
&
Y
(
this
->
Y
());
scalarField
X
(
composition
.
liquids
().
X
(
Y
));
this
->
Cp
()
=
composition
.
liquids
().
Cp
(
this
->
pc_
,
T0
,
X
);
sigma_
=
composition
.
liquids
().
sigma
(
this
->
pc_
,
T0
,
X
);
scalar
rho0
=
composition
.
liquids
().
rho
(
this
->
pc_
,
T0
,
X
);
this
->
Cp
()
=
composition
.
liquids
().
Cp
(
pc0
,
T0
,
X0
);
sigma_
=
composition
.
liquids
().
sigma
(
pc0
,
T0
,
X0
);
scalar
rho0
=
composition
.
liquids
().
rho
(
pc0
,
T0
,
X0
);
this
->
rho
()
=
rho0
;
mu_
=
composition
.
liquids
().
mu
(
pc0
,
T0
,
X0
);
ParcelType
::
calc
(
td
,
dt
,
cellI
);
...
...
@@ -113,11 +111,13 @@ void Foam::SprayParcel<ParcelType>::calc
this
->
Cp
()
=
composition
.
liquids
().
Cp
(
this
->
pc_
,
T1
,
X1
);
sigma_
=
composition
.
liquids
().
sigma
(
this
->
pc_
,
T1
,
X
);
sigma_
=
composition
.
liquids
().
sigma
(
this
->
pc_
,
T1
,
X
1
);
scalar
rho1
=
composition
.
liquids
().
rho
(
this
->
pc_
,
T1
,
X1
);
this
->
rho
()
=
rho1
;
mu_
=
composition
.
liquids
().
mu
(
this
->
pc_
,
T1
,
X1
);
scalar
d1
=
this
->
d
()
*
cbrt
(
rho0
/
rho1
);
this
->
d
()
=
d1
;
...
...
@@ -360,10 +360,6 @@ void Foam::SprayParcel<ParcelType>::solveTABEq
const
scalar
dt
)
{
typedef
typename
TrackData
::
cloudType
::
reactingCloudType
reactingCloudType
;
const
CompositionModel
<
reactingCloudType
>&
composition
=
td
.
cloud
().
composition
();
const
scalar
&
TABCmu
=
td
.
cloud
().
breakup
().
TABCmu
();
const
scalar
&
TABWeCrit
=
td
.
cloud
().
breakup
().
TABWeCrit
();
const
scalar
&
TABComega
=
td
.
cloud
().
breakup
().
TABComega
();
...
...
@@ -372,26 +368,19 @@ void Foam::SprayParcel<ParcelType>::solveTABEq
scalar
r2
=
r
*
r
;
scalar
r3
=
r
*
r2
;
const
scalarField
&
Y
(
this
->
Y
());
scalarField
X
(
composition
.
liquids
().
X
(
Y
));
scalar
rho
=
composition
.
liquids
().
rho
(
this
->
pc
(),
this
->
T
(),
X
);
scalar
mu
=
composition
.
liquids
().
mu
(
this
->
pc
(),
this
->
T
(),
X
);
scalar
sigma
=
composition
.
liquids
().
sigma
(
this
->
pc
(),
this
->
T
(),
X
);
// inverse of characteristic viscous damping time
scalar
rtd
=
0
.
5
*
TABCmu
*
mu
/
(
rho
*
r2
);
scalar
rtd
=
0
.
5
*
TABCmu
*
mu
_
/
(
this
->
rho
()
*
r2
);
// oscillation frequency (squared)
scalar
omega2
=
TABComega
*
sigma
/
(
rho
*
r3
)
-
rtd
*
rtd
;
scalar
omega2
=
TABComega
*
sigma
_
/
(
this
->
rho
()
*
r3
)
-
rtd
*
rtd
;
if
(
omega2
>
0
)
{
scalar
omega
=
sqrt
(
omega2
);
scalar
rhoc
=
this
->
rhoc
();
scalar
We
tmp
=
this
->
We
(
this
->
U
(),
r
,
rhoc
,
sigma
)
/
TABWeCrit
;
scalar
We
=
this
->
We
(
this
->
U
(),
r
,
rhoc
,
sigma
_
)
/
TABWeCrit
;
scalar
y1
=
this
->
y
()
-
We
tmp
;
scalar
y1
=
this
->
y
()
-
We
;
scalar
y2
=
this
->
yDot
()
/
omega
;
// update distortion parameters
...
...
@@ -400,7 +389,7 @@ void Foam::SprayParcel<ParcelType>::solveTABEq
scalar
e
=
exp
(
-
rtd
*
dt
);
y2
=
(
this
->
yDot
()
+
y1
*
rtd
)
/
omega
;
this
->
y
()
=
We
tmp
+
e
*
(
y1
*
c
+
y2
*
s
);
this
->
y
()
=
We
+
e
*
(
y1
*
c
+
y2
*
s
);
if
(
this
->
y
()
<
0
)
{
this
->
y
()
=
0
.
0
;
...
...
@@ -408,7 +397,7 @@ void Foam::SprayParcel<ParcelType>::solveTABEq
}
else
{
this
->
yDot
()
=
(
We
tmp
-
this
->
y
())
*
rtd
+
e
*
omega
*
(
y2
*
c
-
y1
*
s
);
this
->
yDot
()
=
(
We
-
this
->
y
())
*
rtd
+
e
*
omega
*
(
y2
*
c
-
y1
*
s
);
}
}
else
...
...
src/lagrangian/spray/submodels/BreakupModel/ReitzDiwakar/ReitzDiwakar.C
View file @
370bccf1
...
...
@@ -101,11 +101,9 @@ bool Foam::ReitzDiwakar<CloudType>::update
scalar
We
=
0
.
5
*
rhoc
*
sqr
(
Urmag
)
*
d
/
sigma
;
scalar
Re
=
Urmag
*
d
/
nuc
;
scalar
sqRey
=
sqrt
(
Re
);
if
(
We
>
Cbag_
)
{
if
(
We
>
Cstrip_
*
sqRe
y
)
if
(
We
>
Cstrip_
*
sq
rt
(
Re
)
)
{
scalar
dStrip
=
sqr
(
2
.
0
*
Cstrip_
*
sigma
)
/
(
rhoc
*
pow3
(
Urmag
)
*
muc
);
scalar
tauStrip
=
Cs_
*
d
*
sqrt
(
rho
/
rhoc
)
/
Urmag
;
...
...
@@ -117,9 +115,7 @@ bool Foam::ReitzDiwakar<CloudType>::update
else
{
scalar
dBag
=
2
.
0
*
Cbag_
*
sigma
/
(
rhoc
*
sqr
(
Urmag
));
scalar
tauBag
=
Cb_
*
d
*
sqrt
(
rho
*
d
/
sigma
);
scalar
fraction
=
dt
/
tauBag
;
// new droplet diameter, implicit calculation
...
...
tutorials/lagrangian/reactingParcelFoam/verticalChannel/constant/chemistryProperties
View file @
370bccf1
...
...
@@ -17,7 +17,7 @@ FoamFile
chemistryType
{
chemistrySolver
ode
;
chemistrySolver
noChemistrySolver
;
chemistryThermo rho;
}
...
...
Write
Preview
Supports
Markdown
0%
Try again
or
attach a new 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