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