cf98db7fbf3d279238561660a23ed6324a0e61cb
[unres.git] / source / cluster / wham / src / read_coords.F
1       subroutine read_coords(ncon,*)
2       implicit none
3       include "DIMENSIONS"
4       include "sizesclu.dat"
5 #ifdef MPI
6       include "mpif.h"
7       integer IERROR,ERRCODE,STATUS(MPI_STATUS_SIZE)
8       include "COMMON.MPI"
9 #endif
10       include "COMMON.CONTROL"
11       include "COMMON.CHAIN"
12       include "COMMON.INTERACT"
13       include "COMMON.IOUNITS"
14       include "COMMON.VAR"
15       include "COMMON.SBRIDGE"
16       include "COMMON.GEO"
17       include "COMMON.CLUSTER"
18       character*3 liczba
19       integer ncon
20       integer i,j,jj,jjj,jj_old,icount,k,kk,l,ii,if,ib,
21      &  nn,nn1,inan
22       integer ixdrf,iret,itmp
23       real*4 prec,reini,refree,rmsdev
24       integer nrec,nlines,iscor,lenrec,lenrec_in
25       double precision energ,t_acq,tcpu
26       integer ilen,iroof
27       external ilen,iroof
28       double precision rjunk
29       integer ntot_all(0:maxprocs-1)
30       logical lerr
31       double precision energia(0:max_ene),etot
32       real*4 csingle(3,maxres2+2)
33       integer Previous,Next
34       character*256 bprotfiles
35 c      print *,"Processor",me," calls read_protein_data"
36 #ifdef MPI
37       if (me.eq.master) then
38         Previous=MPI_PROC_NULL
39       else
40         Previous=me-1
41       endif
42       if (me.eq.nprocs-1) then
43         Next=MPI_PROC_NULL
44       else
45         Next=me+1
46       endif
47 c Set the scratchfile names
48       write (liczba,'(bz,i3.3)') me
49 #endif
50 c 1/27/05 AL Change stored coordinates to single precision and don't store 
51 c         energy components in the binary databases.
52       lenrec=12*(nres+nct-nnt+1)+4*(2*nss+2)+16
53       lenrec_in=12*(nres+nct-nnt+1)+4*(2*nss+2)+24
54 #ifdef DEBUG
55       write (iout,*) "nres",nres," nnt",nnt," nct",nct," nss", nss
56       write (iout,*) "lenrec_in",lenrec_in
57 #endif
58       bprotfiles=scratchdir(:ilen(scratchdir))//
59      &       "/"//prefix(:ilen(prefix))//liczba//".xbin"
60
61 #ifdef CHUJ
62       ICON=1
63   123 continue
64       if (from_cart .and. .not. from_bx .and. .not. from_cx) then
65         if (efree) then
66         read (intin,*,end=13,err=11) energy(icon),totfree(icon),
67      &    rmstb(icon),
68      &    nss_all(icon),(ihpb_all(ii,icon),jhpb_all(i,icon),
69      &    i=1,nss_all(icon)),iscore(icon)
70         else
71         read (intin,*,end=13,err=11) energy(icon),rmstb(icon),
72      &    nss_all(icon),(ihpb_all(ii,icon),jhpb_all(i,icon),
73      &    i=1,nss_all(icon)),iscore(icon)
74         endif
75         read (intin,'(8f10.5)',end=13,err=10) 
76      &    ((allcart(j,i,icon),j=1,3),i=1,nres),
77      &    ((allcart(j,i+nres,icon),j=1,3),i=nnt,nct)
78         print *,icon,energy(icon),nss_all(icon),rmstb(icon)
79       else 
80         read(intin,'(a80)',end=13,err=12) lineh
81         read(lineh(:5),*,err=8) ic
82         if (efree) then
83         read(lineh(6:),*,err=8) energy(icon)
84         else
85         read(lineh(6:),*,err=8) energy(icon)
86         endif
87         goto 9
88     8   ic=1
89         print *,'error, assuming e=1d10',lineh
90         energy(icon)=1d10
91         nss=0
92     9   continue
93 cold        read(lineh(18:),*,end=13,err=11) nss_all(icon)
94         ii = index(lineh(15:)," ")+15
95         read(lineh(ii:),*,end=13,err=11) nss_all(icon)
96         IF (NSS_all(icon).LT.9) THEN
97           read (lineh(20:),*,end=102)
98      &    (IHPB_all(I,icon),JHPB_all(I,icon),I=1,NSS_all(icon)),
99      &    iscore(icon)
100         ELSE
101           read (lineh(20:),*,end=102) 
102      &           (IHPB_all(I,icon),JHPB_all(I,icon),I=1,8)
103           read (intin,*) (IHPB_all(I,icon),JHPB_all(I,icon),
104      &      I=9,NSS_all(icon)),iscore(icon)
105         ENDIF
106
107   102   continue  
108
109         PRINT *,'IC:',IC,' ENERGY:',ENERGY(ICON)
110         call read_angles(intin,*13)
111         do i=1,nres
112           phiall(i,icon)=phi(i)
113           thetall(i,icon)=theta(i)
114           alphall(i,icon)=alph(i)
115           omall(i,icon)=omeg(i)
116         enddo
117       endif
118       ICON=ICON+1
119       GOTO 123
120 C
121 C CALCULATE DISTANCES
122 C
123    10 print *,'something wrong with angles'
124       goto 13
125    11 print *,'something wrong with NSS',nss
126       goto 13
127    12 print *,'something wrong with header'
128
129    13 NCON=ICON-1
130
131 #endif
132       call flush(iout)
133       jj_old=1
134       open (icbase,file=bprotfiles,status="unknown",
135      &   form="unformatted",access="direct",recl=lenrec)
136 c Read conformations from binary DA files (one per batch) and write them to 
137 c a binary DA scratchfile.
138       jj=0
139       jjj=0
140 #ifdef MPI
141       write (liczba,'(bz,i3.3)') me
142       IF (ME.EQ.MASTER) THEN
143 c Only the master reads the database; it'll send it to the other procs
144 c through a ring.
145 #endif
146         t_acq = tcpu()
147         icount=0
148
149         if (from_bx) then
150
151           open (intin,file=intinname,status="old",form="unformatted",
152      &            access="direct",recl=lenrec_in)
153
154         else if (from_cx) then
155 #if (defined(AIX) && !defined(JUBL))
156           call xdrfopen_(ixdrf,intinname, "r", iret)
157 #else
158           call xdrfopen(ixdrf,intinname, "r", iret)
159 #endif
160           prec=10000.0
161           write (iout,*) "xdrfopen: iret",iret
162           if (iret.eq.0) then
163             write (iout,*) "Error: coordinate file ",
164      &       intinname(:ilen(intinname))," does not exist."
165             call flush(iout)
166 #ifdef MPI
167             call MPI_ABORT(MPI_COMM_WORLD,IERROR,ERRCODE)
168 #endif
169             stop
170           endif
171         else
172           write (iout,*) "Error: coordinate format not specified"
173           call flush(iout)
174 #ifdef MPI
175           call MPI_ABORT(MPI_COMM_WORLD,IERROR,ERRCODE)
176 #else
177           stop
178 #endif
179         endif
180
181 #define DEBUG
182 #ifdef DEBUG
183         write (iout,*) "Opening file ",intinname(:ilen(intinname))
184         write (iout,*) "lenrec",lenrec_in
185         call flush(iout)
186 #endif
187 #undef DEBUG
188 c        write (iout,*) "maxconf",maxconf
189         i=0
190         do while (.true.)
191            i=i+1
192            if (i.gt.maxconf) then
193              write (iout,*) "Error: too many conformations ",
194      &        "(",maxconf,") maximum."
195 #ifdef MPI
196              call MPI_Abort(MPI_COMM_WORLD,errcode,ierror)
197 #endif
198              stop
199            endif
200 c          write (iout,*) "i",i
201 c          call flush(iout)
202           if (from_bx) then
203             read(intin,err=101,end=101) 
204      &       ((csingle(l,k),l=1,3),k=1,nres),
205      &       ((csingle(l,k+nres),l=1,3),k=nnt,nct),
206      &       nss,(ihpb(k),jhpb(k),k=1,nss),
207      &       energy(jj+1),
208      &       entfac(jj+1),rmstb(jj+1),iscor
209              do j=1,2*nres
210                do k=1,3
211                  c(k,j)=csingle(k,j)
212                enddo
213              enddo
214           else
215 #if (defined(AIX) && !defined(JUBL))
216             call xdrf3dfcoord_(ixdrf, csingle, itmp, prec, iret)
217             if (iret.eq.0) goto 101
218             call xdrfint_(ixdrf, nss, iret)
219             if (iret.eq.0) goto 101
220             do j=1,nss
221 cc              if (dyn_ss) then
222 cc                call xdrfint_(ixdrf, idssb(j), iret)
223 cc                if (iret.eq.0) goto 101
224 cc                call xdrfint_(ixdrf, jdssb(j), iret)
225 cc                if (iret.eq.0) goto 101
226 cc             idssb(j)=idssb(j)-nres
227 cc             jdssb(j)=jdssb(j)-nres
228 cc              else
229                 call xdrfint_(ixdrf, ihpb(j), iret)
230                 if (iret.eq.0) goto 101
231                 call xdrfint_(ixdrf, jhpb(j), iret)
232                 if (iret.eq.0) goto 101
233 cc              endif
234             enddo
235             call xdrffloat_(ixdrf,reini,iret)
236             if (iret.eq.0) goto 101
237             call xdrffloat_(ixdrf,refree,iret)
238             if (iret.eq.0) goto 101
239             call xdrffloat_(ixdrf,rmsdev,iret)
240             if (iret.eq.0) goto 101
241             call xdrfint_(ixdrf,iscor,iret)
242             if (iret.eq.0) goto 101
243 #else
244 c            write (iout,*) "calling xdrf3dfcoord"
245             call xdrf3dfcoord(ixdrf, csingle, itmp, prec, iret)
246 c            write (iout,*) "iret",iret
247 c            call flush(iout)
248             if (iret.eq.0) goto 101
249             call xdrfint(ixdrf, nss, iret)
250 c            write (iout,*) "iret",iret
251 c            write (iout,*) "nss",nss
252             call flush(iout)
253             if (iret.eq.0) goto 101
254             do k=1,nss
255 cc              if (dyn_ss) then
256 cc                call xdrfint(ixdrf, idssb(k), iret)
257 cc                if (iret.eq.0) goto 101
258 cc               call xdrfint(ixdrf, jdssb(k), iret)
259 cc                if (iret.eq.0) goto 101
260 cc                idssb(k)=idssb(k)-nres
261 cc                jdssb(k)=jdssb(k)-nres
262 cc              write(iout,*) "TUTU", idssb(k),jdssb(k)
263 cc              else
264                 call xdrfint(ixdrf, ihpb(k), iret)
265                 if (iret.eq.0) goto 101
266                 call xdrfint(ixdrf, jhpb(k), iret)
267                 if (iret.eq.0) goto 101
268 cc              endif
269             enddo
270             call xdrffloat(ixdrf,reini,iret)
271             if (iret.eq.0) goto 101
272             call xdrffloat(ixdrf,refree,iret)
273             if (iret.eq.0) goto 101
274             call xdrffloat(ixdrf,rmsdev,iret)
275             if (iret.eq.0) goto 101
276             call xdrfint(ixdrf,iscor,iret)
277             if (iret.eq.0) goto 101
278 #endif
279             energy(jj+1)=reini
280 cc         write(iout,*) 'reini=', reini, jj+1
281             entfac(jj+1)=dble(refree)
282 cc         write(iout,*) 'refree=', refree,jj+1
283             rmstb(jj+1)=rmsdev
284             do k=1,nres
285               do l=1,3
286                 c(l,k)=csingle(l,k)
287               enddo
288             enddo
289             do k=nnt,nct
290               do l=1,3
291                 c(l,nres+k)=csingle(l,nres+k-nnt+1)
292               enddo
293             enddo
294           endif
295 #ifdef DEBUG
296           write (iout,'(5hREAD ,i5,3f15.4,i10)') 
297      &     jj+1,energy(jj+1),entfac(jj+1),
298      &     rmstb(jj+1),iscor
299           write (iout,*) "Conformation",jjj+1,jj+1
300           write (iout,'(8f10.5)') ((c(j,i),j=1,3),i=1,nres)
301           write (iout,'(8f10.5)') ((c(j,i+nres),j=1,3),i=nnt,nct)
302           call flush(iout)
303 #endif
304           call add_new_cconf(jjj,jj,jj_old,icount,Next)
305         enddo
306   101   continue
307         write (iout,*) i-1," conformations read from DA file ",
308      &    intinname(:ilen(intinname))
309         write (iout,*) jj," conformations read so far"
310         if (from_bx) then
311           close(intin)
312         else
313 #if (defined(AIX) && !defined(JUBL))
314           call xdrfclose_(ixdrf, iret)
315 #else
316           call xdrfclose(ixdrf, iret)
317 #endif
318         endif
319 #ifdef MPI
320 #ifdef DEBUG    
321         write (iout,*) "jj_old",jj_old," jj",jj
322 #endif
323         call write_and_send_cconf(icount,jj_old,jj,Next)
324         call MPI_Send(0,1,MPI_INTEGER,Next,570,
325      &             MPI_COMM_WORLD,IERROR)
326         jj_old=jj+1
327 #else
328         call write_and_send_cconf(icount,jj_old,jj,Next)
329 #endif
330         t_acq = tcpu() - t_acq
331 #ifdef MPI
332         write (iout,*) "Processor",me,
333      &    " time for conformation read/send",t_acq
334       ELSE
335 c A worker gets the confs from the master and sends them to its neighbor
336         t_acq = tcpu()
337         call receive_and_pass_cconf(icount,jj_old,jj,
338      &    Previous,Next)
339         t_acq = tcpu() - t_acq
340       ENDIF
341 #endif
342       ncon=jj
343 c      close(icbase)
344       close(intin)
345
346       write(iout,*)"A total of",ncon," conformations read."
347
348 #ifdef MPI
349 c Check if everyone has the same number of conformations
350       call MPI_Allgather(ncon,1,MPI_INTEGER,
351      &  ntot_all(0),1,MPI_INTEGER,MPI_Comm_World,IERROR)
352       lerr=.false.
353       do i=0,nprocs-1
354         if (i.ne.me) then
355             if (ncon.ne.ntot_all(i)) then
356               write (iout,*) "Number of conformations at processor",i,
357      &         " differs from that at processor",me,
358      &         ncon,ntot_all(i)
359               lerr = .true.
360             endif
361         endif
362       enddo
363       if (lerr) then
364         write (iout,*)
365         write (iout,*) "Number of conformations read by processors"
366         write (iout,*)
367         do i=0,nprocs-1
368           write (iout,'(8i10)') i,ntot_all(i)
369         enddo
370         write (iout,*) "Calculation terminated."
371         call flush(iout)
372         return1
373       endif
374       return
375 #endif
376  1111 write(iout,*) "Error opening coordinate file ",
377      & intinname(:ilen(intinname))
378       call flush(iout)
379       return1
380       end
381 c------------------------------------------------------------------------------
382       subroutine add_new_cconf(jjj,jj,jj_old,icount,Next)
383       implicit none
384       include "DIMENSIONS"
385       include "sizesclu.dat"
386       include "COMMON.CLUSTER"
387       include "COMMON.CHAIN"
388       include "COMMON.INTERACT"
389       include "COMMON.LOCAL"
390       include "COMMON.IOUNITS"
391       include "COMMON.NAMES"
392       include "COMMON.VAR"
393       include "COMMON.SBRIDGE"
394       include "COMMON.GEO"
395       integer i,j,jj,jjj,jj_old,icount,k,kk,l,ii,ib
396      &  nn,nn1,inan,Next,itj
397       double precision etot,energia(0:max_ene)
398       jjj=jjj+1
399       call int_from_cart1(.false.)
400       do j=nnt+1,nct
401         if (vbld(j).lt.2.0d0 .or. vbld(j).gt.5.0d0) then
402           write (iout,*) "Conformation",jjj,jj+1
403           write (iout,*) "Bad CA-CA bond length",j," ",vbld(j)
404           write (iout,*) "The Cartesian geometry is:"
405           write (iout,'(8f10.5)') ((c(l,k),l=1,3),k=1,nres)
406           write (iout,'(8f10.5)') ((c(l,k+nres),l=1,3),k=nnt,nct)
407           write (iout,*) "The internal geometry is:"
408           write (iout,'(8f10.4)') (vbld(k),k=nnt+1,nct)
409           write (iout,'(8f10.4)') (vbld(k),k=nres+nnt,nres+nct)
410           write (iout,'(8f10.4)') (rad2deg*theta(k),k=3,nres)
411           write (iout,'(8f10.4)') (rad2deg*phi(k),k=4,nres)
412           write (iout,'(8f10.4)') (rad2deg*alph(k),k=2,nres-1)
413           write (iout,'(8f10.4)') (rad2deg*omeg(k),k=2,nres-1)
414           write (iout,*) 
415      &      "This conformation WILL NOT be added to the database."
416           return
417         endif
418       enddo
419       do j=nnt,nct
420         itj=itype(j)
421         if (itype(j).ne.10 .and. (vbld(nres+j)-dsc(itj)).gt.2.0d0) then
422           write (iout,*) "Conformation",jjj,jj+1
423           write (iout,*) "Bad CA-SC bond length",j," ",vbld(nres+j)
424           write (iout,*) "The Cartesian geometry is:"
425           write (iout,'(8f10.5)') ((c(l,k),l=1,3),k=1,nres)
426           write (iout,'(8f10.5)') ((c(l,k+nres),l=1,3),k=nnt,nct)
427           write (iout,*) "The internal geometry is:"
428           write (iout,'(8f10.4)') (vbld(k),k=nnt+1,nct)
429           write (iout,'(8f10.4)') (vbld(k),k=nres+nnt,nres+nct)
430           write (iout,'(8f10.4)') (rad2deg*theta(k),k=3,nres)
431           write (iout,'(8f10.4)') (rad2deg*phi(k),k=4,nres)
432           write (iout,'(8f10.4)') (rad2deg*alph(k),k=2,nres-1)
433           write (iout,'(8f10.4)') (rad2deg*omeg(k),k=2,nres-1)
434           write (iout,*) 
435      &      "This conformation WILL NOT be added to the database."
436           return
437         endif
438       enddo
439       do j=3,nres
440         if (theta(j).le.0.0d0) then
441           write (iout,*) 
442      &      "Zero theta angle(s) in conformation",jjj,jj+1
443           write (iout,*) "The Cartesian geometry is:"
444           write (iout,'(8f10.5)') ((c(l,k),l=1,3),k=1,nres)
445           write (iout,'(8f10.5)') ((c(l,k+nres),l=1,3),k=nnt,nct)
446           write (iout,*) "The internal geometry is:"
447           write (iout,'(8f10.4)') (vbld(k),k=nnt+1,nct)
448           write (iout,'(8f10.4)') (vbld(k),k=nres+nnt,nres+nct)
449           write (iout,'(8f10.4)') (rad2deg*theta(k),k=3,nres)
450           write (iout,'(8f10.4)') (rad2deg*phi(k),k=4,nres)
451           write (iout,'(8f10.4)') (rad2deg*alph(k),k=2,nres-1)
452           write (iout,'(8f10.4)') (rad2deg*omeg(k),k=2,nres-1)
453           write (iout,*)
454      &      "This conformation WILL NOT be added to the database."
455           return
456         endif
457         if (theta(j).gt.179.97*deg2rad) theta(j)=179.97*deg2rad
458       enddo
459       jj=jj+1
460 #ifdef DEBUG
461       write (iout,*) "Conformation",jjj,jj
462       write (iout,'(8f10.5)') ((c(j,i),j=1,3),i=1,nres)
463       write (iout,'(8f10.5)') ((c(j,i+nres),j=1,3),i=nnt,nct)
464       write (iout,'(8f10.4)') (vbld(k),k=nnt+1,nct)
465       write (iout,'(8f10.4)') (vbld(k),k=nres+nnt,nres+nct)
466       write (iout,'(8f10.4)') (rad2deg*theta(k),k=3,nres)
467       write (iout,'(8f10.4)') (rad2deg*phi(k),k=4,nres)
468       write (iout,'(8f10.4)') (vbld(k+nres),k=nnt,nct)
469       write (iout,'(8f10.4)') (rad2deg*alph(k),k=2,nres-1)
470       write (iout,'(8f10.4)') (rad2deg*omeg(k),k=2,nres-1)
471       write (iout,'(16i5)') nss,(ihpb(k),jhpb(k),k=1,nss)
472       write (iout,'(e15.5,16i5)') entfac(icount+1),
473      &        iscore(icount+1,0)
474 #endif
475       icount=icount+1
476       call store_cconf_from_file(jj,icount)
477       if (icount.eq.maxstr_proc) then
478 #ifdef DEBUG
479         write (iout,* ) "jj_old",jj_old," jj",jj
480 #endif
481         call write_and_send_cconf(icount,jj_old,jj,Next)
482         jj_old=jj+1
483         icount=0
484       endif
485       return
486       end
487 c------------------------------------------------------------------------------
488       subroutine store_cconf_from_file(jj,icount)
489       implicit none
490       include "DIMENSIONS"
491       include "sizesclu.dat"
492       include "COMMON.CLUSTER"
493       include "COMMON.CHAIN"
494       include "COMMON.SBRIDGE"
495       include "COMMON.INTERACT"
496       include "COMMON.IOUNITS"
497       include "COMMON.VAR"
498       integer i,j,jj,icount
499 c Store the conformation that has been read in
500       do i=1,2*nres
501         do j=1,3
502           allcart(j,i,icount)=c(j,i)
503         enddo
504       enddo
505       nss_all(icount)=nss
506       do i=1,nss
507         ihpb_all(i,icount)=ihpb(i)
508         jhpb_all(i,icount)=jhpb(i)
509       enddo
510       return
511       end
512 c------------------------------------------------------------------------------
513       subroutine write_and_send_cconf(icount,jj_old,jj,Next)
514       implicit none
515       include "DIMENSIONS"
516       include "sizesclu.dat"
517 #ifdef MPI
518       include "mpif.h"
519       integer IERROR
520       include "COMMON.MPI"
521 #endif
522       include "COMMON.CHAIN"
523       include "COMMON.SBRIDGE"
524       include "COMMON.INTERACT"
525       include "COMMON.IOUNITS"
526       include "COMMON.CLUSTER"
527       include "COMMON.VAR"
528       integer icount,jj_old,jj,Next
529 c Write the structures to a scratch file
530 #ifdef MPI
531 c Master sends the portion of conformations that have been read in to the neighbor
532 #ifdef DEBUG
533       write (iout,*) "Processor",me," entered WRITE_AND_SEND_CONF"
534       call flush(iout)
535 #endif
536       call MPI_Send(icount,1,MPI_INTEGER,Next,570,MPI_COMM_WORLD,IERROR)
537       call MPI_Send(nss_all(1),icount,MPI_INTEGER,
538      &    Next,571,MPI_COMM_WORLD,IERROR)
539       call MPI_Send(ihpb_all(1,1),icount,MPI_INTEGER,
540      &    Next,572,MPI_COMM_WORLD,IERROR)
541       call MPI_Send(jhpb_all(1,1),icount,MPI_INTEGER,
542      &    Next,573,MPI_COMM_WORLD,IERROR)
543       call MPI_Send(rmstb(jj_old),icount,MPI_DOUBLE_PRECISION,
544      &    Next,577,MPI_COMM_WORLD,IERROR)
545       call MPI_Send(entfac(jj_old),icount,MPI_DOUBLE_PRECISION,
546      &    Next,579,MPI_COMM_WORLD,IERROR)
547       call MPI_Send(allcart(1,1,1),3*icount*maxres2,
548      &    MPI_REAL,Next,580,MPI_COMM_WORLD,IERROR)
549 #endif
550       call dawrite_ccoords(jj_old,jj,icbase)
551       return
552       end
553 c------------------------------------------------------------------------------
554 #ifdef MPI
555       subroutine receive_and_pass_cconf(icount,jj_old,jj,Previous,
556      &  Next)
557       implicit none
558       include "DIMENSIONS"
559       include "sizesclu.dat"
560       include "mpif.h"
561       integer IERROR,STATUS(MPI_STATUS_SIZE)
562       include "COMMON.MPI"
563       include "COMMON.CHAIN"
564       include "COMMON.SBRIDGE"
565       include "COMMON.INTERACT"
566       include "COMMON.IOUNITS"
567       include "COMMON.VAR"
568       include "COMMON.GEO"
569       include "COMMON.CLUSTER"
570       integer i,j,k,icount,jj_old,jj,Previous,Next
571       icount=1
572 #ifdef DEBUG
573       write (iout,*) "Processor",me," entered RECEIVE_AND_PASS_CONF"
574       call flush(iout)
575 #endif
576       do while (icount.gt.0) 
577       call MPI_Recv(icount,1,MPI_INTEGER,Previous,570,MPI_COMM_WORLD,
578      &     STATUS,IERROR)
579       call MPI_Send(icount,1,MPI_INTEGER,Next,570,MPI_COMM_WORLD,
580      &     IERROR)
581 #ifdef DEBUG
582       write (iout,*) "Processor",me," icount",icount
583 #endif
584       if (icount.eq.0) return
585       call MPI_Recv(nss_all(1),icount,MPI_INTEGER,
586      &    Previous,571,MPI_COMM_WORLD,STATUS,IERROR)
587       call MPI_Send(nss_all(1),icount,MPI_INTEGER,
588      &  Next,571,MPI_COMM_WORLD,IERROR)
589       call MPI_Recv(ihpb_all(1,1),icount,MPI_INTEGER,
590      &    Previous,572,MPI_COMM_WORLD,STATUS,IERROR)
591       call MPI_Send(ihpb_all(1,1),icount,MPI_INTEGER,
592      &  Next,572,MPI_COMM_WORLD,IERROR)
593       call MPI_Recv(jhpb_all(1,1),icount,MPI_INTEGER,
594      &    Previous,573,MPI_COMM_WORLD,STATUS,IERROR)
595       call MPI_Send(jhpb_all(1,1),icount,MPI_INTEGER,
596      &  Next,573,MPI_COMM_WORLD,IERROR)
597       call MPI_Recv(rmstb(jj_old),icount,MPI_DOUBLE_PRECISION,
598      &  Previous,577,MPI_COMM_WORLD,STATUS,IERROR)
599       call MPI_Send(rmstb(jj_old),icount,MPI_DOUBLE_PRECISION,
600      &  Next,577,MPI_COMM_WORLD,IERROR)
601       call MPI_Recv(entfac(jj_old),icount,MPI_DOUBLE_PRECISION,
602      &  Previous,579,MPI_COMM_WORLD,STATUS,IERROR)
603       call MPI_Send(entfac(jj_old),icount,MPI_DOUBLE_PRECISION,
604      &  Next,579,MPI_COMM_WORLD,IERROR)
605       call MPI_Recv(allcart(1,1,1),3*icount*maxres2,
606      &  MPI_REAL,Previous,580,MPI_COMM_WORLD,STATUS,IERROR)
607       call MPI_Send(allcart(1,1,1),3*icount*maxres2,
608      &  MPI_REAL,Next,580,MPI_COMM_WORLD,IERROR)
609       jj=jj_old+icount-1
610       call dawrite_ccoords(jj_old,jj,icbase)
611       jj_old=jj+1
612 #ifdef DEBUG
613       write (iout,*) "Processor",me," received",icount," conformations"
614       do i=1,icount
615         write (iout,'(8f10.4)') (allcart(l,k,i),l=1,3,k=1,nres)
616         write (iout,'(8f10.4)')((allcart(l,k,i+nres),l=1,3,k=nnt,nct)
617         write (iout,'(e15.5,16i5)') entfac(i)
618       enddo
619 #endif
620       enddo
621       return
622       end
623 #endif
624 c------------------------------------------------------------------------------
625       subroutine daread_ccoords(istart_conf,iend_conf)
626       implicit none
627       include "DIMENSIONS"
628       include "sizesclu.dat"
629 #ifdef MPI
630       include "mpif.h"
631       include "COMMON.MPI"
632 #endif
633       include "COMMON.CHAIN"
634       include "COMMON.CLUSTER"
635       include "COMMON.IOUNITS"
636       include "COMMON.INTERACT"
637       include "COMMON.VAR"
638       include "COMMON.SBRIDGE"
639       include "COMMON.GEO"
640       integer istart_conf,iend_conf
641       integer i,j,ij,ii,iii
642       integer len
643       character*16 form,acc
644       character*32 nam
645 c
646 c Read conformations off a DA scratchfile.
647 c
648 #ifdef DEBUG
649       write (iout,*) "DAREAD_COORDS"
650       write (iout,*) "istart_conf",istart_conf," iend_conf",iend_conf
651       inquire(unit=icbase,name=nam,recl=len,form=form,access=acc)
652       write (iout,*) "len=",len," form=",form," acc=",acc
653       write (iout,*) "nam=",nam
654       call flush(iout)
655 #endif
656       do ii=istart_conf,iend_conf
657         ij = ii - istart_conf + 1
658         iii=list_conf(ii)
659 #ifdef DEBUG
660         write (iout,*) "Reading binary file, record",iii," ii",ii
661         call flush(iout)
662 #endif
663         if (dyn_ss) then
664         read(icbase,rec=iii) ((allcart(j,i,ij),j=1,3),i=1,nres),
665      &    ((allcart(j,i,ij),j=1,3),i=nnt+nres,nct+nres),
666 c     &    nss_all(ij),(ihpb_all(i,ij),jhpb_all(i,ij),i=1,nss),
667      &    entfac(ii),rmstb(ii)
668         else
669         read(icbase,rec=iii) ((allcart(j,i,ij),j=1,3),i=1,nres),
670      &    ((allcart(j,i,ij),j=1,3),i=nnt+nres,nct+nres),
671      &    nss_all(ij),(ihpb_all(i,ij),jhpb_all(i,ij),i=1,nss),
672      &    entfac(ii),rmstb(ii)
673         endif
674 #ifdef DEBUG
675         write (iout,*) ii,iii,ij,entfac(ii)
676         write (iout,'(8f10.5)') ((allcart(j,i,ij),j=1,3),i=1,nres)
677         write (iout,'(8f10.4)') ((allcart(j,i,ij),j=1,3),
678      &    i=nnt+nres,nct+nres)
679         write (iout,'(2e15.5)') entfac(ij)
680         write (iout,'(16i5)') nss_all(ij),(ihpb_all(i,ij),
681      &    jhpb_all(i,ij),i=1,nss)
682         call flush(iout)
683 #endif
684       enddo
685       return
686       end
687 c------------------------------------------------------------------------------
688       subroutine dawrite_ccoords(istart_conf,iend_conf,unit_out)
689       implicit none
690       include "DIMENSIONS"
691       include "sizesclu.dat"
692 #ifdef MPI
693       include "mpif.h"
694       include "COMMON.MPI"
695 #endif
696       include "COMMON.CHAIN"
697       include "COMMON.INTERACT"
698       include "COMMON.IOUNITS"
699       include "COMMON.VAR"
700       include "COMMON.SBRIDGE"
701       include "COMMON.GEO"
702       include "COMMON.CLUSTER"
703       integer istart_conf,iend_conf
704       integer i,j,ii,ij,iii,unit_out
705       integer len
706       character*16 form,acc
707       character*32 nam
708 c
709 c Write conformations to a DA scratchfile.
710 c
711 #ifdef DEBUG
712       write (iout,*) "DAWRITE_COORDS"
713       write (iout,*) "istart_conf",istart_conf," iend_conf",iend_conf
714       write (iout,*) "lenrec",lenrec
715       inquire(unit=unit_out,name=nam,recl=len,form=form,access=acc)
716       write (iout,*) "len=",len," form=",form," acc=",acc
717       write (iout,*) "nam=",nam
718       call flush(iout)
719 #endif
720       do ii=istart_conf,iend_conf
721         iii=list_conf(ii)
722         ij = ii - istart_conf + 1
723 #ifdef DEBUG
724         write (iout,*) "Writing binary file, record",iii," ii",ii
725         call flush(iout)
726 #endif
727        if (dyn_ss) then
728         write(unit_out,rec=iii) ((allcart(j,i,ij),j=1,3),i=1,nres),
729      &    ((allcart(j,i,ij),j=1,3),i=nnt+nres,nct+nres),
730 c     &    nss_all(ij),(ihpb_all(i,ij),jhpb_all(i,ij),i=1,nss_all(ij))
731      &    entfac(ii),rmstb(ii)
732         else
733         write(unit_out,rec=iii) ((allcart(j,i,ij),j=1,3),i=1,nres),
734      &    ((allcart(j,i,ij),j=1,3),i=nnt+nres,nct+nres),
735      &    nss_all(ij),(ihpb_all(i,ij),jhpb_all(i,ij),i=1,nss_all(ij)),
736      &    entfac(ii),rmstb(ii)
737         endif
738 #ifdef DEBUG
739         write (iout,'(8f10.5)') ((allcart(j,i,ij),j=1,3),i=1,nres)
740         write (iout,'(8f10.4)') ((allcart(j,i,ij),j=1,3),i=nnt+nres,
741      &   nct+nres)
742         write (iout,'(2e15.5)') entfac(ij)
743         write (iout,'(16i5)') nss_all(ij),(ihpb(i,ij),jhpb(i,ij),i=1,
744      &   nss_all(ij))
745         call flush(iout)
746 #endif
747       enddo
748       return
749       end