File: project_fc.f

package info (click to toggle)
aces3 3.0.6-7
  • links: PTS, VCS
  • area: main
  • in suites: wheezy
  • size: 82,460 kB
  • sloc: fortran: 225,647; ansic: 20,413; cpp: 4,349; makefile: 953; sh: 137
file content (228 lines) | stat: -rw-r--r-- 6,907 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
      Subroutine projec_FC(Coords, Hess, AtmMass, Grad, Hess_project, 
     &                     Work, Threshold, Nreals, Move_CMass, 
     &                     Proj_rots, Proj_grads)
C
      Implicit Double Precision (A-H, O-Z)
C
      Double Precision MO_Inertia(3,3)
      Logical Proj_rots, Proj_grads, Move_CMass
      Integer L_Index, R_Index
C
      Dimension Coords(3*Nreals), Grad(3*Nreals), AtmMass(Nreals), 
     &          Hess(3*Nreals, 3*Nreals), Hess_Project(3*Nreals, 
     &          3*Nreals), Work(3*Nreals, 3*Nreals), Asym_ten(3,3,3),
     &          CM(3), Center_Mass(3)
 
      Data Asym_ten/ 0.0D+00,  0.0D+00,  0.0D+00,
     &               0.0D+00,  0.0D+00, -1.0D+00,
     &               0.0D+00,  1.0D+00,  0.0D+00,
     &               0.0D+00,  0.0D+00,  1.0D+00,
     &               0.0D+00,  0.0D+00,  0.0D+00,
     &              -1.0D+00,  0.0D+00,  0.0D+00,
     &               0.0D+00, -1.0D+00,  0.0D+00,
     &               1.0D+00,  0.0D+00,  0.0D+00, 
     &               0.0D+00,  0.0D+00,  0.0D+00  /
C
C Normalize the gradinet vector. 
C
      If (Proj_grads) Then
         Grad_sqr = Ddot(3*Nreals, Grad, 1, Grad, 1)   
         If (Grad_sqr .GT. Threshold)  Call Normal(Grad, 3*Nreals)
      Endif
C

      Write(6,*)
      Write(6,"(a)") "@-Projec_FC At enetry the Coords and grads"
      Write(6, "(3F17.13)") (Coords(i),i=1,3*Nreals)
      Write(6,*)
      Write(6, "(3F17.13)") (Grad(i),i=1,3*Nreals)

C
C Move to the center of Mass coordinate system.
C
      If (Move_CMass) Then
         Call Dzero(CM, 3)
         TotMass = 0.0D0
         Do Iatoms = 1, Nreals
            Ioff = 3*(Iatoms - 1)
            TotMass = TotMass + AtmMass(Iatoms)
            Do Ixyz = 1, 3
               CM(Ixyz)   = Coords(Ioff + Ixyz)*AtmMass(Iatoms) +
     &                      CM(Ixyz)
            Enddo
         Enddo
C
         Do Ixyz = 1, 3
            If (TotMass .GT. Threshold) Then
                Center_Mass(IxYz) = CM(Ixyz)/TotMass
            Else
                Write(6, "(a)") "@-Project_FC Zero total mass"
                Call Errex
            Endif
         Enddo
C         
         Do Iatoms = 1, Nreals
            Ioff  = 3*(Iatoms - 1)
            Do Ixyz = 1, 3
               Coords(Ioff + Ixyz) = Coords(Ioff + Ixyz) - 
     &                               Center_Mass(Ixyz)
            Enddo
         Enddo
      Endif
C

      Write(6,*)
      Write(6,"(a)") "@-Projec_FC The center of mass Coords"
      Write(6, "(3F17.13)") (Coords(i),i=1,3*Nreals)

C
C compute the inverse of inertia matrix.
c 
      Call CalMOI(NReals, Coords, AtmMass, MO_Inertia)

      Write(6,*)
      Write(6,"(a)") "@-Projec_FC the moments of inertia matrix"
      Call output(MO_Inertia, 1, 3, 1, 3, 3, 3,1)


      Call Minv(MO_Inertia, 3, 3, Work, Det, 1.0D-8, 0, 1)
C
C Mass weigh the incomming Cartesians 
C
      Ioff = 1
      Do Iatoms = 1, Nreals
          Sqrtmass = Dsqrt(AtmMass(Iatoms))
          Do Ixyz = 1, 3
             Coords(Ioff) = Coords(Ioff)*Sqrtmass
             Ioff = Ioff + 1 
          Enddo
      Enddo    
C

      Write(6,*)
      Write(6,"(a)") "@-Projec_FC Mass weighted center of mass Coords"
      Write(6, "(3F17.13)") (Coords(i),i=1,3*Nreals)
      Write(6,*)
      Write(6,"(a)") "@-Projec_FC Inv. of the moms. of inertia matrix"
      Call output(MO_Inertia, 1, 3, 1, 3, 3, 3,1)

C
C Build the Hessian projector See, Miller, Handy and Adams, JCP, 
C 72, 99, (1980).
C
      Do Iatms = 1, Nreals
         Ioff = 3*(Iatms - 1)
CSSS         Itmp = Max(3*(Iatms - 1), 6*(Iatms -1)-3*Nreals)

         Do Jatms = 1, Iatms
          Joff = 3*(Jatms - 1)
CSSS          Jtmp = Max(3*(Jatms - 1), 6*(Jatms -1)-3*Nreals)





            Do Iz = 1, 3
               L_index = Ioff + Iz
CSSS               Itmp_L = Itmp + Iz
               Jz_end = 3
               If (Iatms .EQ. Jatms)  Jz_end = Iz
C
                 Do Jz = 1, Jz_end
                    R_index = Joff + Jz
CSSS                    Jtmp_R =  Jtmp + JZ
                    Sum = 0.0D0
C
                    If (Proj_rots) Then 

                        Do Ix = 1, 3
 
                           Do Iy = 1, 3
C
                              If (Asym_ten(Ix, Iy, Iz) .NE. 0.0D0) 
     &                            Then
                                   
                                   Do Jx = 1, 3

                                      Do Jy = 1, 3
                                         If (Asym_ten(Jx, Jy, Jz) .NE. 
     &                                       0.0D0) Then
                                             Sum = Sum + 
     &                                             Coords(Ioff + IY)*
     &                                             Coords(Joff + Jy)*
     &                                            MO_inertia(Ix, Jx)*
     &                                          Asym_ten(Ix, Iy, Iz)*
     &                                          Asym_ten(Jx, Jy, Jz)
                                         Endif
                                      Enddo
                                   Enddo 
                              Endif
                           Enddo
                        Enddo
                    Endif
C 
                    Hess_project(L_index, R_index) = Sum
C
                    If (Proj_grads) Hess_project(L_index, R_index) = 
     &                              Hess_project(L_index, R_index) + 
     &                              Grad(R_Index)*Grad(L_index)
C
                    If (Iz .EQ. Jz) Then



                        Hess_project(L_index, R_index) =
     &                  Hess_project(L_index, R_index) + 
     &                  Dsqrt(AtmMass(Iatms)*AtmMass(Jatms))/
     &                  Totmass
                    Endif
C
                 Enddo 
            Enddo 
C
         Enddo
      Enddo
C






C
C Build the projector (I - P)
C
      Do Jdeg = 1, 3*Nreals
          Do Ideg = 1, Jdeg
              Hess_project(Jdeg, Ideg) = -Hess_project(Jdeg, Ideg)
             If (Ideg .Eq. Jdeg)  Hess_project(Jdeg, Ideg) = 1.0D0 +
     &                            Hess_project(Jdeg, Ideg)   
             If (Dabs(Hess_project(Jdeg, Ideg)) .LT. Threshold) 
     &                            Hess_project(Jdeg, Ideg) = 0.0D0
              Hess_project(Ideg, Jdeg) = Hess_project(Jdeg, Ideg)
          Enddo
      Enddo
C 
C

      NX = 3*Nreals
      Write(6,*)
      Write(6,"(a)") "The Hessian projector:(I-P)"
      CALL OUTPUT(Hess_project, 1, NX, 1, Nx, Nx, Nx, 1)
      

C
C
C Project the Hessian (I - P)H(I - P)
C
      Call Xgemm("N", "N", 3*Nreals, 3*Nreals, 3*Nreals, 1.0D0, 
     &            Hess_project, 3*Nreals, Hess, 3*Nreals, 0.0D0, 
     &            Work, 3*Nreals)
C
      Call Xgemm("N", "N", 3*Nreals, 3*Nreals, 3*Nreals, 1.0D0, 
     &            Work, 3*Nreals, Hess_project, 3*Nreals, 0.0D0, 
     &            Hess, 3*Nreals)
C
      Return
      End