VTK
|
Given a 3D domain space represented by an unstructured grid composed of tetrahedral cells with bivariate fields, this filter tessellates each cell in the domain to polyhedral fragments by intersecting the projection of the cell into 2-D range space against two sets of cutting planes, one set is defined along the first field, the second set is defined along the second field. More...
#include <vtkContinuousScatterplot.h>
Public Types | |
typedef vtkImageAlgorithm | Superclass |
Public Types inherited from vtkImageAlgorithm | |
typedef vtkAlgorithm | Superclass |
Public Types inherited from vtkAlgorithm | |
enum | DesiredOutputPrecision { SINGLE_PRECISION, DOUBLE_PRECISION, DEFAULT_PRECISION } |
Values used for setting the desired output precision for various algorithms. More... | |
typedef vtkObject | Superclass |
Public Member Functions | |
virtual vtkTypeBool | IsA (const char *type) |
Return 1 if this class is the same type of (or a subclass of) the named class. More... | |
vtkContinuousScatterplot * | NewInstance () const |
void | PrintSelf (ostream &os, vtkIndent indent) override |
Methods invoked by print to print information about the object including superclasses. More... | |
virtual double | GetEpsilon () |
Get the tolerance used when comparing floating point numbers for equality. More... | |
virtual void | SetEpsilon (double) |
Set the tolerance used when comparing floating point numbers for equality. More... | |
void | SetField1 (const char *fieldName, vtkIdType ResX) |
Specify the name of the first field to be used in subdividing the dataset. More... | |
void | SetField2 (const char *fieldName, vtkIdType ResY) |
Specify the name of the second field to be used in subdividing the dataset. More... | |
Public Member Functions inherited from vtkImageAlgorithm | |
vtkImageAlgorithm * | NewInstance () const |
int | ProcessRequest (vtkInformation *, vtkInformationVector **, vtkInformationVector *) override |
Process a request from the executive. More... | |
vtkImageData * | GetOutput () |
Get the output data object for a port on this algorithm. More... | |
vtkImageData * | GetOutput (int) |
Get the output data object for a port on this algorithm. More... | |
virtual void | SetOutput (vtkDataObject *d) |
Get the output data object for a port on this algorithm. More... | |
void | SetInputData (vtkDataObject *) |
Assign a data object as input. More... | |
void | SetInputData (int, vtkDataObject *) |
Assign a data object as input. More... | |
vtkDataObject * | GetInput (int port) |
Get a data object for one of the input port connections. More... | |
vtkDataObject * | GetInput () |
Get a data object for one of the input port connections. More... | |
vtkImageData * | GetImageDataInput (int port) |
Get a data object for one of the input port connections. More... | |
virtual void | AddInputData (vtkDataObject *) |
Assign a data object as input. More... | |
virtual void | AddInputData (int, vtkDataObject *) |
Assign a data object as input. More... | |
Public Member Functions inherited from vtkAlgorithm | |
vtkAlgorithm * | NewInstance () const |
int | HasExecutive () |
Check whether this algorithm has an assigned executive. More... | |
vtkExecutive * | GetExecutive () |
Get this algorithm's executive. More... | |
virtual void | SetExecutive (vtkExecutive *executive) |
Set this algorithm's executive. More... | |
int | ProcessRequest (vtkInformation *request, vtkCollection *inInfo, vtkInformationVector *outInfo) |
Version of ProcessRequest() that is wrapped. More... | |
virtual int | ComputePipelineMTime (vtkInformation *request, vtkInformationVector **inInfoVec, vtkInformationVector *outInfoVec, int requestFromOutputPort, vtkMTimeType *mtime) |
A special version of ProcessRequest meant specifically for the pipeline modified time request. More... | |
virtual int | ModifyRequest (vtkInformation *request, int when) |
This method gives the algorithm a chance to modify the contents of a request before or after (specified in the when argument) it is forwarded. More... | |
vtkInformation * | GetInputPortInformation (int port) |
Get the information object associated with an input port. More... | |
vtkInformation * | GetOutputPortInformation (int port) |
Get the information object associated with an output port. More... | |
int | GetNumberOfInputPorts () |
Get the number of input ports used by the algorithm. More... | |
int | GetNumberOfOutputPorts () |
Get the number of output ports provided by the algorithm. More... | |
void | UpdateProgress (double amount) |
Update the progress of the process object. More... | |
virtual void | SetInputArrayToProcess (int idx, int port, int connection, const char *fieldAssociation, const char *attributeTypeorName) |
String based versions of SetInputArrayToProcess(). More... | |
vtkInformation * | GetInputArrayInformation (int idx) |
Get the info object for the specified input array to this algorithm. More... | |
void | RemoveAllInputs () |
Remove all the input data. More... | |
vtkDataObject * | GetOutputDataObject (int port) |
Get the data object that will contain the algorithm output for the given port. More... | |
vtkDataObject * | GetInputDataObject (int port, int connection) |
Get the data object that will contain the algorithm input for the given port and given connection. More... | |
virtual void | RemoveInputConnection (int port, vtkAlgorithmOutput *input) |
Remove a connection from the given input port index. More... | |
virtual void | RemoveInputConnection (int port, int idx) |
Remove a connection given by index idx. More... | |
virtual void | RemoveAllInputConnections (int port) |
Removes all input connections. More... | |
virtual void | SetInputDataObject (int port, vtkDataObject *data) |
Sets the data-object as an input on the given port index. More... | |
virtual void | SetInputDataObject (vtkDataObject *data) |
virtual void | AddInputDataObject (int port, vtkDataObject *data) |
Add the data-object as an input to this given port. More... | |
virtual void | AddInputDataObject (vtkDataObject *data) |
vtkAlgorithmOutput * | GetOutputPort (int index) |
Get a proxy object corresponding to the given output port of this algorithm. More... | |
vtkAlgorithmOutput * | GetOutputPort () |
int | GetNumberOfInputConnections (int port) |
Get the number of inputs currently connected to a port. More... | |
int | GetTotalNumberOfInputConnections () |
Get the total number of inputs for this algorithm. More... | |
vtkAlgorithmOutput * | GetInputConnection (int port, int index) |
Get the algorithm output port connected to an input port. More... | |
vtkAlgorithm * | GetInputAlgorithm (int port, int index, int &algPort) |
Returns the algorithm and the output port index of that algorithm connected to a port-index pair. More... | |
vtkAlgorithm * | GetInputAlgorithm (int port, int index) |
Returns the algorithm connected to a port-index pair. More... | |
vtkAlgorithm * | GetInputAlgorithm () |
Equivalent to GetInputAlgorithm(0, 0). More... | |
vtkExecutive * | GetInputExecutive (int port, int index) |
Returns the executive associated with a particular input connection. More... | |
vtkExecutive * | GetInputExecutive () |
Equivalent to GetInputExecutive(0, 0) More... | |
vtkInformation * | GetInputInformation (int port, int index) |
Return the information object that is associated with a particular input connection. More... | |
vtkInformation * | GetInputInformation () |
Equivalent to GetInputInformation(0, 0) More... | |
vtkInformation * | GetOutputInformation (int port) |
Return the information object that is associated with a particular output port. More... | |
virtual int | Update (int port, vtkInformationVector *requests) |
This method enables the passing of data requests to the algorithm to be used during execution (in addition to bringing a particular port up-to-date). More... | |
virtual int | Update (vtkInformation *requests) |
Convenience method to update an algorithm after passing requests to its first output port. More... | |
virtual int | UpdatePiece (int piece, int numPieces, int ghostLevels, const int extents[6]=nullptr) |
Convenience method to update an algorithm after passing requests to its first output port. More... | |
virtual int | UpdateExtent (const int extents[6]) |
Convenience method to update an algorithm after passing requests to its first output port. More... | |
virtual int | UpdateTimeStep (double time, int piece=-1, int numPieces=1, int ghostLevels=0, const int extents[6]=nullptr) |
Convenience method to update an algorithm after passing requests to its first output port. More... | |
virtual void | UpdateInformation () |
Bring the algorithm's information up-to-date. More... | |
virtual void | UpdateDataObject () |
Create output object(s). More... | |
virtual void | PropagateUpdateExtent () |
Propagate meta-data upstream. More... | |
virtual void | UpdateWholeExtent () |
Bring this algorithm's outputs up-to-date. More... | |
void | ConvertTotalInputToPortConnection (int ind, int &port, int &conn) |
Convenience routine to convert from a linear ordering of input connections to a port/connection pair. More... | |
virtual vtkInformation * | GetInformation () |
Set/Get the information object associated with this algorithm. More... | |
virtual void | SetInformation (vtkInformation *) |
Set/Get the information object associated with this algorithm. More... | |
void | Register (vtkObjectBase *o) override |
Participate in garbage collection. More... | |
void | UnRegister (vtkObjectBase *o) override |
Participate in garbage collection. More... | |
virtual void | SetAbortExecute (vtkTypeBool) |
Set/Get the AbortExecute flag for the process object. More... | |
virtual vtkTypeBool | GetAbortExecute () |
Set/Get the AbortExecute flag for the process object. More... | |
virtual void | AbortExecuteOn () |
Set/Get the AbortExecute flag for the process object. More... | |
virtual void | AbortExecuteOff () |
Set/Get the AbortExecute flag for the process object. More... | |
virtual void | SetProgress (double) |
Set/Get the execution progress of a process object. More... | |
virtual double | GetProgress () |
Set/Get the execution progress of a process object. More... | |
void | SetProgressText (const char *ptext) |
Set the current text message associated with the progress state. More... | |
virtual char * | GetProgressText () |
Set the current text message associated with the progress state. More... | |
virtual unsigned long | GetErrorCode () |
The error code contains a possible error that occurred while reading or writing the file. More... | |
virtual void | SetInputArrayToProcess (int idx, int port, int connection, int fieldAssociation, const char *name) |
Set the input data arrays that this algorithm will process. More... | |
virtual void | SetInputArrayToProcess (int idx, int port, int connection, int fieldAssociation, int fieldAttributeType) |
Set the input data arrays that this algorithm will process. More... | |
virtual void | SetInputArrayToProcess (int idx, vtkInformation *info) |
Set the input data arrays that this algorithm will process. More... | |
virtual void | SetInputConnection (int port, vtkAlgorithmOutput *input) |
Set the connection for the given input port index. More... | |
virtual void | SetInputConnection (vtkAlgorithmOutput *input) |
Set the connection for the given input port index. More... | |
virtual void | AddInputConnection (int port, vtkAlgorithmOutput *input) |
Add a connection to the given input port index. More... | |
virtual void | AddInputConnection (vtkAlgorithmOutput *input) |
Add a connection to the given input port index. More... | |
virtual void | Update (int port) |
Bring this algorithm's outputs up-to-date. More... | |
virtual void | Update () |
Bring this algorithm's outputs up-to-date. More... | |
virtual void | SetReleaseDataFlag (int) |
Turn release data flag on or off for all output ports. More... | |
virtual int | GetReleaseDataFlag () |
Turn release data flag on or off for all output ports. More... | |
void | ReleaseDataFlagOn () |
Turn release data flag on or off for all output ports. More... | |
void | ReleaseDataFlagOff () |
Turn release data flag on or off for all output ports. More... | |
int | UpdateExtentIsEmpty (vtkInformation *pinfo, vtkDataObject *output) |
This detects when the UpdateExtent will generate no data This condition is satisfied when the UpdateExtent has zero volume (0,-1,...) or the UpdateNumberOfPieces is 0. More... | |
int | UpdateExtentIsEmpty (vtkInformation *pinfo, int extentType) |
This detects when the UpdateExtent will generate no data This condition is satisfied when the UpdateExtent has zero volume (0,-1,...) or the UpdateNumberOfPieces is 0. More... | |
int * | GetUpdateExtent () |
These functions return the update extent for output ports that use 3D extents. More... | |
int * | GetUpdateExtent (int port) |
These functions return the update extent for output ports that use 3D extents. More... | |
void | GetUpdateExtent (int &x0, int &x1, int &y0, int &y1, int &z0, int &z1) |
These functions return the update extent for output ports that use 3D extents. More... | |
void | GetUpdateExtent (int port, int &x0, int &x1, int &y0, int &y1, int &z0, int &z1) |
These functions return the update extent for output ports that use 3D extents. More... | |
void | GetUpdateExtent (int extent[6]) |
These functions return the update extent for output ports that use 3D extents. More... | |
void | GetUpdateExtent (int port, int extent[6]) |
These functions return the update extent for output ports that use 3D extents. More... | |
int | GetUpdatePiece () |
These functions return the update extent for output ports that use piece extents. More... | |
int | GetUpdatePiece (int port) |
These functions return the update extent for output ports that use piece extents. More... | |
int | GetUpdateNumberOfPieces () |
These functions return the update extent for output ports that use piece extents. More... | |
int | GetUpdateNumberOfPieces (int port) |
These functions return the update extent for output ports that use piece extents. More... | |
int | GetUpdateGhostLevel () |
These functions return the update extent for output ports that use piece extents. More... | |
int | GetUpdateGhostLevel (int port) |
These functions return the update extent for output ports that use piece extents. More... | |
void | SetProgressObserver (vtkProgressObserver *) |
If an ProgressObserver is set, the algorithm will report progress through it rather than directly. More... | |
virtual vtkProgressObserver * | GetProgressObserver () |
If an ProgressObserver is set, the algorithm will report progress through it rather than directly. More... | |
Public Member Functions inherited from vtkObject | |
vtkBaseTypeMacro (vtkObject, vtkObjectBase) | |
virtual void | DebugOn () |
Turn debugging output on. More... | |
virtual void | DebugOff () |
Turn debugging output off. More... | |
bool | GetDebug () |
Get the value of the debug flag. More... | |
void | SetDebug (bool debugFlag) |
Set the value of the debug flag. More... | |
virtual void | Modified () |
Update the modification time for this object. More... | |
virtual vtkMTimeType | GetMTime () |
Return this object's modified time. More... | |
void | RemoveObserver (unsigned long tag) |
void | RemoveObservers (unsigned long event) |
void | RemoveObservers (const char *event) |
void | RemoveAllObservers () |
vtkTypeBool | HasObserver (unsigned long event) |
vtkTypeBool | HasObserver (const char *event) |
int | InvokeEvent (unsigned long event) |
int | InvokeEvent (const char *event) |
unsigned long | AddObserver (unsigned long event, vtkCommand *, float priority=0.0f) |
Allow people to add/remove/invoke observers (callbacks) to any VTK object. More... | |
unsigned long | AddObserver (const char *event, vtkCommand *, float priority=0.0f) |
Allow people to add/remove/invoke observers (callbacks) to any VTK object. More... | |
vtkCommand * | GetCommand (unsigned long tag) |
Allow people to add/remove/invoke observers (callbacks) to any VTK object. More... | |
void | RemoveObserver (vtkCommand *) |
Allow people to add/remove/invoke observers (callbacks) to any VTK object. More... | |
void | RemoveObservers (unsigned long event, vtkCommand *) |
Allow people to add/remove/invoke observers (callbacks) to any VTK object. More... | |
void | RemoveObservers (const char *event, vtkCommand *) |
Allow people to add/remove/invoke observers (callbacks) to any VTK object. More... | |
vtkTypeBool | HasObserver (unsigned long event, vtkCommand *) |
Allow people to add/remove/invoke observers (callbacks) to any VTK object. More... | |
vtkTypeBool | HasObserver (const char *event, vtkCommand *) |
Allow people to add/remove/invoke observers (callbacks) to any VTK object. More... | |
template<class U , class T > | |
unsigned long | AddObserver (unsigned long event, U observer, void(T::*callback)(), float priority=0.0f) |
Overloads to AddObserver that allow developers to add class member functions as callbacks for events. More... | |
template<class U , class T > | |
unsigned long | AddObserver (unsigned long event, U observer, void(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f) |
Overloads to AddObserver that allow developers to add class member functions as callbacks for events. More... | |
template<class U , class T > | |
unsigned long | AddObserver (unsigned long event, U observer, bool(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f) |
Allow user to set the AbortFlagOn() with the return value of the callback method. More... | |
int | InvokeEvent (unsigned long event, void *callData) |
This method invokes an event and return whether the event was aborted or not. More... | |
int | InvokeEvent (const char *event, void *callData) |
This method invokes an event and return whether the event was aborted or not. More... | |
Public Member Functions inherited from vtkObjectBase | |
const char * | GetClassName () const |
Return the class name as a string. More... | |
virtual void | Delete () |
Delete a VTK object. More... | |
virtual void | FastDelete () |
Delete a reference to this object. More... | |
void | InitializeObjectBase () |
void | Print (ostream &os) |
Print an object to an ostream. More... | |
int | GetReferenceCount () |
Return the current reference count of this object. More... | |
void | SetReferenceCount (int) |
Sets the reference count. More... | |
void | PrintRevisions (ostream &) |
Legacy. More... | |
virtual void | PrintHeader (ostream &os, vtkIndent indent) |
Methods invoked by print to print information about the object including superclasses. More... | |
virtual void | PrintTrailer (ostream &os, vtkIndent indent) |
Methods invoked by print to print information about the object including superclasses. More... | |
Static Public Member Functions | |
static vtkContinuousScatterplot * | New () |
static vtkTypeBool | IsTypeOf (const char *type) |
static vtkContinuousScatterplot * | SafeDownCast (vtkObjectBase *o) |
Static Public Member Functions inherited from vtkImageAlgorithm | |
static vtkTypeBool | IsTypeOf (const char *type) |
static vtkImageAlgorithm * | SafeDownCast (vtkObjectBase *o) |
Static Public Member Functions inherited from vtkAlgorithm | |
static vtkAlgorithm * | New () |
static vtkTypeBool | IsTypeOf (const char *type) |
static vtkAlgorithm * | SafeDownCast (vtkObjectBase *o) |
static vtkInformationIntegerKey * | INPUT_IS_OPTIONAL () |
Keys used to specify input port requirements. More... | |
static vtkInformationIntegerKey * | INPUT_IS_REPEATABLE () |
static vtkInformationInformationVectorKey * | INPUT_REQUIRED_FIELDS () |
static vtkInformationStringVectorKey * | INPUT_REQUIRED_DATA_TYPE () |
static vtkInformationInformationVectorKey * | INPUT_ARRAYS_TO_PROCESS () |
static vtkInformationIntegerKey * | INPUT_PORT () |
static vtkInformationIntegerKey * | INPUT_CONNECTION () |
static vtkInformationIntegerKey * | CAN_PRODUCE_SUB_EXTENT () |
This key tells the executive that a particular output port is capable of producing an arbitrary subextent of the whole extent. More... | |
static vtkInformationIntegerKey * | CAN_HANDLE_PIECE_REQUEST () |
Key that tells the pipeline that a particular algorithm can or cannot handle piece request. More... | |
static void | SetDefaultExecutivePrototype (vtkExecutive *proto) |
If the DefaultExecutivePrototype is set, a copy of it is created in CreateDefaultExecutive() using NewInstance(). More... | |
Static Public Member Functions inherited from vtkObject | |
static vtkObject * | New () |
Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More... | |
static void | BreakOnError () |
This method is called when vtkErrorMacro executes. More... | |
static void | SetGlobalWarningDisplay (int val) |
This is a global flag that controls whether any debug, warning or error messages are displayed. More... | |
static void | GlobalWarningDisplayOn () |
This is a global flag that controls whether any debug, warning or error messages are displayed. More... | |
static void | GlobalWarningDisplayOff () |
This is a global flag that controls whether any debug, warning or error messages are displayed. More... | |
static int | GetGlobalWarningDisplay () |
This is a global flag that controls whether any debug, warning or error messages are displayed. More... | |
Static Public Member Functions inherited from vtkObjectBase | |
static vtkTypeBool | IsTypeOf (const char *name) |
Return 1 if this class type is the same type of (or a subclass of) the named class. More... | |
static vtkObjectBase * | New () |
Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More... | |
Protected Member Functions | |
virtual vtkObjectBase * | NewInstanceInternal () const |
vtkContinuousScatterplot () | |
int | FillInputPortInformation (int port, vtkInformation *info) override |
These method should be reimplemented by subclasses that have more than a single input or single output. More... | |
int | FillOutputPortInformation (int port, vtkInformation *info) override |
These method should be reimplemented by subclasses that have more than a single input or single output. More... | |
int | RequestData (vtkInformation *, vtkInformationVector **, vtkInformationVector *) override |
This is called in response to a REQUEST_DATA request from the executive. More... | |
Protected Member Functions inherited from vtkImageAlgorithm | |
vtkImageAlgorithm () | |
~vtkImageAlgorithm () override | |
virtual int | RequestInformation (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) |
Subclasses can reimplement this method to collect information from their inputs and set information for their outputs. More... | |
virtual int | RequestUpdateExtent (vtkInformation *, vtkInformationVector **, vtkInformationVector *) |
Subclasses can reimplement this method to translate the update extent requests from each output port into update extent requests for the input connections. More... | |
virtual void | CopyInputArrayAttributesToOutput (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) |
Convenience method to copy the scalar type and number of components from the input data to the output data. More... | |
virtual void | ExecuteDataWithInformation (vtkDataObject *output, vtkInformation *outInfo) |
This is a convenience method that is implemented in many subclasses instead of RequestData. More... | |
virtual void | CopyAttributeData (vtkImageData *in, vtkImageData *out, vtkInformationVector **inputVector) |
Copy the other point and cell data. More... | |
virtual void | ExecuteData (vtkDataObject *output) |
This method is the old style execute method, provided for the sake of backwards compatibility with older filters and readers. More... | |
virtual void | Execute () |
This method is the old style execute method, provided for the sake of backwards compatibility with older filters and readers. More... | |
virtual void | AllocateOutputData (vtkImageData *out, vtkInformation *outInfo, int *uExtent) |
Allocate the output data. More... | |
virtual vtkImageData * | AllocateOutputData (vtkDataObject *out, vtkInformation *outInfo) |
Allocate the output data. More... | |
Protected Member Functions inherited from vtkAlgorithm | |
vtkAlgorithm () | |
~vtkAlgorithm () override | |
virtual void | SetNumberOfInputPorts (int n) |
Set the number of input ports used by the algorithm. More... | |
virtual void | SetNumberOfOutputPorts (int n) |
Set the number of output ports provided by the algorithm. More... | |
int | InputPortIndexInRange (int index, const char *action) |
int | OutputPortIndexInRange (int index, const char *action) |
int | GetInputArrayAssociation (int idx, vtkInformationVector **inputVector) |
Get the assocition of the actual data array for the input array specified by idx, this is only reasonable during the REQUEST_DATA pass. More... | |
vtkInformation * | GetInputArrayFieldInformation (int idx, vtkInformationVector **inputVector) |
This method takes in an index (as specified in SetInputArrayToProcess) and a pipeline information vector. More... | |
virtual vtkExecutive * | CreateDefaultExecutive () |
Create a default executive. More... | |
void | ReportReferences (vtkGarbageCollector *) override |
virtual void | SetNthInputConnection (int port, int index, vtkAlgorithmOutput *input) |
Replace the Nth connection on the given input port. More... | |
virtual void | SetNumberOfInputConnections (int port, int n) |
Set the number of input connections on the given input port. More... | |
void | SetInputDataInternal (int port, vtkDataObject *input) |
These methods are used by subclasses to implement methods to set data objects directly as input. More... | |
void | AddInputDataInternal (int port, vtkDataObject *input) |
int | GetInputArrayAssociation (int idx, int connection, vtkInformationVector **inputVector) |
Filters that have multiple connections on one port can use this signature. More... | |
int | GetInputArrayAssociation (int idx, vtkDataObject *input) |
Filters that have multiple connections on one port can use this signature. More... | |
vtkDataArray * | GetInputArrayToProcess (int idx, vtkInformationVector **inputVector) |
Get the actual data array for the input array specified by idx, this is only reasonable during the REQUEST_DATA pass. More... | |
vtkDataArray * | GetInputArrayToProcess (int idx, vtkInformationVector **inputVector, int &association) |
Get the actual data array for the input array specified by idx, this is only reasonable during the REQUEST_DATA pass. More... | |
vtkDataArray * | GetInputArrayToProcess (int idx, int connection, vtkInformationVector **inputVector) |
Filters that have multiple connections on one port can use this signature. More... | |
vtkDataArray * | GetInputArrayToProcess (int idx, int connection, vtkInformationVector **inputVector, int &association) |
Filters that have multiple connections on one port can use this signature. More... | |
vtkDataArray * | GetInputArrayToProcess (int idx, vtkDataObject *input) |
Filters that have multiple connections on one port can use this signature. More... | |
vtkDataArray * | GetInputArrayToProcess (int idx, vtkDataObject *input, int &association) |
Filters that have multiple connections on one port can use this signature. More... | |
vtkAbstractArray * | GetInputAbstractArrayToProcess (int idx, vtkInformationVector **inputVector) |
Get the actual data array for the input array specified by idx, this is only reasonable during the REQUEST_DATA pass. More... | |
vtkAbstractArray * | GetInputAbstractArrayToProcess (int idx, vtkInformationVector **inputVector, int &association) |
Get the actual data array for the input array specified by idx, this is only reasonable during the REQUEST_DATA pass. More... | |
vtkAbstractArray * | GetInputAbstractArrayToProcess (int idx, int connection, vtkInformationVector **inputVector) |
Filters that have multiple connections on one port can use this signature. More... | |
vtkAbstractArray * | GetInputAbstractArrayToProcess (int idx, int connection, vtkInformationVector **inputVector, int &association) |
Filters that have multiple connections on one port can use this signature. More... | |
vtkAbstractArray * | GetInputAbstractArrayToProcess (int idx, vtkDataObject *input) |
Filters that have multiple connections on one port can use this signature. More... | |
vtkAbstractArray * | GetInputAbstractArrayToProcess (int idx, vtkDataObject *input, int &association) |
Filters that have multiple connections on one port can use this signature. More... | |
virtual void | SetErrorCode (unsigned long) |
The error code contains a possible error that occurred while reading or writing the file. More... | |
Protected Member Functions inherited from vtkObject | |
vtkObject () | |
~vtkObject () override | |
void | RegisterInternal (vtkObjectBase *, vtkTypeBool check) override |
void | UnRegisterInternal (vtkObjectBase *, vtkTypeBool check) override |
void | InternalGrabFocus (vtkCommand *mouseEvents, vtkCommand *keypressEvents=nullptr) |
These methods allow a command to exclusively grab all events. More... | |
void | InternalReleaseFocus () |
These methods allow a command to exclusively grab all events. More... | |
Protected Member Functions inherited from vtkObjectBase | |
vtkObjectBase () | |
virtual | ~vtkObjectBase () |
virtual void | CollectRevisions (ostream &) |
vtkObjectBase (const vtkObjectBase &) | |
void | operator= (const vtkObjectBase &) |
Protected Attributes | |
double | Epsilon |
const char * | Fields [2] |
vtkIdType | ResX |
vtkIdType | ResY |
Protected Attributes inherited from vtkAlgorithm | |
vtkInformation * | Information |
double | Progress |
char * | ProgressText |
vtkProgressObserver * | ProgressObserver |
unsigned long | ErrorCode |
The error code contains a possible error that occurred while reading or writing the file. More... | |
Protected Attributes inherited from vtkObject | |
bool | Debug |
vtkTimeStamp | MTime |
vtkSubjectHelper * | SubjectHelper |
Protected Attributes inherited from vtkObjectBase | |
vtkAtomicInt32 | ReferenceCount |
vtkWeakPointerBase ** | WeakPointers |
Additional Inherited Members | |
Public Attributes inherited from vtkAlgorithm | |
vtkTypeBool | AbortExecute |
Static Protected Member Functions inherited from vtkAlgorithm | |
static vtkInformationIntegerKey * | PORT_REQUIREMENTS_FILLED () |
Static Protected Attributes inherited from vtkAlgorithm | |
static vtkExecutive * | DefaultExecutivePrototype |
Given a 3D domain space represented by an unstructured grid composed of tetrahedral cells with bivariate fields, this filter tessellates each cell in the domain to polyhedral fragments by intersecting the projection of the cell into 2-D range space against two sets of cutting planes, one set is defined along the first field, the second set is defined along the second field.
The volume of these subdivided polyhedral fragments can be computed and aggregated over cells to depict the density distribution of the data projection in the bivariate range space.
Given a bivariate field (f1,f2) defined on an unstructured grid which is composed of tetrahedral cells, we can initially subdivide each cell based on its projection in the range into a number of fragments along the first field f1, we refer to these polyhedral fragments as Frag(f1) = {frag(f1)_1, frag(f1)_2, ... , frag(f1)_n}, where frag(f1)_n refers to the nth fragment along the first field subdivision. Each fragment has a range value and the value difference between the neighbouring fragments is represented as fragment width fw_f1, which is uniformly distributed across the range. Based on the structure of Frag(f1), for each of its cell "frag(f1)_n", we can further subdivide this cell based on the second field f2 using fragment width fw_f2. The tessellation along the second field results in an even finer fragment collection which we refer to as Frag(f1,f2) = {frag(f1,f2)_1, frag(f1,f2)_2, ... , frag(f1,f2)_m}. We can observe that Frag(f1,f2) is a finer tessellation of the domain than Frag(f1) and will be used to compute the density distribution in the bivariate range space. The algorithm for fragment computation is similar to the first stage of the work in [0]. Each fragment "s" in Frag(f1,f2) has range values (f1(s), f2(s)) in the bivariate fields. These values can be further mapped to a 2-D bin with a resolution rexX * resY. The mapped bin index (binIndexX, binIndexY) of the fragment can be computed by linear interpolation on its range values : binIndexX = (int) resX * (f1(s) - f1_min) / (f1_max - f1_min) binIndexY = (int) resY * (f2(s) - f2_min) / (f2_max - f2_min), where (f1_min, f1_max) is the range in first field. Once we know which bin a fragment coincides, the density value in each bin equals to the total geometric volume of the fragments in this bin. This volume distribution over the bins will be exported as a point data array in the output data structure. If we map this 2-D bin to a 2-D image with each bin corresponding to a pixel and bin density to pixel transparency, then the image can be displayed as a continuous scatterplot.
The algorithm of this filter can be described as: Require: R.1 The domain space is an unstructured grid data set composed of tetrahedral cells; R.2 The range space contains two scalar fields, say f1 and f2.
The most important step is to compute the fragments. The implementation processes the input grid one cell at a time, explicitly computing the intersection of the cell with the cutting planes defined by the fragment boundaries in each scalar field. In order to subdivide the cell, we need to define a list of cutting planes in each field. The interval between neighbouring cutting planes is related to the output 2-D bin resolution (resX, resY) and can be computed as : fw_f1 = (f1_max - f1_min) / resX fw_f2 = (f2_max - f2_min) / resY, where (f1_max,f1_min) is the scalar range of first field.
For each tetrahedron T in the input grid:
1.1 Subdivide the cell T based on the first field f1, we will obtain a list of fragments: Frag(f1) = {frag(f1)_1, frag(f1)_2, ... , frag(f1)_n}. The steps for subdivision can be described as:
1.1.1 For each cutting plane s with respect to the first field f1, its field value f1(s) = f1_min + n * fw_f1, where n refers to the n-th cutting plane:
1.1.2. Traverse each edge e starting from point a to b in the cell, we will maintain three data classes, namely fragmentFace, residualFace and cutSet: A. fragmentFace contains vertices in the current fragment. B. cutSet contains vertices whose range values equal to f1(s). This set contains the current cutting plane. C. residualFace contains the rest of the vertices in the cell. In order to classify edge vertices into these classes, the following case table is used for each vertex "a" : case 0 : f1(a)---— f1(s) ---—f1(b) condition: f1(a) < f1(s) , f1(b) > f1(s) class: p(s,e), a -> fragmentFace p(s,e) -> cutSet p(s,e) -> residualFace
case 1 : f1(b)---— f1(s) ---—f1(a) condition: f1(a) > f1(s) , f1(b) < f1(s) class: p(s,e) -> fragmentFace p(s,e) -> cutSet a -> residualFace
case 2 : f1(s),f1(a)----------------—f1(b) condition: f1(s) == f1(a), f1(s) <= f1(b) class: a -> fragmentFace a -> residualFace a -> cutSet
case 3 : f1(a)----------------—f1(b), f1(s) condition: f1(s) > f1(a), f1(s) == f1(b) class: a -> fragmentFace
case 4 : f1(s),f1(b)----------------—f1(a) condition: f1(s) < f1(a), f1(s) == f1(b) class: a -> residualFace Remark: 1. we use "->" to indicate "belongs to" relation.
1.1.3. After we have traversed every edge in a cell for the cutting plane s, three classes for storing fragment, cutting plane and residual faces are updated. The faces of the current fragment frag(f1) are the union of all elements in fragmentFace and cutSet.
1.2 Take the output of step 1.1, traverse each fragment in Frag(f1), define a list of cutting planes with respect to field f2, further subdivide the fragments in Frag(f1) following steps from 1.1.2 to 1.1.3. The output of this step will be the fragment collection Frag(f1,f2). Each fragment in Frag(f1,f2) can be further mapped to a 2-D bin based on its range values. The density value in each bin equals to the total geometric volume of the fragments in this bin. This volume distribution over the bins will be exported as a point data array in the output data structure.
The input and output ports of the filter: Input port : the input data set should be a vtkUnstructuredGrid, with each of its cell defined as a tetrahedron. At least two scalar fields are associated with the data. The user needs to specify the name of the two scalar arrays beforehand. Output port: the output data set is a 2D image stored as a vtkImageData. The resolution of the output image can be set by the user. The volume distribution of fragments in each pixel or bin is stored in an point data array named "volume" in the output vtkImageData.
Suppose we have a tetrahedral mesh stored in a vtkUnstructuredGrid, we call this data set "inputData". This data set has two scalar arrays whose names are "f1" and "f2" respectively. We would like the resolution of output image set to (resX,resY). Given these input, this filter can be called as follows in c++ sample code:
vtkSmartPointer<vtkContinuousScatterplot> csp = vtkSmartPointer<vtkContinuousScatterplot>::New(); csp->SetInputData(inputData); csp->SetField1("f1",resX); csp->SetField2("f2",resY); csp->Update();
Then the output, "csp->GetOutput()", will be a vtkImageData containing a scalar array whose name is "volume". This array contains the volume distribution of the fragments.
[0] H.Carr and D.Duke, Joint contour nets: Topological analysis of multivariate data. IEEE Transactions on Visualization and Computer Graphics, volume 20, issue 08, pages 1100-1113, 2014
Definition at line 178 of file vtkContinuousScatterplot.h.
Definition at line 182 of file vtkContinuousScatterplot.h.
|
protected |
|
static |
|
static |
|
virtual |
Return 1 if this class is the same type of (or a subclass of) the named class.
Returns 0 otherwise. This method works in combination with vtkTypeMacro found in vtkSetGet.h.
Reimplemented from vtkImageAlgorithm.
|
static |
|
protectedvirtual |
Reimplemented from vtkImageAlgorithm.
vtkContinuousScatterplot* vtkContinuousScatterplot::NewInstance | ( | ) | const |
|
overridevirtual |
Methods invoked by print to print information about the object including superclasses.
Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.
Reimplemented from vtkImageAlgorithm.
|
virtual |
Get the tolerance used when comparing floating point numbers for equality.
|
virtual |
Set the tolerance used when comparing floating point numbers for equality.
void vtkContinuousScatterplot::SetField1 | ( | const char * | fieldName, |
vtkIdType | ResX | ||
) |
Specify the name of the first field to be used in subdividing the dataset.
Specify the resolution along x axis of the output image.
void vtkContinuousScatterplot::SetField2 | ( | const char * | fieldName, |
vtkIdType | ResY | ||
) |
Specify the name of the second field to be used in subdividing the dataset.
Specify the resolution along y axis of the output image.
|
overrideprotectedvirtual |
These method should be reimplemented by subclasses that have more than a single input or single output.
See vtkAlgorithm for more information.
Reimplemented from vtkImageAlgorithm.
|
overrideprotectedvirtual |
These method should be reimplemented by subclasses that have more than a single input or single output.
See vtkAlgorithm for more information.
Reimplemented from vtkImageAlgorithm.
|
overrideprotectedvirtual |
This is called in response to a REQUEST_DATA request from the executive.
Subclasses should override either this method or the ExecuteDataWithInformation method in order to generate data for their outputs. For images, the output arrays will already be allocated, so all that is necessary is to fill in the voxel values.
Reimplemented from vtkImageAlgorithm.
|
protected |
Definition at line 219 of file vtkContinuousScatterplot.h.
|
protected |
Definition at line 222 of file vtkContinuousScatterplot.h.
|
protected |
Definition at line 225 of file vtkContinuousScatterplot.h.
|
protected |
Definition at line 225 of file vtkContinuousScatterplot.h.