update new files
[unres.git] / source / maxlik / src_FPy / bilans.F
1 c------------------------------------------------------------------------------
2       subroutine bilans
3       implicit none
4       integer ntasks 
5       parameter (ntasks=11)
6       include "DIMENSIONS.ZSCOPT"
7 #ifdef MPI
8       include "mpif.h"
9       include "COMMON.MPI"
10       integer IERROR
11       double precision ttask_all(ntasks)
12       integer nctask_all(ntasks)
13 #endif
14       include "COMMON.IOUNITS"
15       integer i
16       double precision ttask
17       integer nctask
18       common /timing/ ttask(ntasks),nctask(ntasks)
19       character*16 tname(ntasks) /"function","gradient",9*''/
20
21       if (me.eq.Master) then
22       write (iout,*)
23       write (iout,'(80(1h-))')
24 #ifdef MPI
25       write (iout,*) "Routine call info from the processor ",me," ..."
26 #else
27       write (iout,*) "Routine call info ..."
28 #endif
29       write (iout,*)
30       write (iout,'(65(1h-))')
31       write (iout,100) "task","   # calls","     total time",
32      &  "  ave.time/call"
33       write (iout,'(65(1h-))')
34       do i=1,ntasks
35         write (iout,200) tname(i),nctask(i),ttask(i),
36      &  ttask(i)/(nctask(i)+1.0d-10)
37       enddo
38       write (iout,*) 
39       endif
40 #ifdef MPI
41       call MPI_Reduce(ttask(1),ttask_all(1),ntasks, 
42      &  MPI_DOUBLE_PRECISION, MPI_SUM, Master, WHAM_COMM,IERROR)
43       call MPI_Reduce(nctask(1),nctask_all(1),ntasks, 
44      &  MPI_INTEGER, MPI_SUM, Master, WHAM_COMM,IERROR)
45       if (Me.eq.Master) then
46       write (iout,'(80(1h-))')
47       write (iout,*) "Routine call info from all processors ..."
48       write (iout,*)
49       write (iout,'(65(1h-))')
50       write (iout,100) "task","   # calls","     total time",
51      &  "  ave.time/call"
52       write (iout,'(65(1h-))')
53       do i=1,ntasks
54         write (iout,200) tname(i),nctask_all(i),ttask_all(i),
55      &    ttask_all(i)/(nctask_all(i)+1.0d-10)
56       enddo
57       write (iout,*) 
58       endif
59 #endif
60       return
61   100 format(a,t21,a,t31,a,t46,a)
62   200 format(a,t21,i10,f15.2,f15.8)
63       end
64 c------------------------------------------------------------------------------
65 #ifdef MPI
66       subroutine proc_groups
67 C Split the processors into the Master and Workers group, if needed.
68       implicit none
69       include "DIMENSIONS"
70       include "DIMENSIONS.ZSCOPT"
71       include "mpif.h"
72       include "COMMON.IOUNITS"
73       include "COMMON.MPI"
74       include "COMMON.VMCPAR"
75       integer n,chunk,iprot,i,j,ii,remainder
76       integer kolor,key,ierror,errcode
77       logical lprint
78       lprint=.true.
79 C No secondary structure constraints. 
80       Comm1 = WHAM_COMM
81       Me1 = Me
82       Nprocs1 = Nprocs
83       return
84       end
85 c-------------------------------------------------------------------------------
86       subroutine work_partition(lprint)
87 c Split the conformations between processors
88       implicit none
89       include "DIMENSIONS"
90       include "DIMENSIONS.ZSCOPT"
91       include "mpif.h"
92       include "COMMON.CLASSES"
93       include "COMMON.IOUNITS"
94       include "COMMON.MPI"
95       include "COMMON.VMCPAR"
96       integer n,chunk,iprot,i,j,ii,remainder
97       integer kolor,key,ierror,errcode
98       logical lprint
99 C
100 C Divide conformations between processors; for each proteins the first and
101 C the last conformation to handle by ith processor is stored in 
102 C indstart(i,iprot) and indend(i,iprot), respectively.
103 C
104 C First try to assign equal number of conformations to each processor.
105 C
106       do iprot=1,nprot
107         n=ntot_work(iprot)
108         if (me.eq.Master) write (iout,*) "Protein",iprot," n=",n
109         indstart(0,iprot)=1
110         chunk = N/nprocs1
111         scount(0,iprot) = chunk
112 c        print *,"i",0," indstart",indstart(0,iprot)," scount",
113 c     &     scount(0,iprot)
114         do i=1,nprocs1-1
115           indstart(i,iprot)=chunk+indstart(i-1,iprot) 
116           scount(i,iprot)=scount(i-1,iprot)
117 c          print *,"i",i," indstart",indstart(i,iprot)," scount",
118 c     &     scount(i,iprot)
119         enddo 
120 C
121 C Determine how many conformations remained yet unassigned.
122 C
123         remainder=N-(indstart(nprocs1-1,iprot)
124      &    +scount(nprocs1-1,iprot)-1)
125 c        print *,"remainder",remainder
126 C
127 C Assign the remainder conformations to consecutive processors, starting
128 C from the lowest rank; this continues until the list is exhausted.
129 C
130         if (remainder .gt. 0) then 
131           do i=1,remainder
132             scount(i-1,iprot) = scount(i-1,iprot) + 1
133             indstart(i,iprot) = indstart(i,iprot) + i
134           enddo
135           do i=remainder+1,nprocs1-1
136             indstart(i,iprot) = indstart(i,iprot) + remainder
137           enddo
138         endif
139
140         indstart(nprocs1,iprot)=N+1
141         scount(nprocs1,iprot)=0
142
143         do i=0,NProcs1
144           indend(i,iprot)=indstart(i,iprot)+scount(i,iprot)-1
145           idispl(i,iprot)=indstart(i,iprot)-1
146         enddo
147
148         N=0
149         do i=0,Nprocs1-1
150           N=N+indend(i,iprot)-indstart(i,iprot)+1
151         enddo
152
153 c        print *,"N",n," NTOT",ntot_work(iprot)
154         if (N.ne.ntot_work(iprot)) then
155           write (*,*) "!!! Checksum error on processor",me
156 c          call flush(iout)
157           call MPI_Abort( WHAM_COMM, Ierror, Errcode )
158         endif
159
160         ii=0
161         do i=1,ntot_work(iprot)
162           if (i.ge.indstart(me1,iprot) .and. i.le.indend(me1,iprot)) 
163      &    then
164             ii=ii+1
165             i2ii(i,iprot)=ii
166           else
167             i2ii(i,iprot)=0
168           endif
169 c          write (iout,*) "i",i," iprot",iprot," i2ii",i2ii(i,iprot)
170         enddo
171
172       enddo ! iprot
173       if (lprint) then
174         write (iout,*) "Partition of work between processors"
175         do iprot=1,nprot
176           write (iout,*) "Protein",iprot
177 #ifdef DEBUG
178           write (iout,*) "The i2ii array"
179           do j=1,ntot_work(iprot)
180             write (iout,*) j,i2ii(j,iprot)
181           enddo
182 #endif
183           do i=0,nprocs1-1
184             write (iout,'(a,i5,a,i7,a,i7,a,i7)')
185      &        "Processor",i," indstart",indstart(i,iprot),
186      &        " indend",indend(i,iprot)," count",scount(i,iprot)
187           enddo
188         enddo
189       endif
190       return
191       end
192 c------------------------------------------------------------------------------
193       subroutine jebadelko(nvarr)
194       implicit none
195       include "DIMENSIONS"
196       include "DIMENSIONS.ZSCOPT"
197       include "mpif.h"
198       include "COMMON.IOUNITS"
199       include "COMMON.MPI"
200       include "COMMON.VMCPAR"
201       include "COMMON.CLASSES"
202       include "COMMON.OPTIM"
203       include "COMMON.WEIGHTS"
204       include "COMMON.WEIGHTDER"
205       include "COMMON.ENERGIES"
206       include "COMMON.TIME1"
207       include "COMMON.PROTNAME"
208       include "COMMON.PROTFILES"
209       include "COMMON.COMPAR"
210       integer What, TAG, IERROR, status(MPI_STATUS_SIZE), istop, iprot, 
211      &  nvarr, nf, errcode
212       integer count
213       double precision x(max_paropt), g(max_paropt), viol
214       integer uiparm,i,j
215       integer ibatch,ib
216       external fdum
217       double precision rdum
218       double precision tcpu,t1,t1w,t1_ini,t1w_ini
219       logical lprn,op
220       integer iun
221       lprn=.false.
222       if (me.eq.Master) then
223       write(iout,*) "Processor",me,me1," called JEBADELKO"
224       call flush(iout)
225       endif 
226       if (me.eq.Master) then
227         istop=0
228         call func1(nvarr,istop,x)
229       else
230 #ifdef DEBUG
231       write (iout,*) "ELOWEST at slave  starting JEBADELKO"
232       do iprot=1,nprot
233         do ibatch=1,natlike(iprot)+2
234           do ib=1,nbeta(ibatch,iprot)
235             write (iout,*) "iprot",iprot," ibatch",ibatch," ib",ib,
236      &         " elowest",elowest(ib,ibatch,iprot)
237           enddo
238         enddo
239       enddo
240 #endif
241         istop=1
242         t1w_ini = MPI_WTIME()
243         t1_ini = tcpu()
244         do while (istop.ne.0)
245 #ifdef DEBUG
246       write (iout,*) "ELOWEST at slave  calling FUNC1 from JBADELKO"
247       do iprot=1,nprot
248         do ibatch=1,natlike(iprot)+2
249           do ib=1,nbeta(ibatch,iprot)
250             write (iout,*) "iprot",iprot," ibatch",ibatch," ib",ib,
251      &         " elowest",elowest(ib,ibatch,iprot)
252           enddo
253         enddo
254       enddo
255 #endif
256           call func1(nvarr,istop,x)
257         enddo
258         t1w = mpi_wtime() - t1w_ini
259         t1 = tcpu() - t1_ini
260 #ifdef DEBUG
261         write (iout,*)
262         write (iout,*) "CPU time",t1," wall clock time",t1w
263         call flush(iout)
264 #endif
265       endif
266 #ifdef DEBUG
267       write (iout,*)
268       write (iout,*) "Energies of processor",me
269       do j=1,nprot
270         write (iout,*)
271         write (iout,*) "Protein ",protname(iprot)
272         do i=1,ntot_work(j)
273           if (i.ge.indstart(me1,j).and.i.le.indend(me1,j)) then
274             write(iout,*)i,e_total(i,j),rmstb(i,j),iscore(i,0,j)
275           endif
276         enddo
277       enddo
278 #endif
279 #ifndef PYTHON
280 c      do iprot=1,nprot
281 c        write (iout,*) "Deleting scratchfile",bprotfiles(iprot)
282 c        inquire (file=bprotfiles(iprot),number=iun,opened=op)
283 c        write (iout,*) "unit",iun," icbase",icbase
284 cc        if (.not. op) then
285 c          open (icbase,file=bprotfiles(iprot),status="old",
286 c     &    form="unformatted",access="direct",recl=lenrec(iprot))
287 c          close(icbase,status="delete")
288 c        else
289 c          close(iun,status="delete")
290 c        endif
291 c        call flush(iout)
292 c        write (iout,*) "Deleting scratchfile",benefiles(iprot)
293 c        inquire (file=benefiles(iprot),number=iun,opened=op)
294 c        write (iout,*) "unit",iun," ientout",icbase
295 c        if (.not.op) then
296 c          open (ientout,file=benefiles(iprot),status="old",
297 c     &    form="unformatted",access="direct",recl=lenrec_ene(iprot))
298 c          close(ientout,status="delete")
299 c        else
300 c          close (iun,status="delete")
301 c        endif
302 c        call flush(iout)
303 c      enddo
304 #endif
305 #ifdef DEBUG
306       write (iout,*) "Processor",me,"leaves JEBADELKO"
307       call flush(iout)
308 #endif
309       return
310       end
311 #endif
312