XBeach
C:/repositories/XBeach/trunk/src/xbeachlibrary/xmpinew.F90
Go to the documentation of this file.
00001 module xmpi_module
00002 #ifdef USEMPI
00003    use mpi
00004    implicit none
00005 #ifndef HAVE_MPI_WTIME
00006    real*8, external                :: MPI_Wtime
00007 #endif
00008    integer, parameter              :: MPIBOUNDARY_Y    = 1
00009    integer, parameter              :: MPIBOUNDARY_X    = 2
00010    integer, parameter              :: MPIBOUNDARY_AUTO = 3
00011    integer, parameter              :: MPIBOUNDARY_MAN  = 4
00012    integer                         :: xmpi_rank    ! mpi rank of this process
00013    integer                         :: xmpi_size    ! number of mpi processes
00014    integer                         :: xmpi_comm    ! mpi communicator to use
00015    integer                         :: xmpi_master  ! rank of master process
00016    integer                         :: xmpi_m       ! 1st dimension of processor
00017    ! grid
00018    integer                         :: xmpi_n       ! 2nd dimension of processor
00019    ! grid
00020    integer                         :: xmpi_pcol    ! my column in processor grid (starting at 1)
00021    integer                         :: xmpi_prow    ! my row    in processor grid (starting at 1)
00022    integer                         :: xmpi_left    ! left neighbour
00023    integer                         :: xmpi_right   ! right neighbour
00024    integer                         :: xmpi_bot     ! bottom neighbour
00025    integer                         :: xmpi_top     ! top neighbour
00026    logical                         :: xmpi_isleft  ! submatrix is at the left side
00027    ! ie: shares first column with
00028    ! global matrix
00029    logical                         :: xmpi_isright ! submatrix is at the right side
00030    ! ie: shares last column with
00031    ! global matrix
00032    logical                         :: xmpi_istop   ! submatrix is at the top side
00033    ! ie: shares first row with
00034    ! global matrix
00035    logical                         :: xmpi_isbot   ! submatrix is at the bottom side
00036    ! ie: shares first row with
00037    ! global matrix
00038    logical                         :: xmaster      ! .true. if this process reads
00039    !  and writes files
00040    !
00041    !         1 2 3 4 5 6 7   y-axis
00042    !  X   1  x x x x x x x
00043    !  a   2  x x x x x x x
00044    !  x   3  x x x x x x x
00045    !  i   4  x x x x x x x
00046    !  s   5  x x x x x x x
00047    !
00048    !
00049    integer, parameter              :: SHIFT_X_U = 1  ! shift in x direction from high to low: from bottom to top
00050    integer, parameter              :: SHIFT_X_D = 2  ! shift in x direction from low to high: from top to bottom
00051    integer, parameter              :: SHIFT_Y_R = 3  ! shift in y direction from low to high: from left to right
00052    integer, parameter              :: SHIFT_Y_L = 4  ! shift in y direction from high to low: from right to left
00053 #ifdef USEMPE
00054    integer                         :: event_output_start
00055    integer                         :: event_output_end
00056    integer                         :: event_write_start
00057    integer                         :: event_write_end
00058    integer                         :: event_coll_start
00059    integer                         :: event_coll_end
00060 #endif
00061 
00062 #else
00063    implicit none
00064    integer, parameter              :: xmpirank     = 0
00065    integer, parameter              :: xmpisize     = 1
00066    logical, parameter              :: xmaster      = .true.
00067    logical, parameter              :: xmpi_isleft  = .true.
00068    logical, parameter              :: xmpi_isright = .true.
00069    logical, parameter              :: xmpi_istop   = .true.
00070    logical, parameter              :: xmpi_isbot   = .true.
00071    integer, parameter              :: xmpi_pcol    = 1
00072    integer, parameter              :: xmpi_prow    = 1
00073 #endif
00074 #ifdef USEMPI
00075    interface xmpi_bcast
00076       module procedure xmpi_bcast_int, xmpi_bcast_real8, xmpi_bcast_int8
00077       module procedure xmpi_bcast_real4
00078       module procedure xmpi_bcast_array_logical
00079       module procedure xmpi_bcast_logical
00080       module procedure xmpi_bcast_complex16
00081       module procedure xmpi_bcast_array_real8
00082       module procedure xmpi_bcast_array_integer
00083       module procedure xmpi_bcast_matrix_integer
00084       module procedure xmpi_bcast_matrix_real8
00085       module procedure xmpi_bcast_char
00086    end interface xmpi_bcast
00087 
00088    interface xmpi_allreduce
00089       module procedure xmpi_allreduce_r0
00090       module procedure xmpi_allreduce_r1
00091       module procedure xmpi_allreduce_i0
00092    end interface xmpi_allreduce
00093 
00094    interface xmpi_reduce
00095       module procedure xmpi_reduce_r0
00096       module procedure xmpi_reduce_i0
00097       module procedure xmpi_reduce_r1
00098       module procedure xmpi_reduce_i1
00099    end interface xmpi_reduce
00100 
00101    interface xmpi_sendrecv
00102       module procedure xmpi_sendrecv_r1
00103       module procedure xmpi_sendrecv_r2
00104       module procedure xmpi_sendrecv_r3
00105       module procedure xmpi_sendrecv_i1
00106       module procedure xmpi_sendrecv_i2
00107    end interface xmpi_sendrecv
00108 
00109    interface xmpi_shift
00110       module procedure xmpi_shift_r2
00111       module procedure xmpi_shift_i2
00112       module procedure xmpi_shift_r3
00113       module procedure xmpi_shift_i3
00114 
00115       module procedure xmpi_shift_r2_l
00116       module procedure xmpi_shift_r3_l
00117    end interface xmpi_shift
00118 
00119    interface xmpi_shift_ee
00120       module procedure xmpi_shift_ee_r2
00121       module procedure xmpi_shift_ee_r3
00122    end interface xmpi_shift_ee
00123 
00124    interface xmpi_shift_uu
00125       module procedure xmpi_shift_uu_r2
00126       module procedure xmpi_shift_uu_r3
00127    end interface xmpi_shift_uu
00128 
00129    interface xmpi_shift_vv
00130       module procedure xmpi_shift_vv_r2
00131       module procedure xmpi_shift_vv_r3
00132    end interface xmpi_shift_vv
00133 
00134    interface xmpi_shift_zs
00135       module procedure xmpi_shift_zs_r2
00136       module procedure xmpi_shift_zs_r3
00137    end interface xmpi_shift_zs
00138 
00139 #endif
00140 contains
00141 #ifdef USEMPI
00142 
00143    subroutine xmpi_initialize
00144       ! initialize mpi environment
00145       implicit none
00146       integer ierr
00147       ierr = 0
00148       ! Message buffers in openmpi are not initialized so this call can give a vallgrind error
00149       ! http://www.open-mpi.org/community/lists/users/2009/06/9566.php
00150       ! http://valgrind.org/docs/manual/manual-core.html#manual-core.suppress
00151       call MPI_Init(ierr)
00152 #ifdef USEMPE
00153       call MPE_Log_get_solo_eventid(event_output_start)
00154       call MPE_Log_get_solo_eventid(event_output_end)
00155       call MPE_Log_get_solo_eventid(event_write_start)
00156       call MPE_Log_get_solo_eventid(event_write_end)
00157       call MPE_Log_get_solo_eventid(event_coll_start)
00158       call MPE_Log_get_solo_eventid(event_coll_end)
00159       call MPE_Describe_event(event_output_start,'out_start','red')
00160       call MPE_Describe_event(event_output_end,'out_end','blue')
00161       call MPE_Describe_event(event_write_start,'write_start','orange')
00162       call MPE_Describe_event(event_write_end,'write_end','green')
00163       call MPE_Describe_event(event_coll_start,'coll_start','white')
00164       call MPE_Describe_event(event_coll_end,'coll_end','white')
00165 #endif
00166       xmpi_comm = MPI_COMM_WORLD
00167       call MPI_Comm_rank(xmpi_comm,xmpi_rank,ierr)
00168       call MPI_Comm_size(xmpi_comm,xmpi_size,ierr)
00169       xmpi_master = 0
00170       xmaster  = (xmpi_rank .eq. xmpi_master)
00171    end subroutine xmpi_initialize
00172 
00173    subroutine xmpi_finalize
00174       ! ends mpi environment, collective subroutine
00175       implicit none
00176       integer ierr
00177       call MPI_Finalize(ierr)
00178    end subroutine xmpi_finalize
00179 
00180    subroutine xmpi_abort
00181       ! to be used if program has to end, and there is no way
00182       ! to tell other processes about this fact
00183       implicit none
00184       integer ierr
00185       call MPI_Abort(xmpi_comm,1,ierr)
00186       stop 1
00187    end subroutine xmpi_abort
00188 
00189    subroutine xmpi_determine_processor_grid(m,n,divtype,mmanual,nmanual,error)
00190       implicit none
00191       integer, intent(in)     :: m,n,mmanual,nmanual  ! the dimensions of the global domain
00192       integer, intent(out)    :: error
00193       integer, intent(in)     :: divtype
00194 
00195       integer mm,nn, borderlength, min_borderlength
00196 
00197       ! determine the processor grid (xmpi_m ampi_n), such that
00198       ! - xmpi_m * xmpi_n = xmpi_size
00199       ! - the total length of the internal borders is minimal
00200 
00201       min_borderlength = 1000000000
00202       error = 0
00203       select case(divtype)
00204        case(MPIBOUNDARY_Y)   ! Force all subdivisions to run along y-lines
00205          xmpi_m=xmpi_size
00206          xmpi_n=1
00207        case(MPIBOUNDARY_X)  ! Force all subdivisions to run along x-lines
00208          xmpi_m=1
00209          xmpi_n=xmpi_size
00210        case (MPIBOUNDARY_AUTO)
00211          do mm = 1,xmpi_size
00212             nn = xmpi_size/mm
00213             if (mm * nn .eq. xmpi_size) then
00214                borderlength = (mm - 1)*n + (nn -1)*m
00215                if (borderlength .lt. min_borderlength) then
00216                   xmpi_m = mm
00217                   xmpi_n = nn
00218                   min_borderlength = borderlength
00219                endif
00220             endif
00221          enddo
00222        case (MPIBOUNDARY_MAN)
00223          if (mmanual * nmanual .ne. xmpi_size) then
00224             error = 2
00225             return
00226          endif
00227          xmpi_m = mmanual
00228          xmpi_n = nmanual
00229        case default
00230          error = 1
00231          return
00232       end select
00233 
00234       ! The layout of the processors is as follows:
00235       ! example 12 processors, xmpi_m=3, xmpi_n=4:
00236 
00237       ! 0 3 6  9
00238       ! 1 4 7 10
00239       ! 2 5 8 11
00240 
00241       !        top
00242       !   left  X  right
00243       !        bot
00244 
00245       ! Left neighbour:
00246 
00247       xmpi_left   = xmpi_rank - xmpi_m
00248       xmpi_isleft = .false.
00249       if (xmpi_left .lt. 0) then
00250          xmpi_left   = MPI_PROC_NULL
00251          xmpi_isleft = .true.
00252       endif
00253 
00254       ! Right neighbour:
00255 
00256       xmpi_right   = xmpi_rank + xmpi_m
00257       xmpi_isright = .false.
00258       if (xmpi_right .ge. xmpi_size) then
00259          xmpi_right   = MPI_PROC_NULL
00260          xmpi_isright = .true.
00261       endif
00262 
00263       ! Upper neighbour:
00264 
00265       if (mod (xmpi_rank,xmpi_m) .eq. 0) then
00266          xmpi_top   = MPI_PROC_NULL
00267          xmpi_istop = .true.
00268       else
00269          xmpi_top = xmpi_rank - 1
00270          xmpi_istop = .false.
00271       endif
00272 
00273       ! Lower neighbour:
00274 
00275       if (mod (xmpi_rank+1,xmpi_m) .eq. 0) then
00276          xmpi_bot   = MPI_PROC_NULL
00277          xmpi_isbot = .true.
00278       else
00279          xmpi_bot   = xmpi_rank + 1
00280          xmpi_isbot = .false.
00281       endif
00282 
00283       ! my row and column (starting with 1,1)
00284 
00285       xmpi_pcol = xmpi_rank/xmpi_m
00286       xmpi_prow = xmpi_rank - xmpi_pcol*xmpi_m
00287       xmpi_pcol = xmpi_pcol+1
00288       xmpi_prow = xmpi_prow+1
00289 
00290    end subroutine xmpi_determine_processor_grid
00291 
00292    subroutine xmpi_bcast_array_logical(x)
00293       implicit none
00294       logical, dimension(:) :: x
00295       integer ierror,l
00296       l = size(x)
00297       call MPI_Bcast(x, l, MPI_LOGICAL, xmpi_master, xmpi_comm, ierror)
00298    end subroutine xmpi_bcast_array_logical
00299 
00300    subroutine xmpi_bcast_logical(x)
00301       implicit none
00302       logical x
00303       integer ierror
00304       call MPI_Bcast(x, 1, MPI_LOGICAL, xmpi_master, xmpi_comm, ierror)
00305    end subroutine xmpi_bcast_logical
00306 
00307    subroutine xmpi_bcast_array_real8(x)
00308       implicit none
00309       real*8, dimension(:) :: x
00310       integer ierror,l
00311       l = size(x)
00312       call MPI_Bcast(x, l, MPI_DOUBLE_PRECISION, xmpi_master, xmpi_comm, ierror)
00313    end subroutine xmpi_bcast_array_real8
00314 
00315    subroutine xmpi_bcast_matrix_real8(x)
00316       implicit none
00317       real*8, dimension(:,:) :: x
00318       integer ierror,l
00319       l = size(x)
00320       call MPI_Bcast(x, l, MPI_DOUBLE_PRECISION, xmpi_master, xmpi_comm, ierror)
00321    end subroutine xmpi_bcast_matrix_real8
00322 
00323    subroutine xmpi_bcast_matrix_integer(x)
00324       implicit none
00325       integer, dimension(:,:) :: x
00326       integer ierror,l
00327       l = size(x)
00328       call MPI_Bcast(x, l, MPI_INTEGER, xmpi_master, xmpi_comm, ierror)
00329    end subroutine xmpi_bcast_matrix_integer
00330 
00331    subroutine xmpi_bcast_array_integer(x)
00332       implicit none
00333       integer, dimension(:) :: x
00334       integer ierror,l
00335       l = size(x)
00336       call MPI_Bcast(x, l, MPI_INTEGER, xmpi_master, xmpi_comm, ierror)
00337    end subroutine xmpi_bcast_array_integer
00338 
00339    subroutine xmpi_bcast_real4(x)
00340       implicit none
00341       real*4 x
00342       integer ierror
00343       call MPI_Bcast(x, 1, MPI_REAL, xmpi_master, xmpi_comm, ierror)
00344    end subroutine xmpi_bcast_real4
00345 
00346    subroutine xmpi_bcast_real8(x)
00347       implicit none
00348       real*8 x
00349       integer ierror
00350       call MPI_Bcast(x, 1, MPI_DOUBLE_PRECISION, xmpi_master, xmpi_comm, ierror)
00351    end subroutine xmpi_bcast_real8
00352 
00353    subroutine xmpi_bcast_int(x)
00354       implicit none
00355       integer x
00356       integer ierror
00357       call MPI_Bcast(x, 1, MPI_INTEGER, xmpi_master, xmpi_comm, ierror)
00358    end subroutine xmpi_bcast_int
00359 
00360    subroutine xmpi_bcast_int8(x)
00361       implicit none
00362       integer*8 x
00363       integer ierror
00364       call MPI_Bcast(x, 1, MPI_INTEGER8, xmpi_master, xmpi_comm, ierror)
00365    end subroutine xmpi_bcast_int8
00366 
00367    subroutine xmpi_bcast_complex16(x)
00368       implicit none
00369       complex*16 x
00370       integer ierror
00371       call MPI_Bcast(x, 1, MPI_DOUBLE_COMPLEX, xmpi_master, xmpi_comm, ierror)
00372    end subroutine xmpi_bcast_complex16
00373 
00374    subroutine xmpi_bcast_char(x)
00375       !
00376       ! wwvv convert string to integer array,
00377       !  broadcast integer array and convert back
00378       !  Not very efficient, but this routine is seldom called
00379       !  and now it works for every taste of fortran90
00380       !
00381       implicit none
00382       character(len=*)                   :: x
00383 
00384       integer                            :: l,i
00385       integer, allocatable, dimension(:) :: sx
00386 
00387       if (xmaster) then
00388          l = len_trim(x)
00389       endif
00390 
00391       call xmpi_bcast(l)
00392       allocate(sx(l))
00393 
00394       if (xmaster) then
00395          do i = 1,l
00396             sx(i) = ichar(x(i:i))
00397          enddo
00398       endif
00399 
00400       call xmpi_bcast(sx)
00401 
00402       if ( .not. xmaster) then
00403          x = ' '
00404          do i = 1,l
00405             x(i:i) = char(sx(i))
00406          enddo
00407       endif
00408 
00409       deallocate(sx)
00410 
00411    end subroutine xmpi_bcast_char
00412 
00413    subroutine xmpi_sendrecv_r1(sendbuf,dest,recvbuf,source)
00414       implicit none
00415       real*8, dimension(:), intent(in)  :: sendbuf
00416       real*8, dimension(:), intent(out) :: recvbuf
00417       integer, intent(in)               :: dest,source
00418 
00419       integer                           :: ierror
00420       integer                           :: n
00421 
00422       n = size(sendbuf)
00423 
00424       call MPI_Sendrecv(sendbuf,n,MPI_DOUBLE_PRECISION,dest,100,   &
00425       recvbuf,n,MPI_DOUBLE_PRECISION,source,100, &
00426       xmpi_comm,MPI_STATUS_IGNORE,ierror)
00427 
00428    end subroutine xmpi_sendrecv_r1
00429 
00430    subroutine xmpi_sendrecv_r2(sendbuf,dest,recvbuf,source)
00431       implicit none
00432       real*8, dimension(:,:), intent(in)  :: sendbuf
00433       real*8, dimension(:,:), intent(out) :: recvbuf
00434       integer, intent(in)                 :: dest,source
00435 
00436       integer                             :: ierror
00437       integer                             :: n
00438 
00439       n = size(sendbuf)
00440 
00441       call MPI_Sendrecv(sendbuf,n,MPI_DOUBLE_PRECISION,dest,101,   &
00442       recvbuf,n,MPI_DOUBLE_PRECISION,source,101, &
00443       xmpi_comm,MPI_STATUS_IGNORE,ierror)
00444 
00445    end subroutine xmpi_sendrecv_r2
00446 
00447    subroutine xmpi_sendrecv_r3(sendbuf,dest,recvbuf,source)
00448       implicit none
00449       real*8, dimension(:,:,:), intent(in)  :: sendbuf
00450       real*8, dimension(:,:,:), intent(out) :: recvbuf
00451       integer, intent(in)                   :: dest,source
00452 
00453       integer                               :: ierror
00454       integer                               :: n
00455 
00456       n = size(sendbuf)
00457 
00458       call MPI_Sendrecv(sendbuf,n,MPI_DOUBLE_PRECISION,dest,101,   &
00459       recvbuf,n,MPI_DOUBLE_PRECISION,source,101, &
00460       xmpi_comm,MPI_STATUS_IGNORE,ierror)
00461 
00462    end subroutine xmpi_sendrecv_r3
00463 
00464    subroutine xmpi_sendrecv_i1(sendbuf,dest,recvbuf,source)
00465       implicit none
00466       integer, dimension(:), intent(in)  :: sendbuf
00467       integer, dimension(:), intent(out) :: recvbuf
00468       integer, intent(in)                :: dest,source
00469 
00470       integer                            :: ierror
00471       integer                            :: n
00472 
00473       n = size(sendbuf)
00474 
00475       call MPI_Sendrecv(sendbuf,n,MPI_INTEGER,dest,102,   &
00476       recvbuf,n,MPI_INTEGER,source,102, &
00477       xmpi_comm,MPI_STATUS_IGNORE,ierror)
00478 
00479    end subroutine xmpi_sendrecv_i1
00480 
00481    subroutine xmpi_sendrecv_i2(sendbuf,dest,recvbuf,source)
00482       implicit none
00483       integer, dimension(:,:), intent(in)  :: sendbuf
00484       integer, dimension(:,:), intent(out) :: recvbuf
00485       integer, intent(in)                  :: dest,source
00486 
00487       integer                              :: ierror
00488       integer                              :: n
00489 
00490       n = size(sendbuf)
00491 
00492       call MPI_Sendrecv(sendbuf,n,MPI_INTEGER,dest,103,   &
00493       recvbuf,n,MPI_INTEGER,source,103, &
00494       xmpi_comm,MPI_STATUS_IGNORE,ierror)
00495 
00496    end subroutine xmpi_sendrecv_i2
00497 
00498    subroutine xmpi_allreduce_r0(x,op)
00499       implicit none
00500       real*8,intent(inout)  :: x
00501       integer,intent(in)    :: op
00502 
00503       real*8  :: y
00504       integer :: ierror
00505       y = x
00506       call MPI_Allreduce(y,x,1,MPI_DOUBLE_PRECISION,op,xmpi_comm,ierror)
00507    end subroutine xmpi_allreduce_r0
00508 
00509    subroutine xmpi_allreduce_r1(x,op)
00510       implicit none
00511       real*8,dimension(:), intent(inout)  :: x
00512       real*8,dimension(:), allocatable    :: y
00513       integer,intent(in)    :: op
00514 
00515       integer :: ierror
00516       allocate(y(size(x)))
00517       y = x
00518       call MPI_Allreduce(y,x,size(x),MPI_DOUBLE_PRECISION,op,xmpi_comm,ierror)
00519       deallocate(y)
00520    end subroutine xmpi_allreduce_r1
00521 
00522    subroutine xmpi_allreduce_i0(x,op)
00523       implicit none
00524       integer,intent(inout)  :: x
00525       integer,intent(in)    :: op
00526 
00527       integer :: y
00528       integer :: ierror
00529       y = x
00530       call MPI_Allreduce(y,x,1,MPI_INTEGER,op,xmpi_comm,ierror)
00531    end subroutine xmpi_allreduce_i0
00532 
00533    subroutine xmpi_reduce_r0(x,y,op)
00534       implicit none
00535       real*8, intent(in)   :: x
00536       real*8, intent(out)  :: y
00537       integer, intent(in)  :: op
00538 
00539       integer :: ierror
00540       call MPI_Reduce(x,y,1,MPI_DOUBLE_PRECISION,op,xmpi_master,xmpi_comm,ierror)
00541    end subroutine xmpi_reduce_r0
00542 
00543    subroutine xmpi_reduce_r1(x,y,op)
00544       implicit none
00545       real*8,dimension(:), intent(in)  :: x
00546       real*8,dimension(:), intent(out) :: y
00547       integer, intent(in)              :: op
00548 
00549       integer :: ierror
00550       call MPI_Reduce(x,y,size(x),MPI_DOUBLE_PRECISION,op,xmpi_master,xmpi_comm,ierror)
00551    end subroutine xmpi_reduce_r1
00552 
00553    subroutine xmpi_reduce_i0(x,y,op)
00554       implicit none
00555       integer, intent(in)   :: x
00556       integer, intent(out)  :: y
00557       integer, intent(in)   :: op
00558 
00559       integer :: ierror
00560       call MPI_Reduce(x,y,1,MPI_INTEGER,op,xmpi_master,xmpi_comm,ierror)
00561    end subroutine xmpi_reduce_i0
00562 
00563    subroutine xmpi_reduce_i1(x,y,op)
00564       implicit none
00565       integer,dimension(:),intent(in)  :: x
00566       integer,dimension(:),intent(out) :: y
00567       integer,intent(in)               :: op
00568 
00569       integer :: ierror
00570       call MPI_Reduce(x,y,size(x),MPI_INTEGER,op,xmpi_master,xmpi_comm,ierror)
00571    end subroutine xmpi_reduce_i1
00572 
00573    !
00574    ! shift routines:  x(m,n) is the matrix in this process
00575    ! direction = 'u': shift up,    send to top   x(2,:) ,  receive from bot   x(m,:)
00576    ! direction = 'd': shift down,  send to bot   x(m-1,:), receive from top   x(1,:)
00577    ! direction = 'l': shift left,  send to left  x(:,2),   receive from right x(:,n)
00578    ! direction = 'r': shift right, send to right x(:,n-1), receive from left  x(:,1)
00579    !
00580    ! also 'm:', '1:', ':n' and ':1' can be used, easier to remember:
00581    ! 'm:' :  x(m,:) will be replaced, except for a far bottom process
00582    ! '1:' :  x(1,:) will be replaced, except for a far top process
00583    ! ':n' :  x(:,n) will be replaced, except for a far right process
00584    ! ':1' :  x(:,1) will be replaced, except for a far left process
00585 
00586    subroutine xmpi_shift_r2(x,direction)
00587       implicit none
00588       real*8,dimension(:,:),intent(inout) :: x
00589       character(len=*),intent(in)         :: direction
00590 
00591       integer :: m,n
00592 
00593       m = size(x,1)
00594       n = size(x,2)
00595 
00596       select case(direction)
00597        case('u','m:')
00598          call xmpi_sendrecv(x(2,:),xmpi_top,    x(m,:),xmpi_bot)
00599        case('d','1:')
00600          call xmpi_sendrecv(x(m-1,:),xmpi_bot,  x(1,:),xmpi_top)
00601        case('l',':n')
00602          call xmpi_sendrecv(x(:,2),xmpi_left,   x(:,n),xmpi_right)
00603        case('r',':1')
00604          call xmpi_sendrecv(x(:,n-1),xmpi_right,x(:,1),xmpi_left)
00605        case default
00606          if(xmaster) then
00607             write (*,*) 'Invalid direction parameter for xmpi_shift_r2: "'// &
00608             direction//'"'
00609             call halt_program
00610          endif
00611       end select
00612 
00613    end subroutine xmpi_shift_r2
00614 
00615    subroutine xmpi_shift_i2(x,direction)
00616       implicit none
00617       integer, dimension (:,:), intent(inout) :: x
00618       character(len=*), intent(in)            :: direction
00619 
00620       integer :: m,n
00621 
00622       m = size(x,1)
00623       n = size(x,2)
00624 
00625       select case(direction)
00626        case('u','m:')
00627          call xmpi_sendrecv(x(2,:),   xmpi_top,  x(m,:),xmpi_bot)
00628        case('d','1:')
00629          call xmpi_sendrecv(x(m-1,:), xmpi_bot,  x(1,:),xmpi_top)
00630        case('l',':n')
00631          call xmpi_sendrecv(x(:,2),   xmpi_left, x(:,n),xmpi_right)
00632        case('r',':1')
00633          call xmpi_sendrecv(x(:,n-1), xmpi_right,x(:,1),xmpi_left)
00634        case default
00635          if(xmaster) then
00636             write (*,*) 'Invalid direction parameter for xmpi_shift_r2: "'// &
00637             direction//'"'
00638             call halt_program
00639          endif
00640       end select
00641 
00642    end subroutine xmpi_shift_i2
00643 
00644    subroutine xmpi_shift_r3(x,direction)
00645       implicit none
00646       real*8, dimension (:,:,:), intent(inout) :: x
00647       character(len=*), intent(in)             :: direction
00648 
00649       integer :: m,n,l
00650 
00651       m = size(x,1)
00652       n = size(x,2)
00653       l = size(x,3)
00654 
00655       select case(direction)
00656        case('u','m:')
00657          call xmpi_sendrecv(x(2,:,:),  xmpi_top,    x(m,:,:),xmpi_bot)
00658        case('d','1:')
00659          call xmpi_sendrecv(x(m-1,:,:),xmpi_bot,  x(1,:,:),xmpi_top)
00660        case('l',':n')
00661          call xmpi_sendrecv(x(:,2,:),  xmpi_left,   x(:,n,:),xmpi_right)
00662        case('r',':1')
00663          call xmpi_sendrecv(x(:,n-1,:),xmpi_right,x(:,1,:),xmpi_left)
00664        case default
00665          if(xmaster) then
00666             write (*,*) 'Invalid direction parameter for xmpi_shift_r3: "'// &
00667             direction//'"'
00668             call halt_program
00669          endif
00670       end select
00671 
00672    end subroutine xmpi_shift_r3
00673 
00674    subroutine xmpi_shift_i3(x,direction)
00675       implicit none
00676       integer, dimension (:,:,:), intent(inout) :: x
00677       character(len=*), intent(in)              :: direction
00678 
00679       integer :: m,n,l
00680 
00681       m = size(x,1)
00682       n = size(x,2)
00683       l = size(x,3)
00684 
00685       select case(direction)
00686        case('u','m:')
00687          call xmpi_sendrecv(x(2,:,:),  xmpi_top,  x(m,:,:),xmpi_bot)
00688        case('d','1:')
00689          call xmpi_sendrecv(x(m-1,:,:),xmpi_bot,  x(1,:,:),xmpi_top)
00690        case('l',':n')
00691          call xmpi_sendrecv(x(:,2,:),  xmpi_left, x(:,n,:),xmpi_right)
00692        case('r',':1')
00693          call xmpi_sendrecv(x(:,n-1,:),xmpi_right,x(:,1,:),xmpi_left)
00694        case default
00695          if(xmaster) then
00696             write (*,*) 'Invalid direction parameter for xmpi_shift_i3: "'// &
00697             direction//'"'
00698             call halt_program
00699          endif
00700       end select
00701 
00702    end subroutine xmpi_shift_i3
00703 
00704    !  translation from nx+1, ny+1 to m and n
00705    !
00706    !  m=nx+1  nx=m-1
00707    !  n=ny+1  ny=n-1
00708    !
00709    !  variabele   l->r      r->l        b->t       t->b
00710    !              r         l           u          d      shift
00711    !          SHIFT_Y_R  SHIFT_Y_L  SHIFT_X_U  SHIFT_X_D
00712 
00713 
00714    !    ee,rr    m-3:m-2,:  3:4,:      :,n-3:n-2  :,3:4
00715    !             1:2,:      m-1:m,:    :,1:2      :,n-1:n
00716    !
00717    !    uu       m-3,:      2:3,:      :,n-3:n-2  :,3:4
00718    !             1,:        m-2,m-1,:  :,1:2      :,n-1:n
00719    !
00720    !    vv       m-3:m-2,:  3:4,:      :,n-3      :,2:3
00721    !             1:2        m-1:m,:    :,1        :,n-2:n-1
00722    !
00723    !  zs,ccg,zb  m-2,:      3,:        :,n-2      :,3
00724    !             2,:        m-1,:      :,2        :,n-1
00725    !
00726    !  Dus:
00727    !        i     p    j     q
00728    !
00729    !        1 <-> m-3  1 <-> n-3
00730    !        2 <-> m-2  2 <-> n-2
00731    !        3 <-> m-1  3 <-> n-1
00732    !        4 <-> m    4 <-> n
00733    !
00734    !        p=m-4+i    q=n-4+j
00735 
00736 
00737    subroutine xmpi_shift_r2_l(x,direction,i1,i2)
00738       implicit none
00739       real*8,dimension(:,:),intent(inout) :: x
00740       integer, intent(in)                 :: direction
00741       integer, intent(in)                 :: i1,i2
00742 
00743       integer            :: m,n,s1,s2,r1,r2
00744       integer, parameter :: nbord = 2, nover = 2*nbord
00745 
00746       ! nover is the number of overlapping rows/columns
00747 
00748       ! s1,s2 will contain the indices of the first and last row/column to send
00749       ! r1,r2 will contain the indices of the first and last row/column to send
00750       ! mn will contain the 1st or 2nd dimension as appropriate:
00751       !   shift in x direction: mn = m (= nx+1)
00752       !   shift in y direction: mn = n (= ny+1)
00753 
00754       m = size(x,1)
00755       n = size(x,2)
00756 
00757       ! sanity check
00758       select case(direction)
00759        case(SHIFT_Y_R,SHIFT_Y_L,SHIFT_X_U,SHIFT_X_D)
00760          continue
00761        case default
00762          if (xmaster) then
00763             write(*,*) 'Invalid value for direction in xmpi_shift_r2_l ',direction
00764             call halt_program
00765          endif
00766       endselect
00767 
00768       select case(direction)
00769        case(SHIFT_Y_R)
00770          s1 = n - nover + i1
00771          s2 = n - nover + i2
00772          r1 = i1
00773          r2 = i2
00774          call xmpi_sendrecv(x(:,s1:s2),xmpi_right,x(:,r1:r2),xmpi_left)
00775        case(SHIFT_X_D)
00776          s1 = m - nover + i1
00777          s2 = m - nover + i2
00778          r1 = i1
00779          r2 = i2
00780          call xmpi_sendrecv(x(s1:s2,:),xmpi_bot,  x(r1:r2,:),xmpi_top)
00781        case(SHIFT_Y_L)
00782          s1 = i1
00783          s2 = i2
00784          r1 = n - nover + i1
00785          r2 = n - nover + i2
00786          call xmpi_sendrecv(x(:,s1:s2),xmpi_left, x(:,r1:r2),xmpi_right)
00787        case(SHIFT_X_U)
00788          s1 = i1
00789          s2 = i2
00790          r1 = m - nover + i1
00791          r2 = m - nover + i2
00792          call xmpi_sendrecv(x(s1:s2,:),xmpi_top,  x(r1:r2,:),xmpi_bot)
00793       endselect
00794 
00795    end subroutine xmpi_shift_r2_l
00796 
00797    subroutine xmpi_shift_r3_l(x,direction,i1,i2)
00798       !
00799       ! TODO   niet goed, zie _r2 boven
00800       implicit none
00801       real*8,dimension(:,:,:),intent(inout) :: x
00802       integer, intent(in)                   :: direction
00803       integer, intent(in)                   :: i1,i2
00804 
00805       integer            :: m,n,s1,s2,r1,r2
00806       integer, parameter :: nbord = 2, nover = 2*nbord
00807 
00808 
00809       ! nover is the number of overlapping rows/columns
00810 
00811       ! s1,s2 will contain the indices of the first and last row/column to send
00812       ! r1,r2 will contain the indices of the first and last row/column to send
00813       ! mn will contain the 1st or 2nd dimension as appropriate:
00814       !   shift in x direction: mn = m (= nx+1)
00815       !   shift in y direction: mn = n (= ny+1)
00816 
00817       m = size(x,1)
00818       n = size(x,2)
00819 
00820       ! sanity check
00821       select case(direction)
00822        case(SHIFT_Y_R,SHIFT_Y_L,SHIFT_X_U,SHIFT_X_D)
00823          continue
00824        case default
00825          if (xmaster) then
00826             write(*,*) 'Invalid value for direction in xmpi_shift_r2_l ',direction
00827             call halt_program
00828          endif
00829       endselect
00830 
00831       select case(direction)
00832        case(SHIFT_Y_R)
00833          s1 = n - nover + i1
00834          s2 = n - nover + i2
00835          r1 = i1
00836          r2 = i2
00837          call xmpi_sendrecv(x(:,s1:s2,:),xmpi_right,x(:,r1:r2,:),xmpi_left)
00838        case(SHIFT_X_D)
00839          s1 = m - nover + i1
00840          s2 = m - nover + i2
00841          r1 = i1
00842          r2 = i2
00843          call xmpi_sendrecv(x(s1:s2,:,:),xmpi_bot,  x(r1:r2,:,:),xmpi_top)
00844        case(SHIFT_Y_L)
00845          s1 = i1
00846          s2 = i2
00847          r1 = n - nover + i1
00848          r2 = n - nover + i2
00849          call xmpi_sendrecv(x(:,s1:s2,:),xmpi_left, x(:,r1:r2,:),xmpi_right)
00850        case(SHIFT_X_U)
00851          s1 = i1
00852          s2 = i2
00853          r1 = m - nover + i1
00854          r2 = m - nover + i2
00855          call xmpi_sendrecv(x(s1:s2,:,:),xmpi_top,  x(r1:r2,:,:),xmpi_bot)
00856       endselect
00857    end subroutine xmpi_shift_r3_l
00858 
00859    subroutine xmpi_shift_ee_r2(x)
00860       implicit none
00861       real*8,dimension(:,:),intent(inout) :: x
00862       call xmpi_shift(x,SHIFT_Y_R,1,2)
00863       call xmpi_shift(x,SHIFT_Y_L,3,4)
00864       call xmpi_shift(x,SHIFT_X_U,3,4)
00865       call xmpi_shift(x,SHIFT_X_D,1,2)
00866    end subroutine xmpi_shift_ee_r2
00867 
00868    subroutine xmpi_shift_ee_r3(x)
00869       implicit none
00870       real*8,dimension(:,:,:),intent(inout) :: x
00871       call xmpi_shift(x,SHIFT_Y_R,1,2)
00872       call xmpi_shift(x,SHIFT_Y_L,3,4)
00873       call xmpi_shift(x,SHIFT_X_U,3,4)
00874       call xmpi_shift(x,SHIFT_X_D,1,2)
00875    end subroutine xmpi_shift_ee_r3
00876 
00877    subroutine xmpi_shift_uu_r2(x)
00878       implicit none
00879       real*8,dimension(:,:),intent(inout) :: x
00880       call xmpi_shift(x,SHIFT_Y_R,1,2)
00881       call xmpi_shift(x,SHIFT_Y_L,3,4)
00882       call xmpi_shift(x,SHIFT_X_U,2,3)
00883       call xmpi_shift(x,SHIFT_X_D,1,1)
00884    end subroutine xmpi_shift_uu_r2
00885 
00886    subroutine xmpi_shift_uu_r3(x)
00887       implicit none
00888       real*8,dimension(:,:,:),intent(inout) :: x
00889       call xmpi_shift(x,SHIFT_Y_R,1,2)
00890       call xmpi_shift(x,SHIFT_Y_L,3,4)
00891       call xmpi_shift(x,SHIFT_X_U,2,3)
00892       call xmpi_shift(x,SHIFT_X_D,1,1)
00893    end subroutine xmpi_shift_uu_r3
00894 
00895    subroutine xmpi_shift_vv_r2(x)
00896       implicit none
00897       real*8,dimension(:,:),intent(inout) :: x
00898       call xmpi_shift(x,SHIFT_Y_R,1,1)
00899       call xmpi_shift(x,SHIFT_Y_L,2,3)
00900       call xmpi_shift(x,SHIFT_X_U,3,4)
00901       call xmpi_shift(x,SHIFT_X_D,1,2)
00902    end subroutine xmpi_shift_vv_r2
00903 
00904    subroutine xmpi_shift_vv_r3(x)
00905       implicit none
00906       real*8,dimension(:,:,:),intent(inout) :: x
00907       call xmpi_shift(x,SHIFT_Y_R,1,1)
00908       call xmpi_shift(x,SHIFT_Y_L,2,3)
00909       call xmpi_shift(x,SHIFT_X_U,3,4)
00910       call xmpi_shift(x,SHIFT_X_D,1,2)
00911    end subroutine xmpi_shift_vv_r3
00912 
00913    subroutine xmpi_shift_zs_r2(x)
00914       implicit none
00915       real*8,dimension(:,:),intent(inout) :: x
00916       call xmpi_shift(x,SHIFT_Y_R,2,2)
00917       call xmpi_shift(x,SHIFT_Y_L,3,3)
00918       call xmpi_shift(x,SHIFT_X_U,3,3)
00919       call xmpi_shift(x,SHIFT_X_D,2,2)
00920    end subroutine xmpi_shift_zs_r2
00921 
00922    subroutine xmpi_shift_zs_r3(x)
00923       implicit none
00924       real*8,dimension(:,:,:),intent(inout) :: x
00925       call xmpi_shift(x,SHIFT_Y_R,2,2)
00926       call xmpi_shift(x,SHIFT_Y_L,3,3)
00927       call xmpi_shift(x,SHIFT_X_U,3,3)
00928       call xmpi_shift(x,SHIFT_X_D,2,2)
00929    end subroutine xmpi_shift_zs_r3
00930 
00931    subroutine testje
00932       implicit none
00933       real*8, dimension(10,10) :: x
00934       real*8, dimension(10,10,10) :: y
00935       call xmpi_shift_ee(x)
00936       call xmpi_shift_uu(x)
00937       call xmpi_shift_vv(x)
00938       call xmpi_shift_zs(x)
00939       call xmpi_shift_ee(y)
00940       call xmpi_shift_uu(y)
00941       call xmpi_shift_vv(y)
00942       call xmpi_shift_zs(y)
00943    end subroutine testje
00944 
00945    subroutine xmpi_barrier
00946       implicit none
00947       integer ierror
00948       call MPI_Barrier(xmpi_comm,ierror)
00949    end subroutine xmpi_barrier
00950 
00951    !
00952    ! get a row from a matrix in the same processor column
00953    !
00954    subroutine xmpi_getrow(a,n,l,prow,b)
00955       implicit none
00956       real*8, dimension(:,:), intent(in) :: a    ! the matrix
00957       integer, intent(in)                :: n    ! number of elements in the row
00958       character, intent(in)              :: l    ! '1': get first row
00959       ! 'm': get last row
00960       integer, intent(in)                :: prow ! row number of process
00961       ! to get the row from
00962       real*8, dimension(:), intent(out)  :: b    ! the row from process prow
00963 
00964       ! Note: a and l are only needed at the sending process
00965 
00966       integer                            :: row,ierror
00967       integer                            :: ll,dest,source,tag,r
00968       real*8, dimension(n)               :: rowdata
00969       integer, allocatable, dimension(:) :: requests
00970 
00971       source = (xmpi_pcol-1)*xmpi_m + prow -1 ! Sending process
00972       ll = 1
00973 
00974       if (source .eq. xmpi_rank) then  ! the sending process
00975          select case(l)
00976           case('1')
00977             ll = 1
00978           case('m')
00979             ll = size(a,1)
00980           case default
00981             write(*,*) 'Error in xmpi_getrow, l="'//l//'"'
00982             call halt_program
00983          end select
00984 
00985          rowdata = a(ll,:)
00986          allocate(requests(xmpi_m-1))
00987          dest = (xmpi_pcol-1)*xmpi_m    ! First receiving process
00988          r = 0
00989          do row = 1,xmpi_m
00990             if (dest .eq. xmpi_rank) then       ! do no send to myself
00991                b = rowdata                       ! but copy
00992             else
00993                r = r + 1
00994                tag = dest
00995                call MPI_Isend(rowdata(1), n, MPI_DOUBLE_PRECISION,  &
00996                dest, tag, xmpi_comm, requests(r), ierror)
00997             endif
00998             dest = dest + 1                     ! next receiving process
00999          enddo
01000          call MPI_Waitall(r,requests,MPI_STATUSES_IGNORE,ierror)
01001          deallocate(requests)
01002       else                                    ! receiving process
01003          tag = xmpi_rank
01004          call MPI_Recv(b, n, MPI_DOUBLE_PRECISION,  &
01005          source, tag, xmpi_comm, MPI_STATUS_IGNORE, ierror)
01006       endif
01007 
01008       ! wwvv if this is needed often, than a neat subroutine, using
01009       !  column communicators and MPI_Bcast would be appropriate
01010 
01011    end subroutine xmpi_getrow
01012 
01013 
01014 #endif
01015    subroutine halt_program
01016 #ifdef USEMPI
01017       call xmpi_abort
01018 #else
01019       stop 1
01020 #endif
01021    end subroutine halt_program
01022 end module xmpi_module
 All Classes Files Functions Variables Defines