Isorropia: Partitioning, Load Balancing and more
Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes | Private Member Functions | Private Attributes | List of all members
Isorropia::Epetra::ZoltanLib::QueryObject Class Reference

QueryObject is a class that contains the query functions required by the Zoltan library. More...

#include <QueryObject.hpp>

Public Member Functions

 QueryObject (Teuchos::RCP< const Epetra_CrsGraph > graph, Teuchos::RCP< const Isorropia::Epetra::CostDescriber > costs, int inputType)
 Constructor. More...
 
 QueryObject (Teuchos::RCP< const Epetra_RowMatrix > matrix, Teuchos::RCP< const Isorropia::Epetra::CostDescriber > costs, int inputType)
 Constructor. More...
 
 QueryObject (Teuchos::RCP< const Epetra_MultiVector > coords, Teuchos::RCP< const Epetra_MultiVector > weights)
 Constructor. More...
 
 QueryObject (Teuchos::RCP< const Epetra_CrsGraph > graph, Teuchos::RCP< const Isorropia::Epetra::CostDescriber > costs, Teuchos::RCP< const Epetra_MultiVector > coords, Teuchos::RCP< const Epetra_MultiVector > weights, int inputType)
 Constructor. More...
 
 QueryObject (Teuchos::RCP< const Epetra_RowMatrix > matrix, Teuchos::RCP< const Isorropia::Epetra::CostDescriber > costs, Teuchos::RCP< const Epetra_MultiVector > coords, Teuchos::RCP< const Epetra_MultiVector > weights, int inputType)
 Constructor. More...
 
 QueryObject (Teuchos::RCP< const Epetra_BlockMap > inputMap, int inputType)
 Constructor. More...
 
virtual ~QueryObject ()
 Destructor. More...
 
const Epetra_BlockMap & RowMap (void)
 Return the map associated with the object to be partitioned. More...
 
bool haveVertexWeights ()
 Return true if any of the processes in the application have defined vertex weights. More...
 
bool haveGraphEdgeWeights ()
 Return true if any of the processes in the application have defined graph edge weights. More...
 
bool haveHypergraphEdgeWeights ()
 Return true if any of the processes in the application have defined hypergraph edge weights. More...
 

Static Public Member Functions

static int Number_Objects (void *data, int *ierr)
 The interface to a particular QueryObject's My_Number_Objects query function. More...
 
static void Object_List (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int weight_dim, float *object_weights, int *ierr)
 The interface to a particular QueryObject's My_Object_List query function. More...
 
static void Number_Edges_Multi (void *data, int num_gid_entries, int num_lid_entries, int num_obj, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int *num_edges, int *ierr)
 The interface to a particular QueryObject's My_Number_Edges_Multi query function. More...
 
static void Edge_List_Multi (void *data, int num_gid_entries, int num_lid_entries, int num_obj, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int *num_edges, ZOLTAN_ID_PTR neighbor_global_ids, int *neighbor_procs, int weight_dim, float *edge_weights, int *ierr)
 The interface to a particular QueryObject's My_Edges_Multi query function. More...
 
static void HG_Size_CS (void *data, int *num_lists, int *num_pins, int *format, int *ierr)
 The interface to a particular QueryObject's My_HG_Size_CS query function. More...
 
static void HG_CS (void *data, int num_gid_entries, int num_row_or_col, int num_pins, int format, ZOLTAN_ID_PTR vtxedge_GID, int *vtxedge_ptr, ZOLTAN_ID_PTR pin_GID, int *ierr)
 The interface to a particular QueryObject's My_HG_CS query function. More...
 
static void HG_Size_Edge_Weights (void *data, int *num_edges, int *ierr)
 The interface to a particular QueryObject's My_HG_Size_Edge_Weights query function. More...
 
static void HG_Edge_Weights (void *data, int num_gid_entries, int num_lid_entries, int num_edges, int edge_weight_dim, ZOLTAN_ID_PTR edge_GID, ZOLTAN_ID_PTR edge_LID, float *edge_weights, int *ierr)
 The interface to a particular QueryObject's My_HG_Edge_Weights query function. More...
 
static int Number_Geom (void *data, int *ierr)
 The interface to a particular QueryObject's My_Number_Geom query function. More...
 
static void Geom_Multi (void *data, int num_gid_entries, int num_lid_entries, int num_obj, ZOLTAN_ID_PTR gids, ZOLTAN_ID_PTR lids, int num_dim, double *geom_vec, int *ierr)
 The interface to a particular QueryObject's My_Geom_Multi query function. More...
 

Public Attributes

int input_type_
 The input_type_ indicates how the object to be partitioned is to be interpreted - as a graph or a hypergraph for row partitioning, as a hypergraph for fine-grain partitioning, or as a list of coordinates for geometric partitioning. More...
 

Static Public Attributes

static const int hgraph_input_ = 1
 input_type_ == hgraph_input_. More...
 
static const int hgraph2d_finegrain_input_ = 2
 input_type_ == hgraph2d_finegrain_input_. More...
 
static const int graph_input_ = 3
 input_type_ == graph_input_. More...
 
static const int geometric_input_ = 4
 input_type_ == geometric_input_. More...
 
static const int hgraph_graph_input_ = 5
 input_type_ == hgraph_graph_input_ This indicates that the Epetra_MultiVector represents a hypergraph and graph (see above). More...
 
static const int hgraph_geometric_input_ = 6
 input_type_ == hgraph_geometric_input_ This indicates that the Epetra_MultiVector represents a hypergraph and graph (see above). More...
 
static const int graph_geometric_input_ = 7
 input_type_ == graph_geometric_input_ This indicates that the Epetra_MultiVector represents graph and has geometric coordinates (see above). More...
 
static const int hgraph_graph_geometric_input_ = 8
 input_type_ == hgraph_graph_geometric_input_ This indicates that the Epetra_MultiVector represents a hypergraph and graph and has geometric coordinates(see above). More...
 
static const int simple_input_ = 9
 input_type_ == simple_input_ This indicates that a simple method (block, cyclic, or random) will be used. More...
 
static const int unspecified_input_ = 10
 input_type_ == unspecified_input_. More...
 

Private Member Functions

void fill_procmap ()
 
int My_Number_Objects (int *ierr)
 My_Number_Objects() returns the number of objects currently assigned to this process. More...
 
void My_Object_List (int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int weight_dim, float *object_weights, int *ierr)
 My_ObjectList() returns to Zoltan the global ID and weight of the objects currently assigned to this process. More...
 
void My_Number_Edges_Multi (int num_gid_entries, int num_lid_entries, int num_obj, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int *num_edges, int *ierr)
 My_Number_Edges_Multi() is a query function used for graph partitioning only. More...
 
void My_Edge_List_Multi (int num_gid_entries, int num_lid_entries, int num_obj, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int *num_edges, ZOLTAN_ID_PTR neighbor_global_ids, int *neighbor_procs, int weight_dim, float *edge_weights, int *ierr)
 My_Edge_List_Multi() is a query function used for graph partitioning only. More...
 
void My_HG_Size_CS (int *num_lists, int *num_pins, int *format, int *ierr)
 My_HG_Size_CS() is a query function used for hypergraph partitioning only. More...
 
void My_HG_CS (int num_gid_entries, int num_row_or_col, int num_pins, int format, ZOLTAN_ID_PTR vtxedge_GID, int *vtxedge_ptr, ZOLTAN_ID_PTR pin_GID, int *ierr)
 My_HG_CS() is a query function used for hypergraph partitioning only. More...
 
void My_FGHG_CS (int num_gid_entries, int num_row_or_col, int num_pins, int format, ZOLTAN_ID_PTR vtxedge_GID, int *vtxedge_ptr, ZOLTAN_ID_PTR pin_GID, int *ierr)
 My_FGHG_CS() is a query function used for fine-grain hypergraph partitioning only. More...
 
void My_HG_Size_Edge_Weights (int *num_edges, int *ierr)
 My_HG_Size_Edge_Weights() is a query function used for hypergraph partitioning only. More...
 
void My_HG_Edge_Weights (int num_gid_entries, int num_lid_entries, int num_edges, int edge_weight_dim, ZOLTAN_ID_PTR edge_GID, ZOLTAN_ID_PTR edge_LID, float *edge_weights, int *ierr)
 My_HG_Edge_Weights() is a query function used for hypergraph partitioning only. More...
 
int My_Number_Geom (int *ierr)
 My_Number_Geom() is a query function used for geometric partitioning only. More...
 
void My_Geom_Multi (int num_gid_entries, int num_lid_entries, int num_obj, ZOLTAN_ID_PTR gids, ZOLTAN_ID_PTR lids, int num_dim, double *geom_vec, int *ierr)
 My_Geom_Multi() is a query function used for geometric partitioning only. More...
 

Private Attributes

const bool haveGraph_
 haveGraph is true if we have CrsGraph, and not a CrsMatrix or a MultiVector. More...
 
Teuchos::RCP< const
Epetra_CrsGraph > 
graph_
 The CrsGraph. More...
 
Teuchos::RCP< const
Epetra_RowMatrix > 
matrix_
 The CrsMatrix. More...
 
Teuchos::RCP< const
Epetra_MultiVector > 
coords_
 The MultiVector containing 1, 2 or 3 dimensional coordinates. More...
 
const Epetra_BlockMap * rowMap_
 The graph or matrix row map, or the MultiVector map. More...
 
const Epetra_BlockMap * colMap_
 The graph or matrix column map. More...
 
Teuchos::RCP< const
Isorropia::Epetra::CostDescriber
costs_
 The CostDescriber contains optional vertex and/or edge weights for graph and hypergraph partitioning. More...
 
Teuchos::RCP< const
Epetra_MultiVector > 
weights_
 The MultiVector contains optional object (point) weights for geometric partitioning. More...
 
std::map< int, int > procmap_
 
std::set< int > graph_self_edges_
 
unsigned int myProc_
 
unsigned int base_
 

Detailed Description

QueryObject is a class that contains the query functions required by the Zoltan library.

These methods are not part of the Isorropia API (except to Zoltan). They are called by Isorropia itself and by Zoltan.

For a better understanding of Zoltan's query functions, see the Zoltan User's Guide at http://www.cs.sandia.gov/zoltan

Constructor & Destructor Documentation

Isorropia::Epetra::ZoltanLib::QueryObject::QueryObject ( Teuchos::RCP< const Epetra_CrsGraph >  graph,
Teuchos::RCP< const Isorropia::Epetra::CostDescriber costs,
int  inputType 
)

Constructor.

Isorropia::Epetra::ZoltanLib::QueryObject::QueryObject ( Teuchos::RCP< const Epetra_RowMatrix >  matrix,
Teuchos::RCP< const Isorropia::Epetra::CostDescriber costs,
int  inputType 
)

Constructor.

Isorropia::Epetra::ZoltanLib::QueryObject::QueryObject ( Teuchos::RCP< const Epetra_MultiVector >  coords,
Teuchos::RCP< const Epetra_MultiVector >  weights 
)

Constructor.

Isorropia::Epetra::ZoltanLib::QueryObject::QueryObject ( Teuchos::RCP< const Epetra_CrsGraph >  graph,
Teuchos::RCP< const Isorropia::Epetra::CostDescriber costs,
Teuchos::RCP< const Epetra_MultiVector >  coords,
Teuchos::RCP< const Epetra_MultiVector >  weights,
int  inputType 
)

Constructor.

Isorropia::Epetra::ZoltanLib::QueryObject::QueryObject ( Teuchos::RCP< const Epetra_RowMatrix >  matrix,
Teuchos::RCP< const Isorropia::Epetra::CostDescriber costs,
Teuchos::RCP< const Epetra_MultiVector >  coords,
Teuchos::RCP< const Epetra_MultiVector >  weights,
int  inputType 
)

Constructor.

Isorropia::Epetra::ZoltanLib::QueryObject::QueryObject ( Teuchos::RCP< const Epetra_BlockMap >  inputMap,
int  inputType 
)

Constructor.

virtual Isorropia::Epetra::ZoltanLib::QueryObject::~QueryObject ( )
virtual

Destructor.

Member Function Documentation

void Isorropia::Epetra::ZoltanLib::QueryObject::fill_procmap ( )
private
int Isorropia::Epetra::ZoltanLib::QueryObject::My_Number_Objects ( int *  ierr)
private

My_Number_Objects() returns the number of objects currently assigned to this process.

(The objects are interpreted as graph vertices for Graph partitioning, as hypergraph vertices for hypergraph partitioning, or as coordinates for geometric partitioning.)

void Isorropia::Epetra::ZoltanLib::QueryObject::My_Object_List ( int  num_gid_entries,
int  num_lid_entries,
ZOLTAN_ID_PTR  global_ids,
ZOLTAN_ID_PTR  local_ids,
int  weight_dim,
float *  object_weights,
int *  ierr 
)
private

My_ObjectList() returns to Zoltan the global ID and weight of the objects currently assigned to this process.

void Isorropia::Epetra::ZoltanLib::QueryObject::My_Number_Edges_Multi ( int  num_gid_entries,
int  num_lid_entries,
int  num_obj,
ZOLTAN_ID_PTR  global_ids,
ZOLTAN_ID_PTR  local_ids,
int *  num_edges,
int *  ierr 
)
private

My_Number_Edges_Multi() is a query function used for graph partitioning only.

It returns to Zoltan the number of edges (non-zeroes) that each vertex (row) has.

void Isorropia::Epetra::ZoltanLib::QueryObject::My_Edge_List_Multi ( int  num_gid_entries,
int  num_lid_entries,
int  num_obj,
ZOLTAN_ID_PTR  global_ids,
ZOLTAN_ID_PTR  local_ids,
int *  num_edges,
ZOLTAN_ID_PTR  neighbor_global_ids,
int *  neighbor_procs,
int  weight_dim,
float *  edge_weights,
int *  ierr 
)
private

My_Edge_List_Multi() is a query function used for graph partitioning only.

For each vertex (row), it returns a list of the global ID of each neighbor (non-zero) and the process owning that neighbor (that row).

void Isorropia::Epetra::ZoltanLib::QueryObject::My_HG_Size_CS ( int *  num_lists,
int *  num_pins,
int *  format,
int *  ierr 
)
private

My_HG_Size_CS() is a query function used for hypergraph partitioning only.

Zoltan calls this query to get size of the non-zeros lists from the QueryObject.

void Isorropia::Epetra::ZoltanLib::QueryObject::My_HG_CS ( int  num_gid_entries,
int  num_row_or_col,
int  num_pins,
int  format,
ZOLTAN_ID_PTR  vtxedge_GID,
int *  vtxedge_ptr,
ZOLTAN_ID_PTR  pin_GID,
int *  ierr 
)
private

My_HG_CS() is a query function used for hypergraph partitioning only.

Zoltan calls this query to get the non-zeros lists from the QueryObject.

void Isorropia::Epetra::ZoltanLib::QueryObject::My_FGHG_CS ( int  num_gid_entries,
int  num_row_or_col,
int  num_pins,
int  format,
ZOLTAN_ID_PTR  vtxedge_GID,
int *  vtxedge_ptr,
ZOLTAN_ID_PTR  pin_GID,
int *  ierr 
)
private

My_FGHG_CS() is a query function used for fine-grain hypergraph partitioning only.

Zoltan calls this query to get the non-zeros lists from the QueryObject.

void Isorropia::Epetra::ZoltanLib::QueryObject::My_HG_Size_Edge_Weights ( int *  num_edges,
int *  ierr 
)
private

My_HG_Size_Edge_Weights() is a query function used for hypergraph partitioning only.

Zoltan calls this query to get number of hyperedge weights that this QueryObject will be providing.

void Isorropia::Epetra::ZoltanLib::QueryObject::My_HG_Edge_Weights ( int  num_gid_entries,
int  num_lid_entries,
int  num_edges,
int  edge_weight_dim,
ZOLTAN_ID_PTR  edge_GID,
ZOLTAN_ID_PTR  edge_LID,
float *  edge_weights,
int *  ierr 
)
private

My_HG_Edge_Weights() is a query function used for hypergraph partitioning only.

Zoltan calls this query to get hyperedge weights from this QueryObject.

int Isorropia::Epetra::ZoltanLib::QueryObject::My_Number_Geom ( int *  ierr)
private

My_Number_Geom() is a query function used for geometric partitioning only.

Zoltan calls this query to get the dimension of the geometric coordinates.

void Isorropia::Epetra::ZoltanLib::QueryObject::My_Geom_Multi ( int  num_gid_entries,
int  num_lid_entries,
int  num_obj,
ZOLTAN_ID_PTR  gids,
ZOLTAN_ID_PTR  lids,
int  num_dim,
double *  geom_vec,
int *  ierr 
)
private

My_Geom_Multi() is a query function used for geometric partitioning only.

Zoltan calls this query to get a list of coordinates from the QueryObject.

const Epetra_BlockMap& Isorropia::Epetra::ZoltanLib::QueryObject::RowMap ( void  )
inline

Return the map associated with the object to be partitioned.

bool Isorropia::Epetra::ZoltanLib::QueryObject::haveVertexWeights ( )

Return true if any of the processes in the application have defined vertex weights.

bool Isorropia::Epetra::ZoltanLib::QueryObject::haveGraphEdgeWeights ( )

Return true if any of the processes in the application have defined graph edge weights.

bool Isorropia::Epetra::ZoltanLib::QueryObject::haveHypergraphEdgeWeights ( )

Return true if any of the processes in the application have defined hypergraph edge weights.

static int Isorropia::Epetra::ZoltanLib::QueryObject::Number_Objects ( void *  data,
int *  ierr 
)
static

The interface to a particular QueryObject's My_Number_Objects query function.

static void Isorropia::Epetra::ZoltanLib::QueryObject::Object_List ( void *  data,
int  num_gid_entries,
int  num_lid_entries,
ZOLTAN_ID_PTR  global_ids,
ZOLTAN_ID_PTR  local_ids,
int  weight_dim,
float *  object_weights,
int *  ierr 
)
static

The interface to a particular QueryObject's My_Object_List query function.

static void Isorropia::Epetra::ZoltanLib::QueryObject::Number_Edges_Multi ( void *  data,
int  num_gid_entries,
int  num_lid_entries,
int  num_obj,
ZOLTAN_ID_PTR  global_ids,
ZOLTAN_ID_PTR  local_ids,
int *  num_edges,
int *  ierr 
)
static

The interface to a particular QueryObject's My_Number_Edges_Multi query function.

static void Isorropia::Epetra::ZoltanLib::QueryObject::Edge_List_Multi ( void *  data,
int  num_gid_entries,
int  num_lid_entries,
int  num_obj,
ZOLTAN_ID_PTR  global_ids,
ZOLTAN_ID_PTR  local_ids,
int *  num_edges,
ZOLTAN_ID_PTR  neighbor_global_ids,
int *  neighbor_procs,
int  weight_dim,
float *  edge_weights,
int *  ierr 
)
static

The interface to a particular QueryObject's My_Edges_Multi query function.

static void Isorropia::Epetra::ZoltanLib::QueryObject::HG_Size_CS ( void *  data,
int *  num_lists,
int *  num_pins,
int *  format,
int *  ierr 
)
static

The interface to a particular QueryObject's My_HG_Size_CS query function.

static void Isorropia::Epetra::ZoltanLib::QueryObject::HG_CS ( void *  data,
int  num_gid_entries,
int  num_row_or_col,
int  num_pins,
int  format,
ZOLTAN_ID_PTR  vtxedge_GID,
int *  vtxedge_ptr,
ZOLTAN_ID_PTR  pin_GID,
int *  ierr 
)
static

The interface to a particular QueryObject's My_HG_CS query function.

static void Isorropia::Epetra::ZoltanLib::QueryObject::HG_Size_Edge_Weights ( void *  data,
int *  num_edges,
int *  ierr 
)
static

The interface to a particular QueryObject's My_HG_Size_Edge_Weights query function.

static void Isorropia::Epetra::ZoltanLib::QueryObject::HG_Edge_Weights ( void *  data,
int  num_gid_entries,
int  num_lid_entries,
int  num_edges,
int  edge_weight_dim,
ZOLTAN_ID_PTR  edge_GID,
ZOLTAN_ID_PTR  edge_LID,
float *  edge_weights,
int *  ierr 
)
static

The interface to a particular QueryObject's My_HG_Edge_Weights query function.

static int Isorropia::Epetra::ZoltanLib::QueryObject::Number_Geom ( void *  data,
int *  ierr 
)
static

The interface to a particular QueryObject's My_Number_Geom query function.

static void Isorropia::Epetra::ZoltanLib::QueryObject::Geom_Multi ( void *  data,
int  num_gid_entries,
int  num_lid_entries,
int  num_obj,
ZOLTAN_ID_PTR  gids,
ZOLTAN_ID_PTR  lids,
int  num_dim,
double *  geom_vec,
int *  ierr 
)
static

The interface to a particular QueryObject's My_Geom_Multi query function.

Member Data Documentation

const bool Isorropia::Epetra::ZoltanLib::QueryObject::haveGraph_
private

haveGraph is true if we have CrsGraph, and not a CrsMatrix or a MultiVector.

Teuchos::RCP<const Epetra_CrsGraph> Isorropia::Epetra::ZoltanLib::QueryObject::graph_
private

The CrsGraph.

The QueryObject must be constructed with one of an Epetra_CrsGraph, an Epetra_RowMatrix or an Epetra_MultiVector.

Teuchos::RCP<const Epetra_RowMatrix> Isorropia::Epetra::ZoltanLib::QueryObject::matrix_
private

The CrsMatrix.

The QueryObject must be constructed with one of an Epetra_CrsGraph, an Epetra_RowMatrix or an Epetra_MultiVector.

Teuchos::RCP<const Epetra_MultiVector> Isorropia::Epetra::ZoltanLib::QueryObject::coords_
private

The MultiVector containing 1, 2 or 3 dimensional coordinates.

If supplied, we will perform geometric partitioning. The QueryObject must be constructed with one of an Epetra_CrsGraph, an Epetra_RowMatrix or an Epetra_MultiVector.

const Epetra_BlockMap* Isorropia::Epetra::ZoltanLib::QueryObject::rowMap_
private

The graph or matrix row map, or the MultiVector map.

const Epetra_BlockMap* Isorropia::Epetra::ZoltanLib::QueryObject::colMap_
private

The graph or matrix column map.

Teuchos::RCP<const Isorropia::Epetra::CostDescriber> Isorropia::Epetra::ZoltanLib::QueryObject::costs_
private

The CostDescriber contains optional vertex and/or edge weights for graph and hypergraph partitioning.

Teuchos::RCP<const Epetra_MultiVector> Isorropia::Epetra::ZoltanLib::QueryObject::weights_
private

The MultiVector contains optional object (point) weights for geometric partitioning.

Zoltan currently will use only the weights in the first vector (1 dimensional point weights).

std::map<int,int> Isorropia::Epetra::ZoltanLib::QueryObject::procmap_
private
std::set<int> Isorropia::Epetra::ZoltanLib::QueryObject::graph_self_edges_
private
unsigned int Isorropia::Epetra::ZoltanLib::QueryObject::myProc_
private
unsigned int Isorropia::Epetra::ZoltanLib::QueryObject::base_
private
const int Isorropia::Epetra::ZoltanLib::QueryObject::hgraph_input_ = 1
static

input_type_ == hgraph_input_.

This indicates that the matrix or graph represents a hypergraph. Columns represent hyperedges, and row (vertex) partitioning is to be performed.

const int Isorropia::Epetra::ZoltanLib::QueryObject::hgraph2d_finegrain_input_ = 2
static

input_type_ == hgraph2d_finegrain_input_.

This indicates that the matrix or graph represents a hypergraph. Columns represent hyperedges, and non-zeroes are to be partitioned.

const int Isorropia::Epetra::ZoltanLib::QueryObject::graph_input_ = 3
static

input_type_ == graph_input_.

This indicates that the square symmetric matrix or graph represents a graph in the sense that row/column IDs are vertices and non-zeroes represent edges. The vertices are to be partitioned.

const int Isorropia::Epetra::ZoltanLib::QueryObject::geometric_input_ = 4
static

input_type_ == geometric_input_.

This indicates that the Epetra_MultiVector represents geometric coordinates. The MultiVector should have 1, 2 or 3 vectors, representing 1, 2 or 3 dimensional coordinates. The coordinates are to be partitioned.

const int Isorropia::Epetra::ZoltanLib::QueryObject::hgraph_graph_input_ = 5
static

input_type_ == hgraph_graph_input_ This indicates that the Epetra_MultiVector represents a hypergraph and graph (see above).

This is necessary for hierarchical partitioning with both hypergraph and graph methods.

const int Isorropia::Epetra::ZoltanLib::QueryObject::hgraph_geometric_input_ = 6
static

input_type_ == hgraph_geometric_input_ This indicates that the Epetra_MultiVector represents a hypergraph and graph (see above).

This is necessary for hierarchical partitioning with both hypergraph and geometric methods.

const int Isorropia::Epetra::ZoltanLib::QueryObject::graph_geometric_input_ = 7
static

input_type_ == graph_geometric_input_ This indicates that the Epetra_MultiVector represents graph and has geometric coordinates (see above).

This is necessary for hierarchical partitioning with both graph and geometric methods.

const int Isorropia::Epetra::ZoltanLib::QueryObject::hgraph_graph_geometric_input_ = 8
static

input_type_ == hgraph_graph_geometric_input_ This indicates that the Epetra_MultiVector represents a hypergraph and graph and has geometric coordinates(see above).

This is necessary for hierarchical partitioning using hypergraph, graph, and geometric methods.

const int Isorropia::Epetra::ZoltanLib::QueryObject::simple_input_ = 9
static

input_type_ == simple_input_ This indicates that a simple method (block, cyclic, or random) will be used.

const int Isorropia::Epetra::ZoltanLib::QueryObject::unspecified_input_ = 10
static

input_type_ == unspecified_input_.

This value is the "unset" state for the input_type_ instance variable.

int Isorropia::Epetra::ZoltanLib::QueryObject::input_type_

The input_type_ indicates how the object to be partitioned is to be interpreted - as a graph or a hypergraph for row partitioning, as a hypergraph for fine-grain partitioning, or as a list of coordinates for geometric partitioning.


The documentation for this class was generated from the following file: