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