[6237] | 1 | module Agrif_seq |
---|
| 2 | ! |
---|
| 3 | use Agrif_Init |
---|
| 4 | use Agrif_Procs |
---|
| 5 | use Agrif_Arrays |
---|
| 6 | ! |
---|
| 7 | implicit none |
---|
| 8 | ! |
---|
| 9 | contains |
---|
| 10 | ! |
---|
| 11 | #if defined AGRIF_MPI |
---|
| 12 | !=================================================================================================== |
---|
| 13 | function Agrif_seq_allocate_list ( nb_seqs ) result( seqlist ) |
---|
| 14 | !--------------------------------------------------------------------------------------------------- |
---|
| 15 | integer, intent(in) :: nb_seqs |
---|
| 16 | ! |
---|
| 17 | type(Agrif_Sequence_List), pointer :: seqlist |
---|
| 18 | ! |
---|
| 19 | allocate(seqlist) |
---|
| 20 | seqlist % nb_seqs = nb_seqs |
---|
| 21 | allocate(seqlist % sequences(1:nb_seqs)) |
---|
| 22 | !--------------------------------------------------------------------------------------------------- |
---|
| 23 | end function Agrif_seq_allocate_list |
---|
| 24 | !=================================================================================================== |
---|
| 25 | ! |
---|
| 26 | !=================================================================================================== |
---|
| 27 | subroutine Agrif_seq_add_grid ( seqlist, seq_num, grid ) |
---|
| 28 | !--------------------------------------------------------------------------------------------------- |
---|
| 29 | type(Agrif_Sequence_List), intent(inout) :: seqlist |
---|
| 30 | integer, intent(in) :: seq_num |
---|
| 31 | type(Agrif_Grid), pointer, intent(in) :: grid |
---|
| 32 | ! |
---|
| 33 | call Agrif_gl_append(seqlist % sequences(seq_num) % gridlist, grid ) |
---|
| 34 | !--------------------------------------------------------------------------------------------------- |
---|
| 35 | end subroutine Agrif_seq_add_grid |
---|
| 36 | !=================================================================================================== |
---|
| 37 | ! |
---|
| 38 | !=================================================================================================== |
---|
| 39 | recursive subroutine Agrif_seq_init_sequences ( grid ) |
---|
| 40 | !--------------------------------------------------------------------------------------------------- |
---|
| 41 | type(Agrif_Grid), intent(inout) :: grid |
---|
| 42 | ! |
---|
| 43 | type(Agrif_PGrid), pointer :: gp |
---|
| 44 | ! |
---|
| 45 | #if defined AGRIF_MPI |
---|
| 46 | ! |
---|
| 47 | ! Build list of required procs for each child grid |
---|
| 48 | gp => grid % child_list % first |
---|
| 49 | do while ( associated( gp ) ) |
---|
| 50 | call Agrif_seq_build_required_proclist( gp % gr ) |
---|
| 51 | gp => gp % next |
---|
| 52 | enddo |
---|
| 53 | ! |
---|
| 54 | ! Create integration sequences for the current grid |
---|
| 55 | call Agrif_seq_create_proc_sequences( grid ) |
---|
| 56 | call Agrif_seq_allocate_procs_to_childs( grid ) |
---|
| 57 | ! |
---|
| 58 | ! Create new communicators for sequences |
---|
| 59 | call Agrif_seq_create_communicators( grid ) |
---|
| 60 | ! |
---|
| 61 | #endif |
---|
| 62 | !--------------------------------------------------------------------------------------------------- |
---|
| 63 | end subroutine Agrif_seq_init_sequences |
---|
| 64 | !=================================================================================================== |
---|
| 65 | ! |
---|
| 66 | !=================================================================================================== |
---|
| 67 | subroutine Agrif_seq_build_required_proclist ( grid ) |
---|
| 68 | !--------------------------------------------------------------------------------------------------- |
---|
| 69 | type(Agrif_Grid), intent(inout) :: grid |
---|
| 70 | ! |
---|
| 71 | type(Agrif_Grid), pointer :: parent_grid |
---|
| 72 | type(Agrif_Rectangle), pointer :: grid_rect |
---|
| 73 | type(Agrif_Proc_p), pointer :: proc_rect |
---|
| 74 | type(Agrif_Proc), pointer :: proc |
---|
| 75 | logical :: proc_is_required |
---|
| 76 | integer :: i |
---|
| 77 | ! |
---|
| 78 | if ( grid % fixedrank == 0 ) then |
---|
| 79 | ! grid is the Root |
---|
| 80 | if ( grid % required_proc_list % nitems == 0 ) then |
---|
| 81 | print*, "### Error Agrif_seq_build_required_proclist: empty proc list." |
---|
| 82 | print*, "# -> You should check if Agrif_Init_ProcList() is actually working." |
---|
| 83 | stop |
---|
| 84 | endif |
---|
| 85 | return |
---|
| 86 | endif |
---|
| 87 | ! |
---|
| 88 | parent_grid => grid % parent |
---|
| 89 | grid_rect => grid % rect_in_parent |
---|
| 90 | proc_rect => parent_grid % proc_def_list % first |
---|
| 91 | ! |
---|
| 92 | |
---|
| 93 | do while ( associated( proc_rect ) ) |
---|
| 94 | |
---|
| 95 | ! |
---|
| 96 | proc => proc_rect % proc |
---|
| 97 | ! |
---|
| 98 | proc_is_required = .true. |
---|
| 99 | do i = 1,Agrif_Probdim |
---|
| 100 | proc_is_required = ( proc_is_required .and. & |
---|
| 101 | ( grid_rect % imin(i) <= proc % imax(i) ) .and. & |
---|
| 102 | ( grid_rect % imax(i) >= proc % imin(i) ) ) |
---|
| 103 | enddo |
---|
| 104 | ! |
---|
| 105 | if ( proc_is_required ) then |
---|
| 106 | call Agrif_pl_append(grid % required_proc_list, proc) |
---|
| 107 | endif |
---|
| 108 | proc_rect => proc_rect % next |
---|
| 109 | ! |
---|
| 110 | enddo |
---|
| 111 | !--------------------------------------------------------------------------------------------------- |
---|
| 112 | end subroutine Agrif_seq_build_required_proclist |
---|
| 113 | !=================================================================================================== |
---|
| 114 | ! |
---|
| 115 | !=================================================================================================== |
---|
| 116 | subroutine Agrif_seq_create_proc_sequences ( grid ) |
---|
| 117 | !--------------------------------------------------------------------------------------------------- |
---|
| 118 | type(Agrif_Grid), intent(inout) :: grid |
---|
| 119 | ! |
---|
| 120 | type(Agrif_Grid_List), pointer :: sorted_child_list |
---|
| 121 | type(Agrif_PGrid), pointer :: child_p |
---|
| 122 | type(Agrif_PGrid), pointer :: g1p, g2p |
---|
| 123 | type(Agrif_Proc_p), pointer :: pp1, pp2 |
---|
| 124 | type(Agrif_Proc), pointer :: proc |
---|
| 125 | integer :: nb_seq_max, nb_seqs, cur_seq |
---|
| 126 | ! |
---|
| 127 | nb_seq_max = 0 |
---|
| 128 | ! |
---|
| 129 | if ( grid % child_list % nitems == 0 ) return |
---|
| 130 | ! |
---|
| 131 | ! For each required proc... |
---|
| 132 | pp1 => grid % required_proc_list % first |
---|
| 133 | do while ( associated(pp1) ) |
---|
| 134 | proc => pp1 % proc |
---|
| 135 | proc % nb_seqs = 0 |
---|
| 136 | ! ..loop over all child grids... |
---|
| 137 | child_p => grid % child_list % first |
---|
| 138 | do while ( associated(child_p) ) |
---|
| 139 | ! ..and look for 'proc' in the list of procs required by 'child' |
---|
| 140 | pp2 => child_p % gr % required_proc_list % first |
---|
| 141 | do while ( associated(pp2) ) |
---|
| 142 | if ( proc % pn == pp2 % proc % pn ) then |
---|
| 143 | ! 'proc' is required by this child grid, so we increment its number of sequences |
---|
| 144 | proc % nb_seqs = proc % nb_seqs + 1 |
---|
| 145 | pp2 => NULL() |
---|
| 146 | else |
---|
| 147 | pp2 => pp2 % next |
---|
| 148 | endif |
---|
| 149 | enddo |
---|
| 150 | child_p => child_p % next |
---|
| 151 | enddo |
---|
| 152 | nb_seq_max = max(nb_seq_max, proc % nb_seqs) |
---|
| 153 | pp1 => pp1 % next |
---|
| 154 | enddo |
---|
| 155 | ! |
---|
| 156 | ! For each grid... |
---|
| 157 | g1p => grid % child_list % first |
---|
| 158 | do while ( associated(g1p) ) |
---|
| 159 | ! compare it with the following ones |
---|
| 160 | g2p => g1p % next |
---|
| 161 | do while ( associated(g2p) ) |
---|
| 162 | if ( Agrif_seq_grids_are_connected( g1p % gr, g2p % gr ) ) then |
---|
| 163 | call Agrif_gl_append( g1p % gr % neigh_list, g2p % gr ) |
---|
| 164 | call Agrif_gl_append( g2p % gr % neigh_list, g1p % gr ) |
---|
| 165 | endif |
---|
| 166 | g2p => g2p % next |
---|
| 167 | enddo |
---|
| 168 | g1p => g1p % next |
---|
| 169 | enddo |
---|
| 170 | ! |
---|
| 171 | ! Colorize graph nodes |
---|
| 172 | nb_seqs = Agrif_seq_colorize_grid_list(grid % child_list) |
---|
| 173 | sorted_child_list => Agrif_gl_merge_sort ( grid % child_list, compare_colors ) |
---|
| 174 | ! |
---|
| 175 | ! Create sequence structure |
---|
| 176 | cur_seq = 0 |
---|
| 177 | grid % child_seq => Agrif_seq_allocate_list(nb_seqs) |
---|
| 178 | child_p => sorted_child_list % first |
---|
| 179 | do while ( associated(child_p) ) |
---|
| 180 | if ( cur_seq /= child_p % gr % seq_num ) then |
---|
| 181 | cur_seq = child_p % gr % seq_num |
---|
| 182 | endif |
---|
| 183 | call Agrif_seq_add_grid(grid % child_seq,cur_seq,child_p% gr) |
---|
| 184 | child_p => child_p % next |
---|
| 185 | enddo |
---|
| 186 | ! |
---|
| 187 | call Agrif_gl_delete(sorted_child_list) |
---|
| 188 | !--------------------------------------------------------------------------------------------------- |
---|
| 189 | end subroutine Agrif_seq_create_proc_sequences |
---|
| 190 | !=================================================================================================== |
---|
| 191 | ! |
---|
| 192 | !=================================================================================================== |
---|
| 193 | function Agrif_seq_grids_are_connected( g1, g2 ) result( connection ) |
---|
| 194 | ! |
---|
| 195 | !< Compare required_proc_list for g1 and g2. These are connected if they share a same proc. |
---|
| 196 | !--------------------------------------------------------------------------------------------------- |
---|
| 197 | type(Agrif_Grid), intent(in) :: g1, g2 |
---|
| 198 | ! |
---|
| 199 | logical :: connection |
---|
| 200 | type(Agrif_Proc_p), pointer :: pp1, pp2 |
---|
| 201 | ! |
---|
| 202 | connection = .false. |
---|
| 203 | ! |
---|
| 204 | pp1 => g1 % required_proc_list % first |
---|
| 205 | ! |
---|
| 206 | do while( associated(pp1) .and. (.not. connection) ) |
---|
| 207 | ! |
---|
| 208 | pp2 => g2 % required_proc_list % first |
---|
| 209 | do while ( associated(pp2) .and. (.not. connection) ) |
---|
| 210 | if ( pp1 % proc % pn == pp2 % proc % pn ) then |
---|
| 211 | ! if pp1 and pp2 are the same proc, it means that g1 and g2 are connected. We stop here. |
---|
| 212 | connection = .true. |
---|
| 213 | endif |
---|
| 214 | pp2 => pp2 % next |
---|
| 215 | enddo |
---|
| 216 | pp1 => pp1 % next |
---|
| 217 | ! |
---|
| 218 | enddo |
---|
| 219 | !--------------------------------------------------------------------------------------------------- |
---|
| 220 | end function Agrif_seq_grids_are_connected |
---|
| 221 | !=================================================================================================== |
---|
| 222 | ! |
---|
| 223 | !=================================================================================================== |
---|
| 224 | function Agrif_seq_colorize_grid_list ( gridlist ) result ( ncolors ) |
---|
| 225 | ! |
---|
| 226 | !< 1. Sort nodes in decreasing order of degree. |
---|
| 227 | !< 2. Color the node with highest degree with color 1. |
---|
| 228 | !< 3. Choose the node with the largest DSAT value. In case of conflict, choose the one with the |
---|
| 229 | !! highest degree. Then the one corresponding to the largest grid. |
---|
| 230 | !< 4. Color this node with the smallest possible color. |
---|
| 231 | !< 5. If all nodes are colored, then stop. Otherwise, go to 3. |
---|
| 232 | !--------------------------------------------------------------------------------------------------- |
---|
| 233 | type(Agrif_Grid_List), intent(in) :: gridlist |
---|
| 234 | ! |
---|
| 235 | type(Agrif_Grid_List), pointer :: X, Y |
---|
| 236 | type(Agrif_PGrid), pointer :: gridp |
---|
| 237 | type(Agrif_Grid_List), pointer :: tmp_gl |
---|
| 238 | integer :: ncolors |
---|
| 239 | integer, dimension(1:gridlist%nitems) :: colors |
---|
| 240 | ! |
---|
| 241 | ! Be carefull... |
---|
| 242 | nullify(Y) |
---|
| 243 | ! |
---|
| 244 | ! First initialize the color of each node |
---|
| 245 | gridp => gridlist % first |
---|
| 246 | do while ( associated(gridp) ) |
---|
| 247 | gridp % gr % seq_num = 0 |
---|
| 248 | gridp => gridp % next |
---|
| 249 | enddo |
---|
| 250 | ! |
---|
| 251 | ! Then sort the grids by decreasing degree |
---|
| 252 | X => Agrif_gl_merge_sort( gridlist, compare_grid_degrees ) |
---|
| 253 | gridp => X % first |
---|
| 254 | ! |
---|
| 255 | ! Colorize the first grid in the list |
---|
| 256 | gridp % gr % seq_num = 1 |
---|
| 257 | gridp => gridp % next |
---|
| 258 | ! |
---|
| 259 | ! Then for each of its neighbors... |
---|
| 260 | do while ( associated(gridp) ) |
---|
| 261 | ! |
---|
| 262 | if ( gridp % gr % neigh_list % nitems == 0 ) then |
---|
| 263 | ! this grid is alone... let.s attach it to an existing sequence |
---|
| 264 | call Agrif_seq_attach_grid( X, gridp % gr ) |
---|
| 265 | gridp => gridp % next |
---|
| 266 | cycle |
---|
| 267 | endif |
---|
| 268 | ! |
---|
| 269 | ! Compute dsat value of all non-colored grids |
---|
| 270 | tmp_gl => Agrif_gl_build_from_gp(gridp) |
---|
| 271 | call Agrif_seq_calc_dsat(tmp_gl) |
---|
| 272 | ! |
---|
| 273 | ! Sort non-colored grids by decreasing dsat value, then by size |
---|
| 274 | call Agrif_gl_delete(Y) |
---|
| 275 | Y => Agrif_gl_merge_sort( tmp_gl, compare_dsat_values, compare_size_values ) |
---|
| 276 | ! |
---|
| 277 | ! Next coloration is for the first grid in this list TODO : maybe we could find a better choice ..? |
---|
| 278 | gridp => Y % first |
---|
| 279 | ! |
---|
| 280 | ! Assign a color to the chosen grid |
---|
| 281 | gridp % gr % seq_num = Agrif_seq_smallest_available_color_in_neighborhood( gridp % gr % neigh_list ) |
---|
| 282 | ! |
---|
| 283 | gridp => gridp % next |
---|
| 284 | call Agrif_gl_delete(tmp_gl) |
---|
| 285 | ! |
---|
| 286 | enddo |
---|
| 287 | ! |
---|
| 288 | call Agrif_gl_delete(X) |
---|
| 289 | call Agrif_seq_colors_in_neighborhood( gridlist, colors ) |
---|
| 290 | ncolors = maxval(colors) |
---|
| 291 | !--------------------------------------------------------------------------------------------------- |
---|
| 292 | end function Agrif_seq_colorize_grid_list |
---|
| 293 | !=================================================================================================== |
---|
| 294 | ! |
---|
| 295 | !=================================================================================================== |
---|
| 296 | subroutine Agrif_seq_attach_grid ( gridlist, grid ) |
---|
| 297 | ! |
---|
| 298 | !< 'grid' is not connected to any neighbor. Therefore, we give an existing and well chosen color. |
---|
| 299 | !--------------------------------------------------------------------------------------------------- |
---|
| 300 | type(Agrif_Grid_List), intent(in) :: gridlist |
---|
| 301 | type(Agrif_Grid), intent(inout) :: grid |
---|
| 302 | ! |
---|
| 303 | integer, dimension(gridlist%nitems) :: colors |
---|
| 304 | integer, dimension(:), allocatable :: ngrids_by_color |
---|
| 305 | integer :: i, color, ncolors |
---|
| 306 | ! |
---|
| 307 | call Agrif_seq_colors_in_neighborhood( gridlist, colors ) |
---|
| 308 | ncolors = maxval(colors) |
---|
| 309 | ! |
---|
| 310 | allocate(ngrids_by_color(ncolors)) |
---|
| 311 | ngrids_by_color = 0 |
---|
| 312 | ! |
---|
| 313 | do i = 1,gridlist % nitems |
---|
| 314 | if (colors(i) > 0) ngrids_by_color(colors(i)) = ngrids_by_color(colors(i)) + 1 |
---|
| 315 | enddo |
---|
| 316 | ! |
---|
| 317 | color = ncolors |
---|
| 318 | do i = 1,ncolors |
---|
| 319 | if ( ngrids_by_color(i) < color ) color = i |
---|
| 320 | enddo |
---|
| 321 | ! |
---|
| 322 | grid % seq_num = color |
---|
| 323 | deallocate(ngrids_by_color) |
---|
| 324 | !--------------------------------------------------------------------------------------------------- |
---|
| 325 | end subroutine Agrif_seq_attach_grid |
---|
| 326 | !=================================================================================================== |
---|
| 327 | ! |
---|
| 328 | !=================================================================================================== |
---|
| 329 | subroutine Agrif_seq_colors_in_neighborhood ( neigh_list, colors ) |
---|
| 330 | !--------------------------------------------------------------------------------------------------- |
---|
| 331 | type(Agrif_Grid_List), intent(in) :: neigh_list |
---|
| 332 | integer, dimension(:), intent(out) :: colors |
---|
| 333 | ! |
---|
| 334 | integer :: i |
---|
| 335 | type(Agrif_PGrid), pointer :: gridp |
---|
| 336 | ! |
---|
| 337 | i = lbound(colors,1) |
---|
| 338 | colors = 0 |
---|
| 339 | gridp => neigh_list % first |
---|
| 340 | ! |
---|
| 341 | do while ( associated(gridp) ) |
---|
| 342 | ! |
---|
| 343 | if ( i > ubound(colors,1) ) then |
---|
| 344 | print*,'Error in Agrif_seq_colors_in_neighborhood : "colors" array is too small' |
---|
| 345 | stop |
---|
| 346 | endif |
---|
| 347 | colors(i) = gridp % gr % seq_num |
---|
| 348 | gridp => gridp % next |
---|
| 349 | i = i+1 |
---|
| 350 | ! |
---|
| 351 | enddo |
---|
| 352 | !--------------------------------------------------------------------------------------------------- |
---|
| 353 | end subroutine Agrif_seq_colors_in_neighborhood |
---|
| 354 | !=================================================================================================== |
---|
| 355 | ! |
---|
| 356 | !=================================================================================================== |
---|
| 357 | function Agrif_seq_smallest_available_color_in_neighborhood ( neigh_list ) result ( smallest ) |
---|
| 358 | !--------------------------------------------------------------------------------------------------- |
---|
| 359 | type(Agrif_Grid_List), intent(in) :: neigh_list |
---|
| 360 | ! |
---|
| 361 | integer, dimension(:), allocatable :: color_is_met |
---|
| 362 | integer :: colors_tab(1:neigh_list%nitems) |
---|
| 363 | integer :: i, smallest, max_color |
---|
| 364 | ! |
---|
| 365 | call Agrif_seq_colors_in_neighborhood( neigh_list, colors_tab ) |
---|
| 366 | max_color = maxval(colors_tab) |
---|
| 367 | ! |
---|
| 368 | allocate(color_is_met(1:max_color)) |
---|
| 369 | color_is_met = 0 |
---|
| 370 | ! |
---|
| 371 | do i = 1,neigh_list % nitems |
---|
| 372 | if ( colors_tab(i) /= 0 ) then |
---|
| 373 | color_is_met(colors_tab(i)) = 1 |
---|
| 374 | endif |
---|
| 375 | enddo |
---|
| 376 | ! |
---|
| 377 | smallest = max_color+1 |
---|
| 378 | do i = 1,max_color |
---|
| 379 | if ( color_is_met(i) == 0 ) then |
---|
| 380 | smallest = i |
---|
| 381 | exit |
---|
| 382 | endif |
---|
| 383 | enddo |
---|
| 384 | ! |
---|
| 385 | deallocate(color_is_met) |
---|
| 386 | !--------------------------------------------------------------------------------------------------- |
---|
| 387 | end function Agrif_seq_smallest_available_color_in_neighborhood |
---|
| 388 | !=================================================================================================== |
---|
| 389 | ! |
---|
| 390 | !=================================================================================================== |
---|
| 391 | subroutine Agrif_seq_calc_dsat ( gridlist ) |
---|
| 392 | !< For each node 'v' : |
---|
| 393 | !< if none of its neighbors is colored then |
---|
| 394 | !< DSAT(v) = degree(v) # degree(v) := number of neighbors |
---|
| 395 | !< else |
---|
| 396 | !< DSAT(v) = number of different colors used in the first neighborhood of v. |
---|
| 397 | !--------------------------------------------------------------------------------------------------- |
---|
| 398 | type(Agrif_Grid_List), intent(in) :: gridlist |
---|
| 399 | ! |
---|
| 400 | type(Agrif_PGrid), pointer :: gridp |
---|
| 401 | type(Agrif_Grid), pointer :: grid |
---|
| 402 | integer, dimension(:), allocatable :: colors, color_is_met |
---|
| 403 | integer :: i, ncolors |
---|
| 404 | ! |
---|
| 405 | gridp => gridlist % first |
---|
| 406 | ! |
---|
| 407 | do while ( associated(gridp) ) |
---|
| 408 | ! |
---|
| 409 | grid => gridp % gr |
---|
| 410 | ! |
---|
| 411 | allocate(colors(grid % neigh_list % nitems)) |
---|
| 412 | call Agrif_seq_colors_in_neighborhood( grid % neigh_list, colors ) |
---|
| 413 | |
---|
| 414 | allocate(color_is_met(1:maxval(colors))) |
---|
| 415 | color_is_met = 0 |
---|
| 416 | ! |
---|
| 417 | do i = 1,grid % neigh_list % nitems |
---|
| 418 | if ( colors(i) /= 0 ) then |
---|
| 419 | color_is_met(colors(i)) = 1 |
---|
| 420 | endif |
---|
| 421 | enddo |
---|
| 422 | ncolors = sum(color_is_met) |
---|
| 423 | ! |
---|
| 424 | if ( ncolors == 0 ) then |
---|
| 425 | grid % dsat = grid % neigh_list % nitems |
---|
| 426 | else |
---|
| 427 | grid % dsat = ncolors |
---|
| 428 | endif |
---|
| 429 | deallocate(colors, color_is_met) |
---|
| 430 | gridp => gridp % next |
---|
| 431 | enddo |
---|
| 432 | !--------------------------------------------------------------------------------------------------- |
---|
| 433 | end subroutine Agrif_seq_calc_dsat |
---|
| 434 | !=================================================================================================== |
---|
| 435 | ! |
---|
| 436 | !=================================================================================================== |
---|
| 437 | subroutine Agrif_seq_allocate_procs_to_childs ( coarse_grid ) |
---|
| 438 | !--------------------------------------------------------------------------------------------------- |
---|
| 439 | type(Agrif_Grid), intent(inout) :: coarse_grid |
---|
| 440 | ! |
---|
| 441 | integer :: is, ip, ig, ngrids |
---|
| 442 | type(Agrif_Grid_List), pointer :: gridlist |
---|
| 443 | type(Agrif_PGrid), pointer :: gp |
---|
| 444 | type(Agrif_Grid), pointer :: grid |
---|
| 445 | type(Agrif_Proc_List), pointer :: proclist |
---|
| 446 | type(Agrif_Proc), pointer :: proc |
---|
| 447 | type(Agrif_Proc_p), pointer :: pp |
---|
| 448 | type(Agrif_Proc), dimension(:), allocatable, target :: procarray |
---|
| 449 | type(Agrif_Grid), dimension(:), allocatable :: gridarray |
---|
| 450 | type(Agrif_Sequence_List), pointer :: seqlist |
---|
| 451 | real,dimension(:),allocatable :: grid_costs |
---|
| 452 | integer,dimension(:), allocatable :: nbprocs_per_grid |
---|
| 453 | integer :: i1, i2, j1, j2 |
---|
| 454 | real :: max_cost |
---|
| 455 | integer :: max_index |
---|
| 456 | ! |
---|
| 457 | seqlist => coarse_grid % child_seq |
---|
| 458 | if ( .not. associated(seqlist) ) return |
---|
| 459 | ! |
---|
| 460 | ! Initialize proc allocation |
---|
| 461 | pp => coarse_grid % proc_def_list % first |
---|
| 462 | do while ( associated(pp) ) |
---|
| 463 | pp % proc % grid_id = 0 |
---|
| 464 | pp => pp % next |
---|
| 465 | enddo |
---|
| 466 | ! |
---|
| 467 | ! For each sequence... |
---|
| 468 | do is = 1,seqlist % nb_seqs |
---|
| 469 | ! |
---|
| 470 | proclist => seqlist % sequences(is) % proclist |
---|
| 471 | gridlist => seqlist % sequences(is) % gridlist |
---|
| 472 | ! |
---|
| 473 | ! Copy coarse grid proc list and convert it to an array |
---|
| 474 | call Agrif_pl_deep_copy( coarse_grid % proc_def_list, proclist ) |
---|
| 475 | call Agrif_pl_to_array ( proclist, procarray ) |
---|
| 476 | ! |
---|
| 477 | ! Allocate a temporary array with concerned grid numbers |
---|
| 478 | ngrids = gridlist % nitems |
---|
| 479 | allocate(gridarray(1:ngrids)) |
---|
| 480 | allocate(grid_costs(1:ngrids)) |
---|
| 481 | allocate(nbprocs_per_grid(1:ngrids)) |
---|
| 482 | |
---|
| 483 | nbprocs_per_grid = 0 |
---|
| 484 | ! |
---|
| 485 | ! Allocate required procs to each grid |
---|
| 486 | gp => gridlist % first |
---|
| 487 | ig = 0 |
---|
| 488 | do while ( associated(gp) ) |
---|
| 489 | grid => gp % gr |
---|
| 490 | ig = ig+1 ; gridarray(ig) = grid |
---|
| 491 | pp => grid % required_proc_list % first |
---|
| 492 | do while ( associated(pp) ) |
---|
| 493 | procarray( pp % proc % pn+1 ) % grid_id = grid % fixedrank |
---|
| 494 | nbprocs_per_grid(ig) = nbprocs_per_grid(ig) + 1 |
---|
| 495 | pp => pp % next |
---|
| 496 | enddo |
---|
| 497 | gp => gp % next |
---|
| 498 | enddo |
---|
| 499 | ! |
---|
| 500 | ! Add unused procs to the grids |
---|
| 501 | ! TODO FIXME: This is just a dummy allocation. You should take into account grid size and more |
---|
| 502 | ! information here... |
---|
| 503 | |
---|
| 504 | ! Estimate current costs |
---|
| 505 | |
---|
| 506 | do ig = 1, ngrids |
---|
| 507 | i1 = gridarray(ig)%ix(1) |
---|
| 508 | i2 = gridarray(ig)%ix(1)+gridarray(ig)%nb(1)/gridarray(ig)%spaceref(1)-1 |
---|
| 509 | j1 = gridarray(ig)%ix(2) |
---|
| 510 | j2 = gridarray(ig)%ix(2)+gridarray(ig)%nb(2)/gridarray(ig)%spaceref(2)-1 |
---|
| 511 | Call Agrif_estimate_parallel_cost(i1,i2,j1,j2,nbprocs_per_grid(ig),grid_costs(ig)) |
---|
| 512 | grid_costs(ig) = grid_costs(ig) * gridarray(ig)%timeref(1) |
---|
| 513 | enddo |
---|
| 514 | |
---|
| 515 | ig = 1 |
---|
| 516 | do ip = 1,proclist%nitems |
---|
| 517 | if ( procarray( ip ) % grid_id == 0 ) then |
---|
| 518 | ! this proc is unused |
---|
| 519 | max_cost = 0. |
---|
| 520 | max_index = 1 |
---|
| 521 | do ig = 1,ngrids |
---|
| 522 | if (grid_costs(ig) > max_cost) then |
---|
| 523 | max_cost = grid_costs(ig) |
---|
| 524 | max_index = ig |
---|
| 525 | endif |
---|
| 526 | enddo |
---|
| 527 | |
---|
| 528 | ig = max_index |
---|
| 529 | procarray( ip ) % grid_id = gridarray(ig) % fixedrank |
---|
| 530 | |
---|
| 531 | nbprocs_per_grid(ig) = nbprocs_per_grid(ig) + 1 |
---|
| 532 | i1 = gridarray(ig)%ix(1) |
---|
| 533 | i2 = gridarray(ig)%ix(1)+gridarray(ig)%nb(1)/gridarray(ig)%spaceref(1)-1 |
---|
| 534 | j1 = gridarray(ig)%ix(2) |
---|
| 535 | j2 = gridarray(ig)%ix(2)+gridarray(ig)%nb(2)/gridarray(ig)%spaceref(2)-1 |
---|
| 536 | Call Agrif_estimate_parallel_cost(i1,i2,j1,j2,nbprocs_per_grid(ig),grid_costs(ig)) |
---|
| 537 | grid_costs(ig) = grid_costs(ig) * gridarray(ig)%timeref(1) |
---|
| 538 | |
---|
| 539 | endif |
---|
| 540 | enddo |
---|
| 541 | ! |
---|
| 542 | ! Allocate proc nums to each grid |
---|
| 543 | gp => gridlist % first |
---|
| 544 | do while ( associated(gp) ) |
---|
| 545 | do ip = 1,proclist%nitems |
---|
| 546 | if ( procarray( ip ) % grid_id == gp % gr % fixedrank ) then |
---|
| 547 | allocate(proc) |
---|
| 548 | proc = procarray( ip ) |
---|
| 549 | call Agrif_pl_append(gp % gr % proc_def_in_parent_list, proc) |
---|
| 550 | endif |
---|
| 551 | enddo |
---|
| 552 | gp => gp % next |
---|
| 553 | enddo |
---|
| 554 | ! |
---|
| 555 | ! Clean up |
---|
| 556 | deallocate(procarray, gridarray, grid_costs, nbprocs_per_grid) |
---|
| 557 | ! |
---|
| 558 | enddo |
---|
| 559 | !--------------------------------------------------------------------------------------------------- |
---|
| 560 | end subroutine Agrif_seq_allocate_procs_to_childs |
---|
| 561 | !=================================================================================================== |
---|
| 562 | ! |
---|
| 563 | !=================================================================================================== |
---|
| 564 | subroutine Agrif_seq_create_communicators ( grid ) |
---|
| 565 | !--------------------------------------------------------------------------------------------------- |
---|
| 566 | type(Agrif_Grid), intent(inout) :: grid |
---|
| 567 | ! |
---|
| 568 | include 'mpif.h' |
---|
| 569 | type(Agrif_Sequence_List), pointer :: seqlist ! List of child sequences |
---|
| 570 | type(Agrif_PGrid), pointer :: gridp |
---|
| 571 | type(Agrif_Proc), pointer :: proc |
---|
| 572 | integer :: i, ierr |
---|
| 573 | integer :: current_comm, comm_seq, color_seq |
---|
| 574 | ! |
---|
| 575 | seqlist => grid % child_seq |
---|
| 576 | if ( .not. associated(seqlist) ) return |
---|
| 577 | ! |
---|
| 578 | current_comm = grid % communicator |
---|
| 579 | color_seq = MPI_COMM_NULL |
---|
| 580 | ! |
---|
| 581 | ! For each sequence, split the current communicator into as many groups as needed. |
---|
| 582 | do i = 1,seqlist % nb_seqs |
---|
| 583 | ! |
---|
| 584 | ! Loop over each proclist to give a color to the current process |
---|
| 585 | gridp => seqlist % sequences(i) % gridlist % first |
---|
| 586 | grid_loop : do while ( associated(gridp) ) |
---|
| 587 | proc => Agrif_pl_search_proc( gridp % gr % proc_def_in_parent_list, Agrif_Procrank ) |
---|
| 588 | if ( associated(proc) ) then |
---|
| 589 | if ( gridp % gr % fixedrank /= proc % grid_id ) then |
---|
| 590 | write(*,'("### Error Agrif_seq_create_communicators : ")') |
---|
| 591 | write(*,'(" inconsitancy on proc ",i2,":")') Agrif_Procrank |
---|
| 592 | write(*,'("gr % fixedrank = ",i0,", where proc % grid_id = ",i0)') & |
---|
| 593 | gridp % gr % fixedrank, proc % grid_id |
---|
| 594 | stop |
---|
| 595 | endif |
---|
| 596 | color_seq = gridp % gr % fixedrank |
---|
| 597 | exit grid_loop |
---|
| 598 | endif |
---|
| 599 | gridp => gridp % next |
---|
| 600 | enddo grid_loop |
---|
| 601 | ! |
---|
| 602 | call MPI_COMM_SPLIT(current_comm, color_seq, Agrif_ProcRank, comm_seq, ierr) |
---|
| 603 | gridp % gr % communicator = comm_seq |
---|
| 604 | ! |
---|
| 605 | enddo |
---|
| 606 | !--------------------------------------------------------------------------------------------------- |
---|
| 607 | end subroutine Agrif_seq_create_communicators |
---|
| 608 | !=================================================================================================== |
---|
| 609 | ! |
---|
| 610 | !=================================================================================================== |
---|
| 611 | function Agrif_seq_select_child ( g, is ) result ( gridp ) |
---|
| 612 | !--------------------------------------------------------------------------------------------------- |
---|
| 613 | type(Agrif_Grid), pointer, intent(in) :: g |
---|
| 614 | integer, intent(in) :: is |
---|
| 615 | ! |
---|
| 616 | type(Agrif_PGrid), pointer :: gridp |
---|
| 617 | type(Agrif_Proc), pointer :: proc |
---|
| 618 | ! |
---|
| 619 | call Agrif_Instance( g ) |
---|
| 620 | gridp => g % child_seq % sequences(is) % gridlist % first |
---|
| 621 | ! |
---|
| 622 | do while ( associated(gridp) ) |
---|
| 623 | proc => Agrif_pl_search_proc( gridp % gr % proc_def_in_parent_list, Agrif_Procrank ) |
---|
| 624 | if ( associated(proc) ) then |
---|
| 625 | return |
---|
| 626 | endif |
---|
| 627 | gridp => gridp % next |
---|
| 628 | enddo |
---|
| 629 | write(*,'("### Error Agrif_seq_select_child : no grid found in sequence ",i0," (mother G",i0,") for P",i0)')& |
---|
| 630 | is, g%fixedrank, Agrif_Procrank |
---|
| 631 | stop |
---|
| 632 | !--------------------------------------------------------------------------------------------------- |
---|
| 633 | end function Agrif_seq_select_child |
---|
| 634 | !=================================================================================================== |
---|
| 635 | #else |
---|
| 636 | subroutine dummy_Agrif_seq () |
---|
| 637 | end subroutine dummy_Agrif_seq |
---|
| 638 | #endif |
---|
| 639 | ! |
---|
| 640 | end module Agrif_seq |
---|