int build_elem_dd(MESH_INFO_PTR mesh) { /* Create a distributed directory of the elements so we can track their * processor assignment after migrations. */ int maxelems; MPI_Allreduce(&(mesh->num_elems), &maxelems, 1, MPI_INT, MPI_MAX, MPI_COMM_WORLD); if (Zoltan_DD_Create(&(mesh->dd), MPI_COMM_WORLD, 1, 0, 0, maxelems, 0) != 0){ Gen_Error(0, "fatal: NULL returned from Zoltan_DD_Create()\n"); return 0; } return update_elem_dd(mesh); }
void migrate_post_process(void *data, int num_gid_entries, int num_lid_entries, int num_import, ZOLTAN_ID_PTR import_global_ids, ZOLTAN_ID_PTR import_local_ids, int *import_procs, int *import_to_part, int num_export, ZOLTAN_ID_PTR export_global_ids, ZOLTAN_ID_PTR export_local_ids, int *export_procs, int *export_to_part, int *ierr) { if (data == NULL) { *ierr = ZOLTAN_FATAL; return; } MESH_INFO_PTR mesh = (MESH_INFO_PTR) data; ELEM_INFO *elements = mesh->elements; int proc = 0, num_proc = 0; MPI_Comm_rank(MPI_COMM_WORLD, &proc); MPI_Comm_size(MPI_COMM_WORLD, &num_proc); /* compact elements array, as the application expects the array to be dense */ for (int i = 0; i < New_Elem_Index_Size; i++) { if (New_Elem_Index[i] != ZOLTAN_ID_INVALID) continue; /* Don't want to shift all elements down one position to fill the */ /* blank spot -- too much work to adjust adjacencies! So find the */ /* last element in the array and move it to the blank spot. */ int last = 0; for (last = New_Elem_Index_Size-1; last >= 0; last--) if (New_Elem_Index[last] != ZOLTAN_ID_INVALID) break; /* If (last < i), array is already dense; i is just in some blank spots */ /* at the end of the array. Quit the compacting. */ if (last < i) break; /* Copy elements[last] to elements[i]. */ elements[i] = elements[last]; /* Adjust adjacencies for local elements. Off-processor adjacencies */ /* don't matter here. */ for (int j = 0; j < elements[i].adj_len; j++) { /* Skip NULL adjacencies (sides that are not adjacent to another elem). */ if (elements[i].adj[j] == ZOLTAN_ID_INVALID) continue; ZOLTAN_ID_TYPE adj_elem = elements[i].adj[j]; /* See whether adjacent element is local; if so, adjust its entry */ /* for local element i. */ if (elements[i].adj_proc[j] == proc) { for (int k = 0; k < elements[adj_elem].adj_len; k++) { if (elements[adj_elem].adj[k] == (ZOLTAN_ID_TYPE)last && elements[adj_elem].adj_proc[k] == proc) { /* found adjacency entry for element last; change it to i */ elements[adj_elem].adj[k] = (ZOLTAN_ID_TYPE)i; break; } } } } /* Update New_Elem_Index */ New_Elem_Index[i] = New_Elem_Index[last]; New_Elem_Index[last] = ZOLTAN_ID_INVALID; /* clear elements[last] */ elements[last].globalID = ZOLTAN_ID_INVALID; elements[last].border = 0; elements[last].my_part = -1; elements[last].perm_value = -1; elements[last].invperm_value = -1; elements[last].nadj = 0; elements[last].adj_len = 0; elements[last].elem_blk = -1; for (int k=0; k<MAX_CPU_WGTS; k++) elements[last].cpu_wgt[k] = 0; elements[last].mem_wgt = 0; elements[last].avg_coord[0] = elements[last].avg_coord[1] = elements[last].avg_coord[2] = 0.; elements[last].coord = NULL; elements[last].connect = NULL; elements[last].adj = NULL; elements[last].adj_proc = NULL; elements[last].edge_wgt = NULL; } if (New_Elem_Index != NULL) { delete [] New_Elem_Index; New_Elem_Index = NULL; } New_Elem_Index_Size = 0; if (!build_elem_comm_maps(proc, mesh)) { Gen_Error(0, "Fatal: error rebuilding elem comm maps"); } if (mesh->data_type == HYPERGRAPH && !update_elem_dd(mesh)) { Gen_Error(0, "Fatal: error updating element dd"); } if (mesh->data_type == HYPERGRAPH && mesh->hvertex_proc && !update_hvertex_proc(mesh)) { Gen_Error(0, "Fatal: error updating hyperedges"); } if (Debug_Driver > 3) print_distributed_mesh(proc, num_proc, mesh); }