Zoltan User's Guide  |  Next  |  Previous

Ordering Functions

Zoltan provides limited capability for ordering a set of objects, typically given as a graph. The following functions are the ordering interface functions in the Zoltan library; their descriptions are included below.
Zoltan_Order


C: int Zoltan_Order (
      struct Zoltan_Struct *zz,
      int *num_gid_entries,
      int *num_lid_entries,
      int num_obj,
      ZOLTAN_ID_PTR global_ids,
      ZOLTAN_ID_PTR local_ids,
      int *rank,
      int *iperm,
      struct Zoltan_Order_Struct *order_info);
FORTRAN: FUNCTION Zoltan_Order(zz, num_gid_entries, num_lid_entries, num_obj, global_ids, local_ids, rank, iperm)
INTEGER(Zoltan_INT) :: Zoltan_Order
TYPE(Zoltan_Struct), INTENT(IN) :: zz
INTEGER(Zoltan_INT), INTENT(OUT) :: num_gid_entries, num_lid_entries
INTEGER(Zoltan_INT), INTENT(IN) :: num_obj
INTEGER(Zoltan_INT) :: global_ids(*), local_ids(*)
INTEGER(Zoltan_INT) :: rank(*), iperm(*)
C++: int Zoltan::Order (
      int &num_gid_entries,
      int &num_lid_entries,
      const int &num_obj,
      ZOLTAN_ID_PTR global_ids,
      ZOLTAN_ID_PTR local_ids,
      int *rank,
      int *iperm);

Zoltan_Order invokes the ordering routine specified by the  ORDER_METHOD parameter. Results of the ordering are returned in  the arrays rank and iperm. rank[i] gives the rank of global_ids[i] in the computed ordering, while iperm is the inverse permutation of rank, that is, iperm[rank[i]] = i.  The ordering may be either global or local, depending on ORDER_TYPE. The arrays global_ids, local_ids, rank, and iperm should all be allocated by the application before Zoltan_Order is called. Each array must have space for (at least) num_obj elements, where num_obj is the number of objects residing on a processor.
 
Arguments:
    zz Pointer to the Zoltan structure, created by Zoltan_Create, to be used in this invocation of the load-balancing routine.
    num_gid_entries Upon return, 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 Upon return, 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 Number of objects to order on this processor. At present, num_obj should be the total number of objects residing on a processor. In future releases, ordering only a subset of the objects may be permitted. 
    global_ids An array of global IDs of objects to be ordered on this processor. (size = num_obj * num_gid_entries)
The array may be uninitialized on input (if REORDER is false), but memory must have been allocated before Zoltan_Order is called.
    local_ids [Optional.] An array of local IDs of objects to be ordered on this processor. (size = num_obj * num_lid_entries)
The array may be uninitialized on input (if REORDER is false), but memory must have been allocated before Zoltan_Order is called.
    rank Upon return, an array of length num_obj containing the rank of each object in the computed ordering. When rank[i] = j, that means that the object corresponding to global_ids[i] is the jth object in the ordering. (This array corresponds directly to the perm array in METIS and the order array in ParMETIS.) Note that the rank may refer to either a local or a global ordering, depending on ORDER_TYPE.  Memory for this array must have been allocated before Zoltan_Order is called.
    iperm Upon return, an  array of length num_obj containing the inverse permutation of rank. That is, iperm[rank[i]] = i. In other words, iperm[j] gives the jth object in the ordering. Memory for this array must have been allocated before Zoltan_Order is called.
    order_info Upon return,  this struct contains additional information about the ordering produced. This parameter is currently not used and should always be set to NULL. It is not included in the FORTRAN or C++ interface.
Returned Value:
    int Error code.


[Table of Contents  | Next:  Coloring Functions  |  Previous:  Migration Functions]