the first working version of multichain homology
[unres.git] / source / wham / src-M / readpdb.f
index fa2d05d..e7343d9 100644 (file)
@@ -3,6 +3,7 @@ C Read the PDB file and convert the peptide geometry into virtual-chain
 C geometry.
       implicit none
       include 'DIMENSIONS'
+      include 'DIMENSIONS.FREE'
       include 'DIMENSIONS.ZSCOPT'
       include 'COMMON.CONTROL'
       include 'COMMON.LOCAL'
@@ -26,8 +27,10 @@ C geometry.
           goto 10
         else if (card(:3).eq.'TER') then
 C End current chain
-          ires_old=ires+1 
-          itype(ires_old)=21
+c          ires_old=ires+1 
+          ires_old=ires+2
+          itype(ires_old-1)=ntyp1 
+          itype(ires_old)=ntyp1
           ibeg=2
 c          write (iout,*) "Chain ended",ires,ishift,ires_old
           call sccenter(ires,iii,sccor)
@@ -48,7 +51,7 @@ c            write (iout,'(a80)') card
               ishift=ires-1
               if (res.ne.'GLY' .and. res.ne. 'ACE') then
                 ishift=ishift-1
-                itype(1)=21
+                itype(1)=ntyp1
               endif
 c              write (iout,*) "ires",ires," ibeg",ibeg," ishift",ishift
               ibeg=0          
@@ -85,14 +88,51 @@ C system
       nres=ires
       do i=2,nres-1
 c        write (iout,*) i,itype(i)
-        if (itype(i).eq.21) then
-c          write (iout,*) "dummy",i,itype(i)
-          do j=1,3
-            c(j,i)=((c(j,i-1)+c(j,i+1))/2+2*c(j,i-1)-c(j,i-2))/2
-c            c(j,i)=(c(j,i-1)+c(j,i+1))/2
-            dc(j,i)=c(j,i)
-          enddo
-        endif
+
+        if (itype(i).eq.ntyp1) then
+         if (itype(i+1).eq.ntyp1) then
+C 16/01/2014 by Adasko: Adding to dummy atoms in the chain
+C first is connected prevous chain (itype(i+1).eq.ntyp1)=true
+C second dummy atom is conected to next chain itype(i+1).eq.ntyp1=false
+C           if (unres_pdb) then
+C 2/15/2013 by Adam: corrected insertion of the last dummy residue
+C            call refsys(i-3,i-2,i-1,e1,e2,e3,fail)
+C            if (fail) then
+C              e2(1)=0.0d0
+C              e2(2)=1.0d0
+C              e2(3)=0.0d0
+C            endif !fail
+C            do j=1,3
+C             c(j,i)=c(j,i-1)-1.9d0*e2(j)
+C            enddo
+C           else   !unres_pdb
+           do j=1,3
+             dcj=(c(j,i-2)-c(j,i-3))/2.0
+             c(j,i)=c(j,i-1)+dcj
+             c(j,nres+i)=c(j,i)
+           enddo     
+C          endif   !unres_pdb
+         else     !itype(i+1).eq.ntyp1
+C          if (unres_pdb) then
+C 2/15/2013 by Adam: corrected insertion of the first dummy residue
+C            call refsys(i+1,i+2,i+3,e1,e2,e3,fail)
+C            if (fail) then
+C              e2(1)=0.0d0
+C              e2(2)=1.0d0
+C              e2(3)=0.0d0
+C            endif
+C            do j=1,3
+C              c(j,i)=c(j,i+1)-1.9d0*e2(j)
+C            enddo
+C          else !unres_pdb
+           do j=1,3
+            dcj=(c(j,i+3)-c(j,i+2))/2.0
+            c(j,i)=c(j,i+1)-dcj
+            c(j,nres+i)=c(j,i)
+           enddo
+C          endif !unres_pdb
+         endif !itype(i+1).eq.ntyp1
+        endif  !itype.eq.ntyp1
       enddo
 C Calculate the CM of the last side chain.
       call sccenter(ires,iii,sccor)
@@ -100,9 +140,9 @@ C Calculate the CM of the last side chain.
       nstart_sup=1
       if (itype(nres).ne.10) then
         nres=nres+1
-        itype(nres)=21
+        itype(nres)=ntyp1
         do j=1,3
-          dcj=c(j,nres-2)-c(j,nres-3)
+          dcj=(c(j,nres-2)-c(j,nres-3))/2.0
           c(j,nres)=c(j,nres-1)+dcj
           c(j,2*nres)=c(j,nres)
         enddo
@@ -116,11 +156,11 @@ C Calculate the CM of the last side chain.
         c(j,nres+1)=c(j,1)
         c(j,2*nres)=c(j,nres)
       enddo
-      if (itype(1).eq.21) then
+      if (itype(1).eq.ntyp1) then
         nsup=nsup-1
         nstart_sup=2
         do j=1,3
-          dcj=c(j,4)-c(j,3)
+          dcj=(c(j,4)-c(j,3))/2.0
           c(j,1)=c(j,2)-dcj
           c(j,nres+1)=c(j,1)
         enddo
@@ -131,7 +171,9 @@ C Calculate internal coordinates.
      &    ires,itype(ires),restyp(itype(ires)),(c(j,ires),j=1,3),
      &    (c(j,nres+ires),j=1,3)
       enddo
+      call int_from_cart1(.false.)
       call int_from_cart(.true.,.false.)
+      call sc_loc_geom(.true.)
       write (iout,*) "After int_from_cart"
       call flush(iout)
       do i=1,nres-1
@@ -148,6 +190,16 @@ C Calculate internal coordinates.
 c        write (iout,*) i,(dc(j,i+nres),j=1,3),(dc_norm(j,i+nres),j=1,3),
 c     &   vbld_inv(i+nres)
       enddo
+      do i=1,nres
+        thetaref(i)=theta(i)
+        phiref(i)=phi(i)
+c
+        phi_ref(i)=phi(i)
+        theta_ref(i)=theta(i)
+        alph_ref(i)=alph(i)
+        omeg_ref(i)=omeg(i)
+      enddo
+      
 c      call chainbuild
 C Copy the coordinates to reference coordinates
 c      do i=1,2*nres
@@ -163,7 +215,7 @@ C Splits to single chain if occurs
       lll=lll+1
 cc      write (iout,*) "spraw lancuchy",(c(j,i),j=1,3)
       if (i.gt.1) then 
-      if (itype(i-1).eq.21) then
+      if ((itype(i-1).eq.ntyp1).and.(i.gt.2).and.(i.ne.nres)) then
       chain_length=lll-1
       kkk=kkk+1
 c       write (iout,*) "spraw lancuchy",(c(j,i),j=1,3)
@@ -173,19 +225,19 @@ c       write (iout,*) "spraw lancuchy",(c(j,i),j=1,3)
         do j=1,3
           cref(j,i,cou)=c(j,i)
           cref(j,i+nres,cou)=c(j,i+nres)
-          if ((i.le.nres).and.(symetr.gt.1)) then
+          if (i.le.nres) then
           chain_rep(j,lll,kkk)=c(j,i)
           chain_rep(j,lll+nres,kkk)=c(j,i+nres)
           endif
          enddo
       enddo
-      if (symetr.gt.1) then
+      if (chain_length.eq.0) chain_length=nres
+      write (iout,*) chain_length
       do j=1,3
       chain_rep(j,chain_length,symetr)=chain_rep(j,chain_length,1)
       chain_rep(j,chain_length+nres,symetr)
      &=chain_rep(j,chain_length+nres,1)
       enddo
-      endif
 c diagnostic
 
 c diagnostic
@@ -280,27 +332,18 @@ c---------------------------------------------------------------------------
      & '       Phi'
         endif
       endif
-      do i=2,nres
+      do i=1,nres-1
         iti=itype(i)
-        write (iout,*) i,i-1,(c(j,i),j=1,3),(c(j,i-1),j=1,3),dist(i,i-1)
-        if (itype(i-1).ne.21 .and. itype(i).ne.21 .and.
-     &    (dist(i,i-1).lt.2.0D0 .or. dist(i,i-1).gt.5.0D0)) then
+        if (iti.ne.ntyp1 .and. itype(i+1).ne.ntyp1 .and. 
+     &      (dist(i,i+1).lt.2.0D0 .or. dist(i,i+1).gt.5.0D0)) then
           write (iout,'(a,i4)') 'Bad Cartesians for residue',i
           stop
         endif
-        theta(i+1)=alpha(i-1,i,i+1)
+        vbld(i+1)=dist(i,i+1)
+        vbld_inv(i+1)=1.0d0/vbld(i+1)
+        if (i.gt.1) theta(i+1)=alpha(i-1,i,i+1)
         if (i.gt.2) phi(i+1)=beta(i-2,i-1,i,i+1)
       enddo
-      if (itype(1).eq.21) then
-        do j=1,3
-          c(j,1)=c(j,2)+(c(j,3)-c(j,4))
-        enddo
-      endif
-      if (itype(nres).eq.21) then
-        do j=1,3
-          c(j,nres)=c(j,nres-1)+(c(j,nres-2)-c(j,nres-3))
-        enddo
-      endif
       if (lside) then
         do i=2,nres-1
           do j=1,3
@@ -326,6 +369,97 @@ c---------------------------------------------------------------------------
       endif
       return
       end
+
+c-------------------------------------------------------------------------------
+      subroutine sc_loc_geom(lprn)
+      implicit real*8 (a-h,o-z)
+      include 'DIMENSIONS'
+      include 'DIMENSIONS.ZSCOPT'
+      include 'DIMENSIONS.FREE'
+      include 'COMMON.LOCAL'
+      include 'COMMON.VAR'
+      include 'COMMON.CHAIN'
+      include 'COMMON.INTERACT'
+      include 'COMMON.IOUNITS'
+      include 'COMMON.GEO'
+      include 'COMMON.NAMES'
+      include 'COMMON.CONTROL'
+      include 'COMMON.SETUP'
+      double precision x_prime(3),y_prime(3),z_prime(3)
+      logical lprn
+      do i=1,nres-1
+        do j=1,3
+          dc_norm(j,i)=vbld_inv(i+1)*(c(j,i+1)-c(j,i))
+        enddo
+      enddo
+      do i=2,nres-1
+        if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then
+          do j=1,3
+            dc_norm(j,i+nres)=vbld_inv(i+nres)*(c(j,i+nres)-c(j,i))
+          enddo
+        else
+          do j=1,3
+            dc_norm(j,i+nres)=0.0d0
+          enddo
+        endif
+      enddo
+      do i=2,nres-1
+        costtab(i+1) =dcos(theta(i+1))
+        sinttab(i+1) =dsqrt(1-costtab(i+1)*costtab(i+1))
+        cost2tab(i+1)=dsqrt(0.5d0*(1.0d0+costtab(i+1)))
+        sint2tab(i+1)=dsqrt(0.5d0*(1.0d0-costtab(i+1)))
+        cosfac2=0.5d0/(1.0d0+costtab(i+1))
+        cosfac=dsqrt(cosfac2)
+        sinfac2=0.5d0/(1.0d0-costtab(i+1))
+        sinfac=dsqrt(sinfac2)
+        it=itype(i)
+        if (it.ne.10 .and. itype(i).ne.ntyp1) then
+c
+C  Compute the axes of tghe local cartesian coordinates system; store in
+c   x_prime, y_prime and z_prime 
+c
+        do j=1,3
+          x_prime(j) = 0.00
+          y_prime(j) = 0.00
+          z_prime(j) = 0.00
+        enddo
+        do j = 1,3
+          x_prime(j) = (dc_norm(j,i) - dc_norm(j,i-1))*cosfac
+          y_prime(j) = (dc_norm(j,i) + dc_norm(j,i-1))*sinfac
+        enddo
+        call vecpr(x_prime,y_prime,z_prime)
+c
+C Transform the unit vector of the ith side-chain centroid, dC_norm(*,i),
+C to local coordinate system. Store in xx, yy, zz.
+c
+        xx=0.0d0
+        yy=0.0d0
+        zz=0.0d0
+        do j = 1,3
+          xx = xx + x_prime(j)*dc_norm(j,i+nres)
+          yy = yy + y_prime(j)*dc_norm(j,i+nres)
+          zz = zz + z_prime(j)*dc_norm(j,i+nres)
+        enddo
+
+        xxref(i)=xx
+        yyref(i)=yy
+        zzref(i)=zz
+        else
+        xxref(i)=0.0d0
+        yyref(i)=0.0d0
+        zzref(i)=0.0d0
+        endif
+      enddo
+      if (lprn) then
+        do i=2,nres
+          iti=itype(i)
+          if(me.eq.king.or..not.out1file)
+     &     write (iout,'(a3,i4,3f10.5)') restyp(iti),i,xxref(i),
+     &      yyref(i),zzref(i)
+        enddo
+      endif
+      return
+      end
 c---------------------------------------------------------------------------
       subroutine sccenter(ires,nscat,sccor)
       implicit none