partition.F [SRC] [CPP] [JOB] [SCAN]
SOURCES / MAIN / MASTER



   1 | include(dom.inc)
   2 | 
   3 |       SUBROUTINE partition
   4 | 
   5 | !       ================================================================!
   6 | !                                                                       !
   7 | !       partition.F   : Calculates de 'begin' and 'end' directions      !
   8 | !                       and bands for each processor. Vactors will be   !
   9 | !                       partitioned using this values.                  !
  10 | !                                                                       !
  11 | !       out           : Vectors containing information on first/last    !
  12 | !                       direction/band to treat for each processor.     !
  13 | !                                                                       !
  14 | !       comments      : When using GRAYCASE there is no band integration!
  15 | !                       so the partitioning is only made over the       !
  16 | !                       directions. A 3rd kind of parallelism could be  !
  17 | !                       done over the domain, leting each processor     !
  18 | !                       calculate only part of the cells.               !
  19 | !                                                                       !
  20 | !       nota          : with this kind of partitioning if the number of !
  21 | !                       directions is < than the number of processors   !
  22 | !                       for a gray case, there will be some UNUSED      !
  23 | !                       processors!!                                    !
  24 | !                                                                       !
  25 | !       nota 2        : epsilon = 1e-5 is needed beacause in some cases !
  26 | !                       it would be possible to have a sum equal to the !
  27 | !                       number of processors, but with a floating sum   !
  28 | !                       slightly superior to this value (p.e. 24.00001) !
  29 | !                       To avoid errors an epsilon value is substracted !
  30 | !                       from the sum.                                   !
  31 | !                                                                       !
  32 | !       author        : J. AMAYA (september 2007)                       !
  33 | !                                                                       !
  34 | !       ================================================================!
  35 | 
  36 |         USE mod_pmm
  37 |         USE mod_prissma
  38 |         USE mod_inout
  39 | 
  40 |         IMPLICIT NONE
  41 | 
  42 |         include 'pmm_constants.h'
  43 | 
  44 | !       LOCAL
  45 |         DOM_INT    :: base_step, step, reste, instep, i
  46 |         DOM_INT    :: node_reste, node_step
  47 |         DOM_INT    :: node_beg, node_end
  48 |         DOM_INT    :: cell_reste, cell_step
  49 |         DOM_INT    :: cell_beg, cell_end
  50 |         DOM_INT    :: bf_reste, bf_step
  51 |         DOM_INT    :: bfbeg, bfend
  52 |         DOM_INT    :: idir_reste, idir_step
  53 |         DOM_INT    :: idir_beg, idir_end
  54 |         DOM_INT    :: iproc, addi, ierr, partitype, nspec
  55 |         DOM_INT, PARAMETER :: buffersize =18
  56 |         DOM_INT    :: buffer(buffersize)
  57 |         DOM_REAL   :: flostep, fdirbeg, fdirend
  58 |         DOM_REAL   :: epsilon
  59 | 
  60 | !       -----------------!
  61 | !       Allocate vectors !
  62 | !       -----------------!
  63 | 
  64 |         IF (ALLOCATED(cd)) DEALLOCATE(cd)
  65 |         IF (ALLOCATED(cf)) DEALLOCATE(cf)
  66 |         IF (ALLOCATED(dir_d)) DEALLOCATE(dir_d)
  67 |         IF (ALLOCATED(dir_f)) DEALLOCATE(dir_f)
  68 | 
  69 |         ALLOCATE(cd(pmm_n_p))
  70 |         ALLOCATE(cf(pmm_n_p))
  71 |         ALLOCATE(dir_d(pmm_n_p))
  72 |         ALLOCATE(dir_f(pmm_n_p))
  73 | 
  74 | !       -------------------------------------!
  75 | !       Calculate domain decomposition steps !
  76 | !       -------------------------------------!
  77 | 
  78 |         node_step  = INT(i_dom_nnodes/pmm_n_p)
  79 |         node_reste = MOD(i_dom_nnodes,pmm_n_p)
  80 |         node_end   = 0
  81 | 
  82 |         cell_step  = INT(i_dom_ncells/pmm_n_p)
  83 |         cell_reste = MOD(i_dom_ncells,pmm_n_p)
  84 |         cell_end   = 0
  85 | 
  86 |         bf_step  = INT(i_dom_nbfaces/pmm_n_p)
  87 |         bf_reste = MOD(i_dom_nbfaces,pmm_n_p)
  88 |         bfend    = 0
  89 | 
  90 | !       ----------------------!
  91 | !       Calculate 'step' size !
  92 | !       ----------------------!
  93 | 
  94 |         idir_step  = INT(ndir/pmm_n_p)
  95 |         idir_reste = MOD(ndir,pmm_n_p)
  96 |         idir_end   = 0
  97 | 
  98 |         epsilon = 1e-5
  99 | 
 100 |         IF (trim(mediumtype).eq.'CK') THEN
 101 |           nspec = nallbandes
 102 |         ELSE
 103 |           nspec = nkabs
 104 |         ENDIF
 105 | 
 106 |         base_step = INT(nspec*ndir/pmm_n_p)
 107 |         reste = MOD(nspec*ndir,pmm_n_p)
 108 | 
 109 |         fdirend = 0.
 110 | 
 111 |         partitype   = 2
 112 |         IF (pmm_n_p.gt.ndir) partitype = 1
 113 | 
 114 |         DO iproc=1,pmm_n_p
 115 | 
 116 | !       --------------------!
 117 | !       Domain partitioning !
 118 | !       --------------------!
 119 | 
 120 |         node_beg  = node_end + 1
 121 |         node_end  = node_beg + node_step - 1
 122 |         IF (node_reste.gt.0) THEN
 123 |           node_end   = node_end + 1
 124 |           node_reste = node_reste - 1
 125 |         ENDIF
 126 | 
 127 |         cell_beg  = cell_end + 1
 128 |         cell_end  = cell_beg + cell_step - 1
 129 |         IF (cell_reste.gt.0) THEN
 130 |           cell_end   = cell_end + 1
 131 |           cell_reste = cell_reste - 1
 132 |         ENDIF
 133 | 
 134 | !       ---------------------------!
 135 | !       Boundary face partitioning !
 136 | !       ---------------------------!
 137 | 
 138 |         bfbeg  = bfend + 1
 139 |         bfend  = bfbeg + bf_step - 1
 140 |         IF (bf_reste.gt.0) THEN
 141 |           bfend    = bfend + 1
 142 |           bf_reste = bf_reste - 1
 143 |         ENDIF
 144 | 
 145 | !         --------------------!
 146 | !         Partitioning type 1 !
 147 | !         --------------------!
 148 | 
 149 |           IF (partitype.eq.1) THEN
 150 | 
 151 | !           ----------------------------!
 152 | !           Capture bands for this proc !
 153 | !           ----------------------------!
 154 |             IF (iproc.eq.1) THEN
 155 |               cd(iproc)=1
 156 |               dir_d(iproc) = 1
 157 |             ELSE
 158 |               cd(iproc) = cf(iproc-1)+1
 159 |               dir_d(iproc)= dir_f(iproc-1)
 160 |               IF(cd(iproc).eq.1) dir_d(iproc)= dir_f(iproc-1)+1
 161 |             ENDIF
 162 |             IF (cd(iproc).gt.nspec) THEN
 163 |               cd(iproc) = 1
 164 |               dir_d(iproc)= dir_f(iproc-1)+1
 165 |             ENDIF
 166 | 
 167 |             step = base_step
 168 |             IF (reste > 0) THEN
 169 |               step=step+1
 170 |               reste = reste-1
 171 |             ENDIF
 172 |             cf(iproc) = cd(iproc)+(step-1)
 173 | 
 174 |             DO WHILE(cf(iproc).gt.nspec)
 175 |               cf(iproc) = cf(iproc) - nspec
 176 |             ENDDO
 177 | 
 178 | !           -----------------------!
 179 | !           Capture the directions !
 180 | !           -----------------------!
 181 | 
 182 |             fdirbeg = fdirend
 183 |             flostep = REAL(step)/REAL(nspec)
 184 |             fdirend = fdirbeg + flostep - epsilon
 185 | 
 186 |             dir_f(iproc) = CEILING(fdirend)
 187 | 
 188 | !         ------------------------------!
 189 | !         Partitioning 2 (nproc < ndir) !
 190 | !         ------------------------------!
 191 | 
 192 |           ELSE
 193 | 
 194 |             cd(iproc) = 1
 195 |             cf(iproc) = nspec
 196 | 
 197 |             idir_beg  = idir_end + 1
 198 |             idir_end  = idir_beg + idir_step - 1
 199 |             IF (idir_reste.gt.0) THEN
 200 |               idir_end   = idir_end + 1
 201 |               idir_reste = idir_reste - 1
 202 |             ENDIF
 203 | 
 204 |             dir_d(iproc) = idir_beg
 205 |             dir_f(iproc) = idir_end
 206 | 
 207 |           ENDIF
 208 | 
 209 | !         --------------------------------------------------------!
 210 | !         Filling buffer with partitioning and global information !
 211 | !         --------------------------------------------------------!
 212 | 
 213 |           buffer(1) = i_dom_nnodes
 214 |           buffer(2) = i_dom_ncells
 215 |           buffer(3) = cd(iproc)
 216 |           buffer(4) = cf(iproc)
 217 |           buffer(5) = ndir
 218 |           buffer(6) = dir_d(iproc)
 219 |           buffer(7) = dir_f(iproc)
 220 |           buffer(8) = i_dom_nfacesmax
 221 |           buffer(9) = n_gaz
 222 |           buffer(10) = nallbandes
 223 |           buffer(11)= node_beg
 224 |           buffer(12)= node_end
 225 |           buffer(13)= i_dom_nbfaces
 226 |           buffer(14)= bfbeg
 227 |           buffer(15)= bfend
 228 |           buffer(16)= cell_beg
 229 |           buffer(17)= cell_end
 230 |           buffer(18)= i_dom_nfaces
 231 | 
 232 |           print*, "    Partitionning on proc ", iproc,":"
 233 | !         print*, buffer(8) , "nfacemax"
 234 | !         print*, buffer(18), "Faces"
 235 | !         print*, buffer(13), "Bfaces : ", buffer(14), "-->",buffer(15)
 236 | !         print*, "    ngas     :", buffer(9)
 237 | 
 238 |           print*, i_dom_nbfaces, "Nbfaces"
 239 |           IF (trim(mediumtype).eq.'CK') THEN
 240 |           print*, buffer(2) , "Cells  : ", buffer(16), "-->",buffer(17)
 241 |           print*, buffer(10), "Nbands : ", buffer(3) , "-->",buffer(4)
 242 |           print*, nkabs     , "Nq pts"
 243 |           ELSE
 244 |           print*, buffer(1) , "Nodes  : ", buffer(11), "-->",buffer(12)
 245 |           print*, nkabs     , "Nq pts : ", buffer(3) , "-->",buffer(4)
 246 |           ENDIF
 247 | 
 248 |           print*, buffer(5) , "Ndirs  : ", buffer(6) , "-->",buffer(7)
 249 |           print*
 250 | 
 251 | 
 252 | 
 253 | 
 254 | 
 255 | 
 256 | 
 257 | 
 258 | 
 259 | 
 260 | 
 261 |           CALL pmm_sendpartition(buffer,buffersize,iproc)
 262 | 
 263 |         ENDDO
 264 | 
 265 |       END SUBROUTINE partition


partition.F could be called by:
Makefile [SOURCES] - 141
master_control.F [SOURCES/MAIN/MASTER] - 53 - 54