Skip to content
GitLab
Menu
Projects
Groups
Snippets
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
Menu
Open sidebar
Community
catalyst
Compare Revisions
566968cae355bd15bdfcc133a04c886af604075f...478eeaedc267f3c9965f53bf22855847b3535f4a
Commits (1)
ENH: updated to use OpenFOAM vtk::vtuAdaptor structure
· 478eeaed
Mark Olesen
authored
Feb 11, 2019
478eeaed
Hide whitespace changes
Inline
Side-by-side
src/catalyst/CMakeLists-Project.txt
View file @
478eeaed
...
...
@@ -68,7 +68,6 @@ file(GLOB SOURCE_FILES
volMesh/catalystFvMesh.C
volMesh/foamVtkFvMeshAdaptor.C
volMesh/foamVtkFvMeshAdaptorGeom.C
volMesh/foamVtkFvMeshAdaptorGeomVtu.C
volMesh/foamVtkFvMeshAdaptorFields.C
)
...
...
src/catalyst/volMesh/foamVtkFvMeshAdaptor.H
View file @
478eeaed
...
...
@@ -58,14 +58,7 @@ SourceFiles
#include "PrimitivePatchInterpolation.H"
#include "volPointInterpolation.H"
#include "foamVtkTools.H"
#include "foamVtkMeshMaps.H"
#include <vtkSmartPointer.h>
#include <vtkPoints.h>
#include <vtkPolyData.h>
#include <vtkUnstructuredGrid.h>
#include <vtkMultiBlockDataSet.h>
#include "foamVtkVtuAdaptor.H"
// * * * * * * * * * * * * * Forward Declarations * * * * * * * * * * * * * //
...
...
@@ -136,30 +129,8 @@ private:
//- Bookkeeping for vtkUnstructuredGrid
struct
foamVtuData
:
public
vtk
::
Caching
<
vtkUnstructuredGrid
>
,
public
foamVtkMeshMaps
{
//- The vtk points for the mesh (and decomposition)
vtkSmartPointer
<
vtkPoints
>
points
(
const
fvMesh
&
mesh
)
const
;
//- The vtk points for the mesh (and decomposition)
//- using the provided pointMap
vtkSmartPointer
<
vtkPoints
>
points
(
const
fvMesh
&
mesh
,
const
labelUList
&
pointMap
)
const
;
//- Internal mesh as vtkUnstructuredGrid
vtkSmartPointer
<
vtkUnstructuredGrid
>
internal
(
const
fvMesh
&
mesh
,
const
bool
decompPoly
);
};
public
vtk
::
vtuAdaptor
{};
// Private Data
...
...
@@ -225,15 +196,6 @@ private:
//- Convert specified volume fields
void
convertVolFields
(
const
wordRes
&
selectFields
);
//- Volume field
template
<
class
Type
>
vtkSmartPointer
<
vtkFloatArray
>
convertVolFieldToVTK
(
const
GeometricField
<
Type
,
fvPatchField
,
volMesh
>&
fld
,
const
foamVtuData
&
vtuData
)
const
;
//- Volume field - all types
template
<
class
Type
>
void
convertVolField
...
...
src/catalyst/volMesh/foamVtkFvMeshAdaptorFieldTemplates.C
View file @
478eeaed
...
...
@@ -216,30 +216,17 @@ void Foam::vtk::fvMeshAdaptor::convertVolFieldInternal
foamVtuData
&
vtuData
=
iter
.
object
();
auto
dataset
=
vtuData
.
dataset
;
vtkSmartPointer
<
vtkFloatArray
>
cdata
=
convertVolFieldToVTK
dataset
->
GetCellData
()
->
AddArray
(
fld
,
vtuData
vtuData
.
convertField
(
fld
)
);
dataset
->
GetCellData
()
->
AddArray
(
cdata
);
// double cbound[2];
// cdata->GetRange(cbound);
// Info<< "range cdata: " << cbound[0] << " - " << cbound[1] << nl;
if
(
ptfPtr
.
valid
())
{
vtkSmartPointer
<
vtkFloatArray
>
pdata
=
convertPointField
dataset
->
GetPointData
()
->
AddArray
(
ptfPtr
(),
fld
,
vtuData
convertPointField
(
ptfPtr
(),
fld
,
vtuData
)
);
dataset
->
GetPointData
()
->
AddArray
(
pdata
);
// double bound[2];
// pdata->GetRange(bound);
// Info<< "range pdata: " << bound[0] << " - " << bound[1] << nl;
}
}
...
...
@@ -331,49 +318,6 @@ vtkSmartPointer<vtkFloatArray> Foam::vtk::fvMeshAdaptor::convertPointField
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
//
// low-level conversions
//
template
<
class
Type
>
vtkSmartPointer
<
vtkFloatArray
>
Foam
::
vtk
::
fvMeshAdaptor
::
convertVolFieldToVTK
(
const
GeometricField
<
Type
,
fvPatchField
,
volMesh
>&
fld
,
const
foamVtuData
&
vtuData
)
const
{
const
int
nComp
(
pTraits
<
Type
>::
nComponents
);
const
labelUList
&
cellMap
=
vtuData
.
cellMap
();
auto
data
=
vtkSmartPointer
<
vtkFloatArray
>::
New
();
data
->
SetName
(
fld
.
name
().
c_str
());
data
->
SetNumberOfComponents
(
nComp
);
data
->
SetNumberOfTuples
(
cellMap
.
size
());
DebugInfo
<<
"convert volField: "
<<
fld
.
name
()
<<
" size="
<<
cellMap
.
size
()
<<
" ("
<<
fld
.
size
()
<<
" + "
<<
(
cellMap
.
size
()
-
fld
.
size
())
<<
") nComp="
<<
nComp
<<
endl
;
float
scratch
[
pTraits
<
Type
>::
nComponents
];
vtkIdType
celli
=
0
;
for
(
const
label
meshCelli
:
cellMap
)
{
vtk
::
Tools
::
foamToVtkTuple
(
scratch
,
fld
[
meshCelli
]);
data
->
SetTuple
(
celli
++
,
scratch
);
}
return
data
;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
src/catalyst/volMesh/foamVtkFvMeshAdaptorGeom.C
View file @
478eeaed
...
...
@@ -27,10 +27,10 @@ License
#include "cellCellStencilObject.H"
// VTK includes
#include
<
vtkMultiBlockDataSet.h
>
#include
<
vtkPolyData.h
>
#include
<
vtkUnstructuredGrid.h
>
#include
<
vtkSmartPointer.h
>
#include
"
vtkMultiBlockDataSet.h
"
#include
"
vtkPolyData.h
"
#include
"
vtkUnstructuredGrid.h
"
#include
"
vtkSmartPointer.h
"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
...
...
src/catalyst/volMesh/foamVtkFvMeshAdaptorGeomVtu.C
deleted
100644 → 0
View file @
566968ca
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016-2018 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "foamVtkFvMeshAdaptor.H"
// OpenFOAM includes
#include "fvMesh.H"
#include "foamVtkTools.H"
#include "foamVtuSizing.H"
// VTK includes
#include <vtkUnstructuredGrid.h>
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
vtkSmartPointer
<
vtkPoints
>
Foam
::
vtk
::
fvMeshAdaptor
::
foamVtuData
::
points
(
const
fvMesh
&
mesh
)
const
{
// Convert OpenFOAM mesh vertices to VTK
auto
vtkpoints
=
vtkSmartPointer
<
vtkPoints
>::
New
();
// Normal points
const
pointField
&
pts
=
mesh
.
points
();
// Additional cell centres
const
labelUList
&
addPoints
=
this
->
additionalIds
();
vtkpoints
->
SetNumberOfPoints
(
pts
.
size
()
+
addPoints
.
size
());
// Normal points
vtkIdType
pointId
=
0
;
for
(
const
point
&
p
:
pts
)
{
vtkpoints
->
SetPoint
(
pointId
++
,
p
.
v_
);
}
// Cell centres
for
(
const
label
meshCelli
:
addPoints
)
{
vtkpoints
->
SetPoint
(
pointId
++
,
mesh
.
cellCentres
()[
meshCelli
].
v_
);
}
return
vtkpoints
;
}
vtkSmartPointer
<
vtkPoints
>
Foam
::
vtk
::
fvMeshAdaptor
::
foamVtuData
::
points
(
const
fvMesh
&
mesh
,
const
labelUList
&
pointMap
)
const
{
// Convert OpenFOAM mesh vertices to VTK
auto
vtkpoints
=
vtkSmartPointer
<
vtkPoints
>::
New
();
// Normal points
const
pointField
&
pts
=
mesh
.
points
();
// Additional cell centres
const
labelUList
&
addPoints
=
this
->
additionalIds
();
vtkpoints
->
SetNumberOfPoints
(
pointMap
.
size
()
+
addPoints
.
size
());
// Normal points
vtkIdType
pointId
=
0
;
for
(
const
label
meshPointi
:
pointMap
)
{
vtkpoints
->
SetPoint
(
pointId
++
,
pts
[
meshPointi
].
v_
);
}
// Cell centres
for
(
const
label
meshCelli
:
addPoints
)
{
vtkpoints
->
SetPoint
(
pointId
++
,
mesh
.
cellCentres
()[
meshCelli
].
v_
);
}
return
vtkpoints
;
}
vtkSmartPointer
<
vtkUnstructuredGrid
>
Foam
::
vtk
::
fvMeshAdaptor
::
foamVtuData
::
internal
(
const
fvMesh
&
mesh
,
const
bool
decompPoly
)
{
vtk
::
vtuSizing
sizing
(
mesh
,
decompPoly
);
auto
cellTypes
=
vtkSmartPointer
<
vtkUnsignedCharArray
>::
New
();
auto
cells
=
vtkSmartPointer
<
vtkCellArray
>::
New
();
auto
faces
=
vtkSmartPointer
<
vtkIdTypeArray
>::
New
();
auto
cellLocations
=
vtkSmartPointer
<
vtkIdTypeArray
>::
New
();
auto
faceLocations
=
vtkSmartPointer
<
vtkIdTypeArray
>::
New
();
UList
<
uint8_t
>
cellTypesUL
=
vtk
::
Tools
::
asUList
(
cellTypes
,
sizing
.
nFieldCells
());
UList
<
vtkIdType
>
cellsUL
=
vtk
::
Tools
::
asUList
(
cells
,
sizing
.
nFieldCells
(),
sizing
.
sizeInternal
(
vtk
::
vtuSizing
::
slotType
::
CELLS
)
);
UList
<
vtkIdType
>
cellLocationsUL
=
vtk
::
Tools
::
asUList
(
cellLocations
,
sizing
.
sizeInternal
(
vtk
::
vtuSizing
::
slotType
::
CELLS_OFFSETS
)
);
UList
<
vtkIdType
>
facesUL
=
vtk
::
Tools
::
asUList
(
faces
,
sizing
.
sizeInternal
(
vtk
::
vtuSizing
::
slotType
::
FACES
)
);
UList
<
vtkIdType
>
faceLocationsUL
=
vtk
::
Tools
::
asUList
(
faceLocations
,
sizing
.
sizeInternal
(
vtk
::
vtuSizing
::
slotType
::
FACES_OFFSETS
)
);
sizing
.
populateInternal
(
mesh
,
cellTypesUL
,
cellsUL
,
cellLocationsUL
,
facesUL
,
faceLocationsUL
,
static_cast
<
foamVtkMeshMaps
&>
(
*
this
)
);
auto
vtkmesh
=
vtkSmartPointer
<
vtkUnstructuredGrid
>::
New
();
// Convert OpenFOAM mesh vertices to VTK
// - can only do this *after* populating the decompInfo with cell-ids
// for any additional points (ie, mesh cell-centres)
vtkmesh
->
SetPoints
(
this
->
points
(
mesh
));
if
(
facesUL
.
size
())
{
vtkmesh
->
SetCells
(
cellTypes
,
cellLocations
,
cells
,
faceLocations
,
faces
);
}
else
{
vtkmesh
->
SetCells
(
cellTypes
,
cellLocations
,
cells
,
nullptr
,
nullptr
);
}
return
vtkmesh
;
}
// ************************************************************************* //