Statistiques
| Révision :

root / src / CalcRmsd.f90

Historique | Voir | Annoter | Télécharger (8,54 ko)

1

    
2
! This subroutine calculates RMSD (Root mean square deviation) using quaternions.
3
! It is baNatsed on the F90 routine bu E. Coutsias
4
! http://www.math.unm.edu/~vageli/homepage.html
5
! I (PFL) have just translated it, and I have changed the diagonalization
6
! subroutine.
7
! I also made some changes to make it suitable for Cart package.
8
!----------------------------------------------------------------------
9
!----------------------------------------------------------------------
10
! Copyright (C) 2004, 2005 Chaok Seok, Evangelos Coutsias and Ken Dill
11
!      UCSF, Univeristy of New Mexico, Seoul National University
12
! Witten by Chaok Seok and Evangelos Coutsias 2004.
13
                                                                                
14
! This library is free software; you can redistribute it and/or
15
! modify it under the terms of the GNU Lesser General Public
16
! License as published by the Free Software Foundation; either
17
! version 2.1 of the License, or (at your option) any later version.
18
!
19
                                                                                
20
! This library is distributed in the hope that it will be useful,
21
! but WITHOUT ANY WARRANTY; without even the implied warranty of
22
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
23
! Lesser General Public License for more details.
24
!
25
                                                                                
26
! You should have received a copy of the GNU Lesser General Public
27
! License along with this library; if not, write to the Free Software
28
! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA
29
!----------------------------------------------------------------------------
30
!----------------------------------------------------------------------
31
!  Copyright 2003-2014 Ecole Normale Supérieure de Lyon, 
32
!  Centre National de la Recherche Scientifique,
33
!  Université Claude Bernard Lyon 1. All rights reserved.
34
!
35
!  This work is registered with the Agency for the Protection of Programs 
36
!  as IDDN.FR.001.100009.000.S.P.2014.000.30625
37
!
38
!  Authors: P. Fleurat-Lessard, P. Dayal
39
!  Contact: optnpath@gmail.com
40
!
41
! This file is part of "Opt'n Path".
42
!
43
!  "Opt'n Path" is free software: you can redistribute it and/or modify
44
!  it under the terms of the GNU Affero General Public License as
45
!  published by the Free Software Foundation, either version 3 of the License,
46
!  or (at your option) any later version.
47
!
48
!  "Opt'n Path" is distributed in the hope that it will be useful,
49
!  but WITHOUT ANY WARRANTY; without even the implied warranty of
50
!
51
!  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
52
!  GNU Affero General Public License for more details.
53
!
54
!  You should have received a copy of the GNU Affero General Public License
55
!  along with "Opt'n Path". If not, see <http://www.gnu.org/licenses/>.
56
!
57
! Contact The Office of Technology Licensing, valorisation@ens-lyon.fr,
58
! for commercial licensing opportunities.
59
!----------------------------------------------------------------------
60

    
61
      subroutine CalcRmsd(Nat,x0,y0,z0, x2,y2,z2,U,rmsd,FRot,FAlign)
62
!-----------------------------------------------------------------------
63
!  This subroutine calculates the least square rmsd of two coordinate
64
!  sets coord1(3,n) and coord2(3,n) using a method based on quaternion.
65
!  It then calculate the rotation matrix U and the centers of coord, and uses
66
! them to align the two molecules.
67
! x0(Nat), y0(Nat), z0(Nat) are the coordinates of the reference geometry.
68
! x2(Nat), y2(Nat), z2(Nat) are the coordinates of the geometry which is to
69
! be aligned with the reference geometry. 
70
! The rotation matrix U has INTENT(OUT) in subroutine rotation_matrix(...),
71
! which is called in this CalcRmsd(...).
72
! rmsd: the root mean square deviation and calculated in this subroutine
73
! CalcRmsd(...).
74
! FRot: if .TRUE. the rotation matrix is calculated. This is the matrix that rotates
75
! the first molecule onto the second one.
76
! FAlgin: if .TRUE. then the second molecule is aligned on the first one.
77
!-----------------------------------------------------------------------
78

    
79
        use VarTypes
80

    
81
      IMPLICIT NONE
82

    
83
      INTEGER(KINT) :: Nat
84
      real(KREAL) ::  x0(Nat),y0(Nat),z0(Nat)
85
      real(KREAL) ::  x2(Nat),y2(Nat),z2(Nat)
86
      real(KREAL) ::  U(3,3), rmsd
87
      LOGICAL  FRot,FAlign,Debug
88

    
89
      REAL(KREAL),ALLOCATABLE :: Coord1(:,:), Coord2(:,:) ! 3,Nat
90
      real(KREAL) ::  x0c1,y0c1,z0c1, xc2,yc2,zc2
91
      
92

    
93
      INTEGER(KINT) :: i, j, ia
94
      real(KREAL) :: x_norm, y_norm, lambda
95
      real(KREAL) :: Rmatrix(3,3)
96
      real(KREAL) :: S(4, 4)
97
      real(KREAL) :: EigVec(4, 4), EigVal(4)
98

    
99
      INTERFACE
100
       function valid(string) result (isValid)
101
         CHARACTER(*), intent(in) :: string
102
         logical                  :: isValid
103
         END function VALID
104
      END INTERFACE
105

    
106

    
107
      debug=valid('CalcRmsd').OR.valid('align').OR.valid('alignpartial')
108

    
109
      if (debug) Call Header("Entering CalcRmsd")
110

    
111
      if (debug) THEN
112
         WRITE(*,*) "NAt=",NAt
113
         WRITE(*,*) x0
114
         WRITE(*,*) y0
115
         WRITE(*,*) z0
116
      END IF
117

    
118
      ALLOCATE(Coord1(3,Nat),Coord2(3,Nat))
119
  ! calculate the barycenters, centroidal coordinates, and the norms
120
      x_norm = 0.0d0
121
      y_norm = 0.0d0
122
      x0c1=0.
123
      y0c1=0.
124
      z0c1=0.
125
      xc2=0.
126
      yc2=0.
127
      zc2=0.
128
      do ia=1,Nat
129
         if (debug) WRITE(*,'(A,I5,4(1X,F10.4))') 'ia...',ia,x0(ia), &
130
                             x2(ia),x0c1,xc2
131
         x0c1=x0c1+x0(ia)
132
         xc2=xc2+x2(ia)
133
         y0c1=y0c1+y0(ia)
134
         yc2=yc2+y2(ia)
135
         z0c1=z0c1+z0(ia)
136
         zc2=zc2+z2(ia)
137
         if (debug) WRITE(*,'(A,I5,4(1X,F10.4))') 'ia...',ia,x0(ia), &
138
                             x2(ia),x0c1,xc2
139
      END DO
140
      x0c1=x0c1/dble(Nat)
141
      y0c1=y0c1/dble(Nat)
142
      z0c1=z0c1/dble(Nat)
143
      xc2=xc2/dble(Nat)
144
      yc2=yc2/dble(Nat)
145
      zc2=zc2/dble(Nat)
146

    
147
      IF (debug) WRITE(*,'(1X,A,3(1X,F10.4))') 'Center1',x0c1,y0c1,z0c1
148
      IF (debug) WRITE(*,'(1X,A,3(1X,F10.4))') 'Center2',xc2,yc2,zc2
149
      do i=1,Nat
150
         Coord1(1,i)=x0(i)-x0c1
151
         Coord1(2,i)=y0(i)-y0c1
152
         Coord1(3,i)=z0(i)-z0c1
153
         Coord2(1,i)=x2(i)-xc2
154
         Coord2(2,i)=y2(i)-yc2
155
         Coord2(3,i)=z2(i)-zc2
156
         x_norm=x_norm+Coord1(1,i)**2+Coord1(2,i)**2+Coord1(3,i)**2
157
         y_norm=y_norm+Coord2(1,i)**2+Coord2(2,i)**2+Coord2(3,i)**2
158
      end do
159

    
160
      IF (debug) THEN
161
         WRITE(*,*) "R matrix"
162
         DO I=1,3
163
            WRITE(*,*) (RMatrix(I,j),j=1,3)
164
         END DO
165
      END IF
166

    
167
  ! calculate the R matrix
168
      do i = 1, 3
169
         do j = 1, 3
170
            Rmatrix(i,j)=0.
171
            do ia=1,Nat
172
               Rmatrix(i,j) = Rmatrix(i,j)+Coord1(i,ia)*Coord2(j,ia)
173
            END DO
174
         end do
175
      end do
176

    
177
      IF (debug) THEN
178
         WRITE(*,*) "R matrix"
179
         DO I=1,3
180
            WRITE(*,*) (RMatrix(I,j),j=1,3)
181
         END DO
182
      END IF
183

    
184

    
185
  ! S matrix
186
      S(1, 1) = Rmatrix(1, 1) + Rmatrix(2, 2) + Rmatrix(3, 3)
187
      S(2, 1) = Rmatrix(2, 3) - Rmatrix(3, 2)
188
      S(3, 1) = Rmatrix(3, 1) - Rmatrix(1, 3)
189
      S(4, 1) = Rmatrix(1, 2) - Rmatrix(2, 1)
190

    
191
      S(1, 2) = S(2, 1)
192
      S(2, 2) = Rmatrix(1, 1) - Rmatrix(2, 2) - Rmatrix(3, 3)
193
      S(3, 2) = Rmatrix(1, 2) + Rmatrix(2, 1)
194
      S(4, 2) = Rmatrix(1, 3) + Rmatrix(3, 1)
195

    
196
      S(1, 3) = S(3, 1)
197
      S(2, 3) = S(3, 2)
198
      S(3, 3) =-Rmatrix(1, 1) + Rmatrix(2, 2) - Rmatrix(3, 3)
199
      S(4, 3) = Rmatrix(2, 3) + Rmatrix(3, 2)
200

    
201
      S(1, 4) = S(4, 1)
202
      S(2, 4) = S(4, 2)
203
      S(3, 4) = S(4, 3)
204
      S(4, 4) =-Rmatrix(1, 1) - Rmatrix(2, 2) + Rmatrix(3, 3) 
205

    
206

    
207
! PFL : I use my usual Jacobi diagonalisation
208
  ! Calculate eigenvalues and eigenvectors, and 
209
  ! take the maximum eigenvalue lambda and the corresponding eigenvector q.
210

    
211
      IF (debug) THEN
212
         WRITE(*,*) "S matrix"
213
         DO I=1,4
214
            WRITE(*,*) (S(I,j),j=1,4)
215
         END DO
216
      END IF
217

    
218
      Call Jacobi(S,4,EigVal,EigVec,4)
219

    
220
      Call Trie(4,EigVal,EigVec,4)
221

    
222
      Lambda=EigVal(4)
223

    
224
! RMS Deviation
225
       rmsd=sqrt(max(0.0d0,((x_norm+y_norm)-2.0d0*lambda))/dble(Nat))
226
! The subroutine rotation_matrix(...) constructs rotation matrix U as output
227
! from the input quaternion EigVec(1,4).
228
      if (FRot.OR.FAlign) Call rotation_matrix(EigVec(1,4),U)
229
      IF (FAlign) THEN
230
         DO I=1,Nat
231
            x2(i)=Coord2(1,i)*U(1,1)+Coord2(2,i)*U(2,1) &
232
                +Coord2(3,i)*U(3,1) +x0c1
233
            y2(i)=Coord2(1,i)*U(1,2)+Coord2(2,i)*U(2,2) &
234
                +Coord2(3,i)*U(3,2) +y0c1
235
            z2(i)=Coord2(1,i)*U(1,3)+Coord2(2,i)*U(2,3) &
236
                +Coord2(3,i)*U(3,3) +z0c1
237
         END DO
238
      END IF
239

    
240
      DEALLOCATE(Coord1,Coord2)
241

    
242
    END subroutine CalcRmsd