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