Skip to content
Projects
Groups
Snippets
Help
Loading...
Help
Support
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
VTK
Project overview
Project overview
Details
Activity
Releases
Cycle Analytics
Repository
Repository
Files
Commits
Branches
Tags
Contributors
Graph
Compare
Charts
Issues
0
Issues
0
List
Boards
Labels
Milestones
Merge Requests
1
Merge Requests
1
Packages
Packages
Container Registry
Wiki
Wiki
Snippets
Snippets
Members
Members
Collapse sidebar
Close sidebar
Activity
Graph
Charts
Create a new issue
Commits
Issue Boards
Open sidebar
Michael Migliore
VTK
Commits
a4891a53
Commit
a4891a53
authored
Jul 10, 2018
by
Haocheng LIU
Browse files
Options
Browse Files
Download
Email Patches
Plain Diff
Expose the vtkm Probe filter
parent
2010e287
Changes
6
Hide whitespace changes
Inline
Side-by-side
Showing
6 changed files
with
616 additions
and
0 deletions
+616
-0
Accelerators/Vtkm/CMakeLists.txt
Accelerators/Vtkm/CMakeLists.txt
+3
-0
Accelerators/Vtkm/Testing/Cxx/CMakeLists.txt
Accelerators/Vtkm/Testing/Cxx/CMakeLists.txt
+1
-0
Accelerators/Vtkm/Testing/Cxx/TestVTKMProbe.cxx
Accelerators/Vtkm/Testing/Cxx/TestVTKMProbe.cxx
+152
-0
Accelerators/Vtkm/vtkmProbe.cu
Accelerators/Vtkm/vtkmProbe.cu
+16
-0
Accelerators/Vtkm/vtkmProbe.cxx
Accelerators/Vtkm/vtkmProbe.cxx
+280
-0
Accelerators/Vtkm/vtkmProbe.h
Accelerators/Vtkm/vtkmProbe.h
+164
-0
No files found.
Accelerators/Vtkm/CMakeLists.txt
View file @
a4891a53
...
...
@@ -45,6 +45,7 @@ set(headers
vtkmGradient.h
vtkmPointElevation.h
vtkmPolyDataNormals.h
vtkmProbe.h
vtkmTriangleMeshPointNormals.h
)
...
...
@@ -66,6 +67,7 @@ set(cpu_accelerator_srcs
vtkmGradient.cxx
vtkmPointElevation.cxx
vtkmPolyDataNormals.cxx
vtkmProbe.cxx
vtkmTriangleMeshPointNormals.cxx
vtkmlib/Portals.cxx
vtkmlib/ImplicitFunctionConverter.cxx
...
...
@@ -89,6 +91,7 @@ set(cuda_accelerator_srcs
vtkmGradient.cu
vtkmPointElevation.cu
vtkmPolyDataNormals.cu
vtkmProbe.cu
vtkmTriangleMeshPointNormals.cu
vtkmlib/Portals.cu
vtkmlib/ImplicitFunctionConverter.cu
...
...
Accelerators/Vtkm/Testing/Cxx/CMakeLists.txt
View file @
a4891a53
...
...
@@ -11,6 +11,7 @@ vtk_add_test_cxx(vtkAcceleratorsVtkmCxxTests tests
TestVTKMMarchingCubes.cxx
TestVTKMMarchingCubes2.cxx
TestVTKMPointElevation.cxx
TestVTKMProbe.cxx,NO_VALID
TestVTKMPolyDataNormals.cxx
TestVTKMThreshold.cxx
TestVTKMThreshold2.cxx
...
...
Accelerators/Vtkm/Testing/Cxx/TestVTKMProbe.cxx
0 → 100644
View file @
a4891a53
//=============================================================================
//
// 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 "vtkmProbe.h"
#include "vtkFloatArray.h"
#include "vtkImageData.h"
#include "vtkNew.h"
#include "vtkDataSet.h"
#include "vtkPointData.h"
#include "vtkCellData.h"
#include "vtkDataArray.h"
namespace
{
static
int
inputDim
=
9
;
static
int
sourceDim
=
4
;
void
populatePointAndCellArray
(
vtkFloatArray
*
pointArray
,
vtkFloatArray
*
cellArray
)
{
pointArray
->
SetNumberOfValues
(
sourceDim
*
sourceDim
);
pointArray
->
SetName
(
"pointdata"
);
for
(
vtkIdType
i
=
0
;
i
<
static_cast
<
vtkIdType
>
(
sourceDim
*
sourceDim
);
i
++
)
{
pointArray
->
SetValue
(
i
,
0.3
f
*
i
);
}
cellArray
->
SetName
(
"celldata"
);
cellArray
->
SetNumberOfValues
((
sourceDim
-
1
)
*
(
sourceDim
-
1
));
for
(
vtkIdType
i
=
0
;
i
<
static_cast
<
vtkIdType
>
((
sourceDim
-
1
)
*
(
sourceDim
-
1
));
i
++
)
{
cellArray
->
SetValue
(
i
,
0.7
f
*
i
);
}
}
const
std
::
vector
<
float
>&
GetExpectedPointData
()
{
static
std
::
vector
<
float
>
expected
=
{
1.05
f
,
1.155
f
,
1.26
f
,
1.365
f
,
1.47
f
,
1.575
f
,
1.68
f
,
0.0
f
,
0.0
f
,
1.47
f
,
1.575
f
,
1.68
f
,
1.785
f
,
1.89
f
,
1.995
f
,
2.1
f
,
0.0
f
,
0.0
f
,
1.89
f
,
1.995
f
,
2.1
f
,
2.205
f
,
2.31
f
,
2.415
f
,
2.52
f
,
0.0
f
,
0.0
f
,
2.31
f
,
2.415
f
,
2.52
f
,
2.625
f
,
2.73
f
,
2.835
f
,
2.94
f
,
0.0
f
,
0.0
f
,
2.73
f
,
2.835
f
,
2.94
f
,
3.045
f
,
3.15
f
,
3.255
f
,
3.36
f
,
0.0
f
,
0.0
f
,
3.15
f
,
3.255
f
,
3.36
f
,
3.465
f
,
3.57
f
,
3.675
f
,
3.78
f
,
0.0
f
,
0.0
f
,
3.57
f
,
3.675
f
,
3.78
f
,
3.885
f
,
3.99
f
,
4.095
f
,
4.2
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
};
return
expected
;
}
const
std
::
vector
<
float
>&
GetExpectedCellData
()
{
static
std
::
vector
<
float
>
expected
=
{
0.0
f
,
0.7
f
,
0.7
f
,
0.7
f
,
1.4
f
,
1.4
f
,
1.4
f
,
0.0
f
,
0.0
f
,
2.1
f
,
2.8
f
,
2.8
f
,
2.8
f
,
3.5
f
,
3.5
f
,
3.5
f
,
0.0
f
,
0.0
f
,
2.1
f
,
2.8
f
,
2.8
f
,
2.8
f
,
3.5
f
,
3.5
f
,
3.5
f
,
0.0
f
,
0.0
f
,
2.1
f
,
2.8
f
,
2.8
f
,
2.8
f
,
3.5
f
,
3.5
f
,
3.5
f
,
0.0
f
,
0.0
f
,
4.2
f
,
4.9
f
,
4.9
f
,
4.9
f
,
5.6
f
,
5.6
f
,
5.6
f
,
0.0
f
,
0.0
f
,
4.2
f
,
4.9
f
,
4.9
f
,
4.9
f
,
5.6
f
,
5.6
f
,
5.6
f
,
0.0
f
,
0.0
f
,
4.2
f
,
4.9
f
,
4.9
f
,
4.9
f
,
5.6
f
,
5.6
f
,
5.6
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
,
0.0
f
};
return
expected
;
}
const
std
::
vector
<
size_t
>&
GetExpectedHiddenPoints
()
{
static
std
::
vector
<
size_t
>
expected
=
{
0
,
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
0
,
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
0
,
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
0
,
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
0
,
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
0
,
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
0
,
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
};
return
expected
;
}
const
std
::
vector
<
size_t
>&
GetExpectedHiddenCells
()
{
static
std
::
vector
<
size_t
>
expected
=
{
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
0
,
0
,
0
,
0
,
0
,
0
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
,
2
};
return
expected
;
}
template
<
typename
T
>
void
TestResultArray
(
vtkDataArray
*
result
,
const
std
::
vector
<
T
>&
expected
)
{
if
(
result
->
GetNumberOfValues
()
!=
static_cast
<
vtkIdType
>
(
expected
.
size
()))
{
std
::
cout
<<
"Array "
<<
result
->
GetName
()
<<
" has wrong size"
<<
std
::
endl
;
}
assert
(
result
->
GetNumberOfValues
()
==
static_cast
<
vtkIdType
>
(
expected
.
size
()));
for
(
vtkIdType
i
=
0
;
i
<
result
->
GetNumberOfValues
();
++
i
)
{
if
((
result
->
GetComponent
(
0
,
i
)
-
expected
[
static_cast
<
size_t
>
(
i
)])
>
1e-5
)
{
std
::
cout
<<
"Array "
<<
result
->
GetName
()
<<
" has wrong value"
<<
" at index "
<<
i
<<
". result value="
<<
result
->
GetComponent
(
0
,
i
)
<<
" expected value="
<<
expected
[
static_cast
<
size_t
>
(
i
)]
<<
std
::
endl
;
}
assert
((
result
->
GetComponent
(
0
,
i
)
-
expected
[
static_cast
<
size_t
>
(
i
)])
<
1e-5
);
}
}
}
// Anonymous namespace
int
TestVTKMProbe
(
int
,
char
*
[])
{
vtkNew
<
vtkImageData
>
input
;
input
->
SetOrigin
(
0.7
,
0.7
,
0.0
);
input
->
SetSpacing
(
0.35
,
0.35
,
1.0
);
input
->
SetExtent
(
0
,
inputDim
-
1
,
0
,
inputDim
-
1
,
0
,
0
);
vtkNew
<
vtkImageData
>
source
;
source
->
SetOrigin
(
0.0
,
0.0
,
0.0
);
source
->
SetSpacing
(
1.0
,
1.0
,
1.0
);
source
->
SetExtent
(
0
,
sourceDim
-
1
,
0
,
sourceDim
-
1
,
0
,
0
);
vtkNew
<
vtkFloatArray
>
pointArray
,
cellArray
;
populatePointAndCellArray
(
pointArray
,
cellArray
);
source
->
GetPointData
()
->
AddArray
(
pointArray
);
source
->
GetCellData
()
->
AddArray
(
cellArray
);
vtkNew
<
vtkmProbe
>
probe
;
probe
->
SetValidPointMaskArrayName
(
"validPoint"
);
probe
->
SetValidCellMaskArrayName
(
"validCell"
);
probe
->
SetInputData
(
input
);
probe
->
SetSourceData
(
source
);
probe
->
Update
();
vtkDataSet
*
result
=
probe
->
GetOutput
();
TestResultArray
(
result
->
GetPointData
()
->
GetArray
(
pointArray
->
GetName
()),
GetExpectedPointData
());
TestResultArray
(
result
->
GetCellData
()
->
GetArray
(
cellArray
->
GetName
()),
GetExpectedCellData
());
TestResultArray
(
result
->
GetPointData
()
->
GetArray
(
"validPoint"
),
GetExpectedHiddenPoints
());
TestResultArray
(
result
->
GetCellData
()
->
GetArray
(
"validCell"
),
GetExpectedHiddenCells
());
return
0
;
}
Accelerators/Vtkm/vtkmProbe.cu
0 → 100644
View file @
a4891a53
//=============================================================================
//
// 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 "vtkmProbe.cxx"
Accelerators/Vtkm/vtkmProbe.cxx
0 → 100644
View file @
a4891a53
//=============================================================================
//
// 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 "vtkmProbe.h"
#include "vtkCellData.h"
#include "vtkDataSet.h"
#include "vtkExecutive.h"
#include "vtkImageData.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkObjectFactory.h"
#include "vtkStreamingDemandDrivenPipeline.h"
#include "vtkPointData.h"
#include "vtkmlib/ArrayConverters.h"
#include "vtkmlib/DataSetConverters.h"
#include "vtkmlib/Storage.h"
#include "vtkmCellSetExplicit.h"
#include "vtkmCellSetSingleType.h"
#include "vtkmFilterPolicy.h"
#include "vtkm/filter/Probe.h"
vtkStandardNewMacro
(
vtkmProbe
)
//------------------------------------------------------------------------------
vtkmProbe
::
vtkmProbe
()
{
this
->
SetNumberOfInputPorts
(
2
);
this
->
PassCellArrays
=
false
;
this
->
PassPointArrays
=
false
;
this
->
PassFieldArrays
=
true
;
this
->
ValidPointMaskArrayName
=
"vtkValidPointMask"
;
this
->
ValidCellMaskArrayName
=
"vtkValidCellMask"
;
}
//------------------------------------------------------------------------------
void
vtkmProbe
::
SetSourceData
(
vtkDataObject
*
input
)
{
this
->
SetInputData
(
1
,
input
);
}
//------------------------------------------------------------------------------
vtkDataObject
*
vtkmProbe
::
GetSource
()
{
if
(
this
->
GetNumberOfInputConnections
(
1
)
<
1
)
{
return
nullptr
;
}
return
this
->
GetExecutive
()
->
GetInputData
(
1
,
0
);
}
//------------------------------------------------------------------------------
void
vtkmProbe
::
SetSourceConnection
(
vtkAlgorithmOutput
*
algOutput
)
{
this
->
SetInputConnection
(
1
,
algOutput
);
}
//------------------------------------------------------------------------------
int
vtkmProbe
::
RequestData
(
vtkInformation
*
vtkNotUsed
(
request
),
vtkInformationVector
**
inputVector
,
vtkInformationVector
*
outputVector
)
{
// Get the info objects
vtkInformation
*
inInfo
=
inputVector
[
0
]
->
GetInformationObject
(
0
);
vtkInformation
*
sourceInfo
=
inputVector
[
1
]
->
GetInformationObject
(
0
);
vtkInformation
*
outInfo
=
outputVector
->
GetInformationObject
(
0
);
// Get the input and output
vtkDataSet
*
input
=
vtkDataSet
::
SafeDownCast
(
inInfo
->
Get
(
vtkDataSet
::
DATA_OBJECT
()));
vtkDataSet
*
source
=
vtkDataSet
::
SafeDownCast
(
sourceInfo
->
Get
(
vtkDataSet
::
DATA_OBJECT
()));
vtkDataSet
*
output
=
vtkDataSet
::
SafeDownCast
(
outInfo
->
Get
(
vtkDataSet
::
DATA_OBJECT
()));
// Copy the input to the output as a starting point
output
->
CopyStructure
(
input
);
try
{
// Convert the input dataset to a vtkm::cont::DataSet
vtkm
::
cont
::
DataSet
in
=
tovtkm
::
Convert
(
input
);
// VTK-m's probe filter requires the source to have at least a cellSet.
vtkm
::
cont
::
DataSet
so
=
tovtkm
::
Convert
(
source
,
tovtkm
::
FieldsFlag
::
PointsAndCells
);
if
(
!
so
.
GetNumberOfCellSets
())
{
vtkErrorMacro
(
<<
"The source geometry does not have any cell set,"
"aborting vtkmProbe filter"
);
return
0
;
}
vtkmInputFilterPolicy
policy
;
vtkm
::
filter
::
Probe
probe
;
// The input in VTK is the geometry in VTKM and the source in VTK is the input
// in VTKM.
probe
.
SetGeometry
(
in
);
auto
result
=
probe
.
Execute
(
so
,
policy
);
for
(
vtkm
::
Id
i
=
0
;
i
<
result
.
GetNumberOfFields
();
i
++
)
{
const
vtkm
::
cont
::
Field
&
field
=
result
.
GetField
(
i
);
vtkDataArray
*
fieldArray
=
fromvtkm
::
Convert
(
field
);
if
(
field
.
GetAssociation
()
==
vtkm
::
cont
::
Field
::
Association
::
POINTS
)
{
if
(
strcmp
(
fieldArray
->
GetName
(),
"HIDDEN"
)
==
0
)
{
fieldArray
->
SetName
(
this
->
ValidPointMaskArrayName
.
c_str
());
}
output
->
GetPointData
()
->
AddArray
(
fieldArray
);
}
else
if
(
field
.
GetAssociation
()
==
vtkm
::
cont
::
Field
::
Association
::
CELL_SET
)
{
if
(
strcmp
(
fieldArray
->
GetName
(),
"HIDDEN"
)
==
0
)
{
fieldArray
->
SetName
(
this
->
ValidCellMaskArrayName
.
c_str
());
}
output
->
GetCellData
()
->
AddArray
(
fieldArray
);
}
fieldArray
->
FastDelete
();
}
}
catch
(
const
vtkm
::
cont
::
Error
&
e
)
{
vtkErrorMacro
(
<<
"VTK-m error: "
<<
e
.
GetMessage
());
return
0
;
}
this
->
PassAttributeData
(
input
,
source
,
output
);
return
1
;
}
//------------------------------------------------------------------------------
int
vtkmProbe
::
RequestInformation
(
vtkInformation
*
vtkNotUsed
(
request
),
vtkInformationVector
**
inputVector
,
vtkInformationVector
*
outputVector
)
{
// Update the whole extent in the output
vtkInformation
*
inInfo
=
inputVector
[
0
]
->
GetInformationObject
(
0
);
vtkInformation
*
sourceInfo
=
inputVector
[
1
]
->
GetInformationObject
(
0
);
vtkInformation
*
outInfo
=
outputVector
->
GetInformationObject
(
0
);
int
wholeExtent
[
6
];
if
(
inInfo
&&
outInfo
)
{
outInfo
->
CopyEntry
(
sourceInfo
,
vtkStreamingDemandDrivenPipeline
::
TIME_STEPS
());
outInfo
->
CopyEntry
(
sourceInfo
,
vtkStreamingDemandDrivenPipeline
::
TIME_RANGE
());
inInfo
->
Get
(
vtkStreamingDemandDrivenPipeline
::
WHOLE_EXTENT
(),
wholeExtent
);
outInfo
->
Set
(
vtkStreamingDemandDrivenPipeline
::
WHOLE_EXTENT
(),
wholeExtent
,
6
);
// Make sure that the scalar type and number of components
// are propagated from the source not the input.
if
(
vtkImageData
::
HasScalarType
(
sourceInfo
))
{
vtkImageData
::
SetScalarType
(
vtkImageData
::
GetScalarType
(
sourceInfo
),
outInfo
);
}
if
(
vtkImageData
::
HasNumberOfScalarComponents
(
sourceInfo
))
{
vtkImageData
::
SetNumberOfScalarComponents
(
vtkImageData
::
GetNumberOfScalarComponents
(
sourceInfo
),
outInfo
);
}
return
1
;
}
vtkErrorMacro
(
"Missing input or output info!"
);
return
0
;
}
//------------------------------------------------------------------------------
int
vtkmProbe
::
RequestUpdateExtent
(
vtkInformation
*
vtkNotUsed
(
request
),
vtkInformationVector
**
inputVector
,
vtkInformationVector
*
outputVector
)
{
vtkInformation
*
inInfo
=
inputVector
[
0
]
->
GetInformationObject
(
0
);
vtkInformation
*
sourceInfo
=
inputVector
[
1
]
->
GetInformationObject
(
0
);
vtkInformation
*
outInfo
=
outputVector
->
GetInformationObject
(
0
);
if
(
inInfo
&&
outInfo
)
{
// Source's update exetent should be independent of the resampling extent
inInfo
->
Set
(
vtkStreamingDemandDrivenPipeline
::
EXACT_EXTENT
(),
1
);
sourceInfo
->
Remove
(
vtkStreamingDemandDrivenPipeline
::
UPDATE_EXTENT
());
if
(
sourceInfo
->
Has
(
vtkStreamingDemandDrivenPipeline
::
WHOLE_EXTENT
()))
{
sourceInfo
->
Set
(
vtkStreamingDemandDrivenPipeline
::
UPDATE_EXTENT
(),
sourceInfo
->
Get
(
vtkStreamingDemandDrivenPipeline
::
WHOLE_EXTENT
()),
6
);
}
return
1
;
}
vtkErrorMacro
(
"Missing input or output info!"
);
return
0
;
}
//------------------------------------------------------------------------------
void
vtkmProbe
::
PassAttributeData
(
vtkDataSet
*
input
,
vtkDataObject
*
vtkNotUsed
(
source
),
vtkDataSet
*
output
)
{
if
(
this
->
PassPointArrays
)
{
// Copy point data arrays
int
numPtArrays
=
input
->
GetPointData
()
->
GetNumberOfArrays
();
for
(
int
i
=
0
;
i
<
numPtArrays
;
i
++
)
{
vtkDataArray
*
da
=
input
->
GetPointData
()
->
GetArray
(
i
);
if
(
da
&&
!
output
->
GetPointData
()
->
HasArray
(
da
->
GetName
()))
{
output
->
GetPointData
()
->
AddArray
(
da
);
}
}
// Set active attributes in the ouput to the active attributes in the input
for
(
int
i
=
0
;
i
<
vtkDataSetAttributes
::
NUM_ATTRIBUTES
;
++
i
)
{
vtkAbstractArray
*
da
=
input
->
GetPointData
()
->
GetAttribute
(
i
);
if
(
da
&&
da
->
GetName
()
&&
!
output
->
GetPointData
()
->
GetAttribute
(
i
))
{
output
->
GetPointData
()
->
SetAttribute
(
da
,
i
);
}
}
}
// copy cell data arrays
if
(
this
->
PassCellArrays
)
{
int
numCellArrays
=
input
->
GetCellData
()
->
GetNumberOfArrays
();
for
(
int
i
=
0
;
i
<
numCellArrays
;
++
i
)
{
vtkDataArray
*
da
=
input
->
GetCellData
()
->
GetArray
(
i
);
if
(
!
output
->
GetCellData
()
->
HasArray
(
da
->
GetName
()))
{
output
->
GetCellData
()
->
AddArray
(
da
);
}
}
// Set active attributes in the output to the active attributes in the input
for
(
int
i
=
0
;
i
<
vtkDataSetAttributes
::
NUM_ATTRIBUTES
;
++
i
)
{
vtkAbstractArray
*
da
=
input
->
GetCellData
()
->
GetAttribute
(
i
);
if
(
da
&&
da
->
GetName
()
&&
!
output
->
GetCellData
()
->
GetAttribute
(
i
))
{
output
->
GetCellData
()
->
SetAttribute
(
da
,
i
);
}
}
}
if
(
this
->
PassFieldArrays
)
{
// nothing to do, vtkDemandDrivenPipeline takes care of that.
}
else
{
output
->
GetFieldData
()
->
Initialize
();
}
}
//------------------------------------------------------------------------------
void
vtkmProbe
::
PrintSelf
(
ostream
&
os
,
vtkIndent
indent
)
{
this
->
Superclass
::
PrintSelf
(
os
,
indent
);
os
<<
indent
<<
"PassPointArrays: "
<<
this
->
PassPointArrays
<<
"
\n
"
;
os
<<
indent
<<
"PassCellArrays: "
<<
this
->
PassCellArrays
<<
"
\n
"
;
os
<<
indent
<<
"PassFieldArray: "
<<
this
->
PassFieldArrays
<<
"
\n
"
;
}
Accelerators/Vtkm/vtkmProbe.h
0 → 100644
View file @
a4891a53
//=============================================================================
//
// 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.
//
//=============================================================================
/**
* @class vtkmProbe
* @brief Sample data at specified point locations
*
* vtkmProbe is a filter that computes point attributes(e.g., scalars, vectors,
* etc.) at specific point positions using the probe filter in VTK-m. The
* filter has two inputs: the Input and Source.
* The Input geometric structure is passed through the filter. The point
* attributes are computed at the Input point positions by interpolating into
* the source data. For example, we can compute data values on a plane(plane
* specified as Input from a volume(Source). The source geometry must have cellSet
* defined otherwise the vtkm filter won't work. The cell data of the source data
* is copied to the output based on in which source cell each input point is. If
* an array of the same name exists both in source's point and cell data, only
* the one from the point data is probed. The valid point result is stored as
* a field array whose default name is "vtkValidPointMask" in the point data and
* the valid cell result(Invalid cells are the cells with at least one invalid
* point) is stored as a field array whose default name is "vtkValidCellMask" in
* the cell data.
*
* This filter can be used to resample data, or convert one dataset form into
* another. For example, an unstructured grid (vtkUnstructuredGrid) can be
* probed with a volume (three-dimensional vtkImageData), and then volume
* rendering techniques can be used to visualize the results. Another example:
* a line or curve can be used to probe data to produce x-y plots along
* that line or curve.
*/
#ifndef vtkmProbe_h
#define vtkmProbe_h
#include <algorithm>
#include <utility>
#include <vector>
#include "vtkAcceleratorsVTKmModule.h" //required for export
#include "vtkDataSetAlgorithm.h"
class
VTKACCELERATORSVTKM_EXPORT
vtkmProbe
:
public
vtkDataSetAlgorithm
{
public:
vtkTypeMacro
(
vtkmProbe
,
vtkDataSetAlgorithm
)
void
PrintSelf
(
ostream
&
os
,
vtkIndent
indent
)
override
;
static
vtkmProbe
*
New
();
//@{
/**
* Specify the data set that will be probed at the input points.
* The Input gives the geometry (the points and cells) for the output,
* while the Source is probed (interpolated) to generate the scalars,
* vectors, etc. for the output points based on the point locations.
*/
void
SetSourceData
(
vtkDataObject
*
source
);
vtkDataObject
*
GetSource
();
//@}
//@}
/**
* Specify the data set that will be probed at the input points.
* The Input gives the geometry (the points and cells) for the output,
* while the Source is probed (interpolated) to generate the scalars,
* vectors, etc. for the output points based on the point locations.
*/
void
SetSourceConnection
(
vtkAlgorithmOutput
*
algOutput
);
//@}
//@{
/**
* Shallow copy the input cell data arrays to the output.
* Off by default.
*/
vtkSetMacro
(
PassCellArrays
,
vtkTypeBool
);
vtkBooleanMacro
(
PassCellArrays
,
vtkTypeBool
);
vtkGetMacro
(
PassCellArrays
,
vtkTypeBool
);
//@}
//@{
/**
* Shallow copy the input point data arrays to the output.
* Off by default.
*/
vtkSetMacro
(
PassPointArrays
,
vtkTypeBool
);
vtkBooleanMacro
(
PassPointArrays
,
vtkTypeBool
);
vtkGetMacro
(
PassPointArrays
,
vtkTypeBool
);
//@}
//@{
/**
* Set whether to pass the field-data arrays from the Input i.e. the input
* providing the geometry to the output. On by default.
*/
vtkSetMacro
(
PassFieldArrays
,
vtkTypeBool
);
vtkBooleanMacro
(
PassFieldArrays
,
vtkTypeBool
);
vtkGetMacro
(
PassFieldArrays
,
vtkTypeBool
);
//@}
//@{
/**
* Returns the name of the valid point array added to the output with values 2 for
* hidden points and 0 for valid points.
* Set to "vtkValidPointMask" by default.
*/
vtkSetMacro
(
ValidPointMaskArrayName
,
std
::
string
)
vtkGetMacro
(
ValidPointMaskArrayName
,
std
::
string
)
//@}
//@{
/**
* Returns the name of the valid cell array added to the output with values 2 for
* hidden points and 0 for valid points.
* Set to "vtkValidCellMask" by default.
*/
vtkSetMacro
(
ValidCellMaskArrayName
,
std
::
string
)
vtkGetMacro
(
ValidCellMaskArrayName
,
std
::
string
)
//@}
protected:
vtkmProbe
();
~
vtkmProbe
()
=
default
;
vtkTypeBool
PassCellArrays
;
vtkTypeBool
PassPointArrays
;
vtkTypeBool
PassFieldArrays
;
std
::
string
ValidPointMaskArrayName
;
std
::
string
ValidCellMaskArrayName
;
virtual
int
RequestData
(
vtkInformation
*
,
vtkInformationVector
**
,
vtkInformationVector
*
)
override
;
virtual
int
RequestUpdateExtent
(
vtkInformation
*
,
vtkInformationVector
**
,
vtkInformationVector
*
)
override
;
virtual
int
RequestInformation
(
vtkInformation
*
,
vtkInformationVector
**
,
vtkInformationVector
*
)
override
;
/**
* Call at the end of RequestData() to pass attribute dat a respecting the
* PassCellArrays, PassPointArrays and PassFieldArrays flag
*/
void
PassAttributeData
(
vtkDataSet
*
input
,
vtkDataObject
*
source
,
vtkDataSet
*
output
);
private:
vtkmProbe
(
const
vtkmProbe
&
)
=
delete
;
void
operator
=
(
const
vtkmProbe
&
)
=
delete
;
};
#endif //vtkmProbe_h
// VTK-HeaderTest-Exclude: vtkmProbe.h
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