Skip to content
GitLab
Menu
Projects
Groups
Snippets
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
Menu
Open sidebar
Michael Migliore
VTK
Commits
e31eec14
Commit
e31eec14
authored
Mar 14, 2017
by
Sujin Philip
Browse files
Add vtkmCleanGrid Filter
parent
5459e6ea
Changes
7
Hide whitespace changes
Inline
Side-by-side
Accelerators/Vtkm/CMakeLists.txt
View file @
e31eec14
...
...
@@ -33,6 +33,7 @@ set(lib_srcs
#needed to properly setup language wrappers
set
(
headers
vtkmAverageToPoints.h
vtkmCleanGrid.h
vtkmClip.h
vtkmContour.h
vtkmExternalFaces.h
...
...
@@ -45,6 +46,7 @@ set(headers
#implementation of the algorithms for cpu accelerators
set
(
cpu_accelerator_srcs
vtkmAverageToPoints.cxx
vtkmCleanGrid.cxx
vtkmClip.cxx
vtkmContour.cxx
vtkmExternalFaces.cxx
...
...
@@ -61,6 +63,7 @@ set(cpu_accelerator_srcs
#implementation of the algorithms for gpu accelerators
set
(
cuda_accelerator_srcs
vtkmAverageToPoints.cu
vtkmCleanGrid.cu
vtkmClip.cu
vtkmContour.cu
vtkmExternalFaces.cu
...
...
Accelerators/Vtkm/Testing/Cxx/CMakeLists.txt
View file @
e31eec14
...
...
@@ -3,6 +3,7 @@ find_package(VTKm REQUIRED COMPONENTS Base)
include_directories
(
${
VTKm_INCLUDE_DIRS
}
)
vtk_add_test_cxx
(
${
vtk-module
}
CxxTests tests
TestVTKMCleanGrid.cxx
TestVTKMClip.cxx
# TestVTKMGradientAndVorticity.cxx,NO_VALID
TestVTKMExternalFaces.cxx
...
...
Accelerators/Vtkm/Testing/Cxx/TestVTKMCleanGrid.cxx
0 → 100644
View file @
e31eec14
//=============================================================================
//
// Copyright (c) Kitware, Inc.
// All rights reserved.
// See LICENSE.txt for details.
//
// This software is distributed WITHOUT ANY WARRANTY; without even
// the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
// PURPOSE. See the above copyright notice for more information.
//
// Copyright 2012 Sandia Corporation.
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
//
//=============================================================================
#include
"vtkmCleanGrid.h"
#include
"vtkActor.h"
#include
"vtkGeometryFilter.h"
#include
"vtkNew.h"
#include
"vtkPolyDataMapper.h"
#include
"vtkRegressionTestImage.h"
#include
"vtkRenderer.h"
#include
"vtkRenderWindow.h"
#include
"vtkRenderWindowInteractor.h"
#include
"vtkRTAnalyticSource.h"
int
TestVTKMCleanGrid
(
int
argc
,
char
*
argv
[])
{
vtkNew
<
vtkRTAnalyticSource
>
wavelet
;
wavelet
->
SetWholeExtent
(
-
10
,
10
,
-
10
,
10
,
-
10
,
10
);
wavelet
->
SetCenter
(
0
,
0
,
0
);
vtkNew
<
vtkmCleanGrid
>
cleanGrid
;
cleanGrid
->
SetInputConnection
(
wavelet
->
GetOutputPort
());
vtkNew
<
vtkGeometryFilter
>
geometry
;
geometry
->
SetInputConnection
(
cleanGrid
->
GetOutputPort
());
vtkNew
<
vtkPolyDataMapper
>
mapper
;
mapper
->
SetInputConnection
(
geometry
->
GetOutputPort
());
mapper
->
SetScalarRange
(
37
,
277
);
vtkNew
<
vtkActor
>
actor
;
actor
->
SetMapper
(
mapper
.
GetPointer
());
vtkNew
<
vtkRenderer
>
renderer
;
renderer
->
AddActor
(
actor
.
GetPointer
());
renderer
->
ResetCamera
();
vtkNew
<
vtkRenderWindow
>
renWin
;
renWin
->
AddRenderer
(
renderer
.
GetPointer
());
vtkNew
<
vtkRenderWindowInteractor
>
iren
;
iren
->
SetRenderWindow
(
renWin
.
GetPointer
());
iren
->
Initialize
();
renWin
->
Render
();
int
retVal
=
vtkRegressionTestImage
(
renWin
.
GetPointer
());
if
(
retVal
==
vtkRegressionTester
::
DO_INTERACTOR
)
{
iren
->
Start
();
}
return
!
retVal
;
}
Accelerators/Vtkm/Testing/Data/Baseline/TestVTKMCleanGrid.png.md5
0 → 100644
View file @
e31eec14
604cd6c98e8c0400348770bf989dd79e
Accelerators/Vtkm/vtkmCleanGrid.cu
0 → 100644
View file @
e31eec14
//=============================================================================
//
// Copyright (c) Kitware, Inc.
// All rights reserved.
// See LICENSE.txt for details.
//
// This software is distributed WITHOUT ANY WARRANTY; without even
// the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
// PURPOSE. See the above copyright notice for more information.
//
// Copyright 2012 Sandia Corporation.
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
//
//=============================================================================
#include
"vtkmCleanGrid.cxx"
Accelerators/Vtkm/vtkmCleanGrid.cxx
0 → 100644
View file @
e31eec14
//=============================================================================
//
// Copyright (c) Kitware, Inc.
// All rights reserved.
// See LICENSE.txt for details.
//
// This software is distributed WITHOUT ANY WARRANTY; without even
// the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
// PURPOSE. See the above copyright notice for more information.
//
// Copyright 2012 Sandia Corporation.
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
//
//=============================================================================
#include
"vtkmCleanGrid.h"
#include
"vtkCellData.h"
#include
"vtkDemandDrivenPipeline.h"
#include
"vtkInformation.h"
#include
"vtkInformationVector.h"
#include
"vtkObjectFactory.h"
#include
"vtkPointData.h"
#include
"vtkUnstructuredGrid.h"
#include
"vtkmlib/ArrayConverters.h"
#include
"vtkmlib/DataSetConverters.h"
#include
"vtkmlib/Storage.h"
#include
"vtkmlib/UnstructuredGridConverter.h"
#include
"vtkmCellSetExplicit.h"
#include
"vtkmCellSetSingleType.h"
#include
"vtkmFilterPolicy.h"
#include
<vtkm/filter/CleanGrid.h>
vtkStandardNewMacro
(
vtkmCleanGrid
)
//------------------------------------------------------------------------------
vtkmCleanGrid
::
vtkmCleanGrid
()
:
CompactPoints
(
false
)
{
}
//------------------------------------------------------------------------------
vtkmCleanGrid
::~
vtkmCleanGrid
()
{
}
//------------------------------------------------------------------------------
void
vtkmCleanGrid
::
PrintSelf
(
ostream
&
os
,
vtkIndent
indent
)
{
this
->
Superclass
::
PrintSelf
(
os
,
indent
);
os
<<
indent
<<
"CompactPoints: "
<<
(
this
->
CompactPoints
?
"On"
:
"Off"
)
<<
"
\n
"
;
}
//------------------------------------------------------------------------------
int
vtkmCleanGrid
::
FillInputPortInformation
(
int
,
vtkInformation
*
info
)
{
info
->
Set
(
vtkAlgorithm
::
INPUT_REQUIRED_DATA_TYPE
(),
"vtkDataSet"
);
return
1
;
}
//------------------------------------------------------------------------------
int
vtkmCleanGrid
::
RequestData
(
vtkInformation
*
request
,
vtkInformationVector
**
inputVector
,
vtkInformationVector
*
outputVector
)
{
vtkInformation
*
inInfo
=
inputVector
[
0
]
->
GetInformationObject
(
0
);
vtkInformation
*
outInfo
=
outputVector
->
GetInformationObject
(
0
);
vtkDataSet
*
input
=
vtkDataSet
::
SafeDownCast
(
inInfo
->
Get
(
vtkDataObject
::
DATA_OBJECT
()));
vtkUnstructuredGrid
*
output
=
vtkUnstructuredGrid
::
SafeDownCast
(
outInfo
->
Get
(
vtkDataObject
::
DATA_OBJECT
()));
// convert the input dataset to a vtkm::cont::DataSet
vtkm
::
cont
::
DataSet
in
=
tovtkm
::
Convert
(
input
);
if
(
in
.
GetNumberOfCoordinateSystems
()
<=
0
||
in
.
GetNumberOfCellSets
()
<=
0
)
{
vtkErrorMacro
(
<<
"Could not convert vtk dataset to vtkm dataset"
);
return
0
;
}
// apply the filter
vtkmInputFilterPolicy
policy
;
vtkm
::
filter
::
CleanGrid
filter
;
filter
.
SetCompactPointFields
(
this
->
CompactPoints
);
vtkm
::
filter
::
ResultDataSet
result
=
filter
.
Execute
(
in
,
policy
);
if
(
!
result
.
IsValid
())
{
vtkErrorMacro
(
<<
"VTKm CleanGrid algorithm failed to run"
);
return
0
;
}
// map each point array
vtkPointData
*
pd
=
input
->
GetPointData
();
for
(
vtkIdType
i
=
0
;
i
<
pd
->
GetNumberOfArrays
();
i
++
)
{
vtkDataArray
*
array
=
pd
->
GetArray
(
i
);
if
(
array
==
NULL
)
{
continue
;
}
vtkm
::
cont
::
Field
pfield
=
tovtkm
::
Convert
(
array
,
vtkDataObject
::
FIELD_ASSOCIATION_POINTS
);
try
{
filter
.
MapFieldOntoOutput
(
result
,
pfield
,
policy
);
}
catch
(
vtkm
::
cont
::
Error
&
)
{
vtkWarningMacro
(
<<
"Unable to use VTKm to convert field ("
<<
array
->
GetName
()
<<
") to the output"
);
}
}
// convert back to vtkDataSet (vtkUnstructuredGrid)
if
(
!
fromvtkm
::
Convert
(
result
.
GetDataSet
(),
output
,
input
))
{
vtkErrorMacro
(
<<
"Unable to convert VTKm DataSet back to VTK"
);
return
0
;
}
// pass cell data
output
->
GetCellData
()
->
PassData
(
input
->
GetCellData
());
return
1
;
}
Accelerators/Vtkm/vtkmCleanGrid.h
0 → 100644
View file @
e31eec14
//=============================================================================
//
// Copyright (c) Kitware, Inc.
// All rights reserved.
// See LICENSE.txt for details.
//
// This software is distributed WITHOUT ANY WARRANTY; without even
// the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
// PURPOSE. See the above copyright notice for more information.
//
// Copyright 2012 Sandia Corporation.
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
//
//=============================================================================
#ifndef vtkmCleanGrid_h
#define vtkmCleanGrid_h
#include
"vtkUnstructuredGridAlgorithm.h"
#include
"vtkAcceleratorsVTKmModule.h"
//required for correct implementation
class
vtkDataSet
;
class
vtkUnstructuredGrid
;
class
VTKACCELERATORSVTKM_EXPORT
vtkmCleanGrid
:
public
vtkUnstructuredGridAlgorithm
{
public:
vtkTypeMacro
(
vtkmCleanGrid
,
vtkUnstructuredGridAlgorithm
)
void
PrintSelf
(
ostream
&
os
,
vtkIndent
indent
)
VTK_OVERRIDE
;
static
vtkmCleanGrid
*
New
();
//@{
/**
* Get/Set if the points from the input that are unused in the output should
* be removed. This will take extra time but the result dataset may use
* less memory. Off by default.
*/
vtkSetMacro
(
CompactPoints
,
bool
);
vtkGetMacro
(
CompactPoints
,
bool
);
vtkBooleanMacro
(
CompactPoints
,
bool
);
//@}
protected:
vtkmCleanGrid
();
~
vtkmCleanGrid
();
int
FillInputPortInformation
(
int
,
vtkInformation
*
)
VTK_OVERRIDE
;
int
RequestData
(
vtkInformation
*
,
vtkInformationVector
**
,
vtkInformationVector
*
)
VTK_OVERRIDE
;
bool
CompactPoints
;
private:
vtkmCleanGrid
(
const
vtkmCleanGrid
&
)
VTK_DELETE_FUNCTION
;
void
operator
=
(
const
vtkmCleanGrid
&
)
VTK_DELETE_FUNCTION
;
};
#endif // vtkmCleanGrid_h
// VTK-HeaderTest-Exclude: vtkmCleanGrid.h
Write
Preview
Supports
Markdown
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