added source code
[unres.git] / source / wham / src / wham_multparm.F
1       program WHAM_multparm
2 c Creation/update of the database of conformations
3       implicit none
4 #ifndef ISNAN
5       external proc_proc
6 #endif
7 #ifdef WINPGI
8 cMS$ATTRIBUTES C ::  proc_proc
9 #endif
10       include "DIMENSIONS"
11       include "DIMENSIONS.ZSCOPT"
12       include "DIMENSIONS.FREE"
13 #ifdef MPI
14       include "mpif.h"
15       integer IERROR,ERRCODE
16       include "COMMON.MPI"
17 #endif
18       include "COMMON.IOUNITS"
19       include "COMMON.FREE"
20       include "COMMON.CONTROL"
21       include "COMMON.ALLPARM"
22       include "COMMON.PROT"
23       double precision rr,x(max_paropt)
24       integer idumm
25       integer i,ipar,islice
26 #ifdef MPI
27       call MPI_Init( IERROR )
28       call MPI_Comm_rank( MPI_COMM_WORLD, me, IERROR )
29       call MPI_Comm_size( MPI_COMM_WORLD, nprocs, IERROR )
30       Master = 0
31       if (ierror.gt.0) then
32         write(iout,*) "SEVERE ERROR - Can't initialize MPI."
33         call mpi_finalize(ierror)
34         stop
35       endif
36       if (nprocs.gt.MaxProcs+1) then
37         write (2,*) "Error - too many processors",
38      &   nprocs,MaxProcs+1
39         write (2,*) "Increase MaxProcs and recompile"
40         call MPI_Finalize(IERROR)
41         stop
42       endif
43 #endif
44 c NaNQ initialization
45 #ifndef ISNAN
46       i=-1
47       rr=dacos(100.0d0)
48 #ifdef WINPGI
49       idumm=proc_proc(rr,i)
50 #else
51       call proc_proc(rr,i)
52 #endif
53 #endif
54       call initialize
55       call openunits
56       call cinfo
57       call read_general_data(*10)
58       call flush(iout)
59       call molread(*10)
60       call flush(iout)
61 #ifdef MPI 
62       write (iout,*) "Calling proc_groups"
63       call proc_groups
64       write (iout,*) "proc_groups exited"
65       call flush(iout)
66 #endif
67       do ipar=1,nParmSet
68         write (iout,*) "Calling parmread",ipar
69         call parmread(ipar,*10)
70         if (.not.separate_parset) then
71           call store_parm(ipar)
72           write (iout,*) "Finished storing parameters",ipar
73         else if (ipar.eq.myparm) then
74           call store_parm(1)
75           write (iout,*) "Finished storing parameters",ipar
76         endif
77         call flush(iout)
78       enddo
79       call read_efree(*10)
80       write (iout,*) "Finished READ_EFREE"
81       call flush(iout)
82       call read_protein_data(*10)
83       write (iout,*) "Finished READ_PROTEIN_DATA"
84       call flush(iout)
85       if (indpdb.gt.0) then
86         call promienie
87         call read_compar
88         call read_ref_structure(*10)
89         call proc_cont
90         call fragment_list
91       endif
92       write (iout,*) "Begin read_database"
93       call flush(iout)
94       call read_database(*10)
95       write (iout,*) "Finished read_database"
96       call flush(iout)
97       if (separate_parset) nparmset=1
98       do islice=1,nslice
99         if (ntot(islice).gt.0) then
100 #ifdef MPI 
101           call work_partition(islice,.true.)
102           write (iout,*) "work_partition OK"
103           call flush(iout)
104 #endif
105           call enecalc(islice,*10)
106           write (iout,*) "enecalc OK"
107           call flush(iout)
108           call WHAM_CALC(islice,*10)
109           write (iout,*) "wham_calc OK"
110           call flush(iout)
111           call write_dbase(islice,*10)
112           write (iout,*) "write_dbase OK"
113           call flush(iout)
114           if (ensembles.gt.0) then
115             call make_ensembles(islice,*10)
116             write (iout,*) "make_ensembles OK"
117             call flush(iout)
118           endif
119         endif
120       enddo
121 #ifdef MPI
122       call MPI_Finalize( IERROR )
123 #endif
124       stop
125    10 write (iout,*) "Error termination of the program"
126       call MPI_Finalize( IERROR )
127       stop
128       end
129 c------------------------------------------------------------------------------
130 #ifdef MPI
131       subroutine proc_groups
132 C Split the processors into the Master and Workers group, if needed.
133       implicit none
134       include "DIMENSIONS"
135       include "DIMENSIONS.ZSCOPT"
136       include "DIMENSIONS.FREE"
137       include "mpif.h"
138       include "COMMON.IOUNITS"
139       include "COMMON.MPI"
140       include "COMMON.FREE"
141       integer n,chunk,i,j,ii,remainder
142       integer kolor,key,ierror,errcode
143       logical lprint
144       lprint=.true.
145
146 C Split the communicator if independent runs for different parameter
147 C sets will be performed.
148 C
149       if (nparmset.eq.1 .or. .not.separate_parset) then
150         WHAM_COMM = MPI_COMM_WORLD
151       else if (separate_parset) then
152         if (nprocs.lt.nparmset) then
153           write (iout,*) 
154      & "*** Cannot split parameter sets for fewer processors than sets",
155      &  nprocs,nparmset
156           call MPI_Finalize(ierror)
157           stop
158         endif 
159         write (iout,*) "nparmset",nparmset
160         nprocs = nprocs/nparmset
161         kolor = me/nprocs
162         key = mod(me,nprocs)
163         write (iout,*) "My old rank",me," kolor",kolor," key",key
164         call MPI_Comm_split(MPI_COMM_WORLD,kolor,key,WHAM_COMM,ierror)
165         call MPI_Comm_size(WHAM_COMM,nprocs,ierror)
166         call MPI_Comm_rank(WHAM_COMM,me,ierror)
167         write (iout,*) "My new rank",me," comm size",nprocs
168         write (iout,*) "MPI_COMM_WORLD",MPI_COMM_WORLD,
169      &   " WHAM_COMM",WHAM_COMM
170         myparm=kolor+1
171         write (iout,*) "My parameter set is",myparm
172         call flush(iout)
173       else
174         myparm=nparmset
175       endif
176       Me1 = Me
177       Nprocs1 = Nprocs
178       return
179       end
180 c------------------------------------------------------------------------------
181       subroutine work_partition(islice,lprint)
182 c Split the conformations between processors
183       implicit none
184       include "DIMENSIONS"
185       include "DIMENSIONS.ZSCOPT"
186       include "DIMENSIONS.FREE"
187       include "mpif.h"
188       include "COMMON.IOUNITS"
189       include "COMMON.MPI"
190       include "COMMON.PROT"
191       integer islice
192       integer n,chunk,i,j,ii,remainder
193       integer kolor,key,ierror,errcode
194       logical lprint
195 C
196 C Divide conformations between processors; the first and
197 C the last conformation to handle by ith processor is stored in 
198 C indstart(i) and indend(i), respectively.
199 C
200 C First try to assign equal number of conformations to each processor.
201 C
202         n=ntot(islice)
203         write (iout,*) "n=",n
204         indstart(0)=1
205         chunk = N/nprocs1
206         scount(0) = chunk
207 c        print *,"i",0," indstart",indstart(0)," scount",
208 c     &     scount(0)
209         do i=1,nprocs1-1
210           indstart(i)=chunk+indstart(i-1) 
211           scount(i)=scount(i-1)
212 c          print *,"i",i," indstart",indstart(i)," scount",
213 c     &     scount(i)
214         enddo 
215 C
216 C Determine how many conformations remained yet unassigned.
217 C
218         remainder=N-(indstart(nprocs1-1)
219      &    +scount(nprocs1-1)-1)
220 c        print *,"remainder",remainder
221 C
222 C Assign the remainder conformations to consecutive processors, starting
223 C from the lowest rank; this continues until the list is exhausted.
224 C
225         if (remainder .gt. 0) then 
226           do i=1,remainder
227             scount(i-1) = scount(i-1) + 1
228             indstart(i) = indstart(i) + i
229           enddo
230           do i=remainder+1,nprocs1-1
231             indstart(i) = indstart(i) + remainder
232           enddo
233         endif
234
235         indstart(nprocs1)=N+1
236         scount(nprocs1)=0
237
238         do i=0,NProcs1
239           indend(i)=indstart(i)+scount(i)-1
240           idispl(i)=indstart(i)-1
241         enddo
242
243         N=0
244         do i=0,Nprocs1-1
245           N=N+indend(i)-indstart(i)+1
246         enddo
247
248 c        print *,"N",n," NTOT",ntot(islice)
249         if (N.ne.ntot(islice)) then
250           write (iout,*) "!!! Checksum error on processor",me,
251      &     " slice",islice
252           call flush(iout)
253           call MPI_Abort( MPI_COMM_WORLD, Ierror, Errcode )
254         endif
255
256       if (lprint) then
257         write (iout,*) "Partition of work between processors"
258           do i=0,nprocs1-1
259             write (iout,'(a,i5,a,i7,a,i7,a,i7)')
260      &        "Processor",i," indstart",indstart(i),
261      &        " indend",indend(i)," count",scount(i)
262           enddo
263       endif
264       return
265       end
266 #endif
267 #ifdef AIX
268       subroutine flush(iu)
269       call flush_(iu)
270       return
271       end
272 #endif