modseq.F90 25 KB

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