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
C
catalyst
Project overview
Project overview
Details
Activity
Releases
Repository
Repository
Files
Commits
Branches
Tags
Contributors
Graph
Compare
Issues
2
Issues
2
List
Boards
Labels
Service Desk
Milestones
Merge Requests
0
Merge Requests
0
CI / CD
CI / CD
Pipelines
Jobs
Schedules
Operations
Operations
Incidents
Environments
Packages & Registries
Packages & Registries
Container Registry
Analytics
Analytics
CI / CD
Repository
Value Stream
Wiki
Wiki
Snippets
Snippets
Members
Members
Collapse sidebar
Close sidebar
Activity
Graph
Create a new issue
Jobs
Commits
Issue Boards
Open sidebar
Community
catalyst
Commits
9b969cd4
Commit
9b969cd4
authored
May 14, 2018
by
Mark Olesen
Browse files
Options
Browse Files
Download
Email Patches
Plain Diff
ENH: overset movement for patches (issue
#2
)
BUG: general patch movement was inconsistent
parent
dfd64a1b
Changes
3
Hide whitespace changes
Inline
Side-by-side
Showing
3 changed files
with
167 additions
and
84 deletions
+167
-84
src/catalyst/volMesh/foamVtkFvMeshAdaptor.H
src/catalyst/volMesh/foamVtkFvMeshAdaptor.H
+6
-0
src/catalyst/volMesh/foamVtkFvMeshAdaptorFieldTemplates.C
src/catalyst/volMesh/foamVtkFvMeshAdaptorFieldTemplates.C
+35
-44
src/catalyst/volMesh/foamVtkFvMeshAdaptorGeom.C
src/catalyst/volMesh/foamVtkFvMeshAdaptorGeom.C
+126
-40
No files found.
src/catalyst/volMesh/foamVtkFvMeshAdaptor.H
View file @
9b969cd4
...
@@ -184,6 +184,12 @@ private:
...
@@ -184,6 +184,12 @@ private:
// There will be several for groups, but only one for regular patches.
// There will be several for groups, but only one for regular patches.
void
convertGeometryPatches
();
void
convertGeometryPatches
();
//- Ghost cells - internal mesh
void
applyGhostingInternal
(
const
labelUList
&
types
);
//- Ghost cells - patches
void
applyGhostingPatches
(
const
labelUList
&
types
);
//- Ghost cells to affect the visibility of the geometry.
//- Ghost cells to affect the visibility of the geometry.
// Currently only used for overset.
// Currently only used for overset.
void
applyGhosting
();
void
applyGhosting
();
...
...
src/catalyst/volMesh/foamVtkFvMeshAdaptorFieldTemplates.C
View file @
9b969cd4
...
@@ -83,18 +83,13 @@ void Foam::vtk::fvMeshAdaptor::convertVolField
...
@@ -83,18 +83,13 @@ void Foam::vtk::fvMeshAdaptor::convertVolField
if
(
!
iter
.
found
()
||
!
iter
.
object
().
dataset
)
if
(
!
iter
.
found
()
||
!
iter
.
object
().
dataset
)
{
{
// Should not happen, but for safety require a vtk geometry
// Should not happen, but for safety require a vtk geometry
Pout
<<
"Cache miss for VTP patch "
<<
longName
<<
endl
;
continue
;
continue
;
}
}
foamVtpData
&
vtpData
=
iter
.
object
();
foamVtpData
&
vtpData
=
iter
.
object
();
auto
dataset
=
vtpData
.
dataset
;
auto
dataset
=
vtpData
.
dataset
;
const
labelList
&
patchIds
=
vtpData
.
additionalIds
();
if
(
patchIds
.
empty
())
{
continue
;
}
// This is slightly roundabout, but we deal with groups and with
// This is slightly roundabout, but we deal with groups and with
// single patches.
// single patches.
// For groups (spanning several patches) it is fairly messy to
// For groups (spanning several patches) it is fairly messy to
...
@@ -111,53 +106,48 @@ void Foam::vtk::fvMeshAdaptor::convertVolField
...
@@ -111,53 +106,48 @@ void Foam::vtk::fvMeshAdaptor::convertVolField
);
);
vtkSmartPointer
<
vtkFloatArray
>
pdata
;
vtkSmartPointer
<
vtkFloatArray
>
pdata
;
const
bool
allowPdata
=
(
interpField
&&
patchIds
.
size
()
==
1
);
label
coffset
=
0
;
// The write offset into cell-data
const
fvPatchField
<
Type
>&
ptf
=
fld
.
boundaryField
()[
patchId
];
for
(
label
patchId
:
patchIds
)
{
const
fvPatchField
<
Type
>&
ptf
=
fld
.
boundaryField
()[
patchId
];
if
if
(
isType
<
emptyFvPatchField
<
Type
>>
(
ptf
)
||
(
(
isType
<
emptyFvPatchField
<
Type
>>
(
ptf
)
extrapPatch
||
&&
!
polyPatch
::
constraintType
(
patches
[
patchId
].
type
())
(
extrapPatch
&&
!
polyPatch
::
constraintType
(
patches
[
patchId
].
type
())
)
)
)
{
)
fvPatch
p
(
ptf
.
patch
().
patch
(),
mesh
.
boundary
());
{
fvPatch
p
(
ptf
.
patch
().
patch
(),
mesh
.
boundary
());
tmp
<
Field
<
Type
>>
tpptf
tmp
<
Field
<
Type
>>
tpptf
(
(
fvPatchField
<
Type
>
(
p
,
fld
).
patchInternalField
()
fvPatchField
<
Type
>
(
p
,
fld
).
patchInternalField
()
);
);
coffset
+=
transcribeFloatData
(
cdata
,
tpptf
(),
coffset
);
transcribeFloatData
(
cdata
,
tpptf
()
);
if
(
allowPdata
&&
patchId
<
patchInterpList
.
size
())
if
(
interpField
&&
patchId
<
patchInterpList
.
size
())
{
{
pdata
=
vtk
::
Tools
::
convertFieldToVTK
pdata
=
vtk
::
Tools
::
convertFieldToVTK
(
(
fld
.
name
(),
fld
.
name
(),
patchInterpList
[
patchId
].
faceToPointInterpolate
(
tpptf
)()
patchInterpList
[
patchId
].
faceToPointInterpolate
(
tpptf
)()
);
);
}
}
}
else
}
else
{
transcribeFloatData
(
cdata
,
ptf
);
if
(
interpField
&&
patchId
<
patchInterpList
.
size
())
{
{
coffset
+=
transcribeFloatData
(
cdata
,
ptf
,
coffset
);
pdata
=
vtk
::
Tools
::
convertFieldToVTK
(
if
(
allowPdata
&&
patchId
<
patchInterpList
.
size
())
fld
.
name
(),
{
patchInterpList
[
patchId
].
faceToPointInterpolate
(
ptf
)()
pdata
=
vtk
::
Tools
::
convertFieldToVTK
);
(
fld
.
name
(),
patchInterpList
[
patchId
].
faceToPointInterpolate
(
ptf
)()
);
}
}
}
}
}
...
@@ -214,6 +204,7 @@ void Foam::vtk::fvMeshAdaptor::convertVolFieldInternal
...
@@ -214,6 +204,7 @@ void Foam::vtk::fvMeshAdaptor::convertVolFieldInternal
if
(
!
iter
.
found
()
||
!
iter
.
object
().
dataset
)
if
(
!
iter
.
found
()
||
!
iter
.
object
().
dataset
)
{
{
// Should not happen, but for safety require a vtk geometry
// Should not happen, but for safety require a vtk geometry
Pout
<<
"Cache miss for VTU "
<<
longName
<<
endl
;
return
;
return
;
}
}
foamVtuData
&
vtuData
=
iter
.
object
();
foamVtuData
&
vtuData
=
iter
.
object
();
...
...
src/catalyst/volMesh/foamVtkFvMeshAdaptorGeom.C
View file @
9b969cd4
...
@@ -56,7 +56,7 @@ void Foam::vtk::fvMeshAdaptor::convertGeometryInternal()
...
@@ -56,7 +56,7 @@ void Foam::vtk::fvMeshAdaptor::convertGeometryInternal()
{
{
Info
<<
"reuse "
<<
longName
<<
nl
;
Info
<<
"reuse "
<<
longName
<<
nl
;
}
}
vtuData
.
reuse
();
//
reuse
vtuData
.
reuse
();
// No movement - simply
reuse
return
;
return
;
}
}
else
if
(
meshState_
==
polyMesh
::
POINTS_MOVED
)
else
if
(
meshState_
==
polyMesh
::
POINTS_MOVED
)
...
@@ -87,6 +87,8 @@ void Foam::vtk::fvMeshAdaptor::convertGeometryInternal()
...
@@ -87,6 +87,8 @@ void Foam::vtk::fvMeshAdaptor::convertGeometryInternal()
void
Foam
::
vtk
::
fvMeshAdaptor
::
convertGeometryPatches
()
void
Foam
::
vtk
::
fvMeshAdaptor
::
convertGeometryPatches
()
{
{
// PATCHES
const
polyBoundaryMesh
&
patches
=
mesh_
.
boundaryMesh
();
const
polyBoundaryMesh
&
patches
=
mesh_
.
boundaryMesh
();
const
label
npatches
=
this
->
nPatches
();
const
label
npatches
=
this
->
nPatches
();
...
@@ -102,28 +104,22 @@ void Foam::vtk::fvMeshAdaptor::convertGeometryPatches()
...
@@ -102,28 +104,22 @@ void Foam::vtk::fvMeshAdaptor::convertGeometryPatches()
{
{
if
(
meshState_
==
polyMesh
::
UNCHANGED
)
if
(
meshState_
==
polyMesh
::
UNCHANGED
)
{
{
// Without movement is easy.
if
(
debug
)
if
(
debug
)
{
{
Info
<<
"reuse "
<<
longName
<<
nl
;
Info
<<
"reuse "
<<
longName
<<
nl
;
}
}
vtpData
.
reuse
();
vtpData
.
reuse
();
// No movement - simply reuse
continue
;
continue
;
}
}
else
if
(
meshState_
==
polyMesh
::
POINTS_MOVED
)
else
if
(
meshState_
==
polyMesh
::
POINTS_MOVED
)
{
{
// Point movement on single patch is OK
// Point movement on single patch
if
(
debug
)
const
labelList
&
patchIds
=
vtpData
.
additionalIds
();
if
(
patchIds
.
size
()
==
1
)
{
{
vtkgeom
=
vtpData
.
getCopy
();
Info
<<
"move points "
<<
longName
<<
nl
;
vtkgeom
->
SetPoints
(
vtk
::
Tools
::
Patch
::
points
(
patches
[
patchIds
[
0
]])
);
continue
;
}
}
vtkgeom
=
vtpData
.
getCopy
();
vtkgeom
->
SetPoints
(
vtk
::
Tools
::
Patch
::
points
(
pp
));
}
}
}
}
...
@@ -135,13 +131,13 @@ void Foam::vtk::fvMeshAdaptor::convertGeometryPatches()
...
@@ -135,13 +131,13 @@ void Foam::vtk::fvMeshAdaptor::convertGeometryPatches()
<<
longName
<<
endl
;
<<
longName
<<
endl
;
}
}
//
Store good patch id as additionalIds
//
Unused information
vtpData
.
additionalIds
()
=
{
patchId
}
;
vtpData
.
additionalIds
()
.
clear
()
;
// This is somewhat inconsistent, since we currently only have
// This is somewhat inconsistent, since we currently only have
// normal (non-grouped) patches but this may change in the future.
// normal (non-grouped) patches but this may change in the future.
vtkgeom
=
vtk
::
Tools
::
Patch
::
mesh
(
p
atches
[
patchId
]
);
vtkgeom
=
vtk
::
Tools
::
Patch
::
mesh
(
p
p
);
if
(
vtkgeom
)
if
(
vtkgeom
)
{
{
...
@@ -156,22 +152,11 @@ void Foam::vtk::fvMeshAdaptor::convertGeometryPatches()
...
@@ -156,22 +152,11 @@ void Foam::vtk::fvMeshAdaptor::convertGeometryPatches()
}
}
// These need to be rebuild here, since the component mesh may just have point
void
Foam
::
vtk
::
fvMeshAdaptor
::
applyGhostingInternal
(
const
labelUList
&
types
)
// motion without topology changes.
void
Foam
::
vtk
::
fvMeshAdaptor
::
applyGhosting
()
{
{
if
(
!
usingVolume
())
// MESH = "internal"
{
return
;
}
const
auto
*
stencilPtr
=
mesh_
.
lookupObjectPtr
<
cellCellStencilObject
>
(
cellCellStencilObject
::
typeName
);
if
(
!
stencilPtr
)
if
(
types
.
empty
()
||
!
usingVolume
()
)
{
{
return
;
return
;
}
}
...
@@ -182,12 +167,12 @@ void Foam::vtk::fvMeshAdaptor::applyGhosting()
...
@@ -182,12 +167,12 @@ void Foam::vtk::fvMeshAdaptor::applyGhosting()
if
(
!
iter
.
found
()
||
!
iter
.
object
().
dataset
)
if
(
!
iter
.
found
()
||
!
iter
.
object
().
dataset
)
{
{
// Should not happen, but for safety require a vtk geometry
// Should not happen, but for safety require a vtk geometry
Pout
<<
"Cache miss for VTU "
<<
longName
<<
endl
;
return
;
return
;
}
}
foamVtuData
&
vtuData
=
iter
.
object
();
foamVtuData
&
vtuData
=
iter
.
object
();
auto
dataset
=
vtuData
.
dataset
;
auto
dataset
=
vtuData
.
dataset
;
const
auto
&
stencil
=
*
stencilPtr
;
const
labelUList
&
cellMap
=
vtuData
.
cellMap
();
const
labelUList
&
cellMap
=
vtuData
.
cellMap
();
auto
vtkgcell
=
dataset
->
GetCellGhostArray
();
auto
vtkgcell
=
dataset
->
GetCellGhostArray
();
...
@@ -202,38 +187,139 @@ void Foam::vtk::fvMeshAdaptor::applyGhosting()
...
@@ -202,38 +187,139 @@ void Foam::vtk::fvMeshAdaptor::applyGhosting()
// vtkgpoint = dataset->AllocatePointGhostArray();
// vtkgpoint = dataset->AllocatePointGhostArray();
// }
// }
UList
<
uint8_t
>
gcell
=
UList
<
uint8_t
>
gcell
=
vtk
::
Tools
::
asUList
(
vtkgcell
,
cellMap
.
size
());
vtk
::
Tools
::
asUList
(
vtkgcell
,
cellMap
.
size
());
vtkgcell
->
FillValue
(
0
);
// Initialize to zero
vtkgcell
->
FillValue
(
0
);
// Initialize to zero
const
labelUList
&
types
=
stencil
.
cellTypes
();
forAll
(
cellMap
,
i
)
forAll
(
cellMap
,
i
)
{
{
const
label
cellType
=
types
[
cellMap
[
i
]];
const
label
cellType
=
types
[
cellMap
[
i
]];
if
(
cellType
==
cellCellStencil
::
INTERPOLATED
)
if
(
cellType
==
cellCellStencil
::
HOLE
)
{
{
// Need duplicate (not just HIDDENCELL) for it to be properly
// recognized
gcell
[
i
]
|=
gcell
[
i
]
|=
(
(
vtkDataSetAttributes
::
DUPLICATECELL
vtkDataSetAttributes
::
DUPLICATECELL
|
vtkDataSetAttributes
::
HIDDENCELL
);
);
}
}
else
if
(
cellType
==
cellCellStencil
::
HOLE
)
#if 0
// No special treatment for INTERPOLATED.
// This causes duplicate/overlapping values, but avoids holes
// in the results
else if (cellType == cellCellStencil::INTERPOLATED)
{
{
// Need duplicate (not just HIDDENCELL) for it to be properly
// recognized
gcell[i] |=
gcell[i] |=
(
(
vtkDataSetAttributes::DUPLICATECELL
vtkDataSetAttributes::DUPLICATECELL
|
vtkDataSetAttributes
::
HIDDENCELL
);
);
}
}
#endif
}
}
dataset
->
GetCellData
()
->
AddArray
(
vtkgcell
);
dataset
->
GetCellData
()
->
AddArray
(
vtkgcell
);
}
}
void
Foam
::
vtk
::
fvMeshAdaptor
::
applyGhostingPatches
(
const
labelUList
&
types
)
{
// PATCHES
if
(
types
.
empty
())
{
return
;
}
const
polyBoundaryMesh
&
patches
=
mesh_
.
boundaryMesh
();
const
label
npatches
=
this
->
nPatches
();
for
(
label
patchId
=
0
;
patchId
<
npatches
;
++
patchId
)
{
const
polyPatch
&
pp
=
patches
[
patchId
];
const
word
&
longName
=
pp
.
name
();
auto
iter
=
cachedVtp_
.
find
(
longName
);
if
(
!
iter
.
found
()
||
!
iter
.
object
().
dataset
)
{
// Should not happen, but for safety require a vtk geometry
Pout
<<
"Cache miss for VTP patch "
<<
longName
<<
endl
;
continue
;
}
foamVtpData
&
vtpData
=
iter
.
object
();
auto
dataset
=
vtpData
.
dataset
;
auto
vtkgcell
=
dataset
->
GetCellGhostArray
();
if
(
!
vtkgcell
)
{
vtkgcell
=
dataset
->
AllocateCellGhostArray
();
}
// Determine face ghosting based on interior cells
const
labelUList
&
bCells
=
pp
.
faceCells
();
const
label
len
=
bCells
.
size
();
UList
<
uint8_t
>
gcell
=
vtk
::
Tools
::
asUList
(
vtkgcell
,
len
);
vtkgcell
->
FillValue
(
0
);
// Initialize to zero
for
(
label
i
=
0
;
i
<
len
;
++
i
)
{
const
label
celli
=
bCells
[
i
];
const
label
cellType
=
types
[
celli
];
if
(
cellType
==
cellCellStencil
::
HOLE
)
{
// Need duplicate (not just HIDDENCELL) for it to be properly
// recognized
gcell
[
i
]
|=
(
vtkDataSetAttributes
::
DUPLICATECELL
|
vtkDataSetAttributes
::
HIDDENCELL
);
}
#if 0
// No special treatment for INTERPOLATED.
// This causes duplicate/overlapping values, but avoids holes
// in the results
else if (cellType == cellCellStencil::INTERPOLATED)
{
gcell[i] |=
(
vtkDataSetAttributes::DUPLICATECELL
);
}
#endif
}
dataset
->
GetCellData
()
->
AddArray
(
vtkgcell
);
}
}
// These need to be rebuild here, since the component mesh may just have point
// motion without topology changes.
void
Foam
::
vtk
::
fvMeshAdaptor
::
applyGhosting
()
{
const
auto
*
stencilPtr
=
mesh_
.
lookupObjectPtr
<
cellCellStencilObject
>
(
cellCellStencilObject
::
typeName
);
if
(
stencilPtr
)
{
const
labelUList
&
types
=
stencilPtr
->
cellTypes
();
applyGhostingInternal
(
types
);
applyGhostingPatches
(
types
);
}
}
// ************************************************************************* //
// ************************************************************************* //
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