introduction of genconstr to wham
[unres.git] / source / wham / src-M / 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       if ((constr_dist.gt.0).and.(genconstr.eq.0)) call read_dist_constr
62 #ifdef MPI 
63       write (iout,*) "Calling proc_groups"
64       call proc_groups
65       write (iout,*) "proc_groups exited"
66       call flush(iout)
67 #endif
68       do ipar=1,nParmSet
69         write (iout,*) "Calling parmread",ipar
70         call parmread(ipar,*10)
71         if (.not.separate_parset) then
72           call store_parm(ipar)
73           write (iout,*) "Finished storing parameters",ipar
74         else if (ipar.eq.myparm) then
75           call store_parm(1)
76           write (iout,*) "Finished storing parameters",ipar
77         endif
78         call flush(iout)
79       enddo
80       call read_efree(*10)
81       write (iout,*) "Finished READ_EFREE"
82       call flush(iout)
83       call read_protein_data(*10)
84       write (iout,*) "Finished READ_PROTEIN_DATA"
85       call flush(iout)
86       if (indpdb.gt.0) then
87         call promienie
88         call read_compar
89         call read_ref_structure(*10)
90       if ((constr_dist.gt.0).and.(genconstr.gt.0)) call gen_dist_constr 
91         call proc_cont
92         call fragment_list
93       endif
94
95 C      if (constr_dist.gt.0) call read_dist_constr
96       write (iout,*) "Begin read_database"
97       call flush(iout)
98       call read_database(*10)
99       write (iout,*) "Finished read_database"
100       call flush(iout)
101       if (separate_parset) nparmset=1
102       do islice=1,nslice
103         if (ntot(islice).gt.0) then
104 #ifdef MPI 
105           call work_partition(islice,.true.)
106           write (iout,*) "work_partition OK"
107           call flush(iout)
108 #endif
109           call enecalc(islice,*10)
110           write (iout,*) "enecalc OK"
111           call flush(iout)
112           call WHAM_CALC(islice,*10)
113           write (iout,*) "wham_calc OK"
114           call flush(iout)
115           call write_dbase(islice,*10)
116           write (iout,*) "write_dbase OK"
117           call flush(iout)
118           if (ensembles.gt.0) then
119             call make_ensembles(islice,*10)
120             write (iout,*) "make_ensembles OK"
121             call flush(iout)
122           endif
123         endif
124       enddo
125 #ifdef MPI
126       call MPI_Finalize( IERROR )
127 #endif
128       stop
129    10 write (iout,*) "Error termination of the program"
130       call MPI_Finalize( IERROR )
131       stop
132       end
133 c------------------------------------------------------------------------------
134 #ifdef MPI
135       subroutine proc_groups
136 C Split the processors into the Master and Workers group, if needed.
137       implicit none
138       include "DIMENSIONS"
139       include "DIMENSIONS.ZSCOPT"
140       include "DIMENSIONS.FREE"
141       include "mpif.h"
142       include "COMMON.IOUNITS"
143       include "COMMON.MPI"
144       include "COMMON.FREE"
145       integer n,chunk,i,j,ii,remainder
146       integer kolor,key,ierror,errcode
147       logical lprint
148       lprint=.true.
149
150 C Split the communicator if independent runs for different parameter
151 C sets will be performed.
152 C
153       if (nparmset.eq.1 .or. .not.separate_parset) then
154         WHAM_COMM = MPI_COMM_WORLD
155       else if (separate_parset) then
156         if (nprocs.lt.nparmset) then
157           write (iout,*) 
158      & "*** Cannot split parameter sets for fewer processors than sets",
159      &  nprocs,nparmset
160           call MPI_Finalize(ierror)
161           stop
162         endif 
163         write (iout,*) "nparmset",nparmset
164         nprocs = nprocs/nparmset
165         kolor = me/nprocs
166         key = mod(me,nprocs)
167         write (iout,*) "My old rank",me," kolor",kolor," key",key
168         call MPI_Comm_split(MPI_COMM_WORLD,kolor,key,WHAM_COMM,ierror)
169         call MPI_Comm_size(WHAM_COMM,nprocs,ierror)
170         call MPI_Comm_rank(WHAM_COMM,me,ierror)
171         write (iout,*) "My new rank",me," comm size",nprocs
172         write (iout,*) "MPI_COMM_WORLD",MPI_COMM_WORLD,
173      &   " WHAM_COMM",WHAM_COMM
174         myparm=kolor+1
175         write (iout,*) "My parameter set is",myparm
176         call flush(iout)
177       else
178         myparm=nparmset
179       endif
180       Me1 = Me
181       Nprocs1 = Nprocs
182       return
183       end
184 c------------------------------------------------------------------------------
185       subroutine work_partition(islice,lprint)
186 c Split the conformations between processors
187       implicit none
188       include "DIMENSIONS"
189       include "DIMENSIONS.ZSCOPT"
190       include "DIMENSIONS.FREE"
191       include "mpif.h"
192       include "COMMON.IOUNITS"
193       include "COMMON.MPI"
194       include "COMMON.PROT"
195       integer islice
196       integer n,chunk,i,j,ii,remainder
197       integer kolor,key,ierror,errcode
198       logical lprint
199 C
200 C Divide conformations between processors; the first and
201 C the last conformation to handle by ith processor is stored in 
202 C indstart(i) and indend(i), respectively.
203 C
204 C First try to assign equal number of conformations to each processor.
205 C
206         n=ntot(islice)
207         write (iout,*) "n=",n
208         indstart(0)=1
209         chunk = N/nprocs1
210         scount(0) = chunk
211 c        print *,"i",0," indstart",indstart(0)," scount",
212 c     &     scount(0)
213         do i=1,nprocs1-1
214           indstart(i)=chunk+indstart(i-1) 
215           scount(i)=scount(i-1)
216 c          print *,"i",i," indstart",indstart(i)," scount",
217 c     &     scount(i)
218         enddo 
219 C
220 C Determine how many conformations remained yet unassigned.
221 C
222         remainder=N-(indstart(nprocs1-1)
223      &    +scount(nprocs1-1)-1)
224 c        print *,"remainder",remainder
225 C
226 C Assign the remainder conformations to consecutive processors, starting
227 C from the lowest rank; this continues until the list is exhausted.
228 C
229         if (remainder .gt. 0) then 
230           do i=1,remainder
231             scount(i-1) = scount(i-1) + 1
232             indstart(i) = indstart(i) + i
233           enddo
234           do i=remainder+1,nprocs1-1
235             indstart(i) = indstart(i) + remainder
236           enddo
237         endif
238
239         indstart(nprocs1)=N+1
240         scount(nprocs1)=0
241
242         do i=0,NProcs1
243           indend(i)=indstart(i)+scount(i)-1
244           idispl(i)=indstart(i)-1
245         enddo
246
247         N=0
248         do i=0,Nprocs1-1
249           N=N+indend(i)-indstart(i)+1
250         enddo
251
252 c        print *,"N",n," NTOT",ntot(islice)
253         if (N.ne.ntot(islice)) then
254           write (iout,*) "!!! Checksum error on processor",me,
255      &     " slice",islice
256           call flush(iout)
257           call MPI_Abort( MPI_COMM_WORLD, Ierror, Errcode )
258         endif
259
260       if (lprint) then
261         write (iout,*) "Partition of work between processors"
262           do i=0,nprocs1-1
263             write (iout,'(a,i5,a,i7,a,i7,a,i7)')
264      &        "Processor",i," indstart",indstart(i),
265      &        " indend",indend(i)," count",scount(i)
266           enddo
267       endif
268       return
269       end
270 #endif
271 #ifdef AIX
272       subroutine flush(iu)
273       call flush_(iu)
274       return
275       end
276 #endif