D i L aminokwasy dzialajace dla src_MD
[unres.git] / source / unres / src_MD / geomout.F
1       subroutine pdbout(etot,tytul,iunit)
2       implicit real*8 (a-h,o-z)
3       include 'DIMENSIONS'
4       include 'COMMON.CHAIN'
5       include 'COMMON.INTERACT'
6       include 'COMMON.NAMES'
7       include 'COMMON.IOUNITS'
8       include 'COMMON.HEADER'
9       include 'COMMON.SBRIDGE'
10       include 'COMMON.DISTFIT'
11       include 'COMMON.MD'
12       character*50 tytul
13       dimension ica(maxres)
14       write (iunit,'(3a,1pe15.5)') 'REMARK ',tytul,' ENERGY ',etot
15 cmodel      write (iunit,'(a5,i6)') 'MODEL',1
16       if (nhfrag.gt.0) then
17        do j=1,nhfrag
18         iti=itype(hfrag(1,j))
19         itj=itype(hfrag(2,j))
20         if (j.lt.10) then
21            write (iunit,'(a5,i5,1x,a1,i1,2x,a3,i7,2x,a3,i7,i3,t76,i5)') 
22      &           'HELIX',j,'H',j,
23      &           restyp(iti),hfrag(1,j)-1,
24      &           restyp(itj),hfrag(2,j)-1,1,hfrag(2,j)-hfrag(1,j)
25         else
26              write (iunit,'(a5,i5,1x,a1,i2,1x,a3,i7,2x,a3,i7,i3)') 
27      &           'HELIX',j,'H',j,
28      &           restyp(iti),hfrag(1,j)-1,
29      &           restyp(itj),hfrag(2,j)-1,1,hfrag(2,j)-hfrag(1,j)
30         endif
31        enddo
32       endif
33
34       if (nbfrag.gt.0) then
35
36        do j=1,nbfrag
37
38         iti=itype(bfrag(1,j))
39         itj=itype(bfrag(2,j)-1)
40
41         write (iunit,'(a5,i5,1x,a1,i1,i3,1x,a3,i6,2x,a3,i6,i3)') 
42      &           'SHEET',1,'B',j,2,
43      &           restyp(iti),bfrag(1,j)-1,
44      &           restyp(itj),bfrag(2,j)-2,0
45
46         if (bfrag(3,j).gt.bfrag(4,j)) then
47
48          itk=itype(bfrag(3,j))
49          itl=itype(bfrag(4,j)+1)
50
51          write (iunit,'(a5,i5,1x,a1,i1,i3,1x,a3,i6,2x,a3,i6,i3,
52      &              2x,a1,2x,a3,i6,3x,a1,2x,a3,i6)') 
53      &           'SHEET',2,'B',j,2,
54      &           restyp(itl),bfrag(4,j),
55      &           restyp(itk),bfrag(3,j)-1,-1,
56      &           "N",restyp(itk),bfrag(3,j)-1,
57      &           "O",restyp(iti),bfrag(1,j)-1
58
59         else
60
61          itk=itype(bfrag(3,j))
62          itl=itype(bfrag(4,j)-1)
63
64
65         write (iunit,'(a5,i5,1x,a1,i1,i3,1x,a3,i6,2x,a3,i6,i3,
66      &              2x,a1,2x,a3,i6,3x,a1,2x,a3,i6)') 
67      &           'SHEET',2,'B',j,2,
68      &           restyp(itk),bfrag(3,j)-1,
69      &           restyp(itl),bfrag(4,j)-2,1,
70      &           "N",restyp(itk),bfrag(3,j)-1,
71      &           "O",restyp(iti),bfrag(1,j)-1
72
73
74
75         endif
76          
77        enddo
78       endif 
79
80       if (nss.gt.0) then
81         do i=1,nss
82           write(iunit,'(a6,i4,1x,a3,i7,4x,a3,i7)') 
83      &         'SSBOND',i,'CYS',ihpb(i)-1-nres,
84      &                    'CYS',jhpb(i)-1-nres
85         enddo
86       endif
87       
88       iatom=0
89       do i=nnt,nct
90         ires=i-nnt+1
91         iatom=iatom+1
92         ica(i)=iatom
93         iti=itype(i)
94         write (iunit,10) iatom,restyp(iti),ires,(c(j,i),j=1,3),vtot(i)
95         if (iti.ne.10) then
96           iatom=iatom+1
97           write (iunit,20) iatom,restyp(iti),ires,(c(j,nres+i),j=1,3),
98      &      vtot(i+nres)
99         endif
100       enddo
101       write (iunit,'(a)') 'TER'
102       do i=nnt,nct-1
103         if (itype(i).eq.10) then
104           write (iunit,30) ica(i),ica(i+1)
105         else
106           write (iunit,30) ica(i),ica(i+1),ica(i)+1
107         endif
108       enddo
109       if (itype(nct).ne.10) then
110         write (iunit,30) ica(nct),ica(nct)+1
111       endif
112       do i=1,nss
113         write (iunit,30) ica(ihpb(i)-nres)+1,ica(jhpb(i)-nres)+1
114       enddo
115       write (iunit,'(a6)') 'ENDMDL'     
116   10  FORMAT ('ATOM',I7,'  CA  ',A3,I6,4X,3F8.3,f15.3)
117   20  FORMAT ('ATOM',I7,'  CB  ',A3,I6,4X,3F8.3,f15.3)
118   30  FORMAT ('CONECT',8I5)
119       return
120       end
121 c------------------------------------------------------------------------------
122       subroutine MOL2out(etot,tytul)
123 C Prints the Cartesian coordinates of the alpha-carbons in the Tripos mol2 
124 C format.
125       implicit real*8 (a-h,o-z)
126       include 'DIMENSIONS'
127       include 'COMMON.CHAIN'
128       include 'COMMON.INTERACT'
129       include 'COMMON.NAMES'
130       include 'COMMON.IOUNITS'
131       include 'COMMON.HEADER'
132       include 'COMMON.SBRIDGE'
133       character*50 tytul,fd
134       character*3 zahl
135       character*6 res_num,pom,ucase
136 #ifdef AIX
137       call fdate_(fd)
138 #elif (defined CRAY)
139       call date(fd)
140 #else
141       call fdate(fd)
142 #endif
143       write (imol2,'(a)') '#'
144       write (imol2,'(a)') 
145      & '#         Creating user name:           unres'
146       write (imol2,'(2a)') '#         Creation time:                ',
147      & fd
148       write (imol2,'(/a)') '\@<TRIPOS>MOLECULE'
149       write (imol2,'(a)') tytul
150       write (imol2,'(5i5)') nct-nnt+1,nct-nnt+nss+1,nct-nnt+nss+1,0,0
151       write (imol2,'(a)') 'SMALL'
152       write (imol2,'(a)') 'USER_CHARGES'
153       write (imol2,'(a)') '\@<TRIPOS>ATOM' 
154       do i=nnt,nct
155         write (zahl,'(i3)') i
156         pom=ucase(restyp(itype(i)))
157         res_num = pom(:3)//zahl(2:)
158         write (imol2,10) i,(c(j,i),j=1,3),i,res_num,0.0
159       enddo
160       write (imol2,'(a)') '\@<TRIPOS>BOND'
161       do i=nnt,nct-1
162         write (imol2,'(i5,2i6,i2)') i-nnt+1,i-nnt+1,i-nnt+2,1
163       enddo
164       do i=1,nss
165         write (imol2,'(i5,2i6,i2)') nct-nnt+i,ihpb(i),jhpb(i),1
166       enddo
167       write (imol2,'(a)') '\@<TRIPOS>SUBSTRUCTURE'
168       do i=nnt,nct
169         write (zahl,'(i3)') i
170         pom = ucase(restyp(itype(i)))
171         res_num = pom(:3)//zahl(2:)
172         write (imol2,30) i-nnt+1,res_num,i-nnt+1,0
173       enddo
174   10  FORMAT (I7,' CA      ',3F10.4,' C.3',I8,1X,A,F11.4,' ****')
175   30  FORMAT (I7,1x,A,I14,' RESIDUE',I13,' ****  ****')
176       return
177       end
178 c------------------------------------------------------------------------
179       subroutine intout
180       implicit real*8 (a-h,o-z)
181       include 'DIMENSIONS'
182       include 'COMMON.IOUNITS'
183       include 'COMMON.CHAIN'
184       include 'COMMON.VAR'
185       include 'COMMON.LOCAL'
186       include 'COMMON.INTERACT'
187       include 'COMMON.NAMES'
188       include 'COMMON.GEO'
189       write (iout,'(/a)') 'Geometry of the virtual chain.'
190       write (iout,'(7a)') '  Res  ','         d','     Theta',
191      & '     Gamma','       Dsc','     Alpha','      Beta '
192       do i=1,nres
193         iti=itype(i)
194         write (iout,'(a3,i4,6f10.3)') restyp(iti),i,vbld(i),
195      &     rad2deg*theta(i),rad2deg*phi(i),vbld(nres+i),rad2deg*alph(i),
196      &     rad2deg*omeg(i)
197       enddo
198       return
199       end
200 c---------------------------------------------------------------------------
201       subroutine briefout(it,ener)
202       implicit real*8 (a-h,o-z)
203       include 'DIMENSIONS'
204       include 'COMMON.IOUNITS'
205       include 'COMMON.CHAIN'
206       include 'COMMON.VAR'
207       include 'COMMON.LOCAL'
208       include 'COMMON.INTERACT'
209       include 'COMMON.NAMES'
210       include 'COMMON.GEO'
211       include 'COMMON.SBRIDGE'
212 c     print '(a,i5)',intname,igeom
213 #if defined(AIX) || defined(PGI)
214       open (igeom,file=intname,position='append')
215 #else
216       open (igeom,file=intname,access='append')
217 #endif
218       IF (NSS.LE.9) THEN
219         WRITE (igeom,180) IT,ENER,NSS,(IHPB(I),JHPB(I),I=1,NSS)
220       ELSE
221         WRITE (igeom,180) IT,ENER,NSS,(IHPB(I),JHPB(I),I=1,9)
222         WRITE (igeom,190) (IHPB(I),JHPB(I),I=10,NSS)
223       ENDIF
224 c     IF (nvar.gt.nphi) WRITE (igeom,200) (RAD2DEG*THETA(I),I=3,NRES)
225       WRITE (igeom,200) (RAD2DEG*THETA(I),I=3,NRES)
226       WRITE (igeom,200) (RAD2DEG*PHI(I),I=4,NRES)
227 c     if (nvar.gt.nphi+ntheta) then
228         write (igeom,200) (rad2deg*alph(i),i=2,nres-1)
229         write (igeom,200) (rad2deg*omeg(i),i=2,nres-1)
230 c     endif
231       close(igeom)
232   180 format (I5,F12.3,I2,9(1X,2I3))
233   190 format (3X,11(1X,2I3))
234   200 format (8F10.4)
235       return
236       end
237 #ifdef WINIFL
238       subroutine fdate(fd)
239       character*32 fd
240       write(fd,'(32x)')
241       return
242       end
243 #endif
244 c----------------------------------------------------------------
245 #ifdef NOXDR
246       subroutine cartout(time)
247 #else
248       subroutine cartoutx(time)
249 #endif
250       implicit real*8 (a-h,o-z)
251       include 'DIMENSIONS'
252       include 'COMMON.CHAIN'
253       include 'COMMON.INTERACT'
254       include 'COMMON.NAMES'
255       include 'COMMON.IOUNITS'
256       include 'COMMON.HEADER'
257       include 'COMMON.SBRIDGE'
258       include 'COMMON.DISTFIT'
259       include 'COMMON.MD'
260       double precision time
261 #if defined(AIX) || defined(PGI)
262       open(icart,file=cartname,position="append")
263 #else
264       open(icart,file=cartname,access="append")
265 #endif
266       write (icart,'(e15.8,2e15.5,f12.5,$)') time,potE,uconst,t_bath
267       write (icart,'(i4,$)')
268      &   nss,(ihpb(j),jhpb(j),j=1,nss)
269        write (icart,'(i4,20f7.4)') nfrag+npair+3*nfrag_back,
270      & (qfrag(i),i=1,nfrag),(qpair(i),i=1,npair),
271      & (utheta(i),ugamma(i),uscdiff(i),i=1,nfrag_back)
272       write (icart,'(8f10.5)')
273      & ((c(k,j),k=1,3),j=1,nres),
274      & ((c(k,j+nres),k=1,3),j=nnt,nct)
275       close(icart)
276       return
277       end
278 c-----------------------------------------------------------------
279 #ifndef NOXDR
280       subroutine cartout(time)
281       implicit real*8 (a-h,o-z)
282       include 'DIMENSIONS'
283 #ifdef MPI
284       include 'mpif.h'
285       include 'COMMON.SETUP'
286 #else
287       parameter (me=0)
288 #endif
289       include 'COMMON.CHAIN'
290       include 'COMMON.INTERACT'
291       include 'COMMON.NAMES'
292       include 'COMMON.IOUNITS'
293       include 'COMMON.HEADER'
294       include 'COMMON.SBRIDGE'
295       include 'COMMON.DISTFIT'
296       include 'COMMON.MD'
297       double precision time
298       integer iret,itmp
299       real xcoord(3,maxres2+2),prec
300
301 #ifdef AIX
302       call xdrfopen_(ixdrf,cartname, "a", iret)
303       call xdrffloat_(ixdrf, real(time), iret)
304       call xdrffloat_(ixdrf, real(potE), iret)
305       call xdrffloat_(ixdrf, real(uconst), iret)
306       call xdrffloat_(ixdrf, real(uconst_back), iret)
307       call xdrffloat_(ixdrf, real(t_bath), iret)
308       call xdrfint_(ixdrf, nss, iret) 
309       do j=1,nss
310         call xdrfint_(ixdrf, ihpb(j), iret)
311         call xdrfint_(ixdrf, jhpb(j), iret)
312       enddo
313       call xdrfint_(ixdrf, nfrag+npair+3*nfrag_back, iret)
314       do i=1,nfrag
315         call xdrffloat_(ixdrf, real(qfrag(i)), iret)
316       enddo
317       do i=1,npair
318         call xdrffloat_(ixdrf, real(qpair(i)), iret)
319       enddo
320       do i=1,nfrag_back
321         call xdrffloat_(ixdrf, real(utheta(i)), iret)
322         call xdrffloat_(ixdrf, real(ugamma(i)), iret)
323         call xdrffloat_(ixdrf, real(uscdiff(i)), iret)
324       enddo
325 #else
326       call xdrfopen(ixdrf,cartname, "a", iret)
327       call xdrffloat(ixdrf, real(time), iret)
328       call xdrffloat(ixdrf, real(potE), iret)
329       call xdrffloat(ixdrf, real(uconst), iret)
330       call xdrffloat(ixdrf, real(uconst_back), iret)
331       call xdrffloat(ixdrf, real(t_bath), iret)
332       call xdrfint(ixdrf, nss, iret) 
333       do j=1,nss
334         call xdrfint(ixdrf, ihpb(j), iret)
335         call xdrfint(ixdrf, jhpb(j), iret)
336       enddo
337       call xdrfint(ixdrf, nfrag+npair+3*nfrag_back, iret)
338       do i=1,nfrag
339         call xdrffloat(ixdrf, real(qfrag(i)), iret)
340       enddo
341       do i=1,npair
342         call xdrffloat(ixdrf, real(qpair(i)), iret)
343       enddo
344       do i=1,nfrag_back
345         call xdrffloat(ixdrf, real(utheta(i)), iret)
346         call xdrffloat(ixdrf, real(ugamma(i)), iret)
347         call xdrffloat(ixdrf, real(uscdiff(i)), iret)
348       enddo
349 #endif
350       prec=10000.0
351       do i=1,nres
352        do j=1,3
353         xcoord(j,i)=c(j,i)
354        enddo
355       enddo
356       do i=nnt,nct
357        do j=1,3
358         xcoord(j,nres+i-nnt+1)=c(j,i+nres)
359        enddo
360       enddo
361
362       itmp=nres+nct-nnt+1
363 #ifdef AIX
364       call xdrf3dfcoord_(ixdrf, xcoord, itmp, prec, iret)
365       call xdrfclose_(ixdrf, iret)
366 #else
367       call xdrf3dfcoord(ixdrf, xcoord, itmp, prec, iret)
368       call xdrfclose(ixdrf, iret)
369 #endif
370       return
371       end
372 #endif
373 c-----------------------------------------------------------------
374       subroutine statout(itime)
375       implicit real*8 (a-h,o-z)
376       include 'DIMENSIONS'
377       include 'COMMON.CONTROL'
378       include 'COMMON.CHAIN'
379       include 'COMMON.INTERACT'
380       include 'COMMON.NAMES'
381       include 'COMMON.IOUNITS'
382       include 'COMMON.HEADER'
383       include 'COMMON.SBRIDGE'
384       include 'COMMON.DISTFIT'
385       include 'COMMON.MD'
386       include 'COMMON.REMD'
387       include 'COMMON.SETUP'
388       integer itime
389       double precision energia(0:n_ene)
390       double precision gyrate
391       external gyrate
392       common /gucio/ cm
393       character*256 line1,line2
394       character*4 format1,format2
395       character*30 format
396 #ifdef AIX
397       if(itime.eq.0) then
398        open(istat,file=statname,position="append")
399       endif
400 #else
401 #ifdef PGI
402       open(istat,file=statname,position="append")
403 #else
404       open(istat,file=statname,access="append")
405 #endif
406 #endif
407        if (refstr) then
408          call rms_nac_nnc(rms,frac,frac_nn,co,.false.)
409         if(tnp .or. tnp1 .or. tnh) then
410         write (line1,'(i10,f15.2,3f12.3,f12.6,f7.2,4f6.3,3f12.3,i5,$)')
411      &          itime,totT,EK,potE,totE,hhh,
412      &          rms,frac,frac_nn,co,amax,kinetic_T,t_bath,gyrate(),me
413           format1="a145"
414         else
415           write (line1,'(i10,f15.2,3f12.3,f7.2,4f6.3,3f12.3,i5,$)')
416      &          itime,totT,EK,potE,totE,
417      &          rms,frac,frac_nn,co,amax,kinetic_T,t_bath,gyrate(),me
418           format1="a133"
419         endif
420        else
421         if(tnp .or. tnp1 .or. tnh) then
422           write (line1,'(i10,f15.2,7f12.3,f12.6,i5,$)')
423      &           itime,totT,EK,potE,totE,hhh,
424      &           amax,kinetic_T,t_bath,gyrate(),me
425           format1="a126"
426         else
427           write (line1,'(i10,f15.2,7f12.3,i5,$)')
428      &           itime,totT,EK,potE,totE,
429      &           amax,kinetic_T,t_bath,gyrate(),me
430           format1="a114"
431         endif
432        endif
433         if(usampl.and.totT.gt.eq_time) then
434            write(line2,'(i5,2f9.4,300f7.4)') iset,uconst,uconst_back,
435      &      (qfrag(ii1),ii1=1,nfrag),(qpair(ii2),ii2=1,npair),
436      &      (utheta(i),ugamma(i),uscdiff(i),i=1,nfrag_back)
437            write(format2,'(a1,i3.3)') "a",23+7*nfrag+7*npair
438      &             +21*nfrag_back
439         elseif(hremd.gt.0) then
440            write(line2,'(i5)') iset
441            format2="a005"
442         else
443            format2="a001"
444            line2=' '
445         endif
446         if (print_compon) then
447           write(format,'(a1,a4,a1,a4,a10)') "(",format1,",",format2,
448      &                                                     ",20f12.3)"
449           write (istat,format) line1,line2,
450      &      (potEcomp(print_order(i)),i=1,nprint_ene)
451         else
452           write(format,'(a1,a4,a1,a4,a1)') "(",format1,",",format2,")"
453           write (istat,format) line1,line2
454         endif
455 #if defined(AIX)
456         call flush(istat)
457 #else
458         close(istat)
459 #endif
460        return
461       end
462 c---------------------------------------------------------------                      
463       double precision function gyrate()
464       implicit real*8 (a-h,o-z)
465       include 'DIMENSIONS'
466       include 'COMMON.INTERACT'
467       include 'COMMON.CHAIN'
468       double precision cen(3),rg
469
470       do j=1,3
471        cen(j)=0.0d0
472       enddo
473
474       do i=nnt,nct
475           do j=1,3
476             cen(j)=cen(j)+c(j,i)
477           enddo
478       enddo
479       do j=1,3
480             cen(j)=cen(j)/dble(nct-nnt+1)
481       enddo
482       rg = 0.0d0
483       do i = nnt, nct
484         do j=1,3
485          rg = rg + (c(j,i)-cen(j))**2 
486         enddo
487       end do
488       gyrate = sqrt(rg/dble(nct-nnt+1))
489       return
490       end
491