rename
[unres4.git] / source / wham / work_partition.F90
1       module work_part
2 !------------------------------------------------------------------------------
3       use io_units
4       use MPI_data
5       use wham_data
6       implicit none
7 #ifdef MPI
8 !------------------------------------------------------------------------------
9 !
10 !
11 !-----------------------------------------------------------------------------
12       contains
13 !-----------------------------------------------------------------------------
14 #ifdef CLUSTER
15       subroutine work_partition(lprint,ncon_work)
16 #else
17       subroutine work_partition(islice,lprint)
18 #endif
19 ! Split the conformations between processors
20 !      implicit none
21 !      include "DIMENSIONS"
22 !      include "DIMENSIONS.ZSCOPT"
23 !      include "DIMENSIONS.FREE"
24       include "mpif.h"
25 !      include "COMMON.IOUNITS"
26 !      include "COMMON.MPI"
27 !      include "COMMON.PROT"
28       integer :: islice,ncon_work
29       integer :: n,chunk,i,j,ii,remainder
30 !el      integer :: kolor
31       integer :: key,ierror,errcode
32       logical :: lprint
33 !
34 ! Divide conformations between processors; the first and
35 ! the last conformation to handle by ith processor is stored in 
36 ! indstart(i) and indend(i), respectively.
37 !
38 !el MPI_data
39       if (.not. allocated(indstart)) allocate(indstart(0:nprocs))
40       if (.not. allocated(indend)) allocate(indend(0:nprocs))
41       if (.not. allocated(idispl)) allocate(idispl(0:nprocs))
42       if (.not. allocated(scount)) allocate(scount(0:nprocs))
43 ! First try to assign equal number of conformations to each processor.
44 !
45 #ifdef CLUSTER
46         n=ncon_work
47         write (iout,*) "n=",n," nprocs=",nprocs
48         nprocs1=nprocs
49 #else
50         n=ntot(islice)
51         write (iout,*) "n=",n
52 #endif
53         indstart(0)=1
54         chunk = N/nprocs1
55         scount(0) = chunk
56 write(iout,*)"chunk",chunk,scount(0)
57 flush(iout)
58 !        print *,"i",0," indstart",indstart(0)," scount",&
59 !          scount(0)
60         do i=1,nprocs1-1
61           indstart(i)=chunk+indstart(i-1) 
62           scount(i)=scount(i-1)
63 !          print *,"i",i," indstart",indstart(i)," scount",
64 !     &     scount(i)
65         enddo 
66 !
67 ! Determine how many conformations remained yet unassigned.
68 !
69         remainder=N-(indstart(nprocs1-1) &
70           +scount(nprocs1-1)-1)
71 !        print *,"remainder",remainder
72 !
73 ! Assign the remainder conformations to consecutive processors, starting
74 ! from the lowest rank; this continues until the list is exhausted.
75 !
76         if (remainder .gt. 0) then 
77           do i=1,remainder
78             scount(i-1) = scount(i-1) + 1
79             indstart(i) = indstart(i) + i
80           enddo
81           do i=remainder+1,nprocs1-1
82             indstart(i) = indstart(i) + remainder
83           enddo
84         endif
85
86         indstart(nprocs1)=N+1
87         scount(nprocs1)=0
88
89         do i=0,NProcs1
90           indend(i)=indstart(i)+scount(i)-1
91           idispl(i)=indstart(i)-1
92         enddo
93
94         N=0
95         do i=0,Nprocs1-1
96           N=N+indend(i)-indstart(i)+1
97         enddo
98
99 !        print *,"N",n," NTOT",ntot(islice)
100 #ifdef CLUSTER
101         if (N.ne.ncon_work) then
102           write (iout,*) "!!! Checksum error on processor",me,&
103             n,ncon_work
104 #else
105         if (N.ne.ntot(islice)) then
106           write (iout,*) "!!! Checksum error on processor",me,&
107            " slice",islice
108 #endif
109           call flush(iout)
110           call MPI_Abort( MPI_COMM_WORLD, Ierror, Errcode )
111         endif
112
113       if (lprint) then
114         write (iout,*) "Partition of work between processors"
115           do i=0,nprocs1-1
116             write (iout,'(a,i5,a,i7,a,i7,a,i7)') &
117               "Processor",i," indstart",indstart(i),&
118               " indend",indend(i)," count",scount(i)
119           enddo
120       endif
121       return
122       end subroutine work_partition
123 #endif
124 !----------------------------------------------------------------------------
125       end module work_part
126 !-----------------------------------------------------------------------------
127 !-----------------------------------------------------------------------------