modified: energy_p_new_barrier.F - closed ascii art frames
[unres.git] / source / unres / src_MD-M / together.F
1       Subroutine together
2 c  feeds tasks for parallel processing
3       implicit real*8 (a-h,o-z)
4       include 'DIMENSIONS'
5       include 'mpif.h'
6       real ran1,ran2
7       include 'COMMON.CSA'
8       include 'COMMON.BANK'
9       include 'COMMON.IOUNITS'
10       include 'COMMON.CHAIN'
11       include 'COMMON.TIME1'
12       include 'COMMON.SETUP'
13       include 'COMMON.VAR'
14       include 'COMMON.GEO'
15       include 'COMMON.CONTROL'
16       include 'COMMON.SBRIDGE'
17       real tcpu
18       double precision time_start,time_start_c,time0f,time0i
19       logical ovrtim,sync_iter,timeout,flag,timeout1
20       dimension muster(mpi_status_size)
21       dimension t100(0:100),indx(mxio)
22       dimension xout(maxvar),eout(mxch*(mxch+1)/2+1),ind(9)
23       dimension cout(2)
24       parameter (rad=1.745329252d-2)
25
26 cccccccccccccccccccccccccccccccccccccccccccccccc
27       IF (ME.EQ.KING) THEN
28
29        time0f=MPI_WTIME()
30        ilastnstep=1
31        sync_iter=.false.
32        numch=1
33        nrmsdb=0
34        nrmsdb1=0
35        rmsdbc1c=rmsdbc1
36        nstep=0
37        call csa_read
38        call make_array
39
40        if(iref.ne.0) call from_int(1,0,idum)
41
42 c To minimize input conformation (bank conformation)
43 c Output to $mol.reminimized
44        if (irestart.lt.0) then
45         call read_bank(0,nft,cutdifr)
46         if (irestart.lt.-10) then
47          p_cut=nres*4.d0
48          call prune_bank(p_cut)
49          return
50         endif
51         call reminimize(jlee)
52         return
53        endif
54
55        if (irestart.eq.0) then
56         call initial_write
57         nbank=nconf
58         ntbank=nconf
59         if (ntbankm.eq.0) ntbank=0
60         nstep=0
61         nft=0
62         do i=1,mxio
63          ibank(i)=0
64          jbank(i)=0
65         enddo
66        else
67         call restart_write
68 c!bankt        call read_bankt(jlee,nft,cutdifr)
69         call read_bank(jlee,nft,cutdifr)
70         call read_rbank(jlee,adif)
71         if(iref.ne.0) call from_int(1,0,idum)
72        endif
73
74        nstmax=nstmax+nstep
75        ntrial=n1+n2+n3+n4+n5+n6+n7+n8
76        ntry=ntrial+1
77        ntry=ntry*nseed
78
79 c ntrial : number of trial conformations per seed.
80 c ntry : total number of trial conformations including seed conformations.
81
82        idum2=-123
83 c       imax=2**31-1
84        imax=huge(0)
85        ENDIF
86
87        call mpi_bcast(jend,1,mpi_integer,0,CG_COMM,ierr)
88 cccccccccccccccccccccccccccccccccccccccc
89        do 300 jlee=1,jend
90 cccccccccccccccccccccccccccccccccccccccc
91   331   continue 
92         IF (ME.EQ.KING) THEN
93         if(sync_iter) goto 333
94         idum=-  ran2(idum2)*imax
95         if(jlee.lt.jstart) goto 300
96
97 C Restart the random number generator for conformation generation
98
99         if(irestart.gt.0) then
100          idum2=idum2+nstep
101          if(idum2.le.0) idum2=-idum2+1
102          idum=-  ran2(idum2)*imax
103         endif
104
105         idumm=idum
106         call vrndst(idumm)
107
108         open(icsa_seed,file=csa_seed,status="old")
109         write(icsa_seed,*) "jlee : ",jlee
110         close(icsa_seed)
111
112       call history_append
113       write(icsa_history,*) "number of procs is ",nodes
114       write(icsa_history,*) jlee,idum,idum2
115       close(icsa_history)
116
117 cccccccccccccccccccccccccccccccccccccccccccccccc
118   333 icycle=0
119
120        call history_append
121         write(icsa_history,*) "nbank is ",nbank
122        close(icsa_history)
123
124       if(irestart.eq.1) goto 111
125       if(irestart.eq.2) then
126        icycle=0
127        do i=1,nbank
128         ibank(i)=1
129        enddo
130        do i=nbank+1,nbank+nconf
131         ibank(i)=0
132        enddo
133       endif
134
135 c  start energy minimization
136       nconfr=max0(nconf+nadd,nodes-1)
137       if (sync_iter) nconf_in=0
138 c  king-emperor - feed input and sort output
139        write (iout,*) "NCONF_IN",nconf_in
140        m=0
141        if (nconf_in.gt.0) then
142 c al 7/2/00 - added possibility to read in some of the initial conformations
143         do m=1,nconf_in
144           read (intin,'(i5)',end=11,err=12) iconf
145    12     continue
146           write (iout,*) "write READ_ANGLES",iconf,m
147           call read_angles(intin,*11)
148           if (iref.eq.0) then
149             mm=m
150           else
151             mm=m+1
152           endif
153           do j=2,nres-1
154             dihang_in(1,j,1,mm)=theta(j+1)
155             dihang_in(2,j,1,mm)=phi(j+2)
156             dihang_in(3,j,1,mm)=alph(j)
157             dihang_in(4,j,1,mm)=omeg(j)
158           enddo
159         enddo ! m
160         goto 13
161    11   write (iout,*) nconf_in," conformations requested, but only",
162      &   m-1," found in the angle file."
163         nconf_in=m-1
164    13   continue
165         m=nconf_in
166         write (iout,*) nconf_in,
167      &    " initial conformations have been read in."
168        endif
169        if (iref.eq.0) then
170         if (nconfr.gt.nconf_in) then
171           call make_ranvar(nconfr,m,idum)
172           write (iout,*) nconfr-nconf_in,
173      &     " conformations have been generated randomly."
174         endif
175        else
176         nconfr=nconfr*2
177         call from_int(nconfr,m,idum)
178 c       call from_pdb(nconfr,idum)
179        endif
180        write (iout,*) 'Exitted from make_ranvar nconfr=',nconfr
181        write (*,*) 'Exitted from make_ranvar nconfr=',nconfr
182        do m=1,nconfr
183           write (iout,*) 'Initial conformation',m
184           write(iout,'(8f10.4)') (rad2deg*dihang_in(1,j,1,m),j=2,nres-1)
185           write(iout,'(8f10.4)') (rad2deg*dihang_in(2,j,1,m),j=2,nres-1)
186           write(iout,'(8f10.4)') (rad2deg*dihang_in(3,j,1,m),j=2,nres-1)
187           write(iout,'(8f10.4)') (rad2deg*dihang_in(4,j,1,m),j=2,nres-1)
188        enddo 
189        write(iout,*)'Calling FEEDIN NCONF',nconfr
190        time1i=MPI_WTIME()
191        call feedin(nconfr,nft)
192        write (iout,*) ' Time for first bank min.',MPI_WTIME()-time1i
193        call  history_append
194         write(icsa_history,*) jlee,nft,nbank
195         write(icsa_history,851) (etot(i),i=1,nconfr)
196         write(icsa_history,850) (rmsn(i),i=1,nconfr)
197         write(icsa_history,850) (pncn(i),i=1,nconfr)
198         write(icsa_history,*)
199        close(icsa_history)
200       ELSE
201 c To minimize input conformation (bank conformation)
202 c Output to $mol.reminimized   
203        if (irestart.lt.0) then 
204         call reminimize(jlee)
205         return
206        endif
207        if (irestart.eq.1) goto 111
208 c  soldier - perform energy minimization
209  334   call minim_jlee
210       ENDIF
211
212 ccccccccccccccccccccccccccccccccccc
213 c need to syncronize all procs
214       call mpi_barrier(CG_COMM,ierr)
215       if (ierr.ne.0) then
216        print *, ' cannot synchronize MPI'
217        stop
218       endif
219 ccccccccccccccccccccccccccccccccccc
220
221       IF (ME.EQ.KING) THEN
222
223 c      print *,"ok after minim"
224       nstep=nstep+nconf
225       if(irestart.eq.2) then
226        nbank=nbank+nconf
227 c      ntbank=ntbank+nconf
228        if(ntbank.gt.ntbankm) ntbank=ntbankm
229       endif
230 c      print *,"ok before indexx"
231       if(iref.eq.0) then
232        call indexx(nconfr,etot,indx)
233       else
234 c cc/al 7/6/00
235        do k=1,nconfr
236          indx(k)=k
237        enddo
238        call indexx(nconfr-nconf_in,rmsn(nconf_in+1),indx(nconf_in+1))
239        do k=nconf_in+1,nconfr
240          indx(k)=indx(k)+nconf_in
241        enddo
242 c cc/al
243 c       call indexx(nconfr,rmsn,indx)
244       endif
245 c      print *,"ok after indexx"
246       do im=1,nconf
247        m=indx(im)
248        if (m.gt.mxio .or. m.lt.1) then
249          write (iout,*) 'Dimension ERROR in TOGEHER: IM',im,' M',m
250          call mpi_abort(mpi_comm_world,ierror,ierrcode)
251        endif
252        jbank(im+nbank-nconf)=0
253        bene(im+nbank-nconf)=etot(m)
254        rene(im+nbank-nconf)=etot(m)
255 c!bankt       btene(im)=etot(m)
256 c
257        brmsn(im+nbank-nconf)=rmsn(m)
258        bpncn(im+nbank-nconf)=pncn(m)
259        rrmsn(im+nbank-nconf)=rmsn(m)
260        rpncn(im+nbank-nconf)=pncn(m)
261        if (im+nbank-nconf.gt.mxio .or. im+nbank-nconf.lt.1) then
262          write (iout,*) 'Dimension ERROR in TOGEHER: IM',im,
263      &   ' NBANK',nbank,' NCONF',nconf,' IM+NBANK-NCONF',im+nbank-nconf
264          call mpi_abort(mpi_comm_world,ierror,ierrcode)
265        endif
266        do k=1,numch
267         do j=2,nres-1
268          do i=1,4
269           bvar(i,j,k,im+nbank-nconf)=dihang(i,j,k,m)
270           rvar(i,j,k,im+nbank-nconf)=dihang(i,j,k,m)
271 c!bankt          btvar(i,j,k,im)=dihang(i,j,k,m)
272 c
273          enddo
274         enddo
275        enddo
276        if(iref.eq.1) then
277         if(brmsn(im+nbank-nconf).gt.rmscut.or.
278      &     bpncn(im+nbank-nconf).lt.pnccut) ibank(im+nbank-nconf)=9
279        endif
280        if(vdisulf) then
281            bvar_ns(im+nbank-nconf)=ns-2*nss
282            k=0
283            do i=1,ns
284              j=1
285              do while( iss(i).ne.ihpb(j)-nres .and. 
286      &                 iss(i).ne.jhpb(j)-nres .and. j.le.nss)
287               j=j+1 
288              enddo
289              if (j.gt.nss) then            
290                k=k+1
291                bvar_s(k,im+nbank-nconf)=iss(i)
292              endif
293            enddo
294        endif
295        bvar_nss(im+nbank-nconf)=nss
296        do i=1,nss
297            bvar_ss(1,i,im+nbank-nconf)=ihpb(i)
298            bvar_ss(2,i,im+nbank-nconf)=jhpb(i)
299        enddo
300       enddo
301       ENDIF
302
303   111 continue
304
305       IF (ME.EQ.KING) THEN
306
307       call find_max
308       call find_min
309
310       call get_diff
311       if(nbank.eq.nconf.and.irestart.eq.0) then
312        adif=avedif
313       endif
314
315       cutdif=adif/cut1
316       ctdif1=adif/cut2
317
318 cd      print *,"adif,xctdif,cutdifr"
319 cd      print *,adif,xctdif,cutdifr
320        nst=ntotal/ntrial/nseed
321        xctdif=(cutdif/ctdif1)**(-1.0/nst)
322        if(irestart.ge.1) call estimate_cutdif(adif,xctdif,cutdifr)
323 c       print *,"ok after estimate"
324
325       irestart=0
326
327        call write_rbank(jlee,adif,nft)
328        call write_bank(jlee,nft)
329 c!bankt       call write_bankt(jlee,nft)
330 c       call write_bank1(jlee)
331        call  history_append
332         write(icsa_history,*) "xctdif: ", xctdif,nst,adif/cut1,ctdif1
333         write(icsa_history,851) (bene(i),i=1,nbank)
334         write(icsa_history,850) (brmsn(i),i=1,nbank)
335         write(icsa_history,850) (bpncn(i),i=1,nbank)
336         close(icsa_history)
337   850 format(10f8.3)
338   851 format(5e15.6)
339
340       ifar=nseed/4*3+1
341       ifar=nseed+1
342       ENDIF
343     
344
345       finished=.false.
346       iter = 0
347       irecv = 0
348       isent =0
349       ifrom= 0
350       time0i=MPI_WTIME()
351       time1i=time0i
352       time_start_c=time0i
353       if (.not.sync_iter) then 
354         time_start=time0i
355         nft00=nft
356       else
357         sync_iter=.false.
358       endif
359       nft00_c=nft
360       nft0i=nft
361 ccccccccccccccccccccccccccccccccccccccc
362       do while (.not. finished)
363 ccccccccccccccccccccccccccccccccccccccc
364 crc      print *,"iter ", iter,' isent=',isent
365
366       IF (ME.EQ.KING) THEN
367 c  start energy minimization
368
369        if (isent.eq.0) then
370 c  king-emperor - select seeds & make var & feed input
371 cd        print *,'generating new conf',ntrial,MPI_WTIME()
372         call select_is(nseed,ifar,idum)
373
374         open(icsa_seed,file=csa_seed,status="old")
375         write(icsa_seed,39) 
376      &    jlee,icycle,nstep,(is(i),bene(is(i)),i=1,nseed)
377         close(icsa_seed)
378         call  history_append
379         write(icsa_history,40) jlee,icycle,nstep,cutdif,ibmin,ibmax,
380      *   ebmin,ebmax,nft,iuse,nbank,ntbank
381         close(icsa_history)
382
383          
384
385         call make_var(ntry,idum,iter)
386 cd        print *,'new trial generated',ntrial,MPI_WTIME()
387            time2i=MPI_WTIME()
388            write (iout,'(a20,i4,f12.2)') 
389      &       'Time for make trial',iter+1,time2i-time1i
390        endif
391
392 crc        write(iout,*)'1:Calling FEEDIN NTRY',NTRY,' ntrial',ntrial
393 crc        call feedin(ntry,nft)
394
395        isent=isent+1
396        if (isent.ge.nodes.or.iter.gt.0)  then
397 ct            print *,'waiting ',MPI_WTIME()
398             irecv=irecv+1
399             call recv(0,ifrom,xout,eout,ind,timeout)
400 ct            print *,'   ',irecv,' received from',ifrom,MPI_WTIME()
401        else
402             ifrom=ifrom+1
403        endif
404
405 ct            print *,'sending to',ifrom,MPI_WTIME()
406        call send(isent,ifrom,iter)
407 ct            print *,isent,' sent ',MPI_WTIME()
408
409 c store results -----------------------------------------------
410        if (isent.ge.nodes.or.iter.gt.0)  then
411          nft=nft+ind(3)
412          movernx(irecv)=iabs(ind(5))
413          call getx(ind,xout,eout,cout,rad,iw_pdb,irecv)
414          if(vdisulf) then
415              nss_out(irecv)=nss
416              do i=1,nss
417                iss_out(i,irecv)=ihpb(i)
418                jss_out(i,irecv)=jhpb(i)  
419              enddo
420          endif
421          if(iw_pdb.gt.0) 
422      &          call write_csa_pdb(xout,eout,nft,irecv,iw_pdb)
423        endif
424 c--------------------------------------------------------------
425        if (isent.eq.ntry) then
426            time1i=MPI_WTIME()
427            write (iout,'(a18,f12.2,a14,f10.2)') 
428      &       'Nonsetup time     ',time1i-time_start_c,
429      &       ' sec, Eval/s =',(nft-nft00_c)/(time1i-time_start_c)
430            write (iout,'(a14,i4,f12.2,a14,f10.2)') 
431      &       'Time for iter ',iter+1,time1i-time0i,
432      &       ' sec, Eval/s =',(nft-nft0i)/(time1i-time0i)
433            time0i=time1i
434            nft0i=nft
435            cutdif=cutdif*xctdif
436            if(cutdif.lt.ctdif1) cutdif=ctdif1
437            if (iter.eq.0) then
438               print *,'UPDATING ',ntry-nodes+1,irecv
439               write(iout,*) 'UPDATING ',ntry-nodes+1
440               iter=iter+1
441 c----------------- call update(ntry-nodes+1) -------------------
442               nstep=nstep+ntry-nseed-(nodes-1)
443               call refresh_bank(ntry-nodes+1)
444 c!bankt              call refresh_bankt(ntry-nodes+1)
445            else
446 c----------------- call update(ntry) ---------------------------
447               iter=iter+1
448               print *,'UPDATING ',ntry,irecv
449               write(iout,*) 'UPDATING ',ntry
450               nstep=nstep+ntry-nseed
451               call refresh_bank(ntry)
452 c!bankt              call refresh_bankt(ntry)
453            endif         
454 c----------------------------------------------------------------- 
455
456            call write_bank(jlee,nft)
457 c!bankt           call write_bankt(jlee,nft)
458            call find_min
459
460            time1i=MPI_WTIME()
461            write (iout,'(a20,i4,f12.2)') 
462      &       'Time for refresh ',iter,time1i-time0i
463
464            if(ebmin.lt.estop) finished=.true.
465            if(icycle.gt.icmax) then
466                call write_bank1(jlee)
467                do i=1,nbank
468                  ibank(i)=2
469                  ibank(i)=1
470                enddo
471                nbank=nbank+nconf
472                if(nbank.gt.1000) then 
473                    finished=.true.
474                else
475 crc                   goto 333
476                    sync_iter=.true.
477                endif
478            endif
479            if(nstep.gt.nstmax) finished=.true.
480
481            if(finished.or.sync_iter) then
482             do ij=1,nodes-1
483               call recv(1,ifrom,xout,eout,ind,timeout)
484               if (timeout) then
485                 nstep=nstep+ij-1
486                 print *,'ERROR worker is not responding'
487                 write(iout,*) 'ERROR worker is not responding'
488                 time1i=MPI_WTIME()-time_start_c
489                 print *,'End of cycle, master time for ',iter,' iters ',
490      &             time1i,'sec, Eval/s ',(nft-nft00_c)/time1i
491                 write (iout,*) 'End of cycle, master time for ',iter,
492      &             ' iters ',time1i,' sec'
493                 write (iout,*) 'Total eval/s ',(nft-nft00_c)/time1i
494                 print *,'UPDATING ',ij-1
495                 write(iout,*) 'UPDATING ',ij-1
496                 call flush(iout)
497                 call refresh_bank(ij-1)
498 c!bankt                call refresh_bankt(ij-1)
499                 goto 1002
500               endif
501 c              print *,'node ',ifrom,' finished ',ij,nft
502               write(iout,*) 'node ',ifrom,' finished ',ij,nft
503               call flush(iout)
504               nft=nft+ind(3)
505               movernx(ij)=iabs(ind(5))
506               call getx(ind,xout,eout,cout,rad,iw_pdb,ij)
507               if(vdisulf) then
508                nss_out(ij)=nss
509                do i=1,nss
510                  iss_out(i,ij)=ihpb(i)
511                  jss_out(i,ij)=jhpb(i)  
512                enddo
513               endif
514               if(iw_pdb.gt.0) 
515      &          call write_csa_pdb(xout,eout,nft,ij,iw_pdb)
516             enddo
517             nstep=nstep+nodes-1
518 crc            print *,'---------bcast finished--------',finished
519             time1i=MPI_WTIME()-time_start_c
520             print *,'End of cycle, master time for ',iter,' iters ',
521      &             time1i,'sec, Eval/s ',(nft-nft00_c)/time1i
522             write (iout,*) 'End of cycle, master time for ',iter,
523      &             ' iters ',time1i,' sec'
524             write (iout,*) 'Total eval/s ',(nft-nft00_c)/time1i
525
526 ctimeout            call mpi_bcast(finished,1,mpi_logical,0,CG_COMM,ierr)
527 ctimeout            call mpi_bcast(sync_iter,1,mpi_logical,0,
528 ctimeout     &                                              CG_COMM,ierr)
529             do ij=1,nodes-1 
530                tstart=MPI_WTIME()
531                call mpi_issend(finished,1,mpi_logical,ij,idchar,
532      &             CG_COMM,ireq,ierr)
533                call mpi_issend(sync_iter,1,mpi_logical,ij,idchar,
534      &             CG_COMM,ireq2,ierr)
535                flag=.false.  
536                timeout1=.false.
537                do while(.not. (flag .or. timeout1))
538                  call MPI_TEST(ireq2,flag,muster,ierr)
539                  tend1=MPI_WTIME()
540                  if(tend1-tstart.gt.60) then 
541                   print *,'ERROR worker ',ij,' is not responding'
542                   write(iout,*) 'ERROR worker ',ij,' is not responding'
543                   timeout1=.true.
544                  endif
545                enddo
546                if(timeout1) then
547                 write(iout,*) 'worker ',ij,' NOT OK ',tend1-tstart
548                 timeout=.true.
549                else
550                 write(iout,*) 'worker ',ij,' OK ',tend1-tstart
551                endif
552             enddo
553             print *,'UPDATING ',nodes-1,ij
554             write(iout,*) 'UPDATING ',nodes-1
555             call refresh_bank(nodes-1)
556 c!bankt            call refresh_bankt(nodes-1)
557  1002       continue
558             call write_bank(jlee,nft)
559 c!bankt            call write_bankt(jlee,nft)
560             call find_min
561
562             do i=0,mxmv  
563               do j=1,3
564                 nstatnx_tot(i,j)=nstatnx_tot(i,j)+nstatnx(i,j)
565                 nstatnx(i,j)=0
566               enddo
567             enddo
568
569             write(iout,*)'### Total stats:'
570             do i=0,mxmv
571              if(nstatnx_tot(i,1).ne.0) then
572               if (i.le.9) then
573               write(iout,'(a5,i1,a7,i6,a7,i4,a5,i4,a5,f5.1)') 
574      &        '### N',i,' total=',nstatnx_tot(i,1),
575      &      ' close=',nstatnx_tot(i,2),' far=',nstatnx_tot(i,3),'%acc',
576      &       (nstatnx_tot(i,2)+nstatnx_tot(i,3))*100.0/nstatnx_tot(i,1)
577               else
578               write(iout,'(a4,i2,a7,i6,a7,i4,a5,i4,a5,f5.1)') 
579      &        '###N',i,' total=',nstatnx_tot(i,1),
580      &      ' close=',nstatnx_tot(i,2),' far=',nstatnx_tot(i,3),'%acc',
581      &       (nstatnx_tot(i,2)+nstatnx_tot(i,3))*100.0/nstatnx_tot(i,1)
582               endif
583              else
584               if (i.le.9) then
585               write(iout,'(a5,i1,a7,i6,a7,i4,a5,i4,a5,f5.1)') 
586      &          '### N',i,' total=',nstatnx_tot(i,1),
587      &          ' close=',nstatnx_tot(i,2),' far=',nstatnx_tot(i,3),
588      &          ' %acc',0.0
589               else
590               write(iout,'(a4,i2,a7,i6,a7,i4,a5,i4,a5,f5.1)') 
591      &          '###N',i,' total=',nstatnx_tot(i,1),
592      &          ' close=',nstatnx_tot(i,2),' far=',nstatnx_tot(i,3),
593      &          ' %acc',0.0
594               endif
595              endif
596             enddo
597
598            endif
599            if(sync_iter) goto 331
600
601    39 format(2i3,i7,5(i4,f8.3,1x),19(/,13x,5(i4,f8.3,1x)))
602    40 format(2i2,i8,f8.1,2i4,2(1pe14.5),i10,3i4)
603    43 format(10i8)
604    44 format('jlee =',i3,':',4f10.1,' E =',f8.3,i7,i10)
605
606            isent=0
607            irecv=0
608        endif
609       ELSE
610 c  soldier - perform energy minimization
611         call minim_jlee
612         print *,'End of minim, proc',me,'time ',MPI_WTIME()-time_start
613         write (iout,*) 'End of minim, proc',me,'time ',
614      &                  MPI_WTIME()-time_start
615         call flush(iout)
616 ctimeout        call mpi_bcast(finished,1,mpi_logical,0,CG_COMM,ierr)
617 ctimeout        call mpi_bcast(sync_iter,1,mpi_logical,0,CG_COMM,ierr)
618          call mpi_recv(finished,1,mpi_logical,0,idchar,
619      *                 CG_COMM,muster,ierr)             
620          call mpi_recv(sync_iter,1,mpi_logical,0,idchar,
621      *                 CG_COMM,muster,ierr)             
622         if(sync_iter) goto 331
623       ENDIF
624
625 ccccccccccccccccccccccccccccccccccccccc
626       enddo
627 ccccccccccccccccccccccccccccccccccccccc
628
629       IF (ME.EQ.KING) THEN
630         call  history_append
631         write(icsa_history,40) jlee,icycle,nstep,cutdif,ibmin,ibmax,
632      *  ebmin,ebmax,nft,iuse,nbank,ntbank
633
634         write(icsa_history,44) jlee,0.0,0.0,0.0,
635      &   0.0,ebmin,nstep,nft
636         write(icsa_history,*)
637        close(icsa_history)
638
639        time1i=MPI_WTIME()-time_start
640        print *,'End of RUN, master time ',
641      &             time1i,'sec, Eval/s ',(nft-nft00)/time1i
642        write (iout,*) 'End of RUN, master time  ',
643      &             time1i,' sec'
644        write (iout,*) 'Total eval/s ',(nft-nft00)/time1i
645
646        if(timeout) then 
647         write(iout,*) '!!!! ERROR worker was not responding'
648         write(iout,*) '!!!! cannot finish work normally'
649         write(iout,*) 'Processor0 is calling MPI_ABORT'
650         print *,'!!!! ERROR worker was not responding'
651         print *,'!!!! cannot finish work normally'
652         print *,'Processor0 is calling MPI_ABORT'
653         call flush(iout)
654         call mpi_abort(mpi_comm_world, 111, ierr)
655        endif
656       ENDIF
657
658 cccccccccccccccccccccccccccccc
659   300 continue
660 cccccccccccccccccccccccccccccc
661
662       return
663       end
664 c-------------------------------------------------
665       subroutine feedin(nconf,nft)
666 c  sends out starting conformations and receives results of energy minimization
667       implicit real*8 (a-h,o-z)
668       include 'DIMENSIONS'
669       include 'COMMON.VAR'
670       include 'COMMON.IOUNITS'
671       include 'COMMON.CONTROL'
672       include 'mpif.h'
673       dimension xin(maxvar),xout(maxvar),eout(mxch*(mxch+1)/2+1),
674      *          cout(2),ind(9),info(12)
675       dimension muster(mpi_status_size)
676       include 'COMMON.SETUP'
677       parameter (rad=1.745329252d-2)
678
679       print *,'FEEDIN: NCONF=',nconf
680       mm=0
681 cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
682       if (nconf .lt. nodes-1) then
683         write (*,*) 'FATAL ERROR in FEEDIN, nconf < nodes -1',
684      &   nconf,nodes-1 
685         write (iout,*) 'FATAL ERROR in FEEDIN, nconf < nodes -1',
686      &   nconf,nodes-1 
687         call mpi_abort(mpi_comm_world,ierror,ierrcode)
688       endif
689       do n=1,nconf
690 c  pull out external and internal variables for next start
691         call putx(xin,n,rad)
692 !        write (iout,*) 'XIN from FEEDIN N=',n
693 !        write(iout,'(8f10.4)') (xin(j),j=1,nvar)
694         mm=mm+1
695         if (mm.lt.nodes) then
696 c  feed task to soldier
697 !       print *, ' sending input for start # ',n
698          info(1)=n
699          info(2)=-1
700          info(3)=0
701          info(4)=0
702          info(5)=0
703          info(6)=0
704          call mpi_send(info,12,mpi_integer,mm,idint,CG_COMM,
705      *                  ierr)
706          call mpi_send(xin,nvar,mpi_double_precision,mm,
707      *                  idreal,CG_COMM,ierr)
708         else
709 c  find an available soldier
710          call mpi_recv(ind,9,mpi_integer,mpi_any_source,idint,
711      *                 CG_COMM,muster,ierr)
712 !        print *, ' receiving output from start # ',ind(1)
713          man=muster(mpi_source)
714 c  receive final energies and variables
715          nft=nft+ind(3)
716          call mpi_recv(eout,1,mpi_double_precision,
717      *                 man,idreal,CG_COMM,muster,ierr)
718 !         print *,eout 
719 #ifdef CO_BIAS
720          call mpi_recv(co,1,mpi_double_precision,
721      *                 man,idreal,CG_COMM,muster,ierr)
722          write (iout,'(a15,f3.2,$)') ' BIAS by contact order*100 ',co
723 #endif
724          call mpi_recv(xout,nvar,mpi_double_precision,
725      *                 man,idreal,CG_COMM,muster,ierr)
726 !         print *,nvar , ierr
727 c  feed next task to soldier
728 !        print *, ' sending input for start # ',n
729          info(1)=n
730          info(2)=-1
731          info(3)=0
732          info(4)=0
733          info(5)=0
734          info(6)=0
735          info(7)=0
736          info(8)=0
737          info(9)=0
738          call mpi_send(info,12,mpi_integer,man,idint,CG_COMM,
739      *                  ierr)
740          call mpi_send(xin,nvar,mpi_double_precision,man,
741      *                  idreal,CG_COMM,ierr)
742 c  retrieve latest results
743          call getx(ind,xout,eout,cout,rad,iw_pdb,ind(1))
744          if(iw_pdb.gt.0) 
745      &        call write_csa_pdb(xout,eout,nft,ind(1),iw_pdb)
746         endif
747       enddo
748 cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
749 c  no more input
750 cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
751       do j=1,nodes-1
752 c  wait for a soldier
753        call mpi_recv(ind,9,mpi_integer,mpi_any_source,idint,
754      *               CG_COMM,muster,ierr)
755 crc       if (ierr.ne.0) go to 30
756 !      print *, ' receiving output from start # ',ind(1)
757        man=muster(mpi_source)
758 c  receive final energies and variables
759        nft=nft+ind(3)
760        call mpi_recv(eout,1,
761      *               mpi_double_precision,man,idreal,
762      *               CG_COMM,muster,ierr)
763 !       print *,eout
764 #ifdef CO_BIAS
765          call mpi_recv(co,1,mpi_double_precision,
766      *                 man,idreal,CG_COMM,muster,ierr)
767          write (iout,'(a15,f3.2,$)') ' BIAS by contact order*100 ',co
768 #endif
769 crc       if (ierr.ne.0) go to 30
770        call mpi_recv(xout,nvar,mpi_double_precision,
771      *               man,idreal,CG_COMM,muster,ierr)
772 !       print *,nvar , ierr
773 crc       if (ierr.ne.0) go to 30
774 c  halt soldier
775        info(1)=0
776        info(2)=-1
777        info(3)=0 
778        info(4)=0
779        info(5)=0
780        info(6)=0
781        call mpi_send(info,12,mpi_integer,man,idint,CG_COMM,
782      *                ierr)
783 c  retrieve results
784        call getx(ind,xout,eout,cout,rad,iw_pdb,ind(1))
785        if(iw_pdb.gt.0) 
786      &          call write_csa_pdb(xout,eout,nft,ind(1),iw_pdb)
787       enddo
788 cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
789       return
790    10 print *, ' dispatching error'
791       call mpi_abort(mpi_comm_world,ierror,ierrcode)
792       return
793    20 print *, ' communication error'
794       call mpi_abort(mpi_comm_world,ierror,ierrcode)
795       return
796    30 print *, ' receiving error'
797       call mpi_abort(mpi_comm_world,ierror,ierrcode)
798       return
799       end
800 cccccccccccccccccccccccccccccccccccccccccccccccccc
801       subroutine getx(ind,xout,eout,cout,rad,iw_pdb,k)
802 c  receives and stores data from soldiers
803       implicit real*8 (a-h,o-z)
804       include 'DIMENSIONS'
805       include 'COMMON.IOUNITS'
806       include 'COMMON.CSA'
807       include 'COMMON.BANK'
808       include 'COMMON.VAR'
809       include 'COMMON.CHAIN'
810       include 'COMMON.CONTACTS'
811       dimension ind(9),xout(maxvar),eout(mxch*(mxch+1)/2+1)
812 cjlee
813       double precision przes(3),obr(3,3)
814       logical non_conv
815 cjlee
816       iw_pdb=2
817       if (k.gt.mxio .or. k.lt.1) then 
818         write (iout,*) 
819      &   'ERROR - dimensions of ANGMIN have been exceeded K=',k
820         call mpi_abort(mpi_comm_world,ierror,ierrcode)
821       endif
822 c  store ind()
823       do j=1,9
824        indb(k,j)=ind(j)
825       enddo
826 c  store energies
827       etot(k)=eout(1)
828 c  retrieve dihedral angles etc
829       call var_to_geom(nvar,xout)
830       do j=2,nres-1
831         dihang(1,j,1,k)=theta(j+1)
832         dihang(2,j,1,k)=phi(j+2)
833         dihang(3,j,1,k)=alph(j)
834         dihang(4,j,1,k)=omeg(j)
835       enddo
836       dihang(2,nres-1,1,k)=0.0d0
837 cjlee
838       if(iref.eq.0) then 
839        iw_pdb=1
840 cd       write(iout,'(i3,a3,i4,i5,a6,1pe12.4,a4,i3,i4)') 
841 cd     &      ind(2),' e ',ind(3),ind(1),' etot ',etot(k),' mv ',
842 cd     &      ind(5),ind(4)
843        return
844       endif
845       call chainbuild
846 c     call dihang_to_c(dihang(1,1,1,k))
847 c     call fitsq(rms,c(1,1),crefjlee(1,1),nres,przes,obr,non_conv)
848 c     call fitsq(rms,c(1,2),crefjlee(1,2),nres-1,przes,obr,non_conv)
849 c           call fitsq(rms,c(1,nstart_seq),crefjlee(1,nstart_sup),
850 c    &                 nsup,przes,obr,non_conv)
851 c      rmsn(k)=dsqrt(rms)
852
853        call rmsd_csa(rmsn(k))
854        call contact(.false.,ncont,icont,co)
855        pncn(k)=contact_fract(ncont,ncont_ref,icont,icont_ref)     
856
857 cd       write(iout,'(i3,a3,i4,i5,a6,1pe12.4,a5
858 cd     &      ,0pf5.2,a5,f5.1,a,f6.3,a4,i3,i4)') 
859 cd     &    ind(2),' e ',ind(3),ind(1),' etot ',etot(k),' rms ',
860 cd     &    rmsn(k),' %NC ',pncn(k)*100,' cont.order',co,' mv ',
861 cd     &    ind(5),ind(4)
862
863  
864       if (rmsn(k).gt.rmscut.or.pncn(k).lt.pnccut) iw_pdb=0
865       return
866       end
867 cccccccccccccccccccccccccccccccccccccccccccccccccc
868       subroutine putx(xin,n,rad)
869 c  gets starting variables
870       implicit real*8 (a-h,o-z)
871       include 'DIMENSIONS'
872       include 'COMMON.CSA'
873       include 'COMMON.BANK'
874       include 'COMMON.VAR'
875       include 'COMMON.CHAIN'
876       include 'COMMON.IOUNITS'
877       dimension xin(maxvar)
878
879 c  pull out starting values for variables
880 !       write (iout,*)'PUTX: N=',n
881       do m=1,numch
882 !        write (iout,'(8f10.4)') (dihang_in(1,j,m,n),j=2,nres-1)
883 !        write (iout,'(8f10.4)') (dihang_in(2,j,m,n),j=2,nres-1)
884 !        write (iout,'(8f10.4)') (dihang_in(3,j,m,n),j=2,nres-1)
885 !        write (iout,'(8f10.4)') (dihang_in(4,j,m,n),j=2,nres-1)
886        do j=2,nres-1
887         theta(j+1)=dihang_in(1,j,m,n)
888         phi(j+2)=dihang_in(2,j,m,n)
889         alph(j)=dihang_in(3,j,m,n)
890         omeg(j)=dihang_in(4,j,m,n)
891        enddo
892       enddo
893 c  set up array of variables
894       call geom_to_var(nvar,xin)
895 !       write (iout,*) 'xin in PUTX N=',n 
896 !       call intout
897 !       write (iout,'(8f10.4)') (xin(i),i=1,nvar) 
898       return
899       end
900 c--------------------------------------------------------
901       subroutine putx2(xin,iff,n)
902 c  gets starting variables
903       implicit real*8 (a-h,o-z)
904       include 'DIMENSIONS'
905       include 'COMMON.CSA'
906       include 'COMMON.BANK'
907       include 'COMMON.VAR'
908       include 'COMMON.CHAIN'
909       include 'COMMON.IOUNITS'
910       dimension xin(maxvar),iff(maxres)
911
912 c  pull out starting values for variables
913       do m=1,numch
914        do j=2,nres-1
915         theta(j+1)=dihang_in2(1,j,m,n)
916         phi(j+2)=dihang_in2(2,j,m,n)
917         alph(j)=dihang_in2(3,j,m,n)
918         omeg(j)=dihang_in2(4,j,m,n)
919        enddo
920       enddo
921 c  set up array of variables
922       call geom_to_var(nvar,xin)
923        
924       do i=1,nres
925         iff(i)=iff_in(i,n)
926       enddo
927       return
928       end
929
930 c-------------------------------------------------------
931       subroutine prune_bank(p_cut)
932       implicit real*8 (a-h,o-z)
933       include 'DIMENSIONS'
934       include 'mpif.h'
935       include 'COMMON.CSA'
936       include 'COMMON.BANK'
937       include 'COMMON.IOUNITS'
938       include 'COMMON.CHAIN'
939       include 'COMMON.TIME1'
940       include 'COMMON.SETUP'
941 c---------------------------
942 c This subroutine prunes bank conformations using p_cut
943 c---------------------------
944
945       nprune=0
946       nprune=nprune+1
947       m=1
948       do k=1,numch
949        do j=2,nres-1
950         do i=1,4
951          dihang(i,j,k,nprune)=bvar(i,j,k,m)
952         enddo
953        enddo
954       enddo
955       bene(nprune)=bene(m)
956       brmsn(nprune)=brmsn(m)
957       bpncn(nprune)=bpncn(m) 
958
959       do m=2,nbank
960        ddmin=9.d190
961        do ip=1,nprune
962         call get_diff12(dihang(1,1,1,ip),bvar(1,1,1,m),diff) 
963         if(diff.lt.p_cut) goto 100
964         if(diff.lt.ddmin) ddmin=diff
965        enddo
966        nprune=nprune+1
967        do k=1,numch
968         do j=2,nres-1
969          do i=1,4
970           dihang(i,j,k,nprune)=bvar(i,j,k,m)
971          enddo
972         enddo
973        enddo
974        bene(nprune)=bene(m)
975        brmsn(nprune)=brmsn(m)
976        bpncn(nprune)=bpncn(m)
977   100  continue
978        write (iout,*) 'Pruning :',m,nprune,p_cut,ddmin
979       enddo
980       nbank=nprune
981       print *, 'Pruning :',m,nprune,p_cut
982       call write_bank(0,0)
983
984       return
985       end
986 c-------------------------------------------------------
987
988       subroutine reminimize(jlee)
989       implicit real*8 (a-h,o-z)
990       include 'DIMENSIONS'
991       include 'mpif.h'
992       include 'COMMON.CSA'
993       include 'COMMON.BANK'
994       include 'COMMON.IOUNITS'
995       include 'COMMON.CHAIN'
996       include 'COMMON.TIME1'
997       include 'COMMON.SETUP'
998 c---------------------------
999 c This subroutine re-minimizes bank conformations:
1000 c---------------------------
1001
1002        ntry=nbank
1003
1004        call find_max
1005        call find_min
1006
1007        if (me.eq.king) then
1008         open(icsa_history,file=csa_history,status="old")
1009          write(icsa_history,*) "Re-minimization",nodes,"nodes"
1010          write(icsa_history,851) (bene(i),i=1,nbank)
1011          write(icsa_history,40) jlee,icycle,nstep,cutdif,ibmin,ibmax,
1012      *   ebmin,ebmax,nft,iuse,nbank,ntbank
1013         close(icsa_history)
1014         do index=1,ntry
1015          do k=1,numch
1016           do j=2,nres-1
1017            do i=1,4
1018             dihang_in(i,j,k,index)=bvar(i,j,k,index)
1019            enddo
1020           enddo
1021          enddo
1022         enddo
1023         nft=0
1024         call feedin(ntry,nft)
1025        else
1026         call minim_jlee
1027        endif
1028
1029        call find_max
1030        call find_min
1031
1032        if (me.eq.king) then
1033         do i=1,ntry
1034          call replace_bvar(i,i)
1035         enddo
1036         open(icsa_history,file=csa_history,status="old")
1037          write(icsa_history,40) jlee,icycle,nstep,cutdif,ibmin,ibmax,
1038      *   ebmin,ebmax,nft,iuse,nbank,ntbank
1039          write(icsa_history,851) (bene(i),i=1,nbank)
1040         close(icsa_history)
1041         call write_bank_reminimized(jlee,nft)
1042        endif
1043
1044    40 format(2i2,i8,f8.1,2i4,2(1pe14.5),i10,3i4)
1045   851 format(5e15.6)
1046   850 format(5e15.10)
1047 c  850 format(10f8.3)
1048
1049       return
1050       end
1051 c-------------------------------------------------------
1052       subroutine send(n,mm,it)
1053 c  sends out starting conformation for minimization
1054       implicit real*8 (a-h,o-z)
1055       include 'DIMENSIONS'
1056       include 'COMMON.VAR'
1057       include 'COMMON.IOUNITS'
1058       include 'COMMON.CONTROL'
1059       include 'COMMON.BANK'
1060       include 'COMMON.CHAIN'
1061       include 'mpif.h'
1062       dimension xin(maxvar),xout(maxvar),eout(mxch*(mxch+1)/2+1),
1063      *          cout(2),ind(8),xin2(maxvar),iff(maxres),info(12)
1064       dimension muster(mpi_status_size)
1065       include 'COMMON.SETUP'
1066       parameter (rad=1.745329252d-2)
1067
1068       if (isend2(n).eq.0) then
1069 c  pull out external and internal variables for next start
1070         call putx(xin,n,rad)
1071         info(1)=n
1072         info(2)=it
1073         info(3)=movenx(n)
1074         info(4)=nss_in(n)
1075         info(5)=parent(1,n)
1076         info(6)=parent(2,n)
1077
1078         if (movenx(n).eq.14.or.movenx(n).eq.17) then
1079           info(7)=idata(1,n)
1080           info(8)=idata(2,n)
1081         else if (movenx(n).eq.16) then
1082           info(7)=idata(1,n)
1083           info(8)=idata(2,n)
1084           info(10)=idata(3,n)
1085           info(11)=idata(4,n)
1086           info(12)=idata(5,n)
1087         else
1088          info(7)=0
1089          info(8)=0
1090          info(10)=0
1091          info(11)=0
1092          info(12)=0
1093         endif
1094
1095         if (movenx(n).eq.15) then
1096          info(9)=parent(3,n)
1097         else
1098          info(9)=0
1099         endif
1100         call mpi_send(info,12,mpi_integer,mm,idint,CG_COMM,
1101      *                  ierr)
1102         call mpi_send(xin,nvar,mpi_double_precision,mm,
1103      *                  idreal,CG_COMM,ierr)
1104       else
1105 c  distfit & minimization for n7 move
1106         info(1)=-n
1107         info(2)=it
1108         info(3)=movenx(n)
1109         info(4)=nss_in(n)
1110         info(5)=parent(1,n)
1111         info(6)=parent(2,n)
1112         info(7)=0
1113         info(8)=0
1114         info(9)=0
1115         call mpi_send(info,12,mpi_integer,mm,idint,CG_COMM,
1116      *                  ierr)
1117         call putx2(xin,iff,isend2(n))
1118         call mpi_send(xin,nvar,mpi_double_precision,mm,
1119      *                  idreal,CG_COMM,ierr)
1120         call mpi_send(iff,nres,mpi_integer,mm,
1121      *                  idint,CG_COMM,ierr)
1122         call putx(xin2,n,rad)
1123         call mpi_send(xin2,nvar,mpi_double_precision,mm,
1124      *                  idreal,CG_COMM,ierr)
1125       endif 
1126       if (vdisulf.and.nss_in(n).ne.0) then
1127         call mpi_send(iss_in(1,n),nss_in(n),mpi_integer,mm,
1128      *                  idint,CG_COMM,ierr)
1129         call mpi_send(jss_in(1,n),nss_in(n),mpi_integer,mm,
1130      *                  idint,CG_COMM,ierr)
1131       endif
1132       return
1133       end
1134 c-------------------------------------------------
1135
1136       subroutine recv(ihalt,man,xout,eout,ind,tout)
1137 c  receives results of energy minimization
1138       implicit real*8 (a-h,o-z)
1139       include 'DIMENSIONS'
1140       include 'COMMON.VAR'
1141       include 'COMMON.IOUNITS'
1142       include 'COMMON.CONTROL'
1143       include 'COMMON.SBRIDGE'
1144       include 'COMMON.BANK'
1145       include 'COMMON.CHAIN'
1146       include 'mpif.h'
1147       dimension xin(maxvar),xout(maxvar),eout(mxch*(mxch+1)/2+1),
1148      *          cout(2),ind(9),info(12)
1149       dimension muster(mpi_status_size)
1150       include 'COMMON.SETUP'
1151       logical tout,flag
1152       double precision twait,tstart,tend1
1153       parameter(twait=600.0d0)
1154
1155 c  find an available soldier
1156        tout=.false.
1157        flag=.false.
1158        tstart=MPI_WTIME()
1159        do while(.not. (flag .or. tout))
1160          call MPI_IPROBE(mpi_any_source,idint,CG_COMM,flag, 
1161      *            muster,ierr)
1162          tend1=MPI_WTIME()
1163          if(tend1-tstart.gt.twait .and. ihalt.eq.1) tout=.true.
1164 c_error         if(tend1-tstart.gt.twait) tout=.true.
1165        enddo
1166        if (tout) then 
1167          write(iout,*) 'ERROR = timeout for recv ',tend1-tstart
1168          call flush(iout)
1169          return
1170        endif
1171        man=muster(mpi_source)        
1172
1173 ctimeout         call mpi_recv(ind,9,mpi_integer,mpi_any_source,idint,
1174 ctimeout     *                 CG_COMM,muster,ierr)
1175 !        print *, ' receiving output from start # ',ind(1)
1176 ct         print *,'receiving ',MPI_WTIME()
1177 ctimeout         man=muster(mpi_source)
1178          call mpi_recv(ind,9,mpi_integer,man,idint,
1179      *                 CG_COMM,muster,ierr)
1180 ctimeout
1181 c  receive final energies and variables
1182          call mpi_recv(eout,1,mpi_double_precision,
1183      *                 man,idreal,CG_COMM,muster,ierr)
1184 !         print *,eout 
1185 #ifdef CO_BIAS
1186          call mpi_recv(co,1,mpi_double_precision,
1187      *                 man,idreal,CG_COMM,muster,ierr)
1188          write (iout,'(a15,f3.2,$)') ' BIAS by contact order*100 ',co
1189 #endif
1190          call mpi_recv(xout,nvar,mpi_double_precision,
1191      *                 man,idreal,CG_COMM,muster,ierr)
1192 !         print *,nvar , ierr
1193          if(vdisulf) nss=ind(6)
1194          if(vdisulf.and.nss.ne.0) then
1195           call mpi_recv(ihpb,nss,mpi_integer,
1196      *                 man,idint,CG_COMM,muster,ierr)         
1197           call mpi_recv(jhpb,nss,mpi_integer,
1198      *                 man,idint,CG_COMM,muster,ierr)
1199          endif
1200 c  halt soldier
1201        if(ihalt.eq.1) then 
1202 c        print *,'sending halt to ',man
1203         write(iout,*) 'sending halt to ',man
1204         info(1)=0
1205         call mpi_send(info,12,mpi_integer,man,idint,CG_COMM,ierr)
1206        endif
1207       return
1208       end
1209
1210 c---------------------------------------------------------- 
1211       subroutine history_append 
1212       implicit real*8 (a-h,o-z)
1213       include 'DIMENSIONS'
1214       include 'COMMON.IOUNITS'
1215
1216 #if defined(AIX) || defined(PGI)
1217        open(icsa_history,file=csa_history,position="append")
1218 #else
1219        open(icsa_history,file=csa_history,access="append")
1220 #endif
1221        return
1222        end