Merge branch 'devel' of mmka:unres into devel
[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               call xdrfint_(ixdrf, ihpb(j), iret)
222               if (iret.eq.0) goto 101
223               call xdrfint_(ixdrf, jhpb(j), iret)
224               if (iret.eq.0) goto 101
225             enddo
226             call xdrffloat_(ixdrf,reini,iret)
227             if (iret.eq.0) goto 101
228             call xdrffloat_(ixdrf,refree,iret)
229             if (iret.eq.0) goto 101
230             call xdrffloat_(ixdrf,rmsdev,iret)
231             if (iret.eq.0) goto 101
232             call xdrfint_(ixdrf,iscor,iret)
233             if (iret.eq.0) goto 101
234 #else
235 c            write (iout,*) "calling xdrf3dfcoord"
236             call xdrf3dfcoord(ixdrf, csingle, itmp, prec, iret)
237 c            write (iout,*) "iret",iret
238 c            call flush(iout)
239             if (iret.eq.0) goto 101
240             call xdrfint(ixdrf, nss, iret)
241 c            write (iout,*) "iret",iret
242 c            write (iout,*) "nss",nss
243             call flush(iout)
244             if (iret.eq.0) goto 101
245             do k=1,nss
246               call xdrfint(ixdrf, ihpb(k), iret)
247               if (iret.eq.0) goto 101
248               call xdrfint(ixdrf, jhpb(k), iret)
249               if (iret.eq.0) goto 101
250             enddo
251             call xdrffloat(ixdrf,reini,iret)
252             if (iret.eq.0) goto 101
253             call xdrffloat(ixdrf,refree,iret)
254             if (iret.eq.0) goto 101
255             call xdrffloat(ixdrf,rmsdev,iret)
256             if (iret.eq.0) goto 101
257             call xdrfint(ixdrf,iscor,iret)
258             if (iret.eq.0) goto 101
259 #endif
260             energy(jj+1)=reini
261             entfac(jj+1)=refree
262             rmstb(jj+1)=rmsdev
263             do k=1,nres
264               do l=1,3
265                 c(l,k)=csingle(l,k)
266               enddo
267             enddo
268             do k=nnt,nct
269               do l=1,3
270                 c(l,nres+k)=csingle(l,nres+k-nnt+1)
271               enddo
272             enddo
273           endif
274 #ifdef DEBUG
275           write (iout,'(5hREAD ,i5,3f15.4,i10)') 
276      &     jj+1,energy(jj+1),entfac(jj+1),
277      &     rmstb(jj+1),iscor
278           write (iout,*) "Conformation",jjj+1,jj+1
279           write (iout,'(8f10.5)') ((c(j,i),j=1,3),i=1,nres)
280           write (iout,'(8f10.5)') ((c(j,i+nres),j=1,3),i=nnt,nct)
281           call flush(iout)
282 #endif
283           call add_new_cconf(jjj,jj,jj_old,icount,Next)
284         enddo
285   101   continue
286         write (iout,*) i-1," conformations read from DA file ",
287      &    intinname(:ilen(intinname))
288         write (iout,*) jj," conformations read so far"
289         if (from_bx) then
290           close(intin)
291         else
292 #if (defined(AIX) && !defined(JUBL))
293           call xdrfclose_(ixdrf, iret)
294 #else
295           call xdrfclose(ixdrf, iret)
296 #endif
297         endif
298 #ifdef MPI
299 #ifdef DEBUG    
300         write (iout,*) "jj_old",jj_old," jj",jj
301 #endif
302         call write_and_send_cconf(icount,jj_old,jj,Next)
303         call MPI_Send(0,1,MPI_INTEGER,Next,570,
304      &             MPI_COMM_WORLD,IERROR)
305         jj_old=jj+1
306 #else
307         call write_and_send_cconf(icount,jj_old,jj,Next)
308 #endif
309         t_acq = tcpu() - t_acq
310 #ifdef MPI
311         write (iout,*) "Processor",me,
312      &    " time for conformation read/send",t_acq
313       ELSE
314 c A worker gets the confs from the master and sends them to its neighbor
315         t_acq = tcpu()
316         call receive_and_pass_cconf(icount,jj_old,jj,
317      &    Previous,Next)
318         t_acq = tcpu() - t_acq
319       ENDIF
320 #endif
321       ncon=jj
322 c      close(icbase)
323       close(intin)
324
325       write(iout,*)"A total of",ncon," conformations read."
326
327 #ifdef MPI
328 c Check if everyone has the same number of conformations
329       call MPI_Allgather(ncon,1,MPI_INTEGER,
330      &  ntot_all(0),1,MPI_INTEGER,MPI_Comm_World,IERROR)
331       lerr=.false.
332       do i=0,nprocs-1
333         if (i.ne.me) then
334             if (ncon.ne.ntot_all(i)) then
335               write (iout,*) "Number of conformations at processor",i,
336      &         " differs from that at processor",me,
337      &         ncon,ntot_all(i)
338               lerr = .true.
339             endif
340         endif
341       enddo
342       if (lerr) then
343         write (iout,*)
344         write (iout,*) "Number of conformations read by processors"
345         write (iout,*)
346         do i=0,nprocs-1
347           write (iout,'(8i10)') i,ntot_all(i)
348         enddo
349         write (iout,*) "Calculation terminated."
350         call flush(iout)
351         return1
352       endif
353       return
354 #endif
355  1111 write(iout,*) "Error opening coordinate file ",
356      & intinname(:ilen(intinname))
357       call flush(iout)
358       return1
359       end
360 c------------------------------------------------------------------------------
361       subroutine add_new_cconf(jjj,jj,jj_old,icount,Next)
362       implicit none
363       include "DIMENSIONS"
364       include "sizesclu.dat"
365       include "COMMON.CLUSTER"
366       include "COMMON.CHAIN"
367       include "COMMON.INTERACT"
368       include "COMMON.LOCAL"
369       include "COMMON.IOUNITS"
370       include "COMMON.NAMES"
371       include "COMMON.VAR"
372       include "COMMON.SBRIDGE"
373       include "COMMON.GEO"
374       integer i,j,jj,jjj,jj_old,icount,k,kk,l,ii,ib
375      &  nn,nn1,inan,Next,itj
376       double precision etot,energia(0:max_ene)
377       jjj=jjj+1
378       call int_from_cart1(.false.)
379       do j=nnt+1,nct
380         if (vbld(j).lt.2.0d0 .or. vbld(j).gt.5.0d0) then
381           write (iout,*) "Conformation",jjj,jj+1
382           write (iout,*) "Bad CA-CA bond length",j," ",vbld(j)
383           write (iout,*) "The Cartesian geometry is:"
384           write (iout,'(8f10.5)') ((c(l,k),l=1,3),k=1,nres)
385           write (iout,'(8f10.5)') ((c(l,k+nres),l=1,3),k=nnt,nct)
386           write (iout,*) "The internal geometry is:"
387           write (iout,'(8f10.4)') (vbld(k),k=nnt+1,nct)
388           write (iout,'(8f10.4)') (vbld(k),k=nres+nnt,nres+nct)
389           write (iout,'(8f10.4)') (rad2deg*theta(k),k=3,nres)
390           write (iout,'(8f10.4)') (rad2deg*phi(k),k=4,nres)
391           write (iout,'(8f10.4)') (rad2deg*alph(k),k=2,nres-1)
392           write (iout,'(8f10.4)') (rad2deg*omeg(k),k=2,nres-1)
393           write (iout,*) 
394      &      "This conformation WILL NOT be added to the database."
395           return
396         endif
397       enddo
398       do j=nnt,nct
399         itj=itype(j)
400         if (itype(j).ne.10 .and. (vbld(nres+j)-dsc(itj)).gt.2.0d0) then
401           write (iout,*) "Conformation",jjj,jj+1
402           write (iout,*) "Bad CA-SC bond length",j," ",vbld(nres+j)
403           write (iout,*) "The Cartesian geometry is:"
404           write (iout,'(8f10.5)') ((c(l,k),l=1,3),k=1,nres)
405           write (iout,'(8f10.5)') ((c(l,k+nres),l=1,3),k=nnt,nct)
406           write (iout,*) "The internal geometry is:"
407           write (iout,'(8f10.4)') (vbld(k),k=nnt+1,nct)
408           write (iout,'(8f10.4)') (vbld(k),k=nres+nnt,nres+nct)
409           write (iout,'(8f10.4)') (rad2deg*theta(k),k=3,nres)
410           write (iout,'(8f10.4)') (rad2deg*phi(k),k=4,nres)
411           write (iout,'(8f10.4)') (rad2deg*alph(k),k=2,nres-1)
412           write (iout,'(8f10.4)') (rad2deg*omeg(k),k=2,nres-1)
413           write (iout,*) 
414      &      "This conformation WILL NOT be added to the database."
415           return
416         endif
417       enddo
418       do j=3,nres
419         if (theta(j).le.0.0d0) then
420           write (iout,*) 
421      &      "Zero theta angle(s) in conformation",jjj,jj+1
422           write (iout,*) "The Cartesian geometry is:"
423           write (iout,'(8f10.5)') ((c(l,k),l=1,3),k=1,nres)
424           write (iout,'(8f10.5)') ((c(l,k+nres),l=1,3),k=nnt,nct)
425           write (iout,*) "The internal geometry is:"
426           write (iout,'(8f10.4)') (vbld(k),k=nnt+1,nct)
427           write (iout,'(8f10.4)') (vbld(k),k=nres+nnt,nres+nct)
428           write (iout,'(8f10.4)') (rad2deg*theta(k),k=3,nres)
429           write (iout,'(8f10.4)') (rad2deg*phi(k),k=4,nres)
430           write (iout,'(8f10.4)') (rad2deg*alph(k),k=2,nres-1)
431           write (iout,'(8f10.4)') (rad2deg*omeg(k),k=2,nres-1)
432           write (iout,*)
433      &      "This conformation WILL NOT be added to the database."
434           return
435         endif
436         if (theta(j).gt.179.97*deg2rad) theta(j)=179.97*deg2rad
437       enddo
438       jj=jj+1
439 #ifdef DEBUG
440       write (iout,*) "Conformation",jjj,jj
441       write (iout,'(8f10.5)') ((c(j,i),j=1,3),i=1,nres)
442       write (iout,'(8f10.5)') ((c(j,i+nres),j=1,3),i=nnt,nct)
443       write (iout,'(8f10.4)') (vbld(k),k=nnt+1,nct)
444       write (iout,'(8f10.4)') (vbld(k),k=nres+nnt,nres+nct)
445       write (iout,'(8f10.4)') (rad2deg*theta(k),k=3,nres)
446       write (iout,'(8f10.4)') (rad2deg*phi(k),k=4,nres)
447       write (iout,'(8f10.4)') (vbld(k+nres),k=nnt,nct)
448       write (iout,'(8f10.4)') (rad2deg*alph(k),k=2,nres-1)
449       write (iout,'(8f10.4)') (rad2deg*omeg(k),k=2,nres-1)
450       write (iout,'(16i5)') nss,(ihpb(k),jhpb(k),k=1,nss)
451       write (iout,'(e15.5,16i5)') entfac(icount+1),
452      &        iscore(icount+1,0)
453 #endif
454       icount=icount+1
455       call store_cconf_from_file(jj,icount)
456       if (icount.eq.maxstr_proc) then
457 #ifdef DEBUG
458         write (iout,* ) "jj_old",jj_old," jj",jj
459 #endif
460         call write_and_send_cconf(icount,jj_old,jj,Next)
461         jj_old=jj+1
462         icount=0
463       endif
464       return
465       end
466 c------------------------------------------------------------------------------
467       subroutine store_cconf_from_file(jj,icount)
468       implicit none
469       include "DIMENSIONS"
470       include "sizesclu.dat"
471       include "COMMON.CLUSTER"
472       include "COMMON.CHAIN"
473       include "COMMON.SBRIDGE"
474       include "COMMON.INTERACT"
475       include "COMMON.IOUNITS"
476       include "COMMON.VAR"
477       integer i,j,jj,icount
478 c Store the conformation that has been read in
479       do i=1,2*nres
480         do j=1,3
481           allcart(j,i,icount)=c(j,i)
482         enddo
483       enddo
484       nss_all(icount)=nss
485       do i=1,nss
486         ihpb_all(i,icount)=ihpb(i)
487         jhpb_all(i,icount)=jhpb(i)
488       enddo
489       return
490       end
491 c------------------------------------------------------------------------------
492       subroutine write_and_send_cconf(icount,jj_old,jj,Next)
493       implicit none
494       include "DIMENSIONS"
495       include "sizesclu.dat"
496 #ifdef MPI
497       include "mpif.h"
498       integer IERROR
499       include "COMMON.MPI"
500 #endif
501       include "COMMON.CHAIN"
502       include "COMMON.SBRIDGE"
503       include "COMMON.INTERACT"
504       include "COMMON.IOUNITS"
505       include "COMMON.CLUSTER"
506       include "COMMON.VAR"
507       integer icount,jj_old,jj,Next
508 c Write the structures to a scratch file
509 #ifdef MPI
510 c Master sends the portion of conformations that have been read in to the neighbor
511 #ifdef DEBUG
512       write (iout,*) "Processor",me," entered WRITE_AND_SEND_CONF"
513       call flush(iout)
514 #endif
515       call MPI_Send(icount,1,MPI_INTEGER,Next,570,MPI_COMM_WORLD,IERROR)
516       call MPI_Send(nss_all(1),icount,MPI_INTEGER,
517      &    Next,571,MPI_COMM_WORLD,IERROR)
518       call MPI_Send(ihpb_all(1,1),icount,MPI_INTEGER,
519      &    Next,572,MPI_COMM_WORLD,IERROR)
520       call MPI_Send(jhpb_all(1,1),icount,MPI_INTEGER,
521      &    Next,573,MPI_COMM_WORLD,IERROR)
522       call MPI_Send(rmstb(jj_old),icount,MPI_DOUBLE_PRECISION,
523      &    Next,577,MPI_COMM_WORLD,IERROR)
524       call MPI_Send(entfac(jj_old),icount,MPI_DOUBLE_PRECISION,
525      &    Next,579,MPI_COMM_WORLD,IERROR)
526       call MPI_Send(allcart(1,1,1),3*icount*maxres2,
527      &    MPI_REAL,Next,580,MPI_COMM_WORLD,IERROR)
528 #endif
529       call dawrite_ccoords(jj_old,jj,icbase)
530       return
531       end
532 c------------------------------------------------------------------------------
533 #ifdef MPI
534       subroutine receive_and_pass_cconf(icount,jj_old,jj,Previous,
535      &  Next)
536       implicit none
537       include "DIMENSIONS"
538       include "sizesclu.dat"
539       include "mpif.h"
540       integer IERROR,STATUS(MPI_STATUS_SIZE)
541       include "COMMON.MPI"
542       include "COMMON.CHAIN"
543       include "COMMON.SBRIDGE"
544       include "COMMON.INTERACT"
545       include "COMMON.IOUNITS"
546       include "COMMON.VAR"
547       include "COMMON.GEO"
548       include "COMMON.CLUSTER"
549       integer i,j,k,icount,jj_old,jj,Previous,Next
550       icount=1
551 #ifdef DEBUG
552       write (iout,*) "Processor",me," entered RECEIVE_AND_PASS_CONF"
553       call flush(iout)
554 #endif
555       do while (icount.gt.0) 
556       call MPI_Recv(icount,1,MPI_INTEGER,Previous,570,MPI_COMM_WORLD,
557      &     STATUS,IERROR)
558       call MPI_Send(icount,1,MPI_INTEGER,Next,570,MPI_COMM_WORLD,
559      &     IERROR)
560 #ifdef DEBUG
561       write (iout,*) "Processor",me," icount",icount
562 #endif
563       if (icount.eq.0) return
564       call MPI_Recv(nss_all(1),icount,MPI_INTEGER,
565      &    Previous,571,MPI_COMM_WORLD,STATUS,IERROR)
566       call MPI_Send(nss_all(1),icount,MPI_INTEGER,
567      &  Next,571,MPI_COMM_WORLD,IERROR)
568       call MPI_Recv(ihpb_all(1,1),icount,MPI_INTEGER,
569      &    Previous,572,MPI_COMM_WORLD,STATUS,IERROR)
570       call MPI_Send(ihpb_all(1,1),icount,MPI_INTEGER,
571      &  Next,572,MPI_COMM_WORLD,IERROR)
572       call MPI_Recv(jhpb_all(1,1),icount,MPI_INTEGER,
573      &    Previous,573,MPI_COMM_WORLD,STATUS,IERROR)
574       call MPI_Send(jhpb_all(1,1),icount,MPI_INTEGER,
575      &  Next,573,MPI_COMM_WORLD,IERROR)
576       call MPI_Recv(rmstb(jj_old),icount,MPI_DOUBLE_PRECISION,
577      &  Previous,577,MPI_COMM_WORLD,STATUS,IERROR)
578       call MPI_Send(rmstb(jj_old),icount,MPI_DOUBLE_PRECISION,
579      &  Next,577,MPI_COMM_WORLD,IERROR)
580       call MPI_Recv(entfac(jj_old),icount,MPI_DOUBLE_PRECISION,
581      &  Previous,579,MPI_COMM_WORLD,STATUS,IERROR)
582       call MPI_Send(entfac(jj_old),icount,MPI_DOUBLE_PRECISION,
583      &  Next,579,MPI_COMM_WORLD,IERROR)
584       call MPI_Recv(allcart(1,1,1),3*icount*maxres2,
585      &  MPI_REAL,Previous,580,MPI_COMM_WORLD,STATUS,IERROR)
586       call MPI_Send(allcart(1,1,1),3*icount*maxres2,
587      &  MPI_REAL,Next,580,MPI_COMM_WORLD,IERROR)
588       jj=jj_old+icount-1
589       call dawrite_ccoords(jj_old,jj,icbase)
590       jj_old=jj+1
591 #ifdef DEBUG
592       write (iout,*) "Processor",me," received",icount," conformations"
593       do i=1,icount
594         write (iout,'(8f10.4)') (allcart(l,k,i),l=1,3,k=1,nres)
595         write (iout,'(8f10.4)')((allcart(l,k,i+nres),l=1,3,k=nnt,nct)
596         write (iout,'(e15.5,16i5)') entfac(i)
597       enddo
598 #endif
599       enddo
600       return
601       end
602 #endif
603 c------------------------------------------------------------------------------
604       subroutine daread_ccoords(istart_conf,iend_conf)
605       implicit none
606       include "DIMENSIONS"
607       include "sizesclu.dat"
608 #ifdef MPI
609       include "mpif.h"
610       include "COMMON.MPI"
611 #endif
612       include "COMMON.CHAIN"
613       include "COMMON.CLUSTER"
614       include "COMMON.IOUNITS"
615       include "COMMON.INTERACT"
616       include "COMMON.VAR"
617       include "COMMON.SBRIDGE"
618       include "COMMON.GEO"
619       integer istart_conf,iend_conf
620       integer i,j,ij,ii,iii
621       integer len
622       character*16 form,acc
623       character*32 nam
624 c
625 c Read conformations off a DA scratchfile.
626 c
627 #ifdef DEBUG
628       write (iout,*) "DAREAD_COORDS"
629       write (iout,*) "istart_conf",istart_conf," iend_conf",iend_conf
630       inquire(unit=icbase,name=nam,recl=len,form=form,access=acc)
631       write (iout,*) "len=",len," form=",form," acc=",acc
632       write (iout,*) "nam=",nam
633       call flush(iout)
634 #endif
635       do ii=istart_conf,iend_conf
636         ij = ii - istart_conf + 1
637         iii=list_conf(ii)
638 #ifdef DEBUG
639         write (iout,*) "Reading binary file, record",iii," ii",ii
640         call flush(iout)
641 #endif
642         read(icbase,rec=iii) ((allcart(j,i,ij),j=1,3),i=1,nres),
643      &    ((allcart(j,i,ij),j=1,3),i=nnt+nres,nct+nres),
644      &    nss_all(ij),(ihpb_all(i,ij),jhpb_all(i,ij),i=1,nss),
645      &    entfac(ii),rmstb(ii)
646 #ifdef DEBUG
647         write (iout,*) ii,iii,ij,entfac(ii)
648         write (iout,'(8f10.5)') ((allcart(j,i,ij),j=1,3),i=1,nres)
649         write (iout,'(8f10.4)') ((allcart(j,i,ij),j=1,3),
650      &    i=nnt+nres,nct+nres)
651         write (iout,'(2e15.5)') entfac(ij)
652         write (iout,'(16i5)') nss_all(ij),(ihpb_all(i,ij),
653      &    jhpb_all(i,ij),i=1,nss)
654         call flush(iout)
655 #endif
656       enddo
657       return
658       end
659 c------------------------------------------------------------------------------
660       subroutine dawrite_ccoords(istart_conf,iend_conf,unit_out)
661       implicit none
662       include "DIMENSIONS"
663       include "sizesclu.dat"
664 #ifdef MPI
665       include "mpif.h"
666       include "COMMON.MPI"
667 #endif
668       include "COMMON.CHAIN"
669       include "COMMON.INTERACT"
670       include "COMMON.IOUNITS"
671       include "COMMON.VAR"
672       include "COMMON.SBRIDGE"
673       include "COMMON.GEO"
674       include "COMMON.CLUSTER"
675       integer istart_conf,iend_conf
676       integer i,j,ii,ij,iii,unit_out
677       integer len
678       character*16 form,acc
679       character*32 nam
680 c
681 c Write conformations to a DA scratchfile.
682 c
683 #ifdef DEBUG
684       write (iout,*) "DAWRITE_COORDS"
685       write (iout,*) "istart_conf",istart_conf," iend_conf",iend_conf
686       write (iout,*) "lenrec",lenrec
687       inquire(unit=unit_out,name=nam,recl=len,form=form,access=acc)
688       write (iout,*) "len=",len," form=",form," acc=",acc
689       write (iout,*) "nam=",nam
690       call flush(iout)
691 #endif
692       do ii=istart_conf,iend_conf
693         iii=list_conf(ii)
694         ij = ii - istart_conf + 1
695 #ifdef DEBUG
696         write (iout,*) "Writing binary file, record",iii," ii",ii
697         call flush(iout)
698 #endif
699         write(unit_out,rec=iii) ((allcart(j,i,ij),j=1,3),i=1,nres),
700      &    ((allcart(j,i,ij),j=1,3),i=nnt+nres,nct+nres),
701      &    nss_all(ij),(ihpb_all(i,ij),jhpb_all(i,ij),i=1,nss_all(ij)),
702      &    entfac(ii),rmstb(ii)
703 #ifdef DEBUG
704         write (iout,'(8f10.5)') ((allcart(j,i,ij),j=1,3),i=1,nres)
705         write (iout,'(8f10.4)') ((allcart(j,i,ij),j=1,3),i=nnt+nres,
706      &   nct+nres)
707         write (iout,'(2e15.5)') entfac(ij)
708         write (iout,'(16i5)') nss_all(ij),(ihpb(i,ij),jhpb(i,ij),i=1,
709      &   nss_all(ij))
710         call flush(iout)
711 #endif
712       enddo
713       return
714       end