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