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