% % This section includes brief examples, including all source % code and gnuplot images where possible, of a geometric problem, % a graph problem, and a hypergraph problem. % \chapter{Examples} \label{cha:ex} This chapter begins with a very simple example. We use the Zoltan library in a parallel application to partition a small set of objects, each of which has only a global ID and a weight. We follow this with three more realistic examples. These examples use the Zoltan library to apply a geometric method, a graph method and a then hypergraph method to a small mesh, graph and hypergraph respectively. If you have the Zoltan source code, you can find these examples in the \textbf{examples/C} directory. The versions found in the source code include code to display the initial partitioning followed by Zoltan's computed partitioning. It may be helpful to compile and run them, and then try changing parameters to see how the partitioning is affected. Consult the Zoltan User's Guide at \url{http://cs.sandia.gov/Zoltan/ug_html/ug.html} for a complete list of all partitioning methods and all of the parameters for each partitioning method. The next four examples are written in C. At the end of this section, we also provide a C++ version of the geometric example. As you read through these examples you will see that, regardless of the load balancing method chosen, every program must: \begin{itemize} \item initialize the Zoltan library \item supply load balancing parameters to Zoltan \item define query functions that Zoltan will use to obtain information from the application about the objects to be partitioned \item call Zoltan\_LB\_Partition to begin the load balancing calculation \item free memory allocated by Zoltan when done \end{itemize} \clearpage \section{A very simple partitioning method} This example is a C program which uses Zoltan to partition in parallel 36 objects, each of which has a weight. The goal is to balance the weights across the parallel application. It uses the Zoltan method called BLOCK, which simply assignes the first $n/p$ data objects to the first process and so on. This method really exists as an example for developers to learn how to write partitioning methods. It can also be used for testing in applications, but was never intended to be a serious partitioning strategy for real applications. However, it is a good choice for our first example. If you are a new user of Zoltan, you may want to initially get your application working with the BLOCK method. Once you have that first step accomplished, you can modify it to use one of the other partitioning methods. In order to use the BLOCK method, we must define two query functions: \begin{itemize} \item A function of type ZOLTAN\_NUM\_OBJ\_FN which provides the number of objects belonging to this process. In our example we call this function \textbf{get\_number\_of\_objects}. \item A function of type ZOLTAN\_OBJ\_LIST\_FN which supplies the global IDs and optional weights for the objects belonging to the process. We name this function \textbf{get\_object\_list}. \end{itemize} The \textbf{ZOLTAN\_*} function prototyes named above are defined in the Zoltan file \textbf{include/zoltan.h}. Zoltan requires that the parallel application supply a unique global ID to identify each of the objects to be partitioned. Zoltan is very flexible about the data type used to represent a global ID, and only requires that it be allowed to represent your global ID internally as an array of integers of some user-defined length. Zoltan also permits the use of an additional set of IDs that are local to a process. If you supply these optional local IDs, Zoltan will also include them when querying your application for data. A version of the code listed below may be found in the Zoltan source in file \textbf{examples/C/simpleBLOCK.c}. \begin{flushleft} \begin{verbatim} #include #include #include #include "zoltan.h" static int numObjects = 36; static int numProcs, myRank, myNumObj; static int myGlobalIDs[36]; static float objWeight(int globalID) { float w; /* simulate an initial imbalance */ if (globalID % numProcs == 0) w = 3; else w = (globalID % 3 + 1); return w; } static int get_number_of_objects(void *data, int *ierr) { *ierr = ZOLTAN_OK; return myNumObj; } static void get_object_list(void *data, int sizeGID, int sizeLID, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int wgt_dim, float *obj_wgts, int *ierr) { int i; *ierr = ZOLTAN_OK; for (i=0; i #include #include #include "zoltan.h" #include "simpleGraph.h" #include "simpleQueries.h" int main(int argc, char *argv[]) { int rc, i, ngids, nextIdx; struct Zoltan_Struct *zz; int changes, numGidEntries, numLidEntries, numImport, numExport; ZOLTAN_ID_PTR importGlobalGids, importLocalGids; ZOLTAN_ID_PTR exportGlobalGids, exportLocalGids; int *importProcs, *importToPart, *exportProcs, *exportToPart; float *wgt_list, ver; int *gid_flags, *gid_list, *lid_list; /****************************************************************** ** Initialize MPI and Zoltan ******************************************************************/ MPI_Init(&argc, &argv); MPI_Comm_rank(MPI_COMM_WORLD, &myRank); MPI_Comm_size(MPI_COMM_WORLD, &numProcs); rc = Zoltan_Initialize(argc, argv, &ver); if (rc != ZOLTAN_OK){ MPI_Finalize(); exit(0); } \end{verbatim} \end{flushleft} \clearpage \begin{flushleft} \begin{verbatim} /****************************************************************** ** Create a Zoltan library structure for this instance of load ** balancing. Set the parameters and query functions that will ** govern the library's calculation. See the Zoltan User's ** Guide for the definition of these and many other parameters. ******************************************************************/ zz = Zoltan_Create(MPI_COMM_WORLD); /* General parameters */ Zoltan_Set_Param(zz, "DEBUG_LEVEL", "0"); Zoltan_Set_Param(zz, "LB_METHOD", "RCB"); Zoltan_Set_Param(zz, "NUM_GID_ENTRIES", "1"); Zoltan_Set_Param(zz, "NUM_LID_ENTRIES", "1"); Zoltan_Set_Param(zz, "OBJ_WEIGHT_DIM", "1"); Zoltan_Set_Param(zz, "RETURN_LISTS", "ALL"); /* RCB parameters */ Zoltan_Set_Param(zz, "KEEP_CUTS", "1"); Zoltan_Set_Param(zz, "RCB_OUTPUT_LEVEL", "0"); Zoltan_Set_Param(zz, "RCB_RECTILINEAR_BLOCKS", "1"); /* Query functions - defined in simpleQueries.h, to return * information about objects defined in simpleGraph.h */ Zoltan_Set_Num_Obj_Fn(zz, get_number_of_objects, NULL); Zoltan_Set_Obj_List_Fn(zz, get_object_list, NULL); Zoltan_Set_Num_Geom_Fn(zz, get_num_geometry, NULL); Zoltan_Set_Geom_Multi_Fn(zz, get_geometry_list, NULL); \end{verbatim} \end{flushleft} \clearpage \begin{flushleft} \begin{verbatim} /****************************************************************** ** Zoltan can now partition the vertices in the simple mesh. ** In this simple example, we assume the number of partitions is ** equal to the number of processes. Process rank 0 will own ** partition 0, process rank 1 will own partition 1, and so on. ******************************************************************/ rc = Zoltan_LB_Partition(zz, /* input (all remaining fields are output) */ &changes, /* 1 if partitioning was changed, 0 otherwise */ &numGidEntries, /* Number of integers used for a global ID */ &numLidEntries, /* Number of integers used for a local ID */ &numImport, /* Number of vertices to be sent to me */ &importGlobalGids, /* Global IDs of vertices to be sent to me */ &importLocalGids, /* Local IDs of vertices to be sent to me */ &importProcs, /* Process rank for source of each incoming vertex */ &importToPart, /* New partition for each incoming vertex */ &numExport, /* Number of vertices I must send to other processes*/ &exportGlobalGids, /* Global IDs of the vertices I must send */ &exportLocalGids, /* Local IDs of the vertices I must send */ &exportProcs, /* Process to which I send each of the vertices */ &exportToPart); /* Partition to which each vertex will belong */ if (rc != ZOLTAN_OK){ printf("sorry...\n"); MPI_Finalize(); Zoltan_Destroy(&zz); exit(0); } /****************************************************************** ** In a real application, you would rebalance the problem now by ** sending the objects to their new partitions. ******************************************************************/ \end{verbatim} \end{flushleft} \clearpage \begin{flushleft} \begin{verbatim} /****************************************************************** ** Free the arrays allocated by Zoltan_LB_Partition, and free ** the storage allocated for the Zoltan structure. ******************************************************************/ Zoltan_LB_Free_Part(&importGlobalGids, &importLocalGids, &importProcs, &importToPart); Zoltan_LB_Free_Part(&exportGlobalGids, &exportLocalGids, &exportProcs, &exportToPart); Zoltan_Destroy(&zz); /********************** ** all done *********** **********************/ MPI_Finalize(); return 0; } \end{verbatim} \end{flushleft} \clearpage \subsection{A graph problem} In this section we provide an example of a parallel C program that uses the Zoltan library to partition the simple graph that was pictured in Figure~\ref{fig:mesh25}. You can compile and run this example if you have the Zoltan source code. It is in the examples directory, in file \textbf{examples/C/simpleGRAPH.c}. \subsubsection{Application defined query functions} We must define four query functions in order to use the particular graph submethod we chose for this example. See the User's Guide at \url{http://cs.sandia.gov/Zoltan/ug_html/ug_alg_parmetis.html} for the requirements of other submethods. The function prototypes referred to below are defined in the Zoltan file \textbf{include/zoltan.h}. The four query functions are: \begin{itemize} \item A function of type ZOLTAN\_NUM\_OBJ\_FN (Figure~\ref{fig:NumObj}) which provides the number of objects belonging to this process \item A function of type ZOLTAN\_OBJ\_LIST\_FN (Figure~\ref{fig:ObjList}) which supplies the global IDs and optional weights for the objects belonging to the process \item A function of type ZOLTAN\_NUM\_EDGES\_MULTI\_FN (Figure~\ref{fig:NumEdges}) which provides the number of edges for each object \item A function of type ZOLTAN\_EDGE\_LIST\_MULTI\_FN (Figure~\ref{fig:EdgeListMulti}) which, for a given object ID, responds with the IDs for each adjacent object, the process owning that adjacent object, and optionally a weight for each edge \end{itemize} An alternative to defining a single ZOLTAN\_OBJ\_LIST\_FN is to define a ZOLTAN\_FIRST\_OBJ\_FN and a ZOLTAN\_NEXT\_OBJ\_FN. The first function would return the global ID for first object. Each call to the second function would return a subsequent global ID. An alternative to defining a ZOLTAN\_NUM\_EDGES\_MULTI\_FN is to define a ZOLTAN\_NUM\_EDGES\_FN, which returns the number of adjacent objects for a single object ID. An alternative to defining a ZOLTAN\_EDGE\_LIST\_MULTI\_FN is to define a ZOLTAN\_EDGE\_LIST\_FN, which returns the same information but for a single object. The variables used in these query functions to describe our simple graph refer to those defined in \textbf{examples/C/simpleGraph.h} which were shown in Figure~\ref{fig:simpleGraphDotH}. \begin{figure} \begin{flushleft} \begin{verbatim} static void get_num_edges_list(void *data, int sizeGID, int sizeLID, int num_obj, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int *numEdges, int *ierr) { int i, idx; if ( (sizeGID != 1) || (sizeLID != 1)){ *ierr = ZOLTAN_FATAL; return; } for (i=0; i < num_obj ; i++){ idx = globalID[i] - 1; numEdges[i] = simpleNumEdges[idx]; } *ierr = ZOLTAN_OK; return; } \end{verbatim} \end{flushleft} \caption{A ZOLTAN\_NUM\_EDGES\_MULTI\_FN query function} \label{fig:NumEdges} \end{figure} \begin{figure} \begin{flushleft} \begin{verbatim} static void get_edge_list(void *data, int sizeGID, int sizeLID, int num_obj, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int *num_edges, ZOLTAN_ID_PTR nborGID, int *nborProc, int wgt_dim, float *ewgts, int *ierr) { int i, j, idx; ZOLTAN_ID_PTR nextID; int *nextProc; if ( (sizeGID != 1) || (sizeLID != 1) || (wgt_dim != 0)){ *ierr = ZOLTAN_FATAL; return; } nextID = nborGID; nextProc = nborProc; for (i=0; i < num_obj ; i++){ idx = globalID[i] - 1; if (num_edges[i] != simpleNumEdges[idx]){ *ierr = ZOLTAN_FATAL; return; } for (j=0; j #include #include #include "zoltan.h" #include "simpleGraph.h" #include "simpleQueries.h" int main(int argc, char *argv[]) { int rc, i, ngids, nextIdx; float ver; struct Zoltan_Struct *zz; int changes, numGidEntries, numLidEntries, numImport, numExport; ZOLTAN_ID_PTR importGlobalGids, importLocalGids; ZOLTAN_ID_PTR exportGlobalGids, exportLocalGids; int *importProcs, *importToPart, *exportProcs, *exportToPart; float *wgt_list; int *gid_flags, *gid_list, *lid_list; /****************************************************************** ** Initialize MPI and Zoltan ******************************************************************/ MPI_Init(&argc, &argv); MPI_Comm_rank(MPI_COMM_WORLD, &myRank); MPI_Comm_size(MPI_COMM_WORLD, &numProcs); rc = Zoltan_Initialize(argc, argv, &ver); if (rc != ZOLTAN_OK){ MPI_Finalize(); exit(0); } \end{verbatim} \end{flushleft} \clearpage \begin{flushleft} \begin{verbatim} /****************************************************************** ** Create a Zoltan library structure for this instance of load ** balancing. Set the parameters and query functions that will ** govern the library's calculation. See the Zoltan User's ** Guide for the definition of these and many other parameters. ******************************************************************/ zz = Zoltan_Create(MPI_COMM_WORLD); /* General parameters */ Zoltan_Set_Param(zz, "DEBUG_LEVEL", "0"); Zoltan_Set_Param(zz, "LB_METHOD", "GRAPH"); Zoltan_Set_Param(zz, "NUM_GID_ENTRIES", "1"); Zoltan_Set_Param(zz, "NUM_LID_ENTRIES", "1"); Zoltan_Set_Param(zz, "OBJ_WEIGHT_DIM", "1"); Zoltan_Set_Param(zz, "RETURN_LISTS", "ALL"); /* Graph parameters */ Zoltan_Set_Param(zz, "PARMETIS_METHOD", "PARTKWAY"); Zoltan_Set_Param(zz, "PARMETIS_COARSE_ALG", "2"); Zoltan_Set_Param(zz, "CHECK_GRAPH", "2"); /* Query functions - defined in simpleQueries.h */ Zoltan_Set_Num_Obj_Fn(zz, get_number_of_objects, NULL); Zoltan_Set_Obj_List_Fn(zz, get_object_list, NULL); Zoltan_Set_Num_Edges_Multi_Fn(zz, get_num_edges_list, NULL); Zoltan_Set_Edge_List_Multi_Fn(zz, get_edge_list, NULL); \end{verbatim} \end{flushleft} \clearpage \begin{flushleft} \begin{verbatim} /****************************************************************** ** Zoltan can now partition the simple graph. ** In this simple example, we assume the number of partitions is ** equal to the number of processes. Process rank 0 will own ** partition 0, process rank 1 will own partition 1, and so on. ******************************************************************/ rc = Zoltan_LB_Partition(zz, /* input (all remaining fields are output) */ &changes, /* 1 if partitioning was changed, 0 otherwise */ &numGidEntries, /* Number of integers used for a global ID */ &numLidEntries, /* Number of integers used for a local ID */ &numImport, /* Number of vertices to be sent to me */ &importGlobalGids, /* Global IDs of vertices to be sent to me */ &importLocalGids, /* Local IDs of vertices to be sent to me */ &importProcs, /* Process rank for source of each incoming vertex */ &importToPart, /* New partition for each incoming vertex */ &numExport, /* Number of vertices I must send to other processes*/ &exportGlobalGids, /* Global IDs of the vertices I must send */ &exportLocalGids, /* Local IDs of the vertices I must send */ &exportProcs, /* Process to which I send each of the vertices */ &exportToPart); /* Partition to which each vertex will belong */ if (rc != ZOLTAN_OK){ printf("sorry...\n"); MPI_Finalize(); Zoltan_Destroy(&zz); exit(0); } /****************************************************************** ** In a real application, you would rebalance the problem now by ** sending the objects to their new partitions. Your query ** functions would need to reflect the new partitioning. ******************************************************************/ \end{verbatim} \end{flushleft} \clearpage \begin{flushleft} \begin{verbatim} /****************************************************************** ** Free the arrays allocated by Zoltan_LB_Partition, and free ** the storage allocated for the Zoltan structure. ******************************************************************/ Zoltan_LB_Free_Part(&importGlobalGids, &importLocalGids, &importProcs, &importToPart); Zoltan_LB_Free_Part(&exportGlobalGids, &exportLocalGids, &exportProcs, &exportToPart); Zoltan_Destroy(&zz); /********************** ** all done *********** **********************/ MPI_Finalize(); return 0; } \end{verbatim} \end{flushleft} \clearpage \subsection{A hypergraph problem} In this example, we use the simple graph in Figure~\ref{fig:mesh25} to create a hypergraph, and then partition it using Zoltan's parallel hypergraph partitioner. In a hypergraph, each edge can be associated with more than two vertices. This edge is called a hyperedge. For each vertex in the simple graph, we create a hyperedge with a vertex list composed of that vertex and all of its graph neighbors. So our hypergraph has 25 hyperedges, one for each vertex in the simple graph. We are creating this hypergraph in the hypergraph query functions. But you can use the Zoltan parameter \emph{PHG\_FROM\_GRAPH\_METHOD=neighbors} described in the User's Guide at \url{http://cs.sandia.gov/Zoltan/ug_html/ug_alg_phg.html}, to accomplish the same result. In this case you use graph query functions instead of hypergraph query funtions to supply the graph to Zoltan. The Zoltan library would then do the conversion to a hypergraph. First we will define the required query functions. \subsubsection{Application defined query functions} We must define six query functions. The function prototyes named below are defined in the Zoltan file \textbf{include/zoltan.h}. \begin{itemize} \item A function of type ZOLTAN\_NUM\_OBJ\_FN (Figure~\ref{fig:NumObj}) which provides the number of objects (vertices) belonging to this process . \item A function of type ZOLTAN\_OBJ\_LIST\_FN (Figure~\ref{fig:ObjList2}) which supplies the object global IDs and optional weights. It differs from the query function created for the last two examples (Figure~\ref{fig:ObjList}) because it omits the optional local IDs and weights. \item A function of type ZOLTAN\_HG\_SIZE\_CS\_FN (Figure~\ref{fig:SizeCS}). Each process will supply a portion of the hypergraph in a compressed storage format. This function tells the Zoltan library the size of the information that will be returned by the application process. \item A function of type ZOLTAN\_HG\_CS\_FN (Figure~\ref{fig:CS}) which supplies a portion of the hypergraph in a compressed storage format. \item A function of type ZOLTAN\_HG\_SIZE\_EDGE\_WTS\_FN (Figure~\ref{fig:SizeEW}) which supplies the number of hyperedges for which the process will provide edge weights. \item A function of type ZOLTAN\_HG\_EDGE\_WTS\_FN (Figure~\ref{fig:EW}) which supplies optional edge weights for the hyperedges. \end{itemize} An alternative to defining a single ZOLTAN\_OBJ\_LIST\_FN is to define a ZOLTAN\_FIRST\_OBJ\_FN and a ZOLTAN\_NEXT\_OBJ\_FN. The first function would return the global ID for first object. Each call to the second function would return a subsequent global ID. In this example the application is maintaining a small structure with information about the size of the hypergraph. The application supplies the address of this structure to the Zoltan library when setting up the query functions. Zoltan supplies the address of the structure to the application when calling the query functions. If you have the Zoltan source, the following query functions are defined in the file \textbf{examples/C/simpleQueries.h}. The variables used in these query functions refer to those defined in \textbf{examples/C/simpleGraph.h} which were shown in Figure~\ref{fig:simpleGraphDotH}. \begin{figure} \begin{flushleft} \begin{verbatim} static void get_hg_object_list(void *data, int sizeGID, int sizeLID, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int wgt_dim, float *obj_wgts, int *ierr) { int i, next; if ((sizeGID != 1) || (sizeLID != 0) || (wgt_dim != 0)) { *ierr = ZOLTAN_FATAL; return; } for (i=0, next=0; inumEdges = 0; hgd->numPins = 0; for (i=0; inumPins += simpleNumEdges[i] + 1; /* my hyperedge */ hgd->numEdges++; } } *num_lists = hgd->numEdges; *num_pins = hgd->numPins; *format = ZOLTAN_COMPRESSED_EDGE; *ierr = ZOLTAN_OK; return; } \end{verbatim} \end{flushleft} \caption{A ZOLTAN\_HG\_SIZE\_CS\_FN query function} \label{fig:SizeCS} \end{figure} \begin{figure} \begin{flushleft} \begin{verbatim} void get_hg(void *data, int sizeGID, int num_rows, int num_pins, int format, ZOLTAN_ID_PTR edge_GID, int *edge_ptr, ZOLTAN_ID_PTR pin_GID, int *ierr) { int i, j, npins; struct _hg_data *hgd; hgd = (struct _hg_data *)data; if ((num_rows != hgd->numEdges) || (num_pins != hgd->numPins) || (format != ZOLTAN_COMPRESSED_EDGE)){ *ierr = ZOLTAN_FATAL; return; } npins = 0; for (i=0; inumEdges; *ierr = ZOLTAN_OK; return; } \end{verbatim} \end{flushleft} \caption{A ZOLTAN\_HG\_SIZE\_EDGE\_WTS\_FN query function} \label{fig:SizeEW} \end{figure} \begin{figure} \begin{flushleft} \begin{verbatim} void get_hyperedge_weights(void *data, int sizeGID, int sizeLID, int num_edges, int edge_weight_dim, ZOLTAN_ID_PTR edge_GID, ZOLTAN_ID_PTR edge_LID, float *edge_weight, int *ierr) { int i; struct _hg_data *hgd; hgd = (struct _hg_data *)data; if ((sizeGID != 1) || (sizeLID != 0) || (num_edges != hgd->numEdges) || (edge_weight_dim != 1)){ *ierr = ZOLTAN_FATAL; return; } for (i=0; i #include #include #include "zoltan.h" #include "simpleGraph.h" #include "simpleQueries.h" int main(int argc, char *argv[]) { int rc, i, ngids, nextIdx; float ver; struct Zoltan_Struct *zz; int changes, numGidEntries, numLidEntries, numImport, numExport; ZOLTAN_ID_PTR importGlobalGids, importLocalGids; ZOLTAN_ID_PTR exportGlobalGids, exportLocalGids; int *importProcs, *importToPart, *exportProcs, *exportToPart; int *gid_flags, *gid_list; /****************************************************************** ** Initialize MPI and Zoltan ******************************************************************/ MPI_Init(&argc, &argv); MPI_Comm_rank(MPI_COMM_WORLD, &myRank); MPI_Comm_size(MPI_COMM_WORLD, &numProcs); rc = Zoltan_Initialize(argc, argv, &ver); if (rc != ZOLTAN_OK){ MPI_Finalize(); exit(0); } \end{verbatim} \end{flushleft} \clearpage \begin{flushleft} \begin{verbatim} /****************************************************************** ** Create a Zoltan library structure for this instance of load ** balancing. Set the parameters and query functions that will ** govern the library's calculation. See the Zoltan User's ** Guide for the definition of these and many other parameters. ******************************************************************/ zz = Zoltan_Create(MPI_COMM_WORLD); /* General parameters */ Zoltan_Set_Param(zz, "DEBUG_LEVEL", "0"); Zoltan_Set_Param(zz, "LB_METHOD", "HYPERGRAPH"); Zoltan_Set_Param(zz, "NUM_GID_ENTRIES", "1"); Zoltan_Set_Param(zz, "NUM_LID_ENTRIES", "0"); Zoltan_Set_Param(zz, "OBJ_WEIGHT_DIM", "0"); Zoltan_Set_Param(zz, "RETURN_LISTS", "ALL"); /* Hypergraph parameters */ Zoltan_Set_Param(zz, "HYPERGRAPH_PACKAGE", "PHG"); Zoltan_Set_Param(zz, "PHG_COARSENING_METHOD", "ipm"); Zoltan_Set_Param(zz, "PHG_COARSEPARTITION_METHOD", "greedy"); Zoltan_Set_Param(zz, "PHG_REFINEMENT_METHOD", "fm"); Zoltan_Set_Param(zz, "PHG_EDGE_WEIGHT_OPERATION", "error"); Zoltan_Set_Param(zz, "PHG_EDGE_SIZE_THRESHOLD", ".7"); Zoltan_Set_Param(zz, "ADD_OBJ_WEIGHT", "unit"); /* Query functions - defined in simpleQueries.h */ Zoltan_Set_Num_Obj_Fn(zz, get_number_of_objects, &hg_data); Zoltan_Set_Obj_List_Fn(zz, get_hg_object_list, &hg_data); Zoltan_Set_HG_Size_CS_Fn(zz, get_hg_size, &hg_data); Zoltan_Set_HG_CS_Fn(zz, get_hg, &hg_data); Zoltan_Set_HG_Size_Edge_Wts_Fn(zz, get_hg_num_edge_weights, &hg_data); Zoltan_Set_HG_Edge_Wts_Fn(zz, get_hyperedge_weights, &hg_data); \end{verbatim} \end{flushleft} \clearpage \begin{flushleft} \begin{verbatim} /****************************************************************** ** Zoltan can now partition the vertices. ** In this simple example, we assume the number of partitions is ** equal to the number of processes. Process rank 0 will own ** partition 0, process rank 1 will own partition 1, and so on. ******************************************************************/ rc = Zoltan_LB_Partition(zz, /* input (all remaining fields are output) */ &changes, /* 1 if partitioning was changed, 0 otherwise */ &numGidEntries, /* Number of integers used for a global ID */ &numLidEntries, /* Number of integers used for a local ID */ &numImport, /* Number of vertices to be sent to me */ &importGlobalGids, /* Global IDs of vertices to be sent to me */ &importLocalGids, /* Local IDs of vertices to be sent to me */ &importProcs, /* Process rank for source of each incoming vertex */ &importToPart, /* New partition for each incoming vertex */ &numExport, /* Number of vertices I must send to other processes*/ &exportGlobalGids, /* Global IDs of the vertices I must send */ &exportLocalGids, /* Local IDs of the vertices I must send */ &exportProcs, /* Process to which I send each of the vertices */ &exportToPart); /* Partition to which each vertex will belong */ if (rc != ZOLTAN_OK){ printf("sorry...\n"); MPI_Finalize(); Zoltan_Destroy(&zz); exit(0); } /****************************************************************** ** In a real application, you would rebalance the problem now by ** sending the objects to their new partitions. Your query ** functions would need to reflect the new partitioning. ******************************************************************/ \end{verbatim} \end{flushleft} \clearpage \begin{flushleft} \begin{verbatim} /****************************************************************** ** Free the arrays allocated by Zoltan_LB_Partition, and free ** the storage allocated for the Zoltan structure. ******************************************************************/ Zoltan_LB_Free_Part(&importGlobalGids, &importLocalGids, &importProcs, &importToPart); Zoltan_LB_Free_Part(&exportGlobalGids, &exportLocalGids, &exportProcs, &exportToPart); Zoltan_Destroy(&zz); /********************** ** all done *********** **********************/ MPI_Finalize(); return 0; } \end{verbatim} \end{flushleft} \clearpage \section{A C++ example} In this section we show a C++ example that solves the same recursive coordinate bisection problem that the C example in Section~\ref{sec:rcb} solves. We created a C++ class called \textbf{rectangularMesh}, which can represent the rectangular mesh shown in Figure~\ref{fig:mesh25}. The query functions required by Zoltan are implemented as static methods of the \textbf{rectangularMesh} class. They take as an argument an object of type \textbf{rectangularMesh}. We will not list the entire class due to space constraints. If you obtain the Zoltan source you will find the class in \textbf{examples/CPP/rectangularMesh.h}. But just to give you an idea, here is the method that returns the number of objects owned by the process: \begin{flushleft} \begin{verbatim} static int get_number_of_objects(void *data, int *ierr) { // Prototype: ZOLTAN_NUM_OBJ_FN // Return the number of objects I own. rectangularMesh *mesh = (rectangularMesh *)data; *ierr = ZOLTAN_OK; return mesh->get_num_my_IDs(); } \end{verbatim} \end{flushleft} We use static class methods, but your query functions could also be global C or C++ functions. The example below can be found in the Zoltan source in \textbf{examples/CPP/simpleRCB.cpp}. \clearpage \begin{flushleft} \begin{verbatim} #include int main(int argc, char *argv[]) { // Initialize MPI and Zoltan int rank, size; float version; MPI::Init(argc, argv); rank = MPI::COMM_WORLD.Get_rank(); size = MPI::COMM_WORLD.Get_size(); Zoltan_Initialize(argc, argv, &version); // Create a Zoltan object Zoltan *zz = new Zoltan(MPI::COMM_WORLD); if (zz == NULL){ MPI::Finalize(); exit(0); } // Create a simple rectangular mesh. It's vertices are the // objects to be partitioned. rectangularMesh *mesh = new rectangularMesh(); mesh->set_x_dim(10); mesh->set_y_dim(8); mesh->set_x_stride(1); mesh->set_y_stride(4); // General parameters: zz->Set_Param("DEBUG_LEVEL", "0"); // amount of debug messages desired zz->Set_Param("LB_METHOD", "RCB"); // recursive coordinate bisection zz->Set_Param("NUM_GID_ENTRIES", "1"); // number of integers in a global ID zz->Set_Param("NUM_LID_ENTRIES", "1"); // number of integers in a local ID zz->Set_Param("OBJ_WEIGHT_DIM", "1"); // dimension of a vertex weight zz->Set_Param("RETURN_LISTS", "ALL"); // return all lists in LB_Partition // RCB parameters: zz->Set_Param("KEEP_CUTS", "1"); // save decomposition zz->Set_Param("RCB_OUTPUT_LEVEL", "0"); // amount of output desired zz->Set_Param("RCB_RECTILINEAR_BLOCKS", "1"); // create rectilinear regions // Query functions: zz->Set_Num_Obj_Fn(rectangularMesh::get_number_of_objects, (void *)mesh); zz->Set_Obj_List_Fn(rectangularMesh::get_object_list, (void *)mesh); zz->Set_Num_Geom_Fn(rectangularMesh::get_num_geometry, NULL); zz->Set_Geom_Multi_Fn(rectangularMesh::get_geometry_list, (void *)mesh); // Zoltan can now partition the vertices in the simple mesh. int changes; int numGidEntries; int numLidEntries; int numImport; ZOLTAN_ID_PTR importGlobalIds; ZOLTAN_ID_PTR importLocalIds; int *importProcs; int *importToPart; int numExport; ZOLTAN_ID_PTR exportGlobalIds; ZOLTAN_ID_PTR exportLocalIds; int *exportProcs; int *exportToPart; int rc = zz->LB_Partition(changes, numGidEntries, numLidEntries, numImport, importGlobalIds, importLocalIds, importProcs, importToPart, numExport, exportGlobalIds, exportLocalIds, exportProcs, exportToPart); if (rc != ZOLTAN_OK){ printf("Partitioning failed on process %d\n",rank); MPI::Finalize(); delete zz; delete mesh; exit(0); } // In a real application, you would rebalance the problem now by // sending the objects to their new partitions. // Free the arrays allocated by LB_Partition, and free // the storage allocated for the Zoltan structure and the mesh. zz->LB_Free_Part(&importGlobalIds, &importLocalIds, &importProcs, &importToPart); zz->LB_Free_Part(&exportGlobalIds, &exportLocalIds, &exportProcs, &exportToPart); delete zz; delete mesh; // all done //////////////////////////////////////////////////// MPI::Finalize(); return 0; } \end{verbatim} \end{flushleft}