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