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
d92d77cc
Commit
d92d77cc
authored
Sep 28, 2018
by
Mark OLESEN
Browse files
Merge remote-tracking branch 'origin/master' into develop
parents
4a61042f
1135f157
Changes
6
Hide whitespace changes
Inline
Side-by-side
src/OpenFOAM/meshes/polyMesh/polyPatches/polyPatch/polyPatch.C
View file @
d92d77cc
...
...
@@ -336,6 +336,25 @@ Foam::tmp<Foam::vectorField> Foam::polyPatch::faceCellCentres() const
}
Foam
::
tmp
<
Foam
::
scalarField
>
Foam
::
polyPatch
::
areaFraction
()
const
{
tmp
<
scalarField
>
tfraction
(
new
scalarField
(
size
()));
scalarField
&
fraction
=
tfraction
.
ref
();
const
vectorField
::
subField
faceAreas
=
this
->
faceAreas
();
const
pointField
&
points
=
this
->
points
();
forAll
(
*
this
,
facei
)
{
const
face
&
curFace
=
this
->
operator
[](
facei
);
fraction
[
facei
]
=
mag
(
faceAreas
[
facei
])
/
(
curFace
.
mag
(
points
)
+
ROOTVSMALL
);
}
return
tfraction
;
}
const
Foam
::
labelUList
&
Foam
::
polyPatch
::
faceCells
()
const
{
if
(
!
faceCellsPtr_
)
...
...
src/OpenFOAM/meshes/polyMesh/polyPatches/polyPatch/polyPatch.H
View file @
d92d77cc
...
...
@@ -373,6 +373,10 @@ public:
//- Return face cell centres
tmp
<
vectorField
>
faceCellCentres
()
const
;
//- Return the area fraction as the ratio of the stored face area
//- and the area given by the face points
tmp
<
scalarField
>
areaFraction
()
const
;
// Addressing into mesh
...
...
src/lagrangian/basic/InteractionLists/InteractionLists.C
View file @
d92d77cc
...
...
@@ -324,10 +324,15 @@ void Foam::InteractionLists<ParticleType>::buildInteractionLists()
if
(
isA
<
wallPolyPatch
>
(
patch
))
{
localWallFaces
.
append
(
identity
(
patch
.
size
(),
patch
.
start
())
);
const
scalarField
areaFraction
(
patch
.
areaFraction
());
forAll
(
areaFraction
,
facei
)
{
if
(
areaFraction
[
facei
]
>
0
.
5
)
{
localWallFaces
.
append
(
facei
+
patch
.
start
());
}
}
}
}
...
...
@@ -335,10 +340,12 @@ void Foam::InteractionLists<ParticleType>::buildInteractionLists()
forAll
(
wallFaceBbs
,
i
)
{
wallFaceBbs
[
i
]
=
treeBoundBox
(
mesh_
.
faces
()[
localWallFaces
[
i
]].
points
(
mesh_
.
points
())
);
wallFaceBbs
[
i
]
=
treeBoundBox
(
mesh_
.
points
(),
mesh_
.
faces
()[
localWallFaces
[
i
]]
);
}
// IAndT: index and transform
...
...
tutorials/incompressible/SRFPimpleFoam/rotor2D/0/epsilon
View file @
d92d77cc
...
...
@@ -24,7 +24,6 @@ boundaryField
rotor
{
type epsilonWallFunction;
U Urel;
value uniform 3.75e-4;
}
...
...
tutorials/incompressible/SRFSimpleFoam/mixer/0/epsilon
View file @
d92d77cc
...
...
@@ -34,13 +34,11 @@ boundaryField
innerWall
{
type epsilonWallFunction;
U Urel;
value uniform 14.855;
}
outerWall
{
type epsilonWallFunction;
U Urel;
value uniform 14.855;
}
cyclic_half0
...
...
tutorials/incompressible/SRFSimpleFoam/mixer/0/nut
View file @
d92d77cc
...
...
@@ -34,13 +34,11 @@ boundaryField
innerWall
{
type nutkWallFunction;
U Urel;
value uniform 0;
}
outerWall
{
type nutkWallFunction;
U Urel;
value uniform 0;
}
cyclic_half0
...
...
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