Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Allow empty partitions in parallel runs #162

Open
wants to merge 9 commits into
base: main
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
8 changes: 7 additions & 1 deletion assemble/Adapt_State.F90
Original file line number Diff line number Diff line change
Expand Up @@ -28,6 +28,7 @@
#include "fdebug.h"

module adapt_state_module
use iso_c_binding
use spud
use fldebug
use global_parameters, only : OPTION_PATH_LEN, periodic_boundary_option_path, adaptivity_mesh_name, adaptivity_mesh_name, domain_bbox, topology_mesh_name, FIELD_NAME_LEN
Expand Down Expand Up @@ -1063,9 +1064,12 @@ subroutine adapt_state_internal(states, metric, initialise_fields)
type(detector_type), pointer :: detector => null()

real :: global_min_quality, quality_tolerance
integer (c_int) :: MPI_COMM_SAVE

ewrite(1, *) "In adapt_state_internal"

MPI_COMM_SAVE=MPI_COMM_NONEMPTY

nullify(node_ownership)

max_adapt_iteration = adapt_iterations()
Expand Down Expand Up @@ -1157,7 +1161,7 @@ subroutine adapt_state_internal(states, metric, initialise_fields)

! Generate a new mesh field based on the current mesh field and the input
! metric
if (.not. vertical_only) then
if (.not. vertical_only .and. node_count(old_positions)>0) then
call adapt_mesh(old_positions, metric, new_positions, node_ownership = node_ownership, &
& force_preserve_regions=initialise_fields)
else
Expand Down Expand Up @@ -1433,6 +1437,8 @@ subroutine adapt_state_internal(states, metric, initialise_fields)
call compute_domain_statistics(states)
end if

call free_communicator(MPI_COMM_SAVE)

ewrite(1, *) "Exiting adapt_state_internal"

end subroutine adapt_state_internal
Expand Down
2 changes: 1 addition & 1 deletion assemble/Surface_Id_Interleaving.F90
Original file line number Diff line number Diff line change
Expand Up @@ -105,7 +105,7 @@ subroutine interleave_surface_ids_vector(mesh, interleaved_surface_ids, max_copl
if(isparallel()) then
#ifdef HAVE_MPI
! Max. coplanar_id must be global to ensure consistent global surface ids
call mpi_allreduce(max_coplanar_id, all_max_coplanar_id, 1, getpinteger(), MPI_MAX, MPI_COMM_FEMTOOLS, ierr)
call mpi_allreduce(max_coplanar_id, all_max_coplanar_id, 1, getpinteger(), MPI_MAX, MPI_COMM_NONEMPTY, ierr)
assert(ierr == MPI_SUCCESS)
max_coplanar_id = all_max_coplanar_id
#endif
Expand Down
22 changes: 11 additions & 11 deletions assemble/Zoltan_callbacks.F90
Original file line number Diff line number Diff line change
Expand Up @@ -11,7 +11,7 @@ module zoltan_callbacks
use fldebug
use data_structures
use mpi_interfaces
use parallel_tools, only: getrank, getnprocs, getprocno, MPI_COMM_FEMTOOLS
use parallel_tools, only: getrank, getnprocs, getprocno, MPI_COMM_FEMTOOLS, MPI_COMM_NONEMPTY
use sparse_tools
use element_numbering
use elements
Expand Down Expand Up @@ -132,7 +132,7 @@ subroutine zoltan_cb_get_num_edges(data, num_gid_entries, num_lid_entries, num_o
assert(num_gid_entries == 1)
assert(num_lid_entries == 1)

count = zoltan_global_zz_halo%nowned_nodes
count = halo_nowned_nodes(zoltan_global_zz_halo)
assert(count == num_obj)

do node=1,count
Expand Down Expand Up @@ -186,8 +186,8 @@ subroutine zoltan_cb_get_edge_list(data, num_gid_entries, num_lid_entries, num_o
assert(num_gid_entries == 1)
assert(num_lid_entries == 1)
assert(wgt_dim == 1)
count = zoltan_global_zz_halo%nowned_nodes

count = halo_nowned_nodes(zoltan_global_zz_halo)
assert(count == num_obj)

my_num_edges = sum(num_edges(1:num_obj))
Expand All @@ -213,14 +213,14 @@ subroutine zoltan_cb_get_edge_list(data, num_gid_entries, num_lid_entries, num_o
! find global ids for each neighbour
nbor_global_id(head:head+size(neighbours)-1) = halo_universal_number(zoltan_global_zz_halo, neighbours)
! find owning proc for each neighbour
nbor_procs(head:head+size(neighbours)-1) = halo_node_owners(zoltan_global_zz_halo, neighbours) - 1
nbor_procs(head:head+size(neighbours)-1) = global_proc_no(halo_node_owners(zoltan_global_zz_halo, neighbours)) - 1
head = head + size(neighbours)
end do
ierr = ZOLTAN_OK
return
else
call MPI_ALLREDUCE(my_num_edges,total_num_edges,1,MPI_INTEGER,MPI_SUM, &
MPI_COMM_FEMTOOLS,err)
MPI_COMM_NONEMPTY,err)
end if

head = 1
Expand All @@ -242,7 +242,7 @@ subroutine zoltan_cb_get_edge_list(data, num_gid_entries, num_lid_entries, num_o
nbor_global_id(head:head+size(neighbours)-1) = halo_universal_number(zoltan_global_zz_halo, neighbours)

! find owning proc for each neighbour
nbor_procs(head:head+size(neighbours)-1) = halo_node_owners(zoltan_global_zz_halo, neighbours) - 1
nbor_procs(head:head+size(neighbours)-1) = global_proc_no(halo_node_owners(zoltan_global_zz_halo, neighbours))-1

! get elements associated with current node
my_nelist => row_m_ptr(zoltan_global_zz_nelist, local_ids(node))
Expand Down Expand Up @@ -299,10 +299,10 @@ subroutine zoltan_cb_get_edge_list(data, num_gid_entries, num_lid_entries, num_o
my_min_weight = minval(ewgts(1:head-1))

! calculate global maximum edge weight
call MPI_ALLREDUCE(my_max_weight,max_weight,1,MPI_REAL,MPI_MAX, MPI_COMM_FEMTOOLS,err)
call MPI_ALLREDUCE(my_max_weight,max_weight,1,MPI_REAL,MPI_MAX, MPI_COMM_NONEMPTY,err)

! calculate global minimum edge weight
call MPI_ALLREDUCE(my_min_weight,min_weight,1,MPI_REAL,MPI_MIN, MPI_COMM_FEMTOOLS,err)
call MPI_ALLREDUCE(my_min_weight,min_weight,1,MPI_REAL,MPI_MIN, MPI_COMM_NONEMPTY,err)

! calculate the local 90th percentile edge weight
ninety_weight = max_weight * 0.90
Expand Down Expand Up @@ -427,7 +427,7 @@ subroutine zoltan_cb_pack_nodes(data, num_gid_entries, num_lid_entries, num_ids,
buf(head:head+row_length(zoltan_global_zz_sparsity_two, node)-1) = halo_universal_number(zoltan_global_zz_halo, row_m_ptr(zoltan_global_zz_sparsity_two, node))
head = head + row_length(zoltan_global_zz_sparsity_two, node)

buf(head:head+row_length(zoltan_global_zz_sparsity_two, node)-1) = halo_node_owners(zoltan_global_zz_halo, row_m_ptr(zoltan_global_zz_sparsity_two, node))
buf(head:head+row_length(zoltan_global_zz_sparsity_two, node)-1) = global_proc_no(halo_node_owners(zoltan_global_zz_halo, row_m_ptr(zoltan_global_zz_sparsity_two, node)))
head = head + row_length(zoltan_global_zz_sparsity_two, node)

buf(head) = row_length(zoltan_global_zz_nelist, node)
Expand Down Expand Up @@ -509,7 +509,7 @@ subroutine zoltan_cb_unpack_nodes(data, num_gid_entries, num_ids, global_ids, si
universal_number = halo_universal_number(zoltan_global_zz_halo, neighbours(j))
call insert(zoltan_global_new_nodes, universal_number, changed=changed)
if (changed) then ! so it is a halo node
old_owner = halo_node_owner(zoltan_global_zz_halo, neighbours(j)) - 1
old_owner = global_proc_no(halo_node_owner(zoltan_global_zz_halo, neighbours(j))) - 1
assert(old_owner < getnprocs())
if (old_owner == rank) then
call insert(halo_nodes_we_currently_own, neighbours(j))
Expand Down
1 change: 1 addition & 0 deletions assemble/Zoltan_global_variables.F90
Original file line number Diff line number Diff line change
Expand Up @@ -39,6 +39,7 @@ module zoltan_global_variables
type(scalar_field), save, pointer, public :: zoltan_global_max_edge_weight_on_node
logical, save, public :: zoltan_global_output_edge_weights = .false.
type(csr_sparsity), save, pointer, public :: zoltan_global_zz_nelist
integer, dimension(:), save, public, allocatable :: global_proc_no

! Needed for zoltan_cb_pack_node_sizes
! - added vector_field to use fields
Expand Down
96 changes: 75 additions & 21 deletions assemble/Zoltan_integration.F90
Original file line number Diff line number Diff line change
Expand Up @@ -5,6 +5,7 @@ module zoltan_integration

#ifdef HAVE_ZOLTAN
use spud
use iso_c_binding
use fldebug
use global_parameters, only: real_size, OPTION_PATH_LEN, topology_mesh_name,&
FIELD_NAME_LEN
Expand Down Expand Up @@ -101,7 +102,7 @@ subroutine zoltan_drive(states, final_adapt_iteration, global_min_quality, metri
real :: load_imbalance_tolerance
logical :: flredecomp
real :: minimum_quality
integer :: flredecomp_input_procs = -1, flredecomp_target_procs = -1
integer :: flredecomp_input_procs = -1, flredecomp_target_procs = -1, ierr

ewrite(1,*) "In zoltan_drive"

Expand Down Expand Up @@ -140,6 +141,13 @@ subroutine zoltan_drive(states, final_adapt_iteration, global_min_quality, metri
zoltan_global_field_weighted_partitions = &
have_option(trim(zoltan_global_base_option_path) // "/field_weighted_partitions")

if (.not. allocated(global_proc_no)) then
allocate(global_proc_no(getnprocs(MPI_COMM_NONEMPTY)))
call MPI_ALLGATHER(getprocno(MPI_COMM_WORLD),1,MPI_INTEGER, &
global_proc_no,1,MPI_INTEGER, &
MPI_COMM_NONEMPTY,ierr)
end if

call setup_module_variables(states, final_adapt_iteration, zz, flredecomp)

call setup_quality_module_variables(metric, minimum_quality) ! this needs to be called after setup_module_variables
Expand Down Expand Up @@ -222,7 +230,8 @@ subroutine zoltan_drive(states, final_adapt_iteration, global_min_quality, metri

call reconstruct_enlist
call reconstruct_senlist
call reconstruct_halo(zz)
call reconstruct_halo(zz, flredecomp, flredecomp_input_procs, &
flredecomp_target_procs)

if(zoltan_global_migrate_extruded_mesh) then
zoltan_global_new_positions_m1d = zoltan_global_new_positions ! save a reference to the horizontal mesh you've just load balanced
Expand Down Expand Up @@ -265,7 +274,8 @@ subroutine zoltan_drive(states, final_adapt_iteration, global_min_quality, metri

call reconstruct_enlist
call reconstruct_senlist
call reconstruct_halo(zz)
call reconstruct_halo(zz, flredecomp, flredecomp_input_procs, &
flredecomp_target_procs)

if (.not. verify_consistent_local_element_numbering(zoltan_global_new_positions%mesh) ) then
ewrite(-1,*) "For the extruded mesh, the local element numbering of elements in the halo region" // &
Expand Down Expand Up @@ -367,7 +377,7 @@ subroutine setup_module_variables(states, final_adapt_iteration, zz, flredecomp,

zoltan_global_zz_nelist => extract_nelist(zoltan_global_zz_mesh)

zz => Zoltan_Create(halo_communicator(zoltan_global_zz_mesh))
zz => Zoltan_Create(MPI_COMM_FEMTOOLS)

nhalos = halo_count(zoltan_global_zz_mesh)
assert(nhalos == 2)
Expand Down Expand Up @@ -395,7 +405,7 @@ subroutine setup_module_variables(states, final_adapt_iteration, zz, flredecomp,
call insert(zoltan_global_old_local_numbering_to_uen, i, halo_universal_number(zoltan_global_zz_ele_halo, i))
end do

allocate(zoltan_global_receives(halo_proc_count(zoltan_global_zz_halo)))
allocate(zoltan_global_receives(getnprocs(MPI_COMM_FEMTOOLS)))
do i=1,size(zoltan_global_receives)
call allocate(zoltan_global_receives(i))
end do
Expand Down Expand Up @@ -858,6 +868,7 @@ subroutine zoltan_load_balance(zz, changes, num_gid_entries, num_lid_entries, &
integer :: min_num_nodes_after_balance, total_num_nodes_before_balance, total_num_nodes_after_balance
integer :: num_empty_partitions, empty_partition
character (len = 10) :: string_load_imbalance_tolerance
logical :: flag

ewrite(1,*) 'in zoltan_load_balance'

Expand Down Expand Up @@ -927,10 +938,16 @@ subroutine zoltan_load_balance(zz, changes, num_gid_entries, num_lid_entries, &
else

min_num_nodes_after_balance = 0
do while (min_num_nodes_after_balance == 0)

flag =.true.
do while (min_num_nodes_after_balance == 0 .and. flag)

ierr = Zoltan_LB_Balance(zz, changes, num_gid_entries, num_lid_entries, p1_num_import, p1_import_global_ids, &
& p1_import_local_ids, p1_import_procs, p1_num_export, p1_export_global_ids, p1_export_local_ids, p1_export_procs)
ierr = Zoltan_LB_Partition(zz, changes, num_gid_entries, num_lid_entries, p1_num_import, p1_import_global_ids, &
& p1_import_local_ids, p1_import_procs, import_to_part, p1_num_export, p1_export_global_ids, &
& p1_export_local_ids, p1_export_procs, export_to_part)
! ierr = Zoltan_LB_Balance(zz, changes, num_gid_entries, num_lid_entries, p1_num_import, p1_import_global_ids, &
! & p1_import_local_ids, p1_import_procs, p1_num_export, p1_export_global_ids, p1_export_local_ids, p1_export_procs)
assert(ierr == ZOLTAN_OK)

! calculate how many owned nodes we'd have after doing the planned load balancing
Expand Down Expand Up @@ -967,7 +984,8 @@ subroutine zoltan_load_balance(zz, changes, num_gid_entries, num_lid_entries, &
assert(ierr == MPI_SUCCESS)

if (min_num_nodes_after_balance == 0) then
FLAbort("Could not stop Zoltan creating empty partitions.")
ewrite(-1,*) "Could not stop Zoltan creating empty partitions."
flag = .false.
else
ewrite(-1,*) 'Load balancing was carried out without edge-weighting being applied. Mesh may not be of expected quality.'
end if
Expand Down Expand Up @@ -1200,7 +1218,7 @@ subroutine deal_with_exporters
universal_number = halo_universal_number(zoltan_global_zz_halo, neighbours(j))
call insert(zoltan_global_new_nodes, universal_number, changed=changed)
if (changed) then
old_owner = halo_node_owner(zoltan_global_zz_halo, neighbours(j)) - 1
old_owner = global_proc_no(halo_node_owner(zoltan_global_zz_halo, neighbours(j))) - 1
if (old_owner == rank) then
call insert(halo_nodes_we_currently_own, neighbours(j))
else
Expand Down Expand Up @@ -1623,7 +1641,7 @@ subroutine reconstruct_senlist
ewrite(1,*) "Exiting reconstruct_senlist"
end subroutine reconstruct_senlist

subroutine reconstruct_halo(zz)
subroutine reconstruct_halo(zz, flredecomp, input_procs, target_procs)
! At this point, the receives sets have been populated with all
! the universal node numbers we need to receive from each process.
! So, we are going to use zoltan to invert this to compute
Expand All @@ -1633,22 +1651,28 @@ subroutine reconstruct_halo(zz)
! l2e halo.
! Supply the peeps with jeeps, brick apiece, capiche?

type(zoltan_struct), pointer, intent(in) :: zz
type(zoltan_struct), pointer, intent(in) :: zz
logical, intent(in) :: flredecomp
integer, intent(in) :: input_procs, target_procs

integer :: num_import, num_export
integer :: num_import, num_export, nprocs
integer, dimension(:), pointer :: import_global_ids, import_local_ids, import_procs
integer, dimension(:), pointer :: export_global_ids, export_local_ids, export_procs, export_to_part
integer :: ierr, i, head
integer :: ierr, i, j, head
type(integer_set), dimension(size(zoltan_global_receives)) :: sends
integer, dimension(size(zoltan_global_receives)) :: nreceives, nsends
integer, dimension(:), allocatable :: nonempty
integer, dimension(ele_count(zoltan_global_new_positions)) :: ele_renumber_permutation
integer, dimension(node_count(zoltan_global_new_positions)) :: node_renumber_permutation
integer :: universal_element_number, old_new_local_element_number, new_new_local_element_number
integer :: universal_node_number, old_new_local_node_number, new_new_local_node_number

integer, dimension(ele_count(zoltan_global_new_positions)) :: old_new_region_ids

integer (c_int) :: comm

ewrite(1,*) "In reconstruct_halo"

allocate(nonempty(getnprocs(MPI_COMM_FEMTOOLS)))

num_import = 0
do i=1,size(zoltan_global_receives)
Expand Down Expand Up @@ -1706,19 +1730,46 @@ subroutine reconstruct_halo(zz)

! Allocate the halo and such
! We had to grow dreads to change our description, two cops is on a milkbox, missin'

call reset_next_mpi_tag()
allocate(zoltan_global_new_positions%mesh%halos(2))
if (allocated(global_proc_no)) deallocate(global_proc_no)

if (flredecomp .and. input_procs>target_procs) then
allocate(global_proc_no(input_procs))
global_proc_no = [(i,i=1,input_procs)]
nprocs = input_procs
comm = MPI_COMM_FEMTOOLS
nonempty = 1
! In this context we don't need to worry about nonempty partitions
call split_communicator(MPI_COMM_FEMTOOLS, MPI_COMM_NONEMPTY, .true.)
else
call split_communicator(MPI_COMM_FEMTOOLS, MPI_COMM_NONEMPTY,&
node_count(zoltan_global_new_positions)>0)
allocate(global_proc_no(getnprocs(MPI_COMM_NONEMPTY)))
call MPI_ALLGATHER(getprocno(MPI_COMM_WORLD),1,MPI_INTEGER, &
global_proc_no,1,MPI_INTEGER, &
MPI_COMM_NONEMPTY,ierr)
nprocs = getnprocs(MPI_COMM_NONEMPTY)
comm = MPI_COMM_NONEMPTY
call MPI_Allgather(merge(1,0,node_count(zoltan_global_new_positions)>0), 1, MPI_INT,&
nonempty, 1, MPI_INT, MPI_COMM_FEMTOOLS, ierr)
end if

call allocate(zoltan_global_new_positions%mesh%halos(2), &
nsends = nsends, &
nreceives = nreceives, &
nsends = pack(nsends, nonempty == nonempty(getprocno(MPI_COMM_FEMTOOLS))), &
nreceives = pack(nreceives, nonempty == nonempty(getprocno(MPI_COMM_FEMTOOLS))), &
nprocs = nprocs, &
name = halo_name(zoltan_global_zz_halo), &
communicator = halo_communicator(zoltan_global_zz_halo), &
communicator = comm, &
nowned_nodes = key_count(zoltan_global_new_nodes) - num_import, &
data_type = halo_data_type(zoltan_global_zz_halo))

j=1
do i=1,size(zoltan_global_receives)
call set_halo_sends(zoltan_global_new_positions%mesh%halos(2), i, fetch(zoltan_global_universal_to_new_local_numbering, set2vector(sends(i))))
call set_halo_receives(zoltan_global_new_positions%mesh%halos(2), i, fetch(zoltan_global_universal_to_new_local_numbering, set2vector(zoltan_global_receives(i))))
if (nonempty(i) == 0) cycle
call set_halo_sends(zoltan_global_new_positions%mesh%halos(2), j, fetch(zoltan_global_universal_to_new_local_numbering, set2vector(sends(i))))
call set_halo_receives(zoltan_global_new_positions%mesh%halos(2), j, fetch(zoltan_global_universal_to_new_local_numbering, set2vector(zoltan_global_receives(i))))
j=j+1
end do

! Now derive all the other halos ...
Expand Down Expand Up @@ -1772,6 +1823,7 @@ subroutine reconstruct_halo(zz)
end do

call reorder_element_numbering(zoltan_global_new_positions)


ewrite(1,*) "Exiting reconstruct_halo"

Expand Down Expand Up @@ -2040,7 +2092,7 @@ subroutine transfer_fields(zz)
integer :: old_ele
integer, dimension(:), pointer :: old_local_nodes, nodes
type(element_type), pointer :: eshape
type(integer_set), dimension(halo_proc_count(zoltan_global_zz_halo)) :: sends
type(integer_set), dimension(:), allocatable :: sends
integer :: i, j, new_owner, universal_element_number
type(integer_set) :: self_sends
integer :: num_import, num_export
Expand All @@ -2061,6 +2113,8 @@ subroutine transfer_fields(zz)
type(detector_type), pointer :: detector => null(), add_detector => null()

ewrite(1,*) 'in transfer_fields'

allocate(sends(getnprocs(MPI_COMM_FEMTOOLS)))

do i=1,size(sends)
call allocate(sends(i))
Expand Down
Loading