[2] | 1 | Module class_DataDecomposition !============================================ |
---|
| 2 | use class_TensorProductGrid |
---|
| 3 | implicit none |
---|
| 4 | public :: NewDataDecomposition |
---|
| 5 | public :: DeleteDataDecomposition |
---|
| 6 | public :: DescribeDecomposition |
---|
| 7 | |
---|
| 8 | type DataDecomposition |
---|
| 9 | integer :: level !!multigrid level 0=fine grid |
---|
| 10 | integer :: nx !!global data size, dim 1 |
---|
| 11 | integer :: ny !!global data size, dim 2 |
---|
| 12 | integer :: nz !!global data size, dim 3 |
---|
| 13 | integer :: numprocs !!number of processors |
---|
| 14 | integer,pointer,dimension(:,:) :: SlabDims,ColDims !!size of local storage arrays |
---|
| 15 | integer,pointer,dimension(:,:) :: SlabVals,ColVals !!num of data vals held |
---|
| 16 | integer,pointer,dimension(:) :: my_global_kvals !!global indexing of local data |
---|
| 17 | integer,pointer,dimension(:) :: my_global_ivals !!global indexing of local data |
---|
| 18 | integer :: locnz |
---|
| 19 | integer :: locnx |
---|
| 20 | character (len=80 ) :: DerivTypes(6) |
---|
| 21 | character (len=80 ) :: label |
---|
| 22 | end type DataDecomposition |
---|
| 23 | |
---|
| 24 | interface assignment (=) |
---|
| 25 | module procedure EqualDataDecompositions |
---|
| 26 | end interface |
---|
| 27 | |
---|
| 28 | Contains |
---|
| 29 | |
---|
| 30 | subroutine DeleteDataDecomposition(DataMap) |
---|
| 31 | type( DataDecomposition ),intent(inout) :: DataMap |
---|
| 32 | deallocate( DataMap%SlabDims ) |
---|
| 33 | deallocate( DataMap%ColDims ) |
---|
| 34 | deallocate( DataMap%SlabVals ) |
---|
| 35 | deallocate( DataMap%ColVals ) |
---|
| 36 | deallocate( DataMap%my_global_kvals ) |
---|
| 37 | deallocate( DataMap%my_global_ivals ) |
---|
| 38 | end subroutine DeleteDataDecomposition |
---|
| 39 | |
---|
| 40 | |
---|
| 41 | subroutine EqualDataDecompositions(new,old) |
---|
| 42 | implicit none |
---|
| 43 | type( DataDecomposition ),intent(inout):: new |
---|
| 44 | type( DataDecomposition ),intent(in) :: old |
---|
| 45 | integer :: numprocs,rcode(6)=0 |
---|
| 46 | |
---|
| 47 | numprocs=old%numprocs |
---|
| 48 | |
---|
| 49 | if( associated(new%SlabDims) ) call DeleteDataDecomposition(new) |
---|
| 50 | |
---|
| 51 | new%level = old%level |
---|
| 52 | new%nx = old%nx |
---|
| 53 | new%ny = old%ny |
---|
| 54 | new%nz = old%nz |
---|
| 55 | new%numprocs = old%numprocs |
---|
| 56 | new%DerivTypes = old%DerivTypes |
---|
| 57 | new%label = old%label |
---|
| 58 | new%locnz = old%locnz |
---|
| 59 | new%locnx = old%locnx |
---|
| 60 | |
---|
| 61 | allocate( new%SlabDims(3,-1:numprocs-1), stat=rcode(1) ) |
---|
| 62 | allocate( new%SlabVals(3,-1:numprocs-1), stat=rcode(2) ) |
---|
| 63 | allocate( new%ColDims(3,-1:numprocs-1), stat=rcode(3) ) |
---|
| 64 | allocate( new%ColVals(3,-1:numprocs-1), stat=rcode(4) ) |
---|
| 65 | allocate( new%my_global_kvals(new%locnz), stat=rcode(5) ) |
---|
| 66 | allocate( new%my_global_ivals(new%locnx), stat=rcode(6) ) |
---|
| 67 | if( any( rcode /= 0 ) ) Stop 'Allocation failed: EqualDataDecompositions' |
---|
| 68 | |
---|
| 69 | new%SlabDims = old%SlabDims |
---|
| 70 | new%ColDims = old%ColDims |
---|
| 71 | new%SlabVals = old%SlabVals |
---|
| 72 | new%ColVals = old%ColVals |
---|
| 73 | new%my_global_kvals = old%my_global_kvals |
---|
| 74 | new%my_global_ivals = old%my_global_ivals |
---|
| 75 | end subroutine EqualDataDecompositions |
---|
| 76 | |
---|
| 77 | subroutine NewDataDecomposition(label,grid,level,DerivTypes,DataMap) |
---|
| 78 | use class_TensorProductGrid |
---|
| 79 | implicit none |
---|
| 80 | |
---|
| 81 | type( DataDecomposition ),intent(inout) :: DataMap |
---|
| 82 | type( TensorProductGrid ),intent(in) :: grid |
---|
| 83 | character (len=80 ) :: DerivTypes(6) |
---|
| 84 | character (len=80 ),intent(in) :: label |
---|
| 85 | integer :: level |
---|
| 86 | integer :: numnodes |
---|
| 87 | integer :: i,k0,rem |
---|
| 88 | integer :: myid,ierr,rcode(6)=0 |
---|
| 89 | include 'mpif.h' |
---|
| 90 | |
---|
| 91 | call mpi_comm_rank (MPI_COMM_WORLD,myid,ierr) |
---|
| 92 | call mpi_comm_size (MPI_COMM_WORLD,numnodes,ierr) |
---|
| 93 | |
---|
| 94 | rem = mod(DataMap%nx,numnodes) |
---|
| 95 | if( rem /= 0 .and. rem /= 1 ) then |
---|
| 96 | if(myid==0) & |
---|
| 97 | write(0,*) ' nx should be evenly divisible by np or have remainder 1' |
---|
| 98 | stop |
---|
| 99 | endif |
---|
| 100 | |
---|
| 101 | rem = mod(DataMap%nz,numnodes) |
---|
| 102 | if( rem /= 0 .and. rem /= 1 ) then |
---|
| 103 | if(myid==0) & |
---|
| 104 | write(0,*) ' nz should be evenly divisible by np or have remainder 1' |
---|
| 105 | stop |
---|
| 106 | endif |
---|
| 107 | |
---|
| 108 | DataMap%level=level |
---|
| 109 | DataMap%label=label |
---|
| 110 | DataMap%DerivTypes(1:3) = DerivTypes(1:3) |
---|
| 111 | DataMap%nx=grid%DataDims(1) |
---|
| 112 | DataMap%ny=grid%DataDims(2) |
---|
| 113 | DataMap%nz=grid%DataDims(3) |
---|
| 114 | DataMap%numprocs=numnodes |
---|
| 115 | |
---|
| 116 | do i=1,3 |
---|
| 117 | if( trim( DerivTypes(i) ) == 'cos' ) then |
---|
| 118 | DataMap%DerivTypes(i+3)='sin' |
---|
| 119 | elseif( trim( DerivTypes(i) ) == 'sin' ) then |
---|
| 120 | DataMap%DerivTypes(i+3)='cos' |
---|
| 121 | else |
---|
| 122 | DataMap%DerivTypes(i+3)=DerivTypes(i) |
---|
| 123 | endif |
---|
| 124 | enddo |
---|
| 125 | |
---|
| 126 | if( .not. associated(DataMap%SlabDims) ) then |
---|
| 127 | allocate( DataMap%SlabDims(3,-1:numnodes-1), stat=rcode(1) ) |
---|
| 128 | allocate( DataMap%SlabVals(3,-1:numnodes-1), stat=rcode(2) ) |
---|
| 129 | allocate( DataMap%ColDims(3,-1:numnodes-1), stat=rcode(3) ) |
---|
| 130 | allocate( DataMap%ColVals(3,-1:numnodes-1), stat=rcode(4) ) |
---|
| 131 | if( any( rcode /= 0 ) ) Stop 'Allocation failed: NewDataDecomposition A' |
---|
| 132 | endif |
---|
| 133 | |
---|
| 134 | i=1 |
---|
| 135 | if( DerivTypes(i) == 'fourier' ) then |
---|
| 136 | DataMap%SlabVals(i,:) = DataMap%nx |
---|
| 137 | DataMap%SlabDims(i,:) = DataMap%nx + 0 !! DON"T accomodate conjugate symmetry |
---|
| 138 | DataMap%ColDims(i,:) = (DataMap%nx)/numnodes !! size counted in real words, but the |
---|
| 139 | DataMap%ColVals(i,:) = (DataMap%nx)/numnodes !! storage is for local part of kx wavenumber space |
---|
| 140 | else |
---|
| 141 | select case ( mod(DataMap%nx,2) ) |
---|
| 142 | case(0) ! nx is even |
---|
| 143 | DataMap%SlabVals(i,:) = DataMap%nx |
---|
| 144 | DataMap%SlabDims(i,:) = DataMap%nx |
---|
| 145 | DataMap%ColDims(i,:) = DataMap%nx/numnodes |
---|
| 146 | DataMap%ColVals(i,:) = DataMap%nx/numnodes |
---|
| 147 | |
---|
| 148 | case(1) ! nx is odd |
---|
| 149 | DataMap%SlabVals(i,:) = DataMap%nx |
---|
| 150 | DataMap%SlabDims(i,:) = DataMap%nx |
---|
| 151 | DataMap%ColDims(i,:) = (DataMap%nx-1)/numnodes |
---|
| 152 | DataMap%ColVals(i,:) = (DataMap%nx-1)/numnodes |
---|
| 153 | DataMap%ColDims(i,numnodes-1) = (DataMap%nx-1)/numnodes + 1 |
---|
| 154 | DataMap%ColVals(i,numnodes-1) = (DataMap%nx-1)/numnodes + 1 |
---|
| 155 | end select |
---|
| 156 | endif |
---|
| 157 | |
---|
| 158 | i=2 |
---|
| 159 | DataMap%SlabVals(i,:) = DataMap%ny |
---|
| 160 | DataMap%SlabDims(i,:) = DataMap%ny |
---|
| 161 | DataMap%ColDims(i,:) = DataMap%ny |
---|
| 162 | DataMap%ColVals(i,:) = DataMap%ny |
---|
| 163 | |
---|
| 164 | i=3 |
---|
| 165 | select case ( mod(DataMap%nz,2) ) |
---|
| 166 | case(0) ! nz is even |
---|
| 167 | DataMap%SlabVals(i,:) = DataMap%nz/numnodes |
---|
| 168 | DataMap%SlabDims(i,:) = DataMap%nz/numnodes |
---|
| 169 | DataMap%ColDims(i,:) = DataMap%nz |
---|
| 170 | DataMap%ColVals(i,:) = DataMap%nz |
---|
| 171 | case(1) ! nz is odd |
---|
| 172 | DataMap%SlabVals(i,:) = (DataMap%nz-1)/numnodes |
---|
| 173 | DataMap%SlabDims(i,:) = (DataMap%nz-1)/numnodes |
---|
| 174 | DataMap%ColDims(i,:) = DataMap%nz |
---|
| 175 | DataMap%ColVals(i,:) = DataMap%nz |
---|
| 176 | DataMap%SlabDims(i,numnodes-1) = (DataMap%nz-1)/numnodes + 1 |
---|
| 177 | DataMap%SlabVals(i,numnodes-1) = (DataMap%nz-1)/numnodes + 1 |
---|
| 178 | end select |
---|
| 179 | |
---|
| 180 | !! determine the global k indices for processor myid |
---|
| 181 | !! first, I already know how many planes are allocated to myid |
---|
| 182 | DataMap%locnz = DataMap%SlabDims(3,myid) |
---|
| 183 | |
---|
| 184 | !! next, add up the k locations below myid |
---|
| 185 | k0=0 |
---|
| 186 | do i=0,myid-1 |
---|
| 187 | k0=k0+DataMap%SlabDims(3,i) !! num of k planes assigned to processor i |
---|
| 188 | enddo |
---|
| 189 | |
---|
| 190 | allocate(DataMap%my_global_kvals( DataMap%locnz ), stat=rcode(1) ) |
---|
| 191 | if( rcode(1) /= 0 ) Stop 'Allocation failed: NewDataDecomposition B' |
---|
| 192 | !!fill the array w/ the k values |
---|
| 193 | do i=1,DataMap%locnz |
---|
| 194 | DataMap%my_global_kvals(i) = k0+i |
---|
| 195 | enddo |
---|
| 196 | |
---|
| 197 | !! determine the global i indices for processor myid (Column orientation) |
---|
| 198 | !! first, I already know how many planes are allocated to myid |
---|
| 199 | DataMap%locnx = DataMap%ColDims(1,myid) |
---|
| 200 | |
---|
| 201 | !! next, add up the i locations left of myid |
---|
| 202 | k0=0 |
---|
| 203 | do i=0,myid-1 |
---|
| 204 | k0=k0+DataMap%ColDims(1,i) !! num of X planes assigned to processor i |
---|
| 205 | enddo |
---|
| 206 | |
---|
| 207 | allocate(DataMap%my_global_ivals( DataMap%locnx ), stat=rcode(1) ) |
---|
| 208 | if( rcode(1) /= 0 ) Stop 'Allocation failed: NewDataDecomposition C' |
---|
| 209 | !!fill the array w/ the i values |
---|
| 210 | do i=1,DataMap%locnx |
---|
| 211 | DataMap%my_global_ivals(i) = k0+i |
---|
| 212 | enddo |
---|
| 213 | |
---|
| 214 | if(myid==0) write(0,*) '......creating MPI derived data types' |
---|
| 215 | call create_MPI_data_types(DataMap) |
---|
| 216 | if(myid==0) write(0,*) '....... data types done' |
---|
| 217 | |
---|
| 218 | end subroutine NewDataDecomposition |
---|
| 219 | |
---|
| 220 | subroutine DescribeDecomposition(DataMap,docfile) |
---|
| 221 | implicit none |
---|
| 222 | |
---|
| 223 | type( DataDecomposition ) :: DataMap |
---|
| 224 | character(len=80) :: docfile |
---|
| 225 | integer :: myid |
---|
| 226 | integer :: numprocs |
---|
| 227 | integer :: ierr |
---|
| 228 | |
---|
| 229 | include 'mpif.h' |
---|
| 230 | |
---|
| 231 | call mpi_comm_rank(MPI_COMM_WORLD,myid,ierr) |
---|
| 232 | call mpi_comm_size(MPI_COMM_WORLD,numprocs,ierr) |
---|
| 233 | |
---|
| 234 | if(myid .ne. 0 ) return |
---|
| 235 | open(1,file=docfile,position='append') |
---|
| 236 | write(1,*) '============================================================' |
---|
| 237 | write(1,*) ' Data Decomposition ',trim(DataMap%label),' for : ',DataMap%numprocs,' processors' |
---|
| 238 | write(1,*) '============================================================' |
---|
| 239 | |
---|
| 240 | write(1,*) ' Grid Level: ',DataMap%level |
---|
| 241 | write(1,*) ' Global number of data points in the ' |
---|
| 242 | write(1,*) ' (x,y,z)<-->(1,2,3)<-->(i,j,k) directions:' |
---|
| 243 | write(1,*) ' nx = ',DataMap%nx |
---|
| 244 | write(1,*) ' ny = ',DataMap%ny |
---|
| 245 | write(1,*) ' nz = ',DataMap%nz |
---|
| 246 | write(1,*) ' Derivative types specified for each dimension: ' |
---|
| 247 | write(1,*) ' in x ', trim(DataMap%DerivTypes(1)) |
---|
| 248 | write(1,*) ' in y ', trim(DataMap%DerivTypes(2)) |
---|
| 249 | write(1,*) ' in z ', trim(DataMap%DerivTypes(3)) |
---|
| 250 | write(1,*) ' In {slab} format, each processor is allocated ' |
---|
| 251 | write(1,*) ' nx*ny*locnz data values per variable, where' |
---|
| 252 | write(1,*) ' nx = ',DataMap%SlabVals(1,0) |
---|
| 253 | write(1,*) ' ny = ',DataMap%SlabVals(2,0) |
---|
| 254 | write(1,*) ' locnz = ',DataMap%SlabVals(3,0) |
---|
| 255 | if(DataMap%SlabVals(3,0) .ne. DataMap%SlabVals(3,DataMap%numprocs-1)) then |
---|
| 256 | write(1,*) ' except for processor ',DataMap%numprocs-1 |
---|
| 257 | write(1,*) ' which is allocated 1 extra plane. ' |
---|
| 258 | write(1,*) ' locnz = ',DataMap%SlabVals(3,DataMap%numprocs-1) |
---|
| 259 | endif |
---|
| 260 | write(1,*) ' ' |
---|
| 261 | write(1,*) ' The local {slab} data arrays should be dimensioned ' |
---|
| 262 | write(1,*) ' nx = ',DataMap%SlabDims(1,0) |
---|
| 263 | write(1,*) ' ny = ',DataMap%SlabDims(2,0) |
---|
| 264 | write(1,*) ' locnz = ',DataMap%SlabDims(3,0) |
---|
| 265 | if(DataMap%SlabDims(3,0) .ne. DataMap%SlabDims(3,DataMap%numprocs-1)) then |
---|
| 266 | write(1,*) ' except for processor ',DataMap%numprocs-1 |
---|
| 267 | write(1,*) ' which is allocated 1 extra plane. ' |
---|
| 268 | write(1,*) ' locnz = ',DataMap%SlabVals(3,DataMap%numprocs-1) |
---|
| 269 | write(1,*) ' ' |
---|
| 270 | write(1,*) ' myid = ',myid |
---|
| 271 | write(1,*) ' DataMap%SlabDims(1,myid) = nx = ',DataMap%SlabDims(1,myid) |
---|
| 272 | write(1,*) ' DataMap%SlabDims(2,myid) = ny = ',DataMap%SlabDims(2,myid) |
---|
| 273 | write(1,*) ' DataMap%SlabDims(3,myid) = locnz = ',DataMap%SlabDims(3,myid) |
---|
| 274 | write(1,*) ' DataMap%locnz = locnz = ',DataMap%locnz |
---|
| 275 | write(1,*) ' DataMap%my_global_kvals = ',DataMap%my_global_kvals |
---|
| 276 | |
---|
| 277 | |
---|
| 278 | endif |
---|
| 279 | write(1,*) ' ' |
---|
| 280 | write(1,*) '=========== ' |
---|
| 281 | write(1,*) ' ' |
---|
| 282 | write(1,*) ' In {column} format, measured in REAL words,' |
---|
| 283 | write(1,*) ' (even though the col arrays are sometimes complex),' |
---|
| 284 | write(1,*) ' each processor is allocated locnx*ny*nz ' |
---|
| 285 | write(1,*) ' REAL data values per variable, where ' |
---|
| 286 | write(1,*) ' locnx = ',DataMap%ColVals(1,0) |
---|
| 287 | write(1,*) ' ny = ',DataMap%ColVals(2,0) |
---|
| 288 | write(1,*) ' nz = ',DataMap%ColVals(3,0) |
---|
| 289 | if(DataMap%ColVals(1,0) .ne. DataMap%ColVals(1,DataMap%numprocs-1)) then |
---|
| 290 | write(1,*) ' except for processor ',DataMap%numprocs-1 |
---|
| 291 | write(1,*) ' which is allocated 1 extra slot. ' |
---|
| 292 | write(1,*) ' locnx = ',DataMap%ColVals(1,DataMap%numprocs-1) |
---|
| 293 | endif |
---|
| 294 | write(1,*) ' ' |
---|
| 295 | write(1,*) ' The local {Col} data arrays should be dimensioned ' |
---|
| 296 | write(1,*) ' locnx = ',DataMap%ColVals(1,0) |
---|
| 297 | write(1,*) ' ny = ',DataMap%ColVals(2,0) |
---|
| 298 | write(1,*) ' nz = ',DataMap%ColVals(3,0) |
---|
| 299 | if(DataMap%ColVals(1,0) .ne. DataMap%ColVals(1,DataMap%numprocs-1)) then |
---|
| 300 | write(1,*) ' except for processor ',DataMap%numprocs-1 |
---|
| 301 | write(1,*) ' which is allocated 1 extra plane. ' |
---|
| 302 | write(1,*) ' locnx = ',DataMap%ColVals(1,DataMap%numprocs-1) |
---|
| 303 | write(1,*) ' ' |
---|
| 304 | write(1,*) ' myid = ',myid |
---|
| 305 | write(1,*) ' DataMap%ColDims(1,myid) = nx = ',DataMap%ColDims(1,myid) |
---|
| 306 | write(1,*) ' DataMap%ColDims(2,myid) = ny = ',DataMap%ColDims(2,myid) |
---|
| 307 | write(1,*) ' DataMap%ColDims(3,myid) = locnx = ',DataMap%ColDims(3,myid) |
---|
| 308 | endif |
---|
| 309 | write(1,*) '============================================================' |
---|
| 310 | write(1,*) ' ' |
---|
| 311 | close(1) |
---|
| 312 | end subroutine DescribeDecomposition |
---|
| 313 | |
---|
| 314 | |
---|
| 315 | end Module class_DataDecomposition |
---|