Zoltan User's Guide  |  Next  |  Previous

General Zoltan Query Functions

The following registered functions are used by various Zoltan algorithms in the Zoltan library. No single algorithm uses all the query functions; the algorithm descriptions indicate which query functions are required by individual algorithms.
Object ID Functions
ZOLTAN_NUM_OBJ_FN
ZOLTAN_OBJ_LIST_FN
ZOLTAN_FIRST_OBJ_FN
ZOLTAN_NEXT_OBJ_FN
ZOLTAN_PARTITION_MULTI_FN or ZOLTAN_PARTITION_FN
Geometry-Based Functions
ZOLTAN_NUM_GEOM_FN
ZOLTAN_GEOM_MULTI_FN or ZOLTAN_GEOM_FN
Graph-Based Functions
ZOLTAN_NUM_EDGES_MULTI_FN or ZOLTAN_NUM_EDGES_FN
ZOLTAN_EDGE_LIST_MULTI_FN or ZOLTAN_EDGE_LIST_FN
Hypergraph-Based Functions
ZOLTAN_HG_SIZE_CS_FN
ZOLTAN_HG_CS_FN

ZOLTAN_HG_SIZE_EDGE_WTS_FN

ZOLTAN_HG_EDGE_WTS_FN
Tree-Based Functions
ZOLTAN_NUM_COARSE_OBJ_FN
ZOLTAN_COARSE_OBJ_LIST_FN
ZOLTAN_FIRST_COARSE_OBJ_FN
ZOLTAN_NEXT_COARSE_OBJ_FN
ZOLTAN_NUM_CHILD_FN
ZOLTAN_CHILD_LIST_FN
ZOLTAN_CHILD_WEIGHT_FN
Border Object Functions (currently unused)
ZOLTAN_NUM_BORDER_OBJ_FN
ZOLTAN_BORDER_OBJ_LIST_FN
ZOLTAN_FIRST_BORDER_OBJ_FN
ZOLTAN_NEXT_BORDER_OBJ_FN


Object ID Functions



C and C++: typedef int ZOLTAN_NUM_OBJ_FN (void *data, int *ierr);
FORTRAN: FUNCTION Get_Num_Obj(data, ierr
INTEGER(Zoltan_INT) :: Get_Num_Obj 
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_NUM_OBJ_FN query function returns the number of objects that are currently assigned to the processor.
 
Function Type: ZOLTAN_NUM_OBJ_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    ierr Error code to be set by function.
Returned Value:
    int The number of objects that are assigned to the processor.



C and C++: typedef void ZOLTAN_OBJ_LIST_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int wgt_dim, float *obj_wgts, int *ierr); 
FORTRAN: SUBROUTINE Get_Obj_List(data, num_gid_entries, num_lid_entries, global_ids, local_ids, wgt_dim, obj_wgts, ierr
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: global_ids 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: local_ids 
INTEGER(Zoltan_INT), INTENT(IN) :: wgt_dim 
REAL(Zoltan_FLOAT), INTENT(OUT), DIMENSION(*) :: obj_wgts 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_OBJ_LIST_FN query function fills two (three if weights are used) arrays with information about the objects currently assigned to the processor. Both arrays are allocated (and subsequently freed) by Zoltan; their size is determined by a call to a ZOLTAN_NUM_OBJ_FN query function to get the array size. For many algorithms, either a ZOLTAN_OBJ_LIST_FN query function or a ZOLTAN_FIRST_OBJ_FN/ZOLTAN_NEXT_OBJ_FN query-function pair must be registered; however, both query options need not be provided.
 
Function Type: ZOLTAN_OBJ_LIST_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    global_ids Upon return, an array of unique global IDs for all objects assigned to the processor.
    local_ids Upon return, an array of local IDs, the meaning of which can be determined by the application, for all objects assigned to the processor.
    wgt_dim The number of weights associated with an object (typically 1), or 0 if weights are not requested. This value is set through the parameter OBJ_WEIGHT_DIM.
    obj_wgts Upon return, an array of object weights. Weights for object i are stored in obj_wgts[(i-1)*wgt_dim:i*wgt_dim-1].  If wgt_dim=0, the return value of obj_wgts is undefined and may be NULL.
    ierr Error code to be set by function.



C and C++: typedef int ZOLTAN_FIRST_OBJ_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR first_global_id, ZOLTAN_ID_PTR first_local_id, int wgt_dim, float *first_obj_wgt, int *ierr); 
FORTRAN: FUNCTION Get_First_Obj(data, num_gid_entries, num_lid_entries, first_global_id, first_local_id, wgt_dim, first_obj_wgt, ierr
INTEGER(Zoltan_INT) :: Get_First_Obj 
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: first_global_id 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: first_local_id 
INTEGER(Zoltan_INT), INTENT(IN) :: wgt_dim 
REAL(Zoltan_FLOAT), INTENT(OUT), DIMENSION(*) :: first_obj_wgt 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_FIRST_OBJ_FN query function initializes an iteration over objects assigned to the processor. It returns the global and local IDs of the first object on the processor. Subsequent calls to a ZOLTAN_NEXT_OBJ_FN query function iterate over and return other objects assigned to the processor. This query-function pair frees the application from having to build an array of objects (as in ZOLTAN_OBJ_LIST_FN) and allows Zoltan's routines to obtain only as much information about objects as they need. For many algorithms, either a ZOLTAN_OBJ_LIST_FN query function or a ZOLTAN_FIRST_OBJ_FN/ZOLTAN_NEXT_OBJ_FN query-function pair must be registered; however, both query options need not be provided.
 
Function Type: ZOLTAN_FIRST_OBJ_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    first_global_id The returned value of the global ID for the first object; the value is ignored if there are no objects.
    first_local_id The returned value of the local ID for the first object; the value is ignored if there are no objects.
    wgt_dim The number of weights associated with an object (typically 1), or 0 if weights are not requested. This value is set through the parameter OBJ_WEIGHT_DIM.
    first_obj_wgt Upon return, the first object's weights; an array of length wgt_dim. Undefined if wgt_dim=0.
    ierr Error code to be set by function.
Returned Value:
    1 If first_global_id and first_local_id contain valid IDs of the first object.
    0 If no objects are available.



C and C++: typedef int ZOLTAN_NEXT_OBJ_FN (void * data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_id, ZOLTAN_ID_PTR local_id, ZOLTAN_ID_PTR next_global_id, ZOLTAN_ID_PTR next_local_id, int wgt_dim, float *next_obj_wgt, int *ierr); 
FORTRAN: FUNCTION Get_Next_Obj(data, num_gid_entries, num_lid_entries, global_id, local_id, next_global_id, next_local_id, wgt_dim, next_obj_wgt, ierr
INTEGER(Zoltan_INT) :: Get_Next_Obj 
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_id 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_id 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: next_global_id 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: next_local_id 
INTEGER(Zoltan_INT), INTENT(IN) :: wgt_dim 
REAL(Zoltan_FLOAT), INTENT(OUT), DIMENSION(*) :: next_obj_wgt 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_NEXT_OBJ_FN query function is an iterator function which, when given an object assigned to the processor, returns the next object assigned to the processor. The first object of the iteration is provided by a ZOLTAN_FIRST_OBJ_FN query function. This query-function pair frees the application from having to build an array of objects (as in ZOLTAN_OBJ_LIST_FN) and allows Zoltan's routines to obtain only as much information about objects as they need. For many algorithms, either a ZOLTAN_OBJ_LIST_FN query function or a ZOLTAN_FIRST_OBJ_FN/ZOLTAN_NEXT_OBJ_FN query-function pair must be registered; however, both query options need not be provided.
 
Function Type: ZOLTAN_NEXT_OBJ_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    global_id The global ID of the previous object.
    local_id The local ID of the previous object.
    next_global_id The returned value of the global ID for the next object; the value is ignored if there are no more objects.
    next_local_id The returned value of the local ID for the next object; the value is ignored if there are no more objects.
    wgt_dim The number of weights associated with an object (typically 1), or 0 if weights are not requested. This value is set through the parameter OBJ_WEIGHT_DIM.
    next_obj_wgt Upon return, the next object's weights; an array of length wgt_dim. Undefined if wgt_dim=0.
    ierr Error code to be set by function.
Returned Value:
    1 If next_global_id and next_local_id contain valid IDs of the next object.
    0 If no more objects are available.



C and C++: typedef void ZOLTAN_PARTITION_MULTI_FN (void *data, int num_gid_entries, int num_lid_entries, int num_obj, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int *parts, int *ierr); 
FORTRAN: SUBROUTINE Get_Partition_Multi(data, num_gid_entries, num_lid_entries, num_obj, global_ids, local_ids, ierr)
<type-data>, INTENT(IN) :: data
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries, num_obj
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_ids
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_ids
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: parts
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation.


A ZOLTAN_PARTITION_MULTI_FN query function returns a list of partitions to which given objects are currently assigned. If a ZOLTAN_PARTITION_MULTI_FN or ZOLTAN_PARTITION_FN is not registered, Zoltan assumes the partition numbers are the processor number of the owning processor. Valid partition numbers are non-negative integers.
 
Function Type: ZOLTAN_PARTITION_MULTI_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    num_obj The number of object IDs in arrays global_ids and local_ids.
    global_ids The global IDs of the objects for which the partition numbers should be returned.
    local_ids The local IDs of the objects for which the partition numbers should be returned.
    parts Upon return, an array of partition numbers corresponding to the global and local IDs.
    ierr Error code to be set by function.





C and C++: typedef int ZOLTAN_PARTITION_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_id, ZOLTAN_ID_PTR local_id, int *ierr); 
FORTRAN: FUNCTION Get_Partition(data, num_gid_entries, num_lid_entries, global_id, local_id, ierr
INTEGER(Zoltan_INT) :: Get_Partition 
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_id 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_id 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_PARTITION_FN query function returns the partition to which a given object is currently assigned. If a ZOLTAN_PARTITION_FN or ZOLTAN_PARTITION_MULTI_FN is not registered, Zoltan assumes the partition numbers are the processor number of the owning processor. Valid partition numbers are non-negative integers.
 
Function Type: ZOLTAN_PARTITION_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    global_id The global ID of the object for which the partition number should be returned.
    local_id The local ID of the object for which the partition number should be returned.
    ierr Error code to be set by function.
Returned Value:
    int The partition number for the object identified by global_id and local_id.



Geometry-based Functions



C and C++: typedef int ZOLTAN_NUM_GEOM_FN (void *data, int *ierr); 
FORTRAN: FUNCTION Get_Num_Geom(data, ierr
INTEGER(Zoltan_INT) :: Get_Num_Geom 
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_NUM_GEOM_FN query function returns the number of values needed to express the geometry of an object. For example, for a two-dimensional mesh-based application, (x,y) coordinates are needed to describe an object's geometry; thus the ZOLTAN_NUM_GEOM_FN query function should return the value of two. For a similar three-dimensional application, the return value should be three.
 
Function Type: ZOLTAN_NUM_GEOM_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    ierr Error code to be set by function.
Returned Value:
    int The number of values needed to express the geometry of an object.



C and C++: typedef void ZOLTAN_GEOM_MULTI_FN (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_dim, double *geom_vec, int *ierr); 
FORTRAN: SUBROUTINE Get_Geom_Multi(data, num_gid_entries, num_lid_entries, num_obj, global_ids, local_ids, num_dim, geom_vec, ierr)
<type-data>, INTENT(IN) :: data
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries
INTEGER(Zoltan_INT), INTENT(IN) :: num_obj, num_dim
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_ids
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_ids
REAL(Zoltan_DOUBLE), INTENT(OUT), DIMENSION(*) :: geom_vec
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_GEOM_MULTI FN query function returns a vector of geometry values for a list of given objects. The geometry vector is allocated by Zoltan to be of size num_obj * num_dim; its format is described below.
 
Function Type: ZOLTAN_GEOM_MULTI_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    num_obj The number of object IDs in arrays global_ids and local_ids.
    global_ids Array of global IDs of objects whose geometry values should be returned.
    local_ids Array of local IDs of objects whose geometry values should be returned.
    num_dim Number of coordinate entries per object (typically 1, 2, or 3).
    geom_vec Upon return, an array containing geometry values. For object i (specified by global_ids[i*num_gid_entries] and local_ids[i*num_lid_entries], i=0,1,...,num_obj-1), coordinate values should be stored in geom_vec[i*num_dim:(i+1)*num_dim-1].
    ierr Error code to be set by function.

 

C and C++: typedef void ZOLTAN_GEOM_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_id, ZOLTAN_ID_PTR local_id, double *geom_vec, int *ierr); 
FORTRAN: SUBROUTINE Get_Geom(data, num_gid_entries, num_lid_entries, global_id, local_id, geom_vec, ierr
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_id 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_id 
REAL(Zoltan_DOUBLE), INTENT(OUT), DIMENSION(*) :: geom_vec 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_GEOM_FN query function returns a vector of geometry values for a given object. The geometry vector is allocated by Zoltan to be of the size returned by a ZOLTAN_NUM_GEOM_FN query function.
 
Function Type: ZOLTAN_GEOM_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    global_id The global ID of the object whose geometry values should be returned.
    local_id The local ID of the object whose geometry values should be returned.
    geom_vec Upon return, an array containing geometry values.
    ierr Error code to be set by function.

 

Graph-based Functions



C and C++: typedef void ZOLTAN_NUM_EDGES_MULTI_FN (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); 
FORTRAN: SUBROUTINE Get_Num_Edges_Multi(data, num_gid_entries, num_lid_entries, num_obj, global_ids, local_ids, num_edges, ierr
INTEGER(Zoltan_INT) :: Get_Num_Edges 
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries, num_obj
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_ids
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_ids
INTEGER(Zoltan_INT), INTENT(OUT),DIMENSION(*) :: num_edges
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_NUM_EDGES_MULTI_FN query function returns the number of edges in the communication graph of the application for each object in a list of objects. That is, for each object in the global_ids/local_ids arrays, the number of objects with which the given object must share information is returned.
 
Function Type: ZOLTAN_NUM_EDGES_MULTI_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    num_obj The number of object IDs in arrays global_ids and local_ids.
    global_ids Array of global IDs of objects whose number of edges should be returned.
    local_ids Array of local IDs of objects whose number of edges should be returned.
    num_edges Upon return, an array containing numbers of edges. For object i (specified by global_ids[i*num_gid_entries] and local_ids[i*num_lid_entries], i=0,1,...,num_obj-1), the number of edges should be stored in num_edges[i].
    ierr Error code to be set by function.



C and C++: typedef int ZOLTAN_NUM_EDGES_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_id, ZOLTAN_ID_PTR local_id, int *ierr); 
FORTRAN: FUNCTION Get_Num_Edges(data, num_gid_entries, num_lid_entries, global_id, local_id, ierr
INTEGER(Zoltan_INT) :: Get_Num_Edges 
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_id 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_id 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_NUM_EDGES_FN query function returns the number of edges for a given object in the communication graph of the application (i.e., the number of objects with which the given object must share information).
 
Function Type: ZOLTAN_NUM_EDGES_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    global_id The global ID of the object for which the number of edges should be returned.
    local_id The local ID of the object for which the number of edges should be returned.
    ierr Error code to be set by function.
Returned Value:
    int The number of edges for the object identified by global_id and local_id.



C and C++: typedef void ZOLTAN_EDGE_LIST_MULTI_FN (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 nbor_global_id, int *nbor_procs, int wgt_dim, float *ewgts, int *ierr); 
FORTRAN: SUBROUTINE Get_Edge_List_Multi(data, num_gid_entries, num_lid_entries, num_obj, global_ids, local_ids, num_edges, nbor_global_id, nbor_procs, wgt_dim, ewgts, ierr
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries, num_obj
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_ids
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_ids
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: num_edges
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: nbor_global_id
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: nbor_procs
INTEGER(Zoltan_INT), INTENT(IN) :: wgt_dim
REAL(Zoltan_FLOAT), INTENT(OUT), DIMENSION(*) :: ewgts
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_EDGE_LIST_MULTI_FN query function returns lists of global IDs, processor IDs, and optionally edge weights for objects sharing edges with objects specified in the global_ids input array; objects share edges when they must share information with other objects. The arrays for the returned neighbor lists are allocated by Zoltan; their size is determined by a calls to ZOLTAN_NUM_EDGES_MULTI_FN or ZOLTAN_NUM_EDGES_FN query functions.
 
Function Type: ZOLTAN_EDGE_LIST_MULTI_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    num_obj The number of object IDs in arrays global_ids and local_ids.
    global_ids Array of global IDs of objects whose edge lists should be returned.
    local_ids Array of local IDs of objects whose edge lists should be returned.
    num_edges An array containing numbers of edges for each object in global_ids. For object i (specified by global_ids[i*num_gid_entries] and local_ids[i*num_lid_entries], i=0,1,...,num_obj-1), the number of edges is stored in num_edges[i].
    nbor_global_id Upon return, an array of global IDs of objects sharing edges with the objects specified in global_ids. For object i (specified by global_ids[i*num_gid_entries] and local_ids[i*num_lid_entries], i=0,1,...,num_obj-1), edges are stored in nbor_global_id[sum*num_gid_entries] to nbor_global_id[(sum+num_edges[i])*num_gid_entries-1], where sum = the sum of num_edges[j] for j=0,1,...,i-1.
    nbor_procs Upon return, an array of processor IDs that identifies where the neighboring objects reside. For neighboring object i (stored in nbor_global_id[i*num_gid_entries]), the processor owning the neighbor is stored in nbor_procs[i].
    wgt_dim The number of weights associated with an edge (typically 1), or 0 if edge weights are not requested. This value is set through the parameter EDGE_WEIGHT_DIM.
    ewgts Upon return, an array of edge weights, where ewgts[i*wgt_dim:(i+1)*wgt_dim-1]
corresponds to the weights for the ith edge. If wgt_dim=0, the return value of ewgts is undefined and may be NULL.
    ierr Error code to be set by function.

 

C and C++: typedef void ZOLTAN_EDGE_LIST_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_id, ZOLTAN_ID_PTR local_id, ZOLTAN_ID_PTR nbor_global_id, int *nbor_procs, int wgt_dim, float *ewgts, int *ierr); 
FORTRAN: SUBROUTINE Get_Edge_List(data, num_gid_entries, num_lid_entries, global_id, local_id, nbor_global_id, nbor_procs, wgt_dim, ewgts, ierr
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_id 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_id 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: nbor_global_id 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: nbor_procs 
INTEGER(Zoltan_INT), INTENT(IN) :: wgt_dim 
REAL(Zoltan_FLOAT), INTENT(OUT), DIMENSION(*) :: ewgts 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_EDGE_LIST_FN query function returns lists of global IDs, processor IDs, and optionally edge weights for objects sharing an edge with a given object (i.e., objects that must share information with the given object). The arrays for the returned neighbor lists are allocated by Zoltan; their size is determined by a call to ZOLTAN_NUM_EDGES_MULTI_FN or ZOLTAN_NUM_EDGES_FN query functions.
 
Function Type: ZOLTAN_EDGE_LIST_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    global_id The global ID of the object for which an edge list should be returned.
    local_id The local ID of the object for which an edge list should be returned.
    nbor_global_id Upon return, an array of global IDs of objects sharing edges with the given object.
    nbor_procs Upon return, an array of processor IDs that identifies where the neighboring objects reside.
    wgt_dim The number of weights associated with an edge (typically 1), or 0 if edge weights are not requested. This value is set through the parameter EDGE_WEIGHT_DIM.
    ewgts Upon return, an array of edge weights, where ewgts[i*wgt_dim:(i+1)*wgt_dim-1]
corresponds to the weights for the ith edge. If wgt_dim=0, the return value of ewgts is undefined and may be NULL.
    ierr Error code to be set by function.

 

Hypergraph-based Functions



C and C++: typedef void ZOLTAN_HG_SIZE_CS_FN (void *data, int *num_lists, int *num_pins, int *format, int *ierr); 
FORTRAN: SUBROUTINE Get_HG_Size_CS(data, num_lists, num_pins, format, ierr
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(OUT) :: num_lists 
INTEGER(Zoltan_INT), INTENT(OUT) :: num_pins 
INTEGER(Zoltan_INT), INTENT(OUT) :: format 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A hypergraph (which may alternatively be viewed as a sparse matrix) can be supplied to the Zoltan library in one of two compressed storage formats. In compressed hyperedge format (ZOLTAN_COMPRESSED_EDGE) a list of global hyperedge IDs is provided. Then a single list of the hypergraph pins, is provided. A pin is the connection between a vertex and a hyperedge (corresponds to a nonzero in a sparse matrix). Pins do not have separate IDs but are rather identified by the global ID of the vertex containing the pin, and implicitly also by the hyperedge ID. An example is provided below.

The other format is compressed vertex (ZOLTAN_COMPRESSED_VERTEX). In this format a list of vertex global IDs is provided. Then a list of pins ordered by vertex and then by hyperedge is provided. The pin ID in this case is the global ID of the row (or hyperedge) in which the pin appears. In both formats, an array must be provided pointing to the start in the list of pins where each row or column begins. Sparse matrix users may think of these two formats as CSR (compressed sparse row) and CSC (compressed sparse column) format, respectively.

The point of this query function is to tell Zoltan in which format the application will supply the hypergraph, how many vertices and hyperedges there will be, and how many pins. The actual hypergraph is supplied with a query function of the type ZOLTAN_HG_CS_FN_TYPE.

This query function is required by all applications using the hypergraph methods of Zoltan (unless they are using the graph-based functions with hypergraph code instead).
 
Function Type: ZOLTAN_HG_SIZE_CS_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_lists Upon return, the number of rows (if using compressed row storage) or columns (if using compressed column storage) that will be supplied to Zoltan by the application process.
    num_pins Upon return, the number of pins (matrix non-zeroes) that will be supplied to Zoltan by the application process.
    format Upon return, the format in which the application process will provide the hypergraph to Zoltan. The options are ZOLTAN_COMPRESSED_EDGE and ZOLTAN_COMPRESSED_VERTEX.
    ierr Error code to be set by function.



C and C++: typedef void ZOLTAN_HG_CS_FN (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); 
FORTRAN: SUBROUTINE Get_HG_CS(data, num_gid_entries, num_row_or_col, num_pins, format, vtxedge_GID, vtxedge_ptr, pin_GID, ierr
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_row_or_col, num_pins, format
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: vtxedge_GID
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: vtxedge_ptr
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: pin_GID
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_HG_CS_FN returns a hypergraph, in a sparse matrix-like style. The size and format of the data to be returned must have been supplied to Zoltan using a ZOLTAN_HG_SIZE_CS_FN_TYPE function.

When a hypergraph is distributed across multiple processes, Zoltan expects that all processes share a consistent global numbering scheme for hyperedges and vertices. Also, no two processes should return the same pin (matrix non-zero) in this query function.

This query function is required by all applications using the hypergraph methods of Zoltan (unless they are using the graph-based functions with hypergraph code instead).
 
Function Type: ZOLTAN_HG_CS_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_vtx_edge The number of global IDs that is expected to appear on return in vtxedge_GID.
    num_pins The number of pins that is expected to appear on return in pin_GID.
    format If format is ZOLTAN_COMPRESSED_EDGE, Zoltan expects that row (hyperedge) global IDs will be returned in vtxedge_GID, and that column (vertex) global IDs will be returned in pin_GIDs. If it is ZOLTAN_COMPRESSED_VERTEX, then column global IDs are expected to be returned in vtxedge_GID and row global IDs are expected to be returned in pin_GIDs.
    vtxedge_GID Upon return, a list of num_row_or_col global IDs.
    vtxedge_ptr Upon return, this array contains num_row_or_col integers. The integer in the i'th array element is the index in array pin_GID where the pins for the i'th row (if format is ZOLTAN_COMPRESSED_EDGE) or i'th column (if format is ZOLTAN_COMPRESSED_VERTEX) begins. Array indices begin at zero.
    pin_GID Upon return, a list of num_pins global IDs. This is the list of the pins (or matrix non-zeros) contained in the rows or columns listed in vtxedge_GID.
    ierr Error code to be set by function.


Example
vertex
hyperedge 10 20 30 40 50
1 0 0 1 1 0
2 0 1 1 0 0
3 1 0 0 0 1
Compressed hyperedge storage:

vtxedge_GID = {1, 2, 3}
vtxedge_ptr = {0, 2, 4}
pin_GID = {30, 40, 20, 30, 10, 50}

Compressed vertex storage:

vtxedge_GID = {10, 20, 30, 40, 50}
vtxedge_ptr = {0, 1, 2, 4, 5}
pin_GID = {3, 2, 1, 2, 1, 3}


C and C++: typedef void ZOLTAN_HG_SIZE_EDGE_WTS_FN (void *data, int *num_edges, int *ierr); 
FORTRAN: SUBROUTINE Get_HG_Size_Edge_Wts(data, num_edges, ierr
<type-data>, INTENT(IN) :: data
INTEGER(Zoltan_INT), INTENT(OUT) :: num_edges
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_HG_SIZE_EDGE_WTS_FN returns the number of hyperedges for which a process will supply edge weights. The number of weights per hyperedge was supplied by the application with the EDGE_WEIGHT_DIM parameter. The actual edge weights will be supplied with a ZOLTAN_HG_EDGE_WTS_FN_TYPE function.

This query function is not required. If no hyperedge weights are supplied, Zoltan will assume every hyperedge has weight 1.0.
 
Function Type: ZOLTAN_HG_SIZE_EDGE_WTS_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_edges Upon return, the number of hyperedges for which edge weights will be supplied.
    ierr Error code to be set by function.


C and C++: typedef void ZOLTAN_HG_EDGE_WTS_FN (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_weight, int *ierr); 
FORTRAN: SUBROUTINE Get_HG_Edge_Wts(data, num_gid_entries, num_lid_entries, num_edges, edge_weight_dim, edge_GID, edge_LID, edge_weight, ierr
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries, num_edges, edge_weight_dim 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: edge_GID 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: edge_LID 
REAL(Zoltan_FLOAT), INTENT(OUT), DIMENSION(*) :: edge_weight 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_HG_EDGE_WTS_FN returns edges weights for a set of hypergraph edges. The number of weights supplied for each hyperedge should equal the value of the EDGE_WEIGHT_DIM parameter. In the case of a hypergraph which is distributed across multiple processes, if more than one process supplies edge weights for the same hyperedge, the different edge weights will be resolved according to the value of the PHG_EDGE_WEIGHT_OPERATION parameter.

This query function is not required. If no hyperedge weights are supplied, Zoltan will assume every hyperedge has weight 1.0.
 
Function Type: ZOLTAN_HG_SIZE_EDGE_WTS_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    num_edges The number of hyperedges for which edge weights should be supplied in the edge_weight array.
    edge_weight_dim The number of weights which should be supplied for each hyperedge. This is also the value of the EDGE_WEIGHT_DIM parameter.
    edge_GID Upon return, this array should contain the global IDs of the num_edges hyperedges for which the application is supplying edge weights.
    edge_LID Upon return, this array can optionally contain the local IDs of the num_edges hyperedges for which the application is supplying edge weights.
    edge_weight Upon return, this array should contain the weights for each edge listed in the edge_GID. If edge_weight_dim is greater than one, all weights for one hyperedge are listed before the weights for the next hyperedge are listed.
    ierr Error code to be set by function.

 

Tree-based Functions



C and C++: typedef int ZOLTAN_NUM_COARSE_OBJ_FN (void *data, int *ierr); 
FORTRAN: FUNCTION Get_Num_Coarse_Obj(data, ierr
INTEGER(Zoltan_INT) :: Get_Num_Coarse_Obj 
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_NUM_COARSE_OBJ_FN query function returns the number of objects (elements) in the initial coarse grid.
 
Function Type: ZOLTAN_NUM_COARSE_OBJ_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    ierr Error code to be set by function.
Returned Value:
    int The number of objects in the coarse grid.



C and C++: typedef void ZOLTAN_COARSE_OBJ_LIST_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int *assigned, int *num_vert, ZOLTAN_ID_PTR vertices, int *in_order, ZOLTAN_ID_PTR in_vertex, ZOLTAN_ID_PTR out_vertex, int *ierr); 
FORTRAN: SUBROUTINE Get_Coarse_Obj_List(data, num_gid_entries, num_lid_entries, global_ids, local_ids, assigned, num_vert, vertices, in_order, in_vertex, out_vertex, ierr
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: global_ids 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: local_ids 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: assigned, num_vert, vertices, in_vertex, out_vertex 
INTEGER(Zoltan_INT), INTENT(OUT) :: in_order, ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_COARSE_OBJ_LIST_FN query function returns lists of global IDs, local IDs, vertices, and order information for all objects (elements) of the initial coarse grid. The vertices are designated by a global ID such that if two elements share a vertex then the same ID designates that vertex in both elements and on all processors. The user may choose to provide the order in which the elements should be traversed or have Zoltan determine the order. If the user provides the order, then entry and exit vertices for a path through the elements may also be provided. The arrays for the returned values are allocated by Zoltan; their size is determined by a call to a ZOLTAN_NUM_COARSE_OBJ_FN query function.
 
Function Type: ZOLTAN_COARSE_OBJ_LIST_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    global_ids Upon return, an array of global IDs of all objects in the coarse grid.
    local_ids Upon return, an array of local IDs of all objects in the coarse grid.
    assigned Upon return, an array of integers indicating whether or not each object is currently assigned to this processor. A value of 1 indicates it is assigned to this processor; a value of 0 indicates it is assigned to some other processor. For elements that have been refined, it is ignored unless weights are assigned to interior nodes of the tree.
    num_vert Upon return, an array containing the number of vertices for each object.
    vertices Upon return, an array of global IDs of the vertices of each object. If the number of vertices for objects 0 through i-1 is N, then the vertices for object i are in vertices[N*num_gid_entries: (N+num_vert[i])*num_gid_entries]
    in_order Upon return, 1 if the user is providing the objects in the order in which they should be traversed, or 0 if Zoltan should determine the order.
    in_vertex Upon return, an array of global IDs of the vertices through which to enter each element in the user provided traversal. It is required only if the user is providing the order for the coarse grid objects (i.e., in_order==1) and allowing Zoltan to select the order of the children in at least one invocation of ZOLTAN_CHILD_LIST_FN.
    out_vertex Upon return, an array of global IDs of the vertex through which to exit each element in the user provided traversal. The same provisions hold as for in_vertex.
    ierr Error code to be set by function.



C and C++: typedef int ZOLTAN_FIRST_COARSE_OBJ_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_id, ZOLTAN_ID_PTR local_id, int *assigned, int *num_vert, ZOLTAN_ID_PTR vertices, int *in_order, ZOLTAN_ID_PTR in_vertex, ZOLTAN_ID_PTR out_vertex, int *ierr); 
FORTRAN: FUNCTION Get_First_Coarse_Obj(data, num_gid_entries, num_lid_entries, global_id, local_id, assigned, num_vert, vertices, in_order, in_vertex, out_vertex, ierr
INTEGER(Zoltan_INT) :: Get_First_Coarse_Obj 
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: global_id 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: local_id 
INTEGER(Zoltan_INT), INTENT(OUT) :: assigned, num_vert, in_order, ierr 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: vertices, in_vertex, out_vertex 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_FIRST_COARSE_OBJ_FN query function initializes an iteration over the objects of the initial coarse grid. It returns the global ID, local ID, vertices, and order information for the first object (element) of the initial coarse grid. Subsequent calls to a ZOLTAN_NEXT_COARSE_OBJ_FN iterate over and return other objects from the coarse grid. The vertices are designated by a global ID such that if two elements share a vertex then the same ID designates that vertex in both elements and on all processors. The user may choose to provide the order in which the elements should be traversed, or have Zoltan determine the order. If the user provides the order, then entry and exit vertices for a path through the elements may also be provided.
 
Function Type: ZOLTAN_FIRST_COARSE_OBJ_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    global_ids Upon return, the global ID of the first object in the coarse grid.
    local_ids Upon return, the local ID of the first object in the coarse grid.
    assigned Upon return, an integer indicating whether or not this object is currently assigned to this processor. A value of 1 indicates it is assigned to this processor; a value of 0 indicates it is assigned to some other processor. For elements that have been refined, it is ignored unless weights are assigned to interior nodes of the tree.
    num_vert Upon return, the number of vertices for this object.
    vertices Upon return, an array of global IDs of the vertices of this object.
    in_order Upon return, 1 if the user is providing the objects in the order in which they should be traversed, or 0 if Zoltan should determine the order.
    in_vertex Upon return, the vertex through which to enter this element in the user provided traversal. It is required only if the user is providing the order for the coarse grid objects (i.e., in_order==1) and allowing Zoltan to select the order of the children in at least one invocation of ZOLTAN_CHILD_LIST_FN.
    out_vertex Upon return, the vertex through which to exit this element in the user provided traversal. The same provisions hold as for in_vertex.
    ierr Error code to be set by function.
Returned Value:
    1 If global_id and local_id contain valid IDs of the first object in the coarse grid.
    0 If no coarse grid is available.



C and C++: typedef int ZOLTAN_NEXT_COARSE_OBJ_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_id, ZOLTAN_ID_PTR local_id, ZOLTAN_ID_PTR next_global_id, ZOLTAN_ID_PTR next_local_id, int *assigned, int *num_vert, ZOLTAN_ID_PTR vertices, ZOLTAN_ID_PTR in_vertex, ZOLTAN_ID_PTR out_vertex, int *ierr); 
FORTRAN: FUNCTION Get_Next_Coarse_Obj(data, num_gid_entries, num_lid_entries, global_id, local_id, next_global_id, next_local_id, assigned, num_vert, vertices, in_vertex, out_vertex, ierr
INTEGER(Zoltan_INT) :: Get_Next_Coarse_Obj 
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_id 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_id 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: next_global_id 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: next_local_id 
INTEGER(Zoltan_INT), INTENT(OUT) :: assigned, num_vertex, ierr 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: vertices, in_vertex, out_vertex 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_NEXT_COARSE_OBJ_FN query function is an iterator function that returns the next object in the initial coarse grid. The first object of the iteration is provided by a ZOLTAN_FIRST_COARSE_OBJ_FN query function.
 
Function Type: ZOLTAN_NEXT_COARSE_OBJ_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    global_id The global ID of the previous object in the coarse grid.
    local_id The local ID of the previous object in the coarse grid.
    next_global_id Upon return, the global ID of the next object in the coarse grid.
    next_local_id Upon return, the local ID of the next object in the coarse grid.
    assigned Upon return, an integer indicating whether or not this object is currently assigned to this processor. A value of 1 indicates it is assigned to this processor; a value of 0 indicates it is assigned to some other processor. For elements that have been refined, it is ignored unless weights are assigned to interior nodes of the tree.
    num_vert Upon return, the number of vertices for this object.
    vertices Upon return, an array of global IDs of the vertices of this object.
    in_vertex Upon return, the vertex through which to enter this element in the user provided traversal. It is required only if the user is providing the order for the coarse grid objects (i.e., in_order==1) and allowing Zoltan to select the order of the children in at least one invocation of ZOLTAN_CHILD_LIST_FN.
    out_vertex Upon return, the vertex through which to exit this element in the user provided traversal. The same provisions hold as for in_vertex.
    ierr Error code to be set by function.
Returned Value:
    1 If global_id and local_id contain valid IDs of the next object in the coarse grid.
    0 If no more objects are available.



C and C++: typedef int ZOLTAN_NUM_CHILD_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_id, ZOLTAN_ID_PTR local_id, int *ierr); 
FORTRAN: FUNCTION Get_Num_Child(data, num_gid_entries, num_lid_entries, global_id, local_id, ierr
INTEGER(Zoltan_INT) :: Get_Num_Child 
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_id 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_id 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_NUM_CHILD_FN query function returns the number of children of the element with the given global and local IDs. If the element has not been refined, the number of children is 0.
 
Function Type: ZOLTAN_NUM_CHILD_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    global_id The global ID of the object for which the number of children is requested.
    local_id The local ID of the object for which the number of children is requested.
    ierr Error code to be set by function.
Returned Value:
    int The number of children.



C and C++: typedef void ZOLTAN_CHILD_LIST_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR parent_gid, ZOLTAN_ID_PTR parent_lid, ZOLTAN_ID_PTR child_gids, ZOLTAN_ID_PTR child_lids, int *assigned, int *num_vert, ZOLTAN_ID_PTR vertices, ZOLTAN_REF_TYPE *ref_type, ZOLTAN_ID_PTR in_vertex, ZOLTAN_ID_PTR out_vertex, int *ierr); 
FORTRAN: SUBROUTINE Get_Child_List(data, num_gid_entries, num_lid_entries, parent_gid, parent_lid, child_gids, child_lids, assigned, num_vert, vertices, ref_type, in_vertex, out_vertex, ierr
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: parent_gid 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: parent_lid 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: child_gids 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: child_lids 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: assigned, num_vert, vertices, in_vertex, out_vertex 
INTEGER(Zoltan_INT), INTENT(OUT) :: ref_type, ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_CHILD_LIST_FN query function returns lists of global IDs, local IDs, vertices, and order information for all children of a refined element. The vertices are designated by a global ID such that if two elements share a vertex then the same ID designates that vertex in both elements and on all processors. The user may choose to provide the order in which the children should be traversed, or have Zoltan determine the order based on the type of element refinement used to create the children. If the user provides the order, then entry and exit vertices for a path through the elements may also be provided. The arrays for the returned values are allocated by Zoltan; their size is determined by a call to a ZOLTAN_NUM_CHILD_FN query function.
 
Function Type: ZOLTAN_CHILD_LIST_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    parent_gid The global ID of the object whose children are requested.
    parent_lid The local ID of the object whose children are requested.
    child_gids Upon return, an array of global IDs of all children of this object.
    child_lids Upon return, an array of local IDs of all children of this object.
    assigned Upon return, an array of integers indicating whether or not each child is currently assigned to this processor. A value of 1 indicates it is assigned to this processor; a value of 0 indicates it is assigned to some other processor. For children that have been further refined, it is ignored unless weights are assigned to interior nodes of the tree.
    num_vert Upon return, an array containing the number of vertices for each object.
    vertices Upon return, an array of global IDs of the vertices of each object. If the number of vertices for objects 0 through i-1 is N, then the vertices for object i are in vertices[N*num_gid_entries: (N+num_vert[i])*num_gid_entries]
    ref_type Upon return, a value indicating what type of refinement was used to create the children. This determines how the children will be ordered. The values currently supported are:
    ZOLTAN_TRI_BISECT Bisection of triangles.
    ZOLTAN_QUAD_QUAD Quadrasection of quadrilaterals.
    ZOLTAN_HEX3D_OCT Octasection of hexahedra.
    ZOLTAN_OTHER_REF All other forms of refinement.
    ZOLTAN_IN_ORDER Traverse the children in the order in which they are provided.
    in_vertex Upon return, an array of global IDs of the vertex through which to enter each element in the user provided traversal. It is required only if the user is providing the order for the children of this element (i.e., ref_type==ZOLTAN_IN_ORDER) but does not provide the order for the children of at least one of those children.
    out_vertex Upon return, an array of global IDs of the vertex through which to exit each element in the user provided traversal. The same provisions hold as for in_vertex.
    ierr Error code to be set by function.



C and C++: typedef void ZOLTAN_CHILD_WEIGHT_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_id, ZOLTAN_ID_PTR local_id, int wgt_dim, float *obj_wgt, int *ierr); 
FORTRAN: SUBROUTINE Get_Child_Weight(data, num_gid_entries, num_lid_entries, global_id, local_id, wgt_dim, obj_wgt, ierr
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_id 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_id 
INTEGER(Zoltan_INT), INTENT(IN) :: wgt_dim 
REAL(Zoltan_FLOAT), INTENT(OUT), DIMENSION(*) :: obj_wgt 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_CHILD_WEIGHT_FN query function returns the weight of an object. Interior nodes of the refinement tree as well as the leaves are allowed to have weights.
 
Function Type: ZOLTAN_CHILD_WEIGHT_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    global_id The global ID of the object whose weight is requested.
    local_id The local ID of the object whose weight is requested.
    wgt_dim The number of weights associated with an object (typically 1), or 0 if weights are not requested. This value is set through the parameter OBJ_WEIGHT_DIM.
    obj_wgt Upon return, an array containing the object's weights. If wgt_dim=0, the return value of obj_wgts is undefined and may be NULL.
    ierr Error code to be set by function.

Border Object Functions (currently not used)



C: typedef int ZOLTAN_NUM_BORDER_OBJ_FN (void *data, int nbor_proc, int *ierr); 
FORTRAN: FUNCTION Get_Num_Border_Obj(data, nbor_proc, ierr
INTEGER(Zoltan_INT) :: Get_Num_Border_Obj
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: nbor_proc 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_NUM_BORDER_OBJ_FN query function returns the number of objects sharing a processor subdomain border (in the communication graph of the application) with a given processor.
 
Function Type: ZOLTAN_NUM_BORDER_OBJ_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    nbor_proc The processor ID of the processor for which the number of border objects should be returned.
    ierr Error code to be set by function.
Returned Value:
    int The number of objects sharing a processor subdomain border with processor nbor_proc.



C: typedef void ZOLTAN_BORDER_OBJ_LIST_FN (void *data, int num_gid_entries, int num_lid_entries, int nbor_proc, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int wgt_dim, float *obj_wgts, int *ierr); 
FORTRAN: SUBROUTINE Get_Border_Obj_List(data, num_gid_entries, num_lid_entries, nbor_proc, global_ids, local_ids, wgt_dim, obj_wgts, ierr
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: nbor_proc 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: global_ids 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: local_ids 
INTEGER(Zoltan_INT), INTENT(IN) :: wgt_dim 
REAL(Zoltan_FLOAT), INTENT(OUT), DIMENSION(*) :: obj_wgts 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_BORDER_OBJ_LIST_FN query function fills two arrays with information about the objects currently assigned to the processor that share a processor subdomain border (in the communication graph of the application) with a given processor. Both arrays are allocated (and subsequently freed) by Zoltan; their size is determined by a call to a ZOLTAN_NUM_BORDER_OBJ_FN query function to get the array size. For certain Zoltan algorithms, either a ZOLTAN_BORDER_OBJ_LIST_FN query function or a ZOLTAN_FIRST_BORDER_OBJ_FN/ZOLTAN_NEXT_BORDER_OBJ_FN query-function pair must be registered; however, both query options need not be provided.
 
Function Type: ZOLTAN_BORDER_OBJ_LIST_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    nbor_proc The processor ID of the processor for which border objects should be returned.
    global_ids Upon return, an array of unique global IDs for all objects assigned to the processor that share a subdomain border with nbor_proc.
    local_ids Upon return, an array of local IDs, the meaning of which can be determined by the application, for all objects assigned to the processor that share a subdomain border with nbor_proc.
    wgt_dim The number of weights associated with an object (typically 1), or 0 if weights are not requested. This value is set through the parameter OBJ_WEIGHT_DIM.
    obj_wgts Upon return, an array of object weights. Weights for object i are stored in obj_wgts[(i-1)*wgt_dim:i*wgt_dim-1].  If wgt_dim=0, obj_wgts is undefined and may be NULL.
    ierr Error code to be set by function.



C: typedef int ZOLTAN_FIRST_BORDER_OBJ_FN (void *data, int num_gid_entries, int num_lid_entries, int nbor_proc, ZOLTAN_ID_PTR first_global_id, ZOLTAN_ID_PTR first_local_id, int wgt_dim, float *first_obj_wgt, int *ierr); 
FORTRAN: FUNCTION Get_First_Border_Obj(data, num_gid_entries, num_lid_entries, nbor_proc, first_global_id, first_local_id, wgt_dim, first_obj_wgt, ierr
INTEGER(Zoltan_INT) :: Get_First_Border_Obj
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(IN) :: nbor_proc 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: first_global_id 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: first_local_id 
INTEGER(Zoltan_INT), INTENT(IN) :: wgt_dim 
REAL(Zoltan_FLOAT), INTENT(OUT), DIMENSION(*) :: first_obj_wgt 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_FIRST_BORDER_OBJ_FN query function initializes an iteration over objects assigned to the processor that share a processor subdomain border with a given processor. It returns the global and local IDs of the first object on the processor along the specified subdomain border. Subsequent calls to a ZOLTAN_NEXT_BORDER_OBJ_FN query function iterate over and return other objects along the requested subdomain border. This query-function pair frees the application from having to build an array of objects (as in ZOLTAN_BORDER_OBJ_LIST_FN) and allows Zoltan to obtain only as much information about objects as it needs. For some algorithms, either a ZOLTAN_BORDER_OBJ_LIST_FN query function or a ZOLTAN_FIRST_BORDER_OBJ_FN/ZOLTAN_NEXT_BORDER_OBJ_FN query-function pair must be registered; however, both query options need not be provided.
 
Function Type: ZOLTAN_FIRST_BORDER_OBJ_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    nbor_proc The processor ID of the processor for which border objects should be returned.
    first_global_id The returned value of the global ID for the first object; the value is ignored if there are no objects along the border.
    first_local_id The returned value of the local ID for the first object; the value is ignored if there are no objects along the border.
    wgt_dim The number of weights associated with an object (typically 1), or 0 if weights are not requested. This value is set through the parameter OBJ_WEIGHT_DIM.
    first_obj_wgt Upon return, the first object's weights; an array of size wgt_dim. Undefined if  wgt_dim=0.
    ierr Error code to be set by function.
Returned Value:
    1 If first_global_id and first_local_id contain valid IDs of the first object along the processor border.
    0 If no objects are available along this processor border.



C: typedef int ZOLTAN_NEXT_BORDER_OBJ_FN (void *data, int num_gid_entries, int num_lid_entries, ZOLTAN_ID_PTR global_id, ZOLTAN_ID_PTR local_id, int nbor_proc, ZOLTAN_ID_PTR next_global_id, ZOLTAN_ID_PTR next_local_id, int wgt_dim, float *next_obj_wgt, int *ierr); 
FORTRAN: FUNCTION Get_Next_Border_Obj(data, num_gid_entries, num_lid_entries, global_id, local_id, nbor_proc, next_global_id, next_local_id, wgt_dim, next_obj_wgt, ierr
INTEGER(Zoltan_INT) :: Get_Next_Border_Obj 
<type-data>, INTENT(IN) :: data 
INTEGER(Zoltan_INT), INTENT(IN) :: num_gid_entries, num_lid_entries 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: global_id 
INTEGER(Zoltan_INT), INTENT(IN), DIMENSION(*) :: local_id 
INTEGER(Zoltan_INT), INTENT(IN) :: nbor_proc 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: next_global_id 
INTEGER(Zoltan_INT), INTENT(OUT), DIMENSION(*) :: next_local_id 
INTEGER(Zoltan_INT), INTENT(IN) :: wgt_dim 
REAL(Zoltan_FLOAT), INTENT(OUT), DIMENSION(*) :: next_obj_wgt 
INTEGER(Zoltan_INT), INTENT(OUT) :: ierr 

<type-data> can be any of INTEGER(Zoltan_INT), DIMENSION(*) or REAL(Zoltan_FLOAT), DIMENSION(*) or REAL(Zoltan_DOUBLE), DIMENSION(*) or TYPE(Zoltan_User_Data_x) where x is 1, 2, 3 or 4. See the section on Fortran query functions for an explanation. 


A ZOLTAN_NEXT_BORDER_OBJ_FN query function is an iterator function which, when given an object assigned to the processor and a neighboring processor ID, returns the next object assigned to the processor that shares a subdomain border with the neighboring processor. The first object of the iteration is provided by a ZOLTAN_FIRST_BORDER_OBJ_FN query function. This query-function pair frees the application from having to build an array of objects (as in ZOLTAN_BORDER_OBJ_LIST_FN) and allows Zoltan to obtain only as much information about objects as it needs. For some algorithms, either a ZOLTAN_BORDER_OBJ_LIST_FN query function or a ZOLTAN_FIRST_BORDER_OBJ_FN/ZOLTAN_NEXT_BORDER_OBJ_FN query-function pair must be registered; however, both query options need not be provided.
 
Function Type: ZOLTAN_NEXT_BORDER_OBJ_FN_TYPE
Arguments:
    data Pointer to user-defined data.
    num_gid_entries The number of array entries used to describe a single global ID.  This value is the maximum value over all processors of the parameter NUM_GID_ENTRIES.
    num_lid_entries The number of array entries used to describe a single local ID.  This value is the maximum value over all processors of the parameter NUM_LID_ENTRIES.
    global_id The global ID of the previous object.
    local_id The local ID of the previous object.
    nbor_proc The processor ID of the processor for which border objects should be returned.
    next_global_id The returned value of the global ID for the next object; the value is ignored if there are no more objects along the border.
    next_local_id The returned value of the local ID for the next object; the value is ignored if there are no more objects along the border.
    wgt_dim The number of weights associated with an object (typically 1), or 0 if weights are not requested. This value is set through the parameter OBJ_WEIGHT_DIM.
    next_obj_wgt Upon return, the weights for the next object; an array of size wgt_dim. Undefined if wgt_dim=0.
    ierr Error code to be set by function.
Returned Value:
    1 If next_global_id and next_local_id contain valid IDs of the next object along the processor border.
    0 If no more objects are available along this processor border.



[Table of Contents  | Next:  Migration Query Functions  |  Previous:  Application-Registered Query Functions]