Statistiques
| Révision :

root / src / CalcRmsd.f90 @ 1

Historique | Voir | Annoter | Télécharger (6,88 ko)

1 1 equemene
2 1 equemene
! This subroutine calculates RMSD (Root mean square deviation) using quaternions.
3 1 equemene
! It is baNatsed on the F90 routine bu E. Coutsias
4 1 equemene
! http://www.math.unm.edu/~vageli/homepage.html
5 1 equemene
! I (PFL) have just translated it, and I have changed the diagonalization
6 1 equemene
! subroutine.
7 1 equemene
! I also made some changes to make it suitable for Cart package.
8 1 equemene
!----------------------------------------------------------------------
9 1 equemene
!----------------------------------------------------------------------
10 1 equemene
! Copyright (C) 2004, 2005 Chaok Seok, Evangelos Coutsias and Ken Dill
11 1 equemene
!      UCSF, Univeristy of New Mexico, Seoul National University
12 1 equemene
! Witten by Chaok Seok and Evangelos Coutsias 2004.
13 1 equemene
14 1 equemene
! This library is free software; you can redistribute it and/or
15 1 equemene
! modify it under the terms of the GNU Lesser General Public
16 1 equemene
! License as published by the Free Software Foundation; either
17 1 equemene
! version 2.1 of the License, or (at your option) any later version.
18 1 equemene
!
19 1 equemene
20 1 equemene
! This library is distributed in the hope that it will be useful,
21 1 equemene
! but WITHOUT ANY WARRANTY; without even the implied warranty of
22 1 equemene
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
23 1 equemene
! Lesser General Public License for more details.
24 1 equemene
!
25 1 equemene
26 1 equemene
! You should have received a copy of the GNU Lesser General Public
27 1 equemene
! License along with this library; if not, write to the Free Software
28 1 equemene
! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA
29 1 equemene
!----------------------------------------------------------------------------
30 1 equemene
31 1 equemene
      subroutine CalcRmsd(Nat,x0,y0,z0, x2,y2,z2,U,rmsd,FRot,FAlign)
32 1 equemene
!-----------------------------------------------------------------------
33 1 equemene
!  This subroutine calculates the least square rmsd of two coordinate
34 1 equemene
!  sets coord1(3,n) and coord2(3,n) using a method based on quaternion.
35 1 equemene
!  It then calculate the rotation matrix U and the centers of coord, and uses
36 1 equemene
! them to align the two molecules.
37 1 equemene
! x0(Nat), y0(Nat), z0(Nat) are the coordinates of the reference geometry.
38 1 equemene
! x2(Nat), y2(Nat), z2(Nat) are the coordinates of the geometry which is to
39 1 equemene
! be aligned with the reference geometry.
40 1 equemene
! The rotation matrix U has INTENT(OUT) in subroutine rotation_matrix(...),
41 1 equemene
! which is called in this CalcRmsd(...).
42 1 equemene
! rmsd: the root mean square deviation and calculated in this subroutine
43 1 equemene
! CalcRmsd(...).
44 1 equemene
! FRot: if .TRUE. the rotation matrix is calculated. This is the matrix that rotates
45 1 equemene
! the first molecule onto the second one.
46 1 equemene
! FAlgin: if .TRUE. then the second molecule is aligned on the first one.
47 1 equemene
!-----------------------------------------------------------------------
48 1 equemene
49 1 equemene
        use VarTypes
50 1 equemene
51 1 equemene
      IMPLICIT NONE
52 1 equemene
53 1 equemene
      INTEGER(KINT) :: Nat
54 1 equemene
      real(KREAL) ::  x0(Nat),y0(Nat),z0(Nat)
55 1 equemene
      real(KREAL) ::  x2(Nat),y2(Nat),z2(Nat)
56 1 equemene
      real(KREAL) ::  U(3,3), rmsd
57 1 equemene
      LOGICAL  FRot,FAlign,Debug
58 1 equemene
59 1 equemene
      REAL(KREAL) :: Coord1(3,Nat), Coord2(3,Nat)
60 1 equemene
      real(KREAL) ::  x0c1,y0c1,z0c1, xc2,yc2,zc2
61 1 equemene
62 1 equemene
63 1 equemene
      INTEGER(KINT) :: i, j, nd, ia
64 1 equemene
      real(KREAL) :: x_norm, y_norm, lambda
65 1 equemene
      real(KREAL) :: Rmatrix(3,3)
66 1 equemene
      real(KREAL) :: S(4,4), dS(4,4), q(4)
67 1 equemene
      real(KREAL) :: tmp(3),tmp1(4),tmp2(4),EigVec(4,4), EigVal(4)
68 1 equemene
69 1 equemene
      INTERFACE
70 1 equemene
       function valid(string) result (isValid)
71 1 equemene
         CHARACTER(*), intent(in) :: string
72 1 equemene
         logical                  :: isValid
73 1 equemene
         END function VALID
74 1 equemene
      END INTERFACE
75 1 equemene
76 1 equemene
77 1 equemene
      debug=valid('CalcRmsd').OR.valid('align').OR.valid('alignpartial')
78 1 equemene
79 1 equemene
  ! calculate the barycenters, centroidal coordinates, and the norms
80 1 equemene
      x_norm = 0.0d0
81 1 equemene
      y_norm = 0.0d0
82 1 equemene
      x0c1=0.
83 1 equemene
      y0c1=0.
84 1 equemene
      z0c1=0.
85 1 equemene
      xc2=0.
86 1 equemene
      yc2=0.
87 1 equemene
      zc2=0.
88 1 equemene
      do ia=1,Nat
89 1 equemene
         x0c1=x0c1+x0(ia)
90 1 equemene
         xc2=xc2+x2(ia)
91 1 equemene
         y0c1=y0c1+y0(ia)
92 1 equemene
         yc2=yc2+y2(ia)
93 1 equemene
         z0c1=z0c1+z0(ia)
94 1 equemene
         zc2=zc2+z2(ia)
95 1 equemene
!         if (debug) WRITE(*,'(A,I5,4(1X,F10.4))') 'ia...',ia,x0(ia),
96 1 equemene
!     &                        x2(ia),x0c1,xc2
97 1 equemene
      END DO
98 1 equemene
      x0c1=x0c1/dble(Nat)
99 1 equemene
      y0c1=y0c1/dble(Nat)
100 1 equemene
      z0c1=z0c1/dble(Nat)
101 1 equemene
      xc2=xc2/dble(Nat)
102 1 equemene
      yc2=yc2/dble(Nat)
103 1 equemene
      zc2=zc2/dble(Nat)
104 1 equemene
105 1 equemene
      IF (debug) WRITE(*,'(1X,A,3(1X,F10.4))') 'Center1',x0c1,y0c1,z0c1
106 1 equemene
      IF (debug) WRITE(*,'(1X,A,3(1X,F10.4))') 'Center2',xc2,yc2,zc2
107 1 equemene
      do i=1,Nat
108 1 equemene
         Coord1(1,i)=x0(i)-x0c1
109 1 equemene
         Coord1(2,i)=y0(i)-y0c1
110 1 equemene
         Coord1(3,i)=z0(i)-z0c1
111 1 equemene
         Coord2(1,i)=x2(i)-xc2
112 1 equemene
         Coord2(2,i)=y2(i)-yc2
113 1 equemene
         Coord2(3,i)=z2(i)-zc2
114 1 equemene
         x_norm=x_norm+Coord1(1,i)**2+Coord1(2,i)**2+Coord1(3,i)**2
115 1 equemene
         y_norm=y_norm+Coord2(1,i)**2+Coord2(2,i)**2+Coord2(3,i)**2
116 1 equemene
      end do
117 1 equemene
118 1 equemene
      IF (debug) THEN
119 1 equemene
         WRITE(*,*) "R matrix"
120 1 equemene
         DO I=1,3
121 1 equemene
            WRITE(*,*) (RMatrix(I,j),j=1,3)
122 1 equemene
         END DO
123 1 equemene
      END IF
124 1 equemene
125 1 equemene
  ! calculate the R matrix
126 1 equemene
      do i = 1, 3
127 1 equemene
         do j = 1, 3
128 1 equemene
            Rmatrix(i,j)=0.
129 1 equemene
            do ia=1,Nat
130 1 equemene
               Rmatrix(i,j) = Rmatrix(i,j)+Coord1(i,ia)*Coord2(j,ia)
131 1 equemene
            END DO
132 1 equemene
         end do
133 1 equemene
      end do
134 1 equemene
135 1 equemene
      IF (debug) THEN
136 1 equemene
         WRITE(*,*) "R matrix"
137 1 equemene
         DO I=1,3
138 1 equemene
            WRITE(*,*) (RMatrix(I,j),j=1,3)
139 1 equemene
         END DO
140 1 equemene
      END IF
141 1 equemene
142 1 equemene
143 1 equemene
  ! S matrix
144 1 equemene
      S(1, 1) = Rmatrix(1, 1) + Rmatrix(2, 2) + Rmatrix(3, 3)
145 1 equemene
      S(2, 1) = Rmatrix(2, 3) - Rmatrix(3, 2)
146 1 equemene
      S(3, 1) = Rmatrix(3, 1) - Rmatrix(1, 3)
147 1 equemene
      S(4, 1) = Rmatrix(1, 2) - Rmatrix(2, 1)
148 1 equemene
149 1 equemene
      S(1, 2) = S(2, 1)
150 1 equemene
      S(2, 2) = Rmatrix(1, 1) - Rmatrix(2, 2) - Rmatrix(3, 3)
151 1 equemene
      S(3, 2) = Rmatrix(1, 2) + Rmatrix(2, 1)
152 1 equemene
      S(4, 2) = Rmatrix(1, 3) + Rmatrix(3, 1)
153 1 equemene
154 1 equemene
      S(1, 3) = S(3, 1)
155 1 equemene
      S(2, 3) = S(3, 2)
156 1 equemene
      S(3, 3) =-Rmatrix(1, 1) + Rmatrix(2, 2) - Rmatrix(3, 3)
157 1 equemene
      S(4, 3) = Rmatrix(2, 3) + Rmatrix(3, 2)
158 1 equemene
159 1 equemene
      S(1, 4) = S(4, 1)
160 1 equemene
      S(2, 4) = S(4, 2)
161 1 equemene
      S(3, 4) = S(4, 3)
162 1 equemene
      S(4, 4) =-Rmatrix(1, 1) - Rmatrix(2, 2) + Rmatrix(3, 3)
163 1 equemene
164 1 equemene
165 1 equemene
! PFL : I use my usual Jacobi diagonalisation
166 1 equemene
  ! Calculate eigenvalues and eigenvectors, and
167 1 equemene
  ! take the maximum eigenvalue lambda and the corresponding eigenvector q.
168 1 equemene
169 1 equemene
      IF (debug) THEN
170 1 equemene
         WRITE(*,*) "S matrix"
171 1 equemene
         DO I=1,4
172 1 equemene
            WRITE(*,*) (S(I,j),j=1,4)
173 1 equemene
         END DO
174 1 equemene
      END IF
175 1 equemene
176 1 equemene
      Call Jacobi(S,4,EigVal,EigVec,4)
177 1 equemene
178 1 equemene
      Call Trie(4,EigVal,EigVec,4)
179 1 equemene
180 1 equemene
      Lambda=EigVal(4)
181 1 equemene
182 1 equemene
! RMS Deviation
183 1 equemene
       rmsd=sqrt(max(0.0d0,((x_norm+y_norm)-2.0d0*lambda))/dble(Nat))
184 1 equemene
! The subroutine rotation_matrix(...) constructs rotation matrix U as output
185 1 equemene
! from the input quaternion EigVec(1,4).
186 1 equemene
      if (FRot.OR.FAlign) Call rotation_matrix(EigVec(1,4),U)
187 1 equemene
      IF (FAlign) THEN
188 1 equemene
         DO I=1,Nat
189 1 equemene
            x2(i)=Coord2(1,i)*U(1,1)+Coord2(2,i)*U(2,1) &
190 1 equemene
                +Coord2(3,i)*U(3,1) +x0c1
191 1 equemene
            y2(i)=Coord2(1,i)*U(1,2)+Coord2(2,i)*U(2,2) &
192 1 equemene
                +Coord2(3,i)*U(3,2) +y0c1
193 1 equemene
            z2(i)=Coord2(1,i)*U(1,3)+Coord2(2,i)*U(2,3) &
194 1 equemene
                +Coord2(3,i)*U(3,3) +z0c1
195 1 equemene
         END DO
196 1 equemene
      END IF
197 1 equemene
198 1 equemene
      END