root / src / Calc_baker.f90
Historique | Voir | Annoter | Télécharger (42,05 ko)
1 |
SUBROUTINE Calc_baker(atome,r_cov,fact,frozen,IGeomRef) |
---|---|
2 |
! |
3 |
! This subroutine analyses a geometry to construct the baker |
4 |
! delocalized internal coordinates |
5 |
! v1.0 |
6 |
! We use only one geometry |
7 |
! |
8 |
|
9 |
!---------------------------------------------------------------------- |
10 |
! Copyright 2003-2014 Ecole Normale Supérieure de Lyon, |
11 |
! Centre National de la Recherche Scientifique, |
12 |
! Université Claude Bernard Lyon 1. All rights reserved. |
13 |
! |
14 |
! This work is registered with the Agency for the Protection of Programs |
15 |
! as IDDN.FR.001.100009.000.S.P.2014.000.30625 |
16 |
! |
17 |
! Authors: P. Fleurat-Lessard, P. Dayal |
18 |
! Contact: optnpath@gmail.com |
19 |
! |
20 |
! This file is part of "Opt'n Path". |
21 |
! |
22 |
! "Opt'n Path" is free software: you can redistribute it and/or modify |
23 |
! it under the terms of the GNU Affero General Public License as |
24 |
! published by the Free Software Foundation, either version 3 of the License, |
25 |
! or (at your option) any later version. |
26 |
! |
27 |
! "Opt'n Path" is distributed in the hope that it will be useful, |
28 |
! but WITHOUT ANY WARRANTY; without even the implied warranty of |
29 |
! |
30 |
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
31 |
! GNU Affero General Public License for more details. |
32 |
! |
33 |
! You should have received a copy of the GNU Affero General Public License |
34 |
! along with "Opt'n Path". If not, see <http://www.gnu.org/licenses/>. |
35 |
! |
36 |
! Contact The Office of Technology Licensing, valorisation@ens-lyon.fr, |
37 |
! for commercial licensing opportunities. |
38 |
!---------------------------------------------------------------------- |
39 |
Use Path_module, only : max_Z,NMaxL,Pi,Massat,BMat_BakerT, & |
40 |
Nat,NCoord,XyzGeomI,NGeomI,NGeomF,UMat,NPrim,BTransInv, & |
41 |
IntCoordI,Coordinate,CurrentCoord,ScanCoord,BprimT,BBT, & |
42 |
BBT_inv,Xprimitive,XPrimitiveF,UMat_local, & |
43 |
BTransInv_local, UMatF,BTransInvF, & |
44 |
indzmat, dzdc,renum,order,OrderInv |
45 |
! BMat_BakerT(3*Nat,NCoord), NCoord=3*Nat or NFree=3*Nat-6-Symmetry_elimination |
46 |
! depending upon the coordinate choice. IntCoordI(NGeomI,NCoord) where |
47 |
! UMat(NPrim,NCoord), NCoord number of vectors in UMat matrix, i.e. NCoord |
48 |
! Baker coordinates. NPrim is the number of primitive internal coordinates. |
49 |
|
50 |
Use Io_module |
51 |
IMPLICIT NONE |
52 |
|
53 |
! IGeom: Geometry index, IGeomRef: Reference Geometry. |
54 |
INTEGER(KINT), INTENT(IN) :: atome(Nat),IGeomRef |
55 |
REAL(KREAL), INTENT(IN) :: fact |
56 |
REAL(KREAL), INTENT(IN) :: r_cov(0:Max_Z) |
57 |
! Frozen contains the indices of frozen atoms |
58 |
INTEGER(KINT), INTENT(IN) :: Frozen(*) |
59 |
|
60 |
INTEGER(KINT), ALLOCATABLE :: LIAISONS(:,:) ! (Nat,0:NMaxL) |
61 |
REAL(KREAL), ALLOCATABLE :: Geom(:,:) !(3,Nat) |
62 |
! NPrim is the number of primitive coordinates and NCoord is the number |
63 |
! of internal coordinates. BMat is actually (NPrim,3*Nat). |
64 |
REAL(KREAL), ALLOCATABLE :: GMat(:,:) !(NPrim,NPrim) |
65 |
! EigVec(..) contains ALL eigevectors of BMat times BprimT, NOT only Baker Coordinate vectors. |
66 |
REAL(KREAL), ALLOCATABLE :: EigVec(:,:), EigVal(:) ! EigVec(NPrim,NPrim) |
67 |
REAL(KREAL), ALLOCATABLE :: x(:), y(:), z(:) |
68 |
REAL(KREAL), ALLOCATABLE :: x_refGeom(:), y_refGeom(:), z_refGeom(:) |
69 |
!REAL(KREAL), ALLOCATABLE :: V(:,:) ! (3*Nat,3*Nat) |
70 |
!REAL(KREAL), ALLOCATABLE :: P(:) ! Used in Baker Coordinates: Matrix Inversion |
71 |
REAL(KREAL), ALLOCATABLE :: UMat_tmp(:,:) |
72 |
INTEGER(KINT) :: NbBonds, NbAngles, NbDihedrals, IGeom |
73 |
|
74 |
real(KREAL) :: vx, vy, vz, Norm |
75 |
real(KREAL) :: vx1,vy1,vz1,norm1 |
76 |
real(KREAL) :: vx2,vy2,vz2,norm2 |
77 |
real(KREAL) :: vx3,vy3,vz3,norm3 |
78 |
real(KREAL) :: vx4,vy4,vz4,norm4 |
79 |
real(KREAL) :: vx5,vy5,vz5,norm5 |
80 |
|
81 |
INTEGER(KINT) :: I, J, IAt, K |
82 |
INTEGER(KINT) :: Kat, Lat, L |
83 |
|
84 |
REAL(KREAL) :: sAngleIatIKat, sAngleIIatLat |
85 |
|
86 |
LOGICAL :: debug, bond, AddPrimitiveCoord |
87 |
|
88 |
INTERFACE |
89 |
function valid(string) result (isValid) |
90 |
CHARACTER(*), intent(in) :: string |
91 |
logical :: isValid |
92 |
END function VALID |
93 |
|
94 |
|
95 |
FUNCTION angle(v1x,v1y,v1z,norm1,v2x,v2y,v2z,norm2) |
96 |
use Path_module, only : Pi,KINT, KREAL |
97 |
real(KREAL) :: v1x,v1y,v1z,norm1 |
98 |
real(KREAL) :: v2x,v2y,v2z,norm2 |
99 |
real(KREAL) :: angle |
100 |
END FUNCTION ANGLE |
101 |
|
102 |
|
103 |
|
104 |
FUNCTION angle_d(v1x,v1y,v1z,norm1,v2x,v2y,v2z,norm2,v3x,v3y,v3z,norm3) |
105 |
use Path_module, only : Pi,KINT, KREAL |
106 |
real(KREAL) :: v1x,v1y,v1z,norm1 |
107 |
real(KREAL) :: v2x,v2y,v2z,norm2 |
108 |
real(KREAL) :: v3x,v3y,v3z,norm3 |
109 |
real(KREAL) :: angle_d,ca,sa |
110 |
END FUNCTION ANGLE_D |
111 |
|
112 |
|
113 |
|
114 |
|
115 |
SUBROUTINE Calc_Xprim(nat,x,y,z,Coordinate,NPrim,XPrimitive,XPrimRef) |
116 |
! |
117 |
! This subroutine uses the description of a list of Coordinates |
118 |
! to compute the values of the coordinates for a given geometry. |
119 |
! |
120 |
!!!!!!!!!! |
121 |
! Input: |
122 |
! Na: INTEGER, Number of atoms |
123 |
! x,y,z(Na): REAL, cartesian coordinates of the considered geometry |
124 |
! Coordinate (Pointer(ListCoord)): description of the wanted coordiantes |
125 |
! NPrim, INTEGER: Number of coordinates to compute |
126 |
! |
127 |
! Optional: XPrimRef(NPrim) REAL: array that contains coordinates values for |
128 |
! a former geometry. Useful for Dihedral angles evolution... |
129 |
|
130 |
!!!!!!!!!!! |
131 |
! Output: |
132 |
! XPrimimite(NPrim) REAL: array that will contain the values of the coordinates |
133 |
! |
134 |
!!!!!!!!! |
135 |
|
136 |
Use VarTypes |
137 |
Use Io_module |
138 |
Use Path_module, only : pi |
139 |
|
140 |
IMPLICIT NONE |
141 |
|
142 |
Type (ListCoord), POINTER :: Coordinate |
143 |
INTEGER(KINT), INTENT(IN) :: Nat,NPrim |
144 |
REAL(KREAL), INTENT(IN) :: x(Nat), y(Nat), z(Nat) |
145 |
REAL(KREAL), INTENT(IN), OPTIONAL :: XPrimRef(NPrim) |
146 |
REAL(KREAL), INTENT(OUT) :: XPrimitive(NPrim) |
147 |
|
148 |
END SUBROUTINE CALC_XPRIM |
149 |
END INTERFACE |
150 |
|
151 |
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
152 |
! |
153 |
!!!!!!!! PFL Debug |
154 |
! |
155 |
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
156 |
INTEGER(KINT) :: I1,I2,I3,I4 |
157 |
LOGICAL :: debugPFL |
158 |
REAL(KREAL), ALLOCATABLE :: val_zmat(:,:) |
159 |
INTEGER(KINT), ALLOCATABLE :: indzmat0(:,:) |
160 |
|
161 |
|
162 |
debug=valid("Calc_baker") |
163 |
debugPFL=valid("BakerPFL") |
164 |
|
165 |
if (debug) Call Header(" Entering Cal_baker ") |
166 |
|
167 |
IF (Nat.le.2) THEN |
168 |
WRITE(*,*) "I do not work for less than 2 atoms :-p" |
169 |
RETURN |
170 |
END IF |
171 |
|
172 |
IF (debug) THEN |
173 |
WRITE(*,*) "DBG Calc_baker, geom" |
174 |
DO I=1,Nat |
175 |
WRITE(*,'(1X,I5,":",I3,3(1X,F15.3))') I,atome(I),XyzGeomI(IGeomRef,1:3,I) |
176 |
END DO |
177 |
END IF |
178 |
|
179 |
|
180 |
|
181 |
|
182 |
! we construct the list of bonds, valence angles and dihedral angles using this connectivity |
183 |
ALLOCATE(Coordinate) |
184 |
NULLIFY(Coordinate%next) |
185 |
NbBonds=0 |
186 |
NbAngles=0 |
187 |
NbDihedrals=0 |
188 |
|
189 |
CurrentCoord => Coordinate |
190 |
|
191 |
ALLOCATE(Liaisons(Nat,0:NMaxL),Geom(3,Nat),x(Nat),y(Nat),z(Nat)) |
192 |
ALLOCATE(x_refGeom(Nat),y_refGeom(Nat),z_refGeom(Nat)) |
193 |
|
194 |
x_refGeom(1:Nat) = XyzGeomI(IGeomRef,1,1:Nat) |
195 |
y_refGeom(1:Nat) = XyzGeomI(IGeomRef,2,1:Nat) |
196 |
z_refGeom(1:Nat) = XyzGeomI(IGeomRef,3,1:Nat) |
197 |
|
198 |
!!!!!!!!!!!!!!!!!!!! |
199 |
! |
200 |
! PFL Debug |
201 |
! |
202 |
!!!!!!!!!!!!!!!!!!!! |
203 |
if (debugPFL) THEN |
204 |
WRITE(*,*) "DBG PFL Calc_BAKER - Calculating Zmat" |
205 |
if (.not.ALLOCATED(indzmat)) THEN |
206 |
ALLOCATE(indzmat(Nat,5)) |
207 |
END IF |
208 |
IF (.NOT.ALLOCATED(DzDc)) THEN |
209 |
ALLOCATE(DzDc(3,nat,3,Nat)) |
210 |
END IF |
211 |
IF (.NOT.ALLOCATED(Val_zmat)) THEN |
212 |
ALLOCATE(val_zmat(nat,3)) |
213 |
END IF |
214 |
IF (.NOT.ALLOCATED(indzmat0)) THEN |
215 |
ALLOCATE(indzmat0(nat,5)) |
216 |
END IF |
217 |
! We construct a Zmat |
218 |
IF (Frozen(1).NE.0) THEN |
219 |
Renum=.TRUE. |
220 |
!!! remplcaer indzmat par indzmat0 puis renumerotation |
221 |
call Calc_zmat_const_frag(Nat,atome,IndZmat0,val_zmat, & |
222 |
x_refGeom,y_refGeom,z_refGeom, & |
223 |
r_cov,fact,frozen) |
224 |
ELSE |
225 |
!no zmat... we create our own |
226 |
call Calc_zmat_frag(Nat,atome,IndZmat0,val_zmat, & |
227 |
x_refGeom,y_refGeom,z_refGeom, & |
228 |
r_cov,fact) |
229 |
END IF |
230 |
|
231 |
DO I=1,Nat |
232 |
Order(IndZmat0(I,1))=I |
233 |
OrderInv(I)=IndZmat0(I,1) |
234 |
END DO |
235 |
IndZmat(1,1)=Order(IndZmat0(1,1)) |
236 |
IndZmat(2,1)=Order(IndZmat0(2,1)) |
237 |
IndZmat(2,2)=Order(IndZmat0(2,2)) |
238 |
IndZmat(3,1)=Order(IndZmat0(3,1)) |
239 |
IndZmat(3,2)=Order(IndZmat0(3,2)) |
240 |
IndZmat(3,3)=Order(IndZmat0(3,3)) |
241 |
DO I=4,Nat |
242 |
DO J=1,4 |
243 |
IndZmat(I,J)=Order(IndZmat0(I,J)) |
244 |
END DO |
245 |
end do |
246 |
|
247 |
|
248 |
WRITE(*,*) "DBG PFL CalcBaker, IndZmat0" |
249 |
DO I=1,Nat |
250 |
WRITE(*,*) I, indZmat0(I,:) |
251 |
end do |
252 |
|
253 |
WRITE(*,*) "DBG PFL CalcBaker, IndZmat" |
254 |
DO I=1,Nat |
255 |
WRITE(*,*) I, indZmat(I,:) |
256 |
end do |
257 |
|
258 |
|
259 |
DO I=1,Nat |
260 |
I1=indzmat0(I,1) |
261 |
I2=indzmat0(I,2) |
262 |
I3=indzmat0(I,3) |
263 |
I4=indzmat0(I,4) |
264 |
! Adding BOND between at1 and at2 |
265 |
WRITE(*,*) "DBG Indzmat",I,I1,I2,I3,I4 |
266 |
if (I2.NE.0) THEN |
267 |
NbBonds=NbBonds+1 |
268 |
! values of the coordinate (bond) is calculated for the reference geometry |
269 |
! whereas the coordinates (bonds) (CurrentCoord%At1, CurrentCoord%At2, etc.) |
270 |
! are determined using all geometries. vecteur is defined in bib_oper.f |
271 |
Call vecteur(I1,I2,x_refGeom,y_refGeom,z_refGeom,vx2,vy2,vz2,Norm2) |
272 |
! CurrentCoord%value=Norm2 |
273 |
WRITE(*,'(1X,A,I5,A,2(1X,F12.6))') "Bond",NbBonds," Norm2,val_zmat",Norm2,val_zmat(I,1) |
274 |
CurrentCoord%value=val_zmat(I,1) |
275 |
CurrentCoord%SignDihedral = 0 |
276 |
CurrentCoord%At1=I1 |
277 |
CurrentCoord%At2=I2 |
278 |
CurrentCoord%Type="BOND" |
279 |
IF (.NOT.ASSOCIATED(CurrentCoord%next)) THEN |
280 |
ALLOCATE(CurrentCoord%next) |
281 |
NULLIFY(CurrentCoord%next%next) |
282 |
END IF |
283 |
CurrentCoord => CurrentCoord%next |
284 |
END IF !IE.NE.0 |
285 |
if (I3.NE.0) THEN |
286 |
|
287 |
NbAngles=NbAngles+1 |
288 |
! values of the coordinate (bond) is calculated for the reference geometry |
289 |
! whereas the coordinates (bonds) (CurrentCoord%At1, CurrentCoord%At2, etc.) are |
290 |
! determined using all geometries. |
291 |
Call vecteur(i2,I3,x_refGeom,y_refGeom,z_refGeom,vx1,vy1,vz1,Norm1) |
292 |
Call vecteur(i2,I1,x_refGeom,y_refGeom,z_refGeom,vx2,vy2,vz2,Norm2) |
293 |
! CurrentCoord%value=angle(vx1,vy1,vz1,Norm1,vx2,vy2,vz2,Norm2)*Pi/180. |
294 |
WRITE(*,'(1X,A,I5,A,2(1X,F12.6))') "Angle",NbAngles," angle,val_zmat", & |
295 |
angle(vx1,vy1,vz1,Norm1,vx2,vy2,vz2,Norm2),val_zmat(I,2) |
296 |
CurrentCoord%value=val_zmat(I,2)*Pi/180. |
297 |
CurrentCoord%SignDihedral = 0 |
298 |
CurrentCoord%At1=I1 |
299 |
CurrentCoord%At2=I2 |
300 |
CurrentCoord%At3=I3 |
301 |
CurrentCoord%Type="ANGLE" |
302 |
IF (.NOT.ASSOCIATED(CurrentCoord%next)) THEN |
303 |
ALLOCATE(CurrentCoord%next) |
304 |
NULLIFY(CurrentCoord%next%next) |
305 |
END IF |
306 |
CurrentCoord => CurrentCoord%next |
307 |
END IF ! matches IF (I3.NE.0) |
308 |
if (I4.NE.0) THEN |
309 |
NbDihedrals=NbDihedrals+1 |
310 |
Call vecteur(I2,I1,x_refGeom,y_refGeom,z_refGeom,vx1,vy1,vz1,Norm1) |
311 |
Call vecteur(I2,I3,x_refGeom,y_refGeom,z_refGeom,vx2,vy2,vz2,Norm2) |
312 |
Call vecteur(I3,I4,x_refGeom,y_refGeom,z_refGeom,vx3,vy3,vz3,Norm3) |
313 |
CALL produit_vect(vx1,vy1,vz1,vx2,vy2,vz2, & |
314 |
vx4,vy4,vz4,norm4) |
315 |
CALL produit_vect(vx3,vy3,vz3,vx2,vy2,vz2, & |
316 |
vx5,vy5,vz5,norm5) |
317 |
! CurrentCoord%value=angle_d(vx4,vy4,vz4,norm4,vx5,vy5,vz5,norm5, & |
318 |
! vx2,vy2,vz2,norm2)*Pi/180. |
319 |
WRITE(*,'(1X,A,I5,A,2(1X,F12.6))') "Dihedral",NbDihedrals, & |
320 |
" angle_d,val_zmat", & |
321 |
angle_d(vx4,vy4,vz4,norm4,vx5,vy5,vz5,norm5, & |
322 |
vx2,vy2,vz2,norm2) & |
323 |
,val_zmat(I,3) |
324 |
|
325 |
CurrentCoord%value = val_zmat(I,3)*Pi/180. |
326 |
CurrentCoord%SignDihedral = int(sign(1.D0,val_zmat(I,3))) |
327 |
CurrentCoord%At1=I1 |
328 |
CurrentCoord%At2=I2 |
329 |
CurrentCoord%At3=I3 |
330 |
CurrentCoord%At4=I4 |
331 |
CurrentCoord%Type="DIHEDRAL" |
332 |
IF (.NOT.ASSOCIATED(CurrentCoord%next)) THEN |
333 |
ALLOCATE(CurrentCoord%next) |
334 |
NULLIFY(CurrentCoord%next%next) |
335 |
END IF |
336 |
CurrentCoord => CurrentCoord%next |
337 |
END IF ! matches IF (I4.NE.0) |
338 |
END DO ! matches DO I=1,Nat |
339 |
deallocate(indzmat0) |
340 |
|
341 |
ELSE !! if (debugPFL) |
342 |
|
343 |
DO IGeom=1, NGeomI |
344 |
!IGeom=1 ! need when using only first geometry. |
345 |
x(1:Nat) = XyzGeomI(IGeom,1,1:Nat) |
346 |
y(1:Nat) = XyzGeomI(IGeom,2,1:Nat) |
347 |
z(1:Nat) = XyzGeomI(IGeom,3,1:Nat) |
348 |
|
349 |
! First step: we calculate the connectivity |
350 |
Liaisons=0 |
351 |
|
352 |
Call CalcCnct(Nat,atome,x,y,z,LIAISONS,r_cov,fact) |
353 |
|
354 |
IF (debug) THEN |
355 |
WRITE(*,*) 'Connectivity done' |
356 |
DO I=1,Nat |
357 |
WRITE(*,'(1X,I5,":",I3,"=",15I4)') I,(Liaisons(I,J),J=0,NMaxL) |
358 |
WRITE(*,'(1X,I5,":",I3,"=",15I4)') I,(Liaisons(I,:)) |
359 |
END DO |
360 |
DO I=1,Nat |
361 |
WRITE(*,'(1X,I5,":",I3," r_cov=",F15.6)') I,atome(I),r_cov(atome(I)) |
362 |
END DO |
363 |
END IF |
364 |
|
365 |
DO I=1,Nat |
366 |
IF (Liaisons(I,0).NE.0) THEN |
367 |
DO J=1,Liaisons(I,0) |
368 |
! We add the bond |
369 |
Iat=Liaisons(I,J) |
370 |
IF (Iat.GT.I) THEN |
371 |
! Checking the uniqueness of bond. |
372 |
! Duplicate bonds should not be added. |
373 |
AddPrimitiveCoord=.True. |
374 |
ScanCoord=>Coordinate |
375 |
DO WHILE (Associated(ScanCoord%next)) |
376 |
IF (ScanCoord%Type .EQ. 'BOND') THEN |
377 |
IF (ScanCoord%At1 .EQ. Iat) THEN |
378 |
IF (ScanCoord%At2 .EQ. I) THEN |
379 |
AddPrimitiveCoord=.False. |
380 |
END IF |
381 |
END IF |
382 |
END IF |
383 |
ScanCoord => ScanCoord%next |
384 |
END DO |
385 |
|
386 |
IF (AddPrimitiveCoord) THEN |
387 |
NbBonds=NbBonds+1 |
388 |
! values of the coordinate (bond) is calculated for the reference geometry |
389 |
! whereas the coordinates (bonds) (CurrentCoord%At1, CurrentCoord%At2, etc.) |
390 |
! are determined using all geometries. vecteur is defined in bib_oper.f |
391 |
Call vecteur(I,Iat,x_refGeom,y_refGeom,z_refGeom,vx2,vy2,vz2,Norm2) |
392 |
CurrentCoord%value=Norm2 |
393 |
CurrentCoord%SignDihedral = 0 |
394 |
CurrentCoord%At1=Iat |
395 |
CurrentCoord%At2=I |
396 |
CurrentCoord%Type="BOND" |
397 |
IF (.NOT.ASSOCIATED(CurrentCoord%next)) THEN |
398 |
ALLOCATE(CurrentCoord%next) |
399 |
NULLIFY(CurrentCoord%next%next) |
400 |
END IF |
401 |
CurrentCoord => CurrentCoord%next |
402 |
END IF ! matches IF (AddPrimitiveCoord) THEN |
403 |
END IF ! matches IF (Iat.GT.I) THEN |
404 |
|
405 |
! Call Annul(Liaisons,IAt,I,Nat,NMaxL) |
406 |
! Liaisons(Iat,0)=Liaisons(Iat,0)-1 |
407 |
|
408 |
! We add the valence angles |
409 |
DO K=J+1,Liaisons(I,0) |
410 |
Kat=Liaisons(I,K) |
411 |
IF (Kat.GT.I) THEN |
412 |
! Checking the uniqueness of valence angle. |
413 |
! Duplicate bonds should not be added. |
414 |
AddPrimitiveCoord=.True. |
415 |
ScanCoord=>Coordinate |
416 |
DO WHILE (Associated(ScanCoord%next)) |
417 |
IF (ScanCoord%Type .EQ. 'ANGLE') THEN |
418 |
IF (ScanCoord%At1 .EQ. Iat) THEN |
419 |
IF (ScanCoord%At2 .EQ. I) THEN |
420 |
IF (ScanCoord%At3 .EQ. Kat) THEN |
421 |
AddPrimitiveCoord=.False. |
422 |
END IF |
423 |
END IF |
424 |
END IF |
425 |
END IF |
426 |
ScanCoord => ScanCoord%next |
427 |
END DO ! matches DO WHILE (Associated(ScanCoord%next)) |
428 |
|
429 |
IF (AddPrimitiveCoord) THEN |
430 |
IF (debug) WRITE(*,*) "Adding angle:",Iat,I,Kat |
431 |
NbAngles=NbAngles+1 |
432 |
! values of the coordinate (bond) is calculated for the reference geometry |
433 |
! whereas the coordinates (bonds) (CurrentCoord%At1, CurrentCoord%At2, etc.) are |
434 |
! determined using all geometries. |
435 |
Call vecteur(i,kat,x_refGeom,y_refGeom,z_refGeom,vx1,vy1,vz1,Norm1) |
436 |
Call vecteur(i,Iat,x_refGeom,y_refGeom,z_refGeom,vx2,vy2,vz2,Norm2) |
437 |
CurrentCoord%value=angle(vx1,vy1,vz1,Norm1,vx2,vy2,vz2,Norm2)*Pi/180. |
438 |
CurrentCoord%SignDihedral = 0 |
439 |
sAngleIatIKat=abs(sin(angle(vx1,vy1,vz1,Norm1,vx2,vy2,vz2,Norm2)*Pi/180.)) |
440 |
CurrentCoord%At1=Iat |
441 |
CurrentCoord%At2=I |
442 |
CurrentCoord%At3=KAt |
443 |
CurrentCoord%Type="ANGLE" |
444 |
IF (.NOT.ASSOCIATED(CurrentCoord%next)) THEN |
445 |
ALLOCATE(CurrentCoord%next) |
446 |
NULLIFY(CurrentCoord%next%next) |
447 |
END IF |
448 |
CurrentCoord => CurrentCoord%next |
449 |
END IF ! matches IF (AddPrimitiveCoord) THEN |
450 |
ELSE ! matches IF (Kat.GT.I) THEN |
451 |
! Kat is less than I. If there is a bond between I and Kat, then this angle |
452 |
! has already been taken into account. |
453 |
if (debug) WRITE(*,*) "Testing angle:",Iat,I,Kat |
454 |
Bond=.FALSE. |
455 |
DO L=1,Liaisons(Iat,0) |
456 |
IF (Liaisons(Iat,L).EQ.Kat) Bond=.TRUE. |
457 |
END DO |
458 |
IF (.NOT.Bond) THEN |
459 |
IF (debug) WRITE(*,*) "Not Bond Adding angle:",Iat,I,Kat |
460 |
|
461 |
! Checking the uniqueness of valence angle. |
462 |
! Duplicate bonds should not be added. |
463 |
AddPrimitiveCoord=.True. |
464 |
ScanCoord=>Coordinate |
465 |
DO WHILE (Associated(ScanCoord%next)) |
466 |
IF (ScanCoord%Type .EQ. 'ANGLE') THEN |
467 |
IF (ScanCoord%At1 .EQ. Iat) THEN |
468 |
IF (ScanCoord%At2 .EQ. I) THEN |
469 |
IF (ScanCoord%At3 .EQ. Kat) THEN |
470 |
AddPrimitiveCoord=.False. |
471 |
END IF |
472 |
END IF |
473 |
END IF |
474 |
END IF |
475 |
ScanCoord => ScanCoord%next |
476 |
END DO |
477 |
|
478 |
IF (AddPrimitiveCoord) THEN |
479 |
NbAngles=NbAngles+1 |
480 |
! values of the coordinate (bond) is calculated for the reference geometry |
481 |
! whereas the coordinates (bonds) (CurrentCoord%At1, CurrentCoord%At2, etc.) |
482 |
! are determined using all geometries. |
483 |
Call vecteur(i,kat,x_refGeom,y_refGeom,z_refGeom,vx1,vy1,vz1,Norm1) |
484 |
Call vecteur(i,Iat,x_refGeom,y_refGeom,z_refGeom,vx2,vy2,vz2,Norm2) |
485 |
CurrentCoord%value=angle(vx1,vy1,vz1,Norm1,vx2,vy2,vz2,Norm2)*Pi/180. |
486 |
CurrentCoord%SignDihedral = 0 |
487 |
CurrentCoord%At1=Iat |
488 |
CurrentCoord%At2=I |
489 |
CurrentCoord%At3=KAt |
490 |
CurrentCoord%Type="ANGLE" |
491 |
IF (.NOT.ASSOCIATED(CurrentCoord%next)) THEN |
492 |
ALLOCATE(CurrentCoord%next) |
493 |
NULLIFY(CurrentCoord%next%next) |
494 |
END IF |
495 |
CurrentCoord => CurrentCoord%next |
496 |
END IF ! matches IF (AddPrimitiveCoord) THEN |
497 |
END IF ! matches IF (.NOT.Bond) THEN |
498 |
END IF ! matches IF (Kat.GT.I) THEN, after else |
499 |
END DO ! matches DO K=J+1,Liaisons(I,0) |
500 |
END DO ! matches DO J=1,Liaisons(I,0) |
501 |
END IF ! matches IF (Liaisons(I,0).NE.0) THEN |
502 |
|
503 |
! We add dihedral angles. We do not concentrate on necessary angles, ie those that are |
504 |
! needed to build the molecule. Instead we collect all of them, and the choice of the best |
505 |
! ones will be done when diagonalizing the Baker matrix. |
506 |
IF (Liaisons(I,0).GE.3) Then |
507 |
DO J=1,Liaisons(I,0) |
508 |
Iat=Liaisons(I,J) |
509 |
! values of the coordinate (bond) is calculated for the reference geometry |
510 |
! whereas the coordinates (bonds) (CurrentCoord%At1, CurrentCoord%At2, etc.) are |
511 |
! determined using all geometries. |
512 |
Call vecteur(Iat,i,x_refGeom,y_refGeom,z_refGeom,vx2,vy2,vz2,Norm2) |
513 |
|
514 |
DO K=J+1,Liaisons(I,0) |
515 |
Kat=Liaisons(I,K) |
516 |
Call vecteur(i,kat,x_refGeom,y_refGeom,z_refGeom,vx1,vy1,vz1,Norm1) |
517 |
sAngleIatIKat=abs(sin(angle(vx1,vy1,vz1,Norm1,-vx2,-vy2,-vz2,Norm2)*Pi/180.)) |
518 |
|
519 |
DO L=K+1,Liaisons(I,0) |
520 |
Lat=Liaisons(I,L) |
521 |
if (debug) WRITE(*,*) "0-Dihe Kat,I,Iat:",Kat,I,Iat,angle(vx1,vy1,vz1,Norm1,vx2,vy2,vz2,Norm2) |
522 |
Call vecteur(iat,lat,x_refGeom,y_refGeom,z_refGeom,vx3,vy3,vz3,Norm3) |
523 |
sAngleIIatLat=abs(sin(angle(vx3,vy3,vz3,Norm3,vx2,vy2,vz2,Norm2)*Pi/180.)) |
524 |
if (debug) WRITE(*,*) "0-Dihe I,Iat,Lat:",angle(vx3,vy3,vz3,Norm3,vx2,vy2,vz2,Norm2) |
525 |
IF ((sAngleIatIKat.ge.0.01d0).AND.(sAngleIIatLat.ge.0.01d0)) THEN |
526 |
if (debug) WRITE(*,*) "Adding dihedral:",Kat,I,Iat,Lat |
527 |
|
528 |
! Checking for the uniqueness of valence angle. |
529 |
! Duplicate bonds should not be added. |
530 |
AddPrimitiveCoord=.True. |
531 |
ScanCoord=>Coordinate |
532 |
DO WHILE (Associated(ScanCoord%next)) |
533 |
IF (ScanCoord%Type .EQ. 'DIHEDRAL') THEN |
534 |
IF (ScanCoord%At1 .EQ. Kat) THEN |
535 |
IF (ScanCoord%At2 .EQ. I) THEN |
536 |
IF (ScanCoord%At3 .EQ. Iat) THEN |
537 |
IF (ScanCoord%At4 .EQ. Lat) THEN |
538 |
AddPrimitiveCoord=.False. |
539 |
END IF |
540 |
END IF |
541 |
END IF |
542 |
END IF |
543 |
END IF |
544 |
ScanCoord => ScanCoord%next |
545 |
END DO ! matches DO WHILE (Associated(ScanCoord%next)) |
546 |
|
547 |
! IF (AddPrimitiveCoord) THEN |
548 |
! NbDihedrals=NbDihedrals+1 |
549 |
! CALL produit_vect(vx3,vy3,vz3,vx2,vy2,vz2, & |
550 |
! vx5,vy5,vz5,norm5) |
551 |
!CALL produit_vect(vx1,vy1,vz1,vx2,vy2,vz2, & |
552 |
! vx4,vy4,vz4,norm4) |
553 |
! CurrentCoord%value=angle_d(vx4,vy4,vz4,norm4,vx5,vy5,vz5,norm5, & |
554 |
! vx2,vy2,vz2,norm2)*Pi/180. |
555 |
! CurrentCoord%SignDihedral = 0 |
556 |
! CurrentCoord%At1=Kat |
557 |
! CurrentCoord%At2=I |
558 |
! CurrentCoord%At3=IAt |
559 |
! CurrentCoord%At4=LAt |
560 |
!WRITE(*,'(1X,":(dihed 1)",I5," - ",I5," - ",I5," - ",I5," = ",F15.6)') & |
561 |
! ScanCoord%At1,ScanCoord%At2,ScanCoord%At3,ScanCoord%At4,ScanCoord%Value*180./Pi |
562 |
! CurrentCoord%Type="DIHEDRAL" |
563 |
! IF (.NOT.ASSOCIATED(CurrentCoord%next)) THEN |
564 |
! ALLOCATE(CurrentCoord%next) |
565 |
! NULLIFY(CurrentCoord%next%next) |
566 |
! END IF |
567 |
! CurrentCoord => CurrentCoord%next |
568 |
!END IF ! matches IF (AddPrimitiveCoord) THEN |
569 |
|
570 |
END IF ! matches IF ((sAngleIatIKat.ge.0.01d0).AND.(sAngleIIatLat.ge.0.01d0)) THEN |
571 |
END DO ! matches DO L=K+1,Liaisons(I,0) |
572 |
END DO ! matches DO K=J+1,Liaisons(I,0) |
573 |
END DO ! matches DO J=1,Liaisons(I,0) |
574 |
END IF !matches IF (Liaisons(I,0).GE.3) Then |
575 |
|
576 |
! we now add other dihedrals |
577 |
DO J=1,Liaisons(I,0) |
578 |
Iat=Liaisons(I,J) |
579 |
If (Iat.LE.I) Cycle |
580 |
! values of the coordinate (bond) is calculated for the reference geometry |
581 |
! whereas the coordinates (bonds) (CurrentCoord%At1, CurrentCoord%At2, etc.) are |
582 |
! determined using all geometries. |
583 |
Call vecteur(Iat,i,x_refGeom,y_refGeom,z_refGeom,vx2,vy2,vz2,Norm2) |
584 |
|
585 |
DO K=1,Liaisons(I,0) |
586 |
IF (Liaisons(I,K).EQ.Iat) Cycle |
587 |
Kat=Liaisons(I,K) |
588 |
Call vecteur(i,kat,x_refGeom,y_refGeom,z_refGeom,vx1,vy1,vz1,Norm1) |
589 |
sAngleIatIKat=abs(sin(angle(vx1,vy1,vz1,Norm1,-vx2,-vy2,-vz2,Norm2)*Pi/180.)) |
590 |
|
591 |
DO L=1,Liaisons(Iat,0) |
592 |
Lat=Liaisons(Iat,L) |
593 |
IF (Lat.eq.I) Cycle |
594 |
|
595 |
if (debug) WRITE(*,*) "Dihe Kat,I,Iat:",angle(vx1,vy1,vz1,Norm1,vx2,vy2,vz2,Norm2),sAngleIatIKat |
596 |
Call vecteur(iat,lat,x_refGeom,y_refGeom,z_refGeom,vx,vy,vz,Norm) |
597 |
|
598 |
sAngleIIatLat=abs(sin(angle(vx,vy,vz,Norm,vx2,vy2,vz2,Norm2)*Pi/180.)) |
599 |
if (debug) WRITE(*,*) "Dihe I,Iat,Lat:",angle(vx,vy,vz,Norm,vx2,vy2,vz2,Norm2) |
600 |
! we add the dihedral angle Kat-I-Iat-Lat |
601 |
if (debug) WRITE(*,*) "Trying dihedral:",Kat,I,Iat,Lat,sAngleIatIKat,sAngleIIatLat |
602 |
IF ((Lat.NE.Kat).AND.(Lat.Ne.I).AND.(sAngleIatIKat.ge.0.01d0) & |
603 |
.AND.(sAngleIIatLat.ge.0.01d0)) THEN |
604 |
if (debug) WRITE(*,*) "Adding dihedral:",Kat,I,Iat,Lat |
605 |
|
606 |
! Checking for the uniqueness of valence angle. |
607 |
! Duplicate bonds should not be added. |
608 |
AddPrimitiveCoord=.True. |
609 |
ScanCoord=>Coordinate |
610 |
DO WHILE (Associated(ScanCoord%next)) |
611 |
IF (ScanCoord%Type .EQ. 'DIHEDRAL') THEN |
612 |
IF (ScanCoord%At1 .EQ. Kat) THEN |
613 |
IF (ScanCoord%At2 .EQ. I) THEN |
614 |
IF (ScanCoord%At3 .EQ. Iat) THEN |
615 |
IF (ScanCoord%At4 .EQ. Lat) THEN |
616 |
AddPrimitiveCoord=.False. |
617 |
END IF |
618 |
END IF |
619 |
END IF |
620 |
END IF |
621 |
END IF |
622 |
ScanCoord => ScanCoord%next |
623 |
END DO |
624 |
|
625 |
IF (AddPrimitiveCoord) THEN |
626 |
NbDihedrals=NbDihedrals+1 |
627 |
Call vecteur(iat,lat,x_refGeom,y_refGeom,z_refGeom,vx3,vy3,vz3,Norm3) |
628 |
CALL produit_vect(vx3,vy3,vz3,vx2,vy2,vz2, & |
629 |
vx5,vy5,vz5,norm5) |
630 |
CALL produit_vect(vx1,vy1,vz1,vx2,vy2,vz2, & |
631 |
vx4,vy4,vz4,norm4) |
632 |
CurrentCoord%value=angle_d(vx4,vy4,vz4,norm4,vx5,vy5,vz5,norm5, & |
633 |
vx2,vy2,vz2,norm2)*Pi/180. |
634 |
CurrentCoord%At1=Kat |
635 |
CurrentCoord%At2=I |
636 |
CurrentCoord%At3=IAt |
637 |
CurrentCoord%At4=LAt |
638 |
CurrentCoord%Type="DIHEDRAL" |
639 |
CurrentCoord%SignDihedral=int(sign(1.d0,angle_d(vx4,vy4,vz4,norm4, & |
640 |
vx5,vy5,vz5,norm5, vx2,vy2,vz2,norm2))) |
641 |
WRITE(*,'(1X,":(dihed 2)",I5," - ",I5," - ",I5," - ",I5," = ",F15.6)') ScanCoord%At1,& |
642 |
ScanCoord%At2, ScanCoord%At3,ScanCoord%At4,ScanCoord%Value*180./Pi |
643 |
IF (.NOT.ASSOCIATED(CurrentCoord%next)) THEN |
644 |
ALLOCATE(CurrentCoord%next) |
645 |
NULLIFY(CurrentCoord%next%next) |
646 |
END IF |
647 |
CurrentCoord => CurrentCoord%next |
648 |
END IF ! matches IF (AddPrimitiveCoord) THEN |
649 |
END IF ! matches IF ((Lat.NE.Kat).AND.(Lat.Ne.I).AND.(sAngleIatIKat.ge.0.01d0) ..... |
650 |
END DO ! matches DO L=1,Liaisons(Iat,0) |
651 |
END DO ! matches DO K=1,Liaisons(I,0) |
652 |
END DO ! matches DO J=1,Liaisons(I,0) |
653 |
END DO ! end of loop over all atoms, DO I=1, Nat |
654 |
END DO ! matches DO IGeom=1, NGeomI |
655 |
|
656 |
|
657 |
|
658 |
END IF ! if (debugPFL) |
659 |
DEALLOCATE(Liaisons) |
660 |
|
661 |
WRITE(*,*) "I have found" |
662 |
WRITE(*,*) NbBonds, " bonds" |
663 |
WRITE(*,*) NbAngles, " angles" |
664 |
WRITE(*,*) NbDihedrals, " dihedrals" |
665 |
|
666 |
WRITE(*,*) 'All in all...' |
667 |
ScanCoord=>Coordinate |
668 |
I=0 |
669 |
DO WHILE (Associated(ScanCoord%next)) |
670 |
I=I+1 |
671 |
SELECT CASE (ScanCoord%Type) |
672 |
CASE ('BOND') |
673 |
!WRITE(IOOUT,'(1X,I3,":",I5," - ",I5," = ",F15.6)') I,ScanCoord%At1,ScanCoord%At2,ScanCoord%Value |
674 |
WRITE(*,'(1X,I3,":(bond )",I5," - ",I5," = ",F15.6)') I,ScanCoord%At1,ScanCoord%At2,ScanCoord%Value |
675 |
CASE ('ANGLE') |
676 |
!WRITE(IOOUT,'(1X,I3,":",I5," - ",I5," - ",I5," = ",F15.6)') I,ScanCoord%At1, & |
677 |
! ScanCoord%At2, ScanCoord%At3,ScanCoord%Value*180./Pi |
678 |
WRITE(*,'(1X,I3,":(angle)",I5," - ",I5," - ",I5," = ",F15.6)') I,ScanCoord%At1, & |
679 |
ScanCoord%At2, ScanCoord%At3,ScanCoord%Value*180./Pi |
680 |
CASE ('DIHEDRAL') |
681 |
!WRITE(IOOUT,'(1X,I3,":",I5," - ",I5," - ",I5," - ",I5," = ",F15.6)') I,ScanCoord%At1,& |
682 |
! ScanCoord%At2, ScanCoord%At3,ScanCoord%At4,ScanCoord%Value*180./Pi |
683 |
WRITE(*,'(1X,I3,":(dihed)",I5," - ",I5," - ",I5," - ",I5," = ",F15.6)') I,ScanCoord%At1,& |
684 |
ScanCoord%At2, ScanCoord%At3,ScanCoord%At4,ScanCoord%Value*180./Pi |
685 |
END SELECT |
686 |
ScanCoord => ScanCoord%next |
687 |
END DO |
688 |
|
689 |
! NPrim is the number of primitive coordinates. |
690 |
NPrim=NbBonds+NbAngles+NbDihedrals |
691 |
|
692 |
! Now calculating values of all primitive bonds for all geometries: |
693 |
ALLOCATE(Xprimitive(NgeomI,NPrim),XPrimitiveF(NGeomF,NPrim)) |
694 |
ALLOCATE(UMatF(NGeomF,NPrim,NCoord)) |
695 |
ALLOCATE(BTransInvF(NGeomF,NCoord,3*Nat)) |
696 |
|
697 |
! We calculate the Primitive values for the first geometry, this plays a special role |
698 |
! as it will serve as a reference for other geometries ! |
699 |
|
700 |
x(1:Nat) = XyzGeomI(1,1,1:Nat) |
701 |
y(1:Nat) = XyzGeomI(1,2,1:Nat) |
702 |
z(1:Nat) = XyzGeomI(1,3,1:Nat) |
703 |
|
704 |
Call Calc_Xprim(nat,x,y,z,Coordinate,NPrim,XPrimitive(1,:)) |
705 |
|
706 |
DO IGeom=2, NGeomI |
707 |
|
708 |
x(1:Nat) = XyzGeomI(IGeom,1,1:Nat) |
709 |
y(1:Nat) = XyzGeomI(IGeom,2,1:Nat) |
710 |
z(1:Nat) = XyzGeomI(IGeom,3,1:Nat) |
711 |
|
712 |
Call Calc_Xprim(nat,x,y,z,Coordinate,NPrim,XPrimitive(IGeom,:),Xprimitive(1,:)) |
713 |
|
714 |
! ScanCoord=>Coordinate |
715 |
! I=0 ! index for the NPrim (NPrim is the number of primitive coordinates). |
716 |
! DO WHILE (Associated(ScanCoord%next)) |
717 |
! I=I+1 |
718 |
! SELECT CASE (ScanCoord%Type) |
719 |
! CASE ('BOND') |
720 |
! Call vecteur(ScanCoord%At2,ScanCoord%At1,x,y,z,vx2,vy2,vz2,Norm2) |
721 |
! Xprimitive(IGeom,I) = Norm2 |
722 |
! CASE ('ANGLE') |
723 |
! Call vecteur(ScanCoord%At2,ScanCoord%At3,x,y,z,vx1,vy1,vz1,Norm1) |
724 |
! Call vecteur(ScanCoord%At2,ScanCoord%At1,x,y,z,vx2,vy2,vz2,Norm2) |
725 |
! Xprimitive(IGeom,I) = angle(vx1,vy1,vz1,Norm1,vx2,vy2,vz2,Norm2)*Pi/180. |
726 |
! CASE ('DIHEDRAL') |
727 |
|
728 |
|
729 |
! Call vecteur(ScanCoord%At2,ScanCoord%At1,x,y,z,vx1,vy1,vz1,Norm1) |
730 |
! Call vecteur(ScanCoord%At2,ScanCoord%At3,x,y,z,vx2,vy2,vz2,Norm2) |
731 |
! Call vecteur(ScanCoord%At3,ScanCoord%At4,x,y,z,vx3,vy3,vz3,Norm3) |
732 |
! Call produit_vect(vx1,vy1,vz1,vx2,vy2,vz2, & |
733 |
! vx4,vy4,vz4,norm4) |
734 |
! Call produit_vect(vx3,vy3,vz3,vx2,vy2,vz2, & |
735 |
! vx5,vy5,vz5,norm5) |
736 |
|
737 |
! DiheTmp= angle_d(vx4,vy4,vz4,norm4,vx5,vy5,vz5,norm5, & |
738 |
! vx2,vy2,vz2,norm2) |
739 |
! Xprimitive(IGeom,I) = DiheTmp*Pi/180. |
740 |
! ! We treat large dihedral angles differently as +180=-180 mathematically and physically |
741 |
! ! but this causes lots of troubles when converting baker to cart. |
742 |
! ! So we ensure that large dihedral angles always have the same sign |
743 |
! if (abs(ScanCoord%SignDihedral).NE.1) THEN |
744 |
! ScanCoord%SignDihedral=Int(Sign(1.D0,DiheTmp)) |
745 |
! ELSE |
746 |
! If ((abs(DiheTmp).GE.170.D0).AND.(Sign(1.,DiheTmp)*ScanCoord%SignDihedral<0)) THEN |
747 |
! Xprimitive(IGeom,I) = DiheTmp*Pi/180.+ ScanCoord%SignDihedral*2.*Pi |
748 |
! END IF |
749 |
! END IF |
750 |
! ! Another test... will all this be consistent ??? |
751 |
! ! We want the shortest path, so we check that the change in dihedral angles is less than Pi: |
752 |
! IF (IGeom.GT.1) THEN |
753 |
! IF (abs(Xprimitive(IGeom,I)-XPrimitive(IGeom-1,I)).GE.Pi) THEN |
754 |
! if ((Xprimitive(IGeom,I)-XPrimitive(IGeom-1,I)).GE.Pi) THEN |
755 |
! Xprimitive(IGeom,I)=Xprimitive(IGeom,I)-2.d0*Pi |
756 |
! ELSE |
757 |
! Xprimitive(IGeom,I)=Xprimitive(IGeom,I)+2.d0*Pi |
758 |
! END IF |
759 |
! END IF |
760 |
! END IF |
761 |
! END SELECT |
762 |
! ScanCoord => ScanCoord%next |
763 |
! END DO ! matches DO WHILE (Associated(ScanCoord%next)) |
764 |
END DO ! matches DO IGeom=1, NGeomI |
765 |
|
766 |
IF (debug) THEN |
767 |
WRITE(*,*) "DBG Calc_Baker Xprimitives for all geometries" |
768 |
DO IGeom=1, NGeomI |
769 |
WRITE(*,*) "Geometry ",IGeom |
770 |
ScanCoord=>Coordinate |
771 |
I=0 ! index for the NPrim (NPrim is the number of primitive coordinates). |
772 |
DO WHILE (Associated(ScanCoord%next)) |
773 |
I=I+1 |
774 |
SELECT CASE (ScanCoord%Type) |
775 |
CASE ('BOND') |
776 |
WRITE(*,'(1X,I3,":",I5," - ",I5," = ",F15.6)') I,ScanCoord%At1,ScanCoord%At2,Xprimitive(IGeom,I) |
777 |
CASE ('ANGLE') |
778 |
WRITE(*,'(1X,I3,":",I5," - ",I5," - ",I5," = ",F15.6)') I,ScanCoord%At1, & |
779 |
ScanCoord%At2, ScanCoord%At3,Xprimitive(IGeom,I)*180./Pi |
780 |
CASE ('DIHEDRAL') |
781 |
WRITE(*,'(1X,I3,":",I5," - ",I5," - ",I5," - ",I5," = ",F15.6)') I,ScanCoord%At1,& |
782 |
ScanCoord%At2, ScanCoord%At3,ScanCoord%At4,Xprimitive(IGeom,I)*180./Pi |
783 |
END SELECT |
784 |
ScanCoord => ScanCoord%next |
785 |
END DO ! matches DO WHILE (Associated(ScanCoord%next)) |
786 |
END DO ! matches DO IGeom=1, NGeomI |
787 |
END IF |
788 |
! Here reference geometry is assigned to Geom. Earlier x,y,z(Nat) were assigned from XyzGeomI(...) |
789 |
! before the call of this (Calc_baker) subroutine and passed to this subroutine as arguments of the |
790 |
! subroutine. |
791 |
! PFL 19 Feb 2008 -> |
792 |
! In order to reproduce the B matrix in the Baker article, one has |
793 |
! to divide Geom by a0. |
794 |
! However, this is inconsitent with our way of storing geometries |
795 |
! so we do not do it ! |
796 |
|
797 |
Geom(1,:)=XyzGeomI(IGeomRef,1,1:Nat) ! XyzGeomI(NGeomI,3,Nat) |
798 |
Geom(2,:)=XyzGeomI(IGeomRef,2,1:Nat) |
799 |
Geom(3,:)=XyzGeomI(IGeomRef,3,1:Nat) |
800 |
|
801 |
! <- PFL 19 Feb 2008 |
802 |
|
803 |
! We now have to compute dq/dx... |
804 |
ALLOCATE(BprimT(3*Nat,NPrim)) |
805 |
BprimT=0.d0 |
806 |
|
807 |
ScanCoord=>Coordinate |
808 |
I=0 |
809 |
DO WHILE (Associated(ScanCoord%next)) |
810 |
I=I+1 |
811 |
SELECT CASE (ScanCoord%Type) |
812 |
CASE ('BOND') |
813 |
CALL CONSTRAINTS_BONDLENGTH_DER(Nat,ScanCoord%at1,ScanCoord%AT2, & |
814 |
Geom,BprimT(1,I)) |
815 |
CASE ('ANGLE') |
816 |
CALL CONSTRAINTS_BONDANGLE_DER(Nat,ScanCoord%At1,ScanCoord%AT2, & |
817 |
ScanCoord%At3,Geom,BprimT(1,I)) |
818 |
CASE ('DIHEDRAL') |
819 |
CALL CONSTRAINTS_TORSION_DER2(Nat,ScanCoord%At1,ScanCoord%AT2, & |
820 |
ScanCoord%At3,ScanCoord%At4,Geom,BprimT(1,I)) |
821 |
! CALL CONSTRAINTS_TORSION_DER(Nat,ScanCoord%At1,ScanCoord%AT2, & |
822 |
! ScanCoord%At3,ScanCoord%At4,Geom,BprimT(1,I)) |
823 |
|
824 |
END SELECT |
825 |
ScanCoord => ScanCoord%next |
826 |
END DO |
827 |
|
828 |
|
829 |
!!!!!!!!!!!!!!!!!!!! |
830 |
! |
831 |
! PFL Debug |
832 |
! |
833 |
!!!!!!!!!!!!!!!!!!!! |
834 |
|
835 |
if (debugPFL) THEN |
836 |
WRITe(*,*) "DBG PFL Calc_Baker =====" |
837 |
DzDc=0.d0 |
838 |
WRITE(*,*) "PFL DBG, DzdC" |
839 |
do I=1,Nat |
840 |
Geom(1,i)=XyzGeomI(IGeomRef,1,orderInv(i)) ! XyzGeomI(NGeomI,3,Nat) |
841 |
Geom(2,i)=XyzGeomI(IGeomRef,2,OrderInv(i)) |
842 |
Geom(3,i)=XyzGeomI(IGeomRef,3,OrderInv(i)) |
843 |
WRITE(*,*) I,OrderInv(I),Geom(:,I) |
844 |
end do |
845 |
CALL CalcBMat_int (nat, Geom, indzmat, dzdc, massat,atome) |
846 |
|
847 |
WRITE(*,*) "PFL DBG, BPrimT" |
848 |
DO I=1,NPrim |
849 |
WRITE(*,'(50(1X,F12.6))') BPrimT(:,I) |
850 |
END DO |
851 |
WRITe(*,*) "DBG PFL Calc_Baker =====" |
852 |
END IF |
853 |
|
854 |
!!!!!!!!!!!!!!!!!!!! |
855 |
! |
856 |
! PFL Debug |
857 |
! |
858 |
!!!!!!!!!!!!!!!!!!!! |
859 |
|
860 |
DEALLOCATE(Geom) |
861 |
|
862 |
! BprimT(3*Nat,NPrim) |
863 |
! We now compute G=B(BT) matrix |
864 |
ALLOCATE(Gmat(NPrim,NPrim)) |
865 |
GMat=0.d0 |
866 |
DO I=1,NPrim |
867 |
DO J=1,3*Nat |
868 |
GMat(:,I)=Gmat(:,I)+BprimT(J,:)*BprimT(J,I) !*1.d0/mass(atome(int(K/3.d0))) |
869 |
END DO |
870 |
END DO |
871 |
|
872 |
!!!!!!!!!!!!!!!!!!!! |
873 |
! |
874 |
! PFL Debug ---> |
875 |
! |
876 |
!!!!!!!!!!!!!!!!!!!! |
877 |
if (debugPFL) THEN |
878 |
GMat=0.d0 |
879 |
DO I=1,NPrim |
880 |
GMat(I,I)=1. |
881 |
END DO |
882 |
END IF |
883 |
|
884 |
!!!!!!!!!!!!!!!!!!!! |
885 |
! |
886 |
! <--- PFL Debug |
887 |
! |
888 |
!!!!!!!!!!!!!!!!!!!! |
889 |
|
890 |
!DO I=1,NPrim |
891 |
! WRITE(*,'(20(1X,F6.3))') GMat(1:min(20,NPrim),I) |
892 |
!END DO |
893 |
|
894 |
! We diagonalize G |
895 |
ALLOCATE(EigVal(NPrim),EigVec(NPrim,NPrim)) |
896 |
Call Jacobi(GMat,NPrim,EigVal,EigVec,NPrim) |
897 |
Call Trie(NPrim,EigVal,EigVec,NPrim) |
898 |
DO I=1,NPrim |
899 |
WRITE(*,'(1X,"Vector ",I3,": e=",F8.3)') I,EigVal(i) |
900 |
WRITE(*,'(20(1X,F8.4))') EigVec(1:min(20,NPrim),I) |
901 |
END DO |
902 |
|
903 |
!DO I=1,NPrim |
904 |
! WRITE(*,'(1X,"Vector ",I3,(20(F8.3)))') I,EigVec(1:min(20,NPrim),I) |
905 |
! END DO |
906 |
|
907 |
!DO I=1,NPrim |
908 |
! WRITE(*,'(1X,"Vector ",I3,F8.3,F8.3,F8.3,F8.3,F8.3,F8.3,F8.3,F8.3,F8.3,F8.3,F8.3,F8.3,& |
909 |
! F8.3,F8.3,F8.3,F8.3,F8.3)') I,EigVec(1,I),EigVec(4,I),EigVec(6,I),EigVec(14,I),EigVec(16,I),& |
910 |
!EigVec(2,I),EigVec(3,I),EigVec(5,I),EigVec(12,I),EigVec(13,I),EigVec(15,I),EigVec(7,I),& |
911 |
!EigVec(8,I),EigVec(9,I),EigVec(10,I),EigVec(11,I),EigVec(17,I) |
912 |
!END DO |
913 |
|
914 |
! We construct the new DzDc(b=baker) matrix |
915 |
! UMat is nonredundant vector set, i.e. set of eigenvectors of BB^T corresponding |
916 |
! to eigenvalues > zero. |
917 |
ALLOCATE(UMat(NPrim,NCoord)) |
918 |
ALLOCATE(UMat_local(NPrim,NCoord)) |
919 |
ALLOCATE(UMat_tmp(NPrim,3*Nat-6)) |
920 |
! BMat_BakerT(3*Nat,NCoord), allocated in Path.f90, |
921 |
! NCoord=3*Nat-6 |
922 |
BMat_BakerT = 0.d0 |
923 |
J=0 |
924 |
UMat=0.d0 |
925 |
UMat_tmp=0.d0 |
926 |
DO I=1,NPrim |
927 |
IF (abs(eigval(I))>=1e-6) THEN |
928 |
J=J+1 |
929 |
DO K=1,NPrim |
930 |
!DzDcb(:,J)=DzDcb(:,J)+Eigvec(K,I)*BprimT(:,K) ! original |
931 |
! BprimT is transpose of B^prim. |
932 |
! B = UMat^T B^prim, B^T = (B^prim)^T UMat |
933 |
BMat_BakerT(:,J)=BMat_BakerT(:,J)+BprimT(:,K)*Eigvec(K,I) |
934 |
!Dzdcb(:,J)=DzDcb(:,J)+Eigvec(K,:)*BprimT(I,K) |
935 |
END DO |
936 |
IF(J .GT. 3*Nat-6) THEN |
937 |
WRITE(*,*) 'Number of vectors in Eigvec with eigval .GT. 1e-6(=UMat) (=' & |
938 |
,J,') exceeded 3*Nat-6=',3*Nat-6, & |
939 |
'Stopping the calculation.' |
940 |
STOP |
941 |
END IF |
942 |
UMat_tmp(:,J) = Eigvec(:,I) |
943 |
END IF |
944 |
END DO |
945 |
|
946 |
if (debug) THEN |
947 |
WRITE(*,*) "DBG Calc_Baker: UMat" |
948 |
DO I=1,3*Nat-6 |
949 |
WRITE(*,'(50(1X,F12.8))') UMat_tmp(:,I) |
950 |
END DO |
951 |
END IF |
952 |
|
953 |
! THIS IS TO BE CORRECTED: |
954 |
UMat = UMat_tmp |
955 |
!J=0 |
956 |
!DO I=1, 3*Nat-6 |
957 |
! IF ((I .NE. 6) .AND. (I .NE. 7) .AND. (I .NE. 9)) Then |
958 |
! J=J+1 |
959 |
! Print *, 'J=', J, ' I=', I |
960 |
! UMat(:,J) = UMat_tmp(:,I) |
961 |
!END IF |
962 |
!END DO |
963 |
! BMat_BakerT=0.d0 |
964 |
!DO I=1,NCoord |
965 |
! DO J=1,NPrim |
966 |
! B = UMat^T B^prim, B^T = (B^prim)^T UMat |
967 |
! BMat_BakerT(:,I)=BMat_BakerT(:,I)+BprimT(:,J)*UMat(J,I) |
968 |
! END DO |
969 |
!END DO |
970 |
! TO BE CORRECTED, ENDS HERE. |
971 |
|
972 |
IntCoordI=0.d0 |
973 |
DO IGeom=1, NGeomI |
974 |
DO I=1, NPrim |
975 |
! Transpose of UMat is needed below, that is why UMat(I,:). |
976 |
IntCoordI(IGeom,:) = IntCoordI(IGeom,:) + UMat(I,:)*Xprimitive(IGeom,I) |
977 |
END DO |
978 |
END DO |
979 |
|
980 |
if (debug) THEN |
981 |
WRITE(*,*) "DBG Calc_Baker IntCoordI" |
982 |
DO IGeom=1, NGeomI |
983 |
WRITE(*,'(1X,I5,":",30(1X,F10.6))') IGeom,IntCoordI(IGeom,:) |
984 |
END DO |
985 |
END IF |
986 |
|
987 |
|
988 |
! the following might be used to get rid of some coordinates that |
989 |
! do not respect the symmetry of the problem. |
990 |
|
991 |
! ! We look at the group symmetry of the molecule. Nothing complicate here |
992 |
! ! we just look for planar symmetry. |
993 |
! ! We first place it into the standard orientation |
994 |
! Call Std_ori(Nat,x,y,z,a,b,c,massat) |
995 |
|
996 |
! ! Is the molecule planar |
997 |
! DO I=1,Nat |
998 |
! WRITE(*,*) nom(atome(i)),x(i),y(i),z(i) |
999 |
! END DO |
1000 |
DEALLOCATE(BprimT,GMat,EigVal,EigVec) ! BprimT is B^prim^T |
1001 |
! Calculation of BTransInv starts here: |
1002 |
! Calculation of BBT(3*Nat-6,3*Nat-6)=BB^T: |
1003 |
! BMat_BakerT(3*Nat,NCoord) is Transpose of B = UMat^TB^prim |
1004 |
ALLOCATE(BBT(NCoord,NCoord)) |
1005 |
BBT = 0.d0 |
1006 |
DO I=1, NCoord |
1007 |
DO J=1, 3*Nat |
1008 |
! BBT(:,I) forms BB^T |
1009 |
BBT(:,I) = BBT(:,I) + BMat_BakerT(J,:)*BMat_BakerT(J,I) |
1010 |
END DO |
1011 |
END DO |
1012 |
|
1013 |
!ALLOCATE(V(3*Nat,3*Nat)) |
1014 |
!Call GINVSE(BBT,3*Nat,3*Nat,BBT_inv,3*Nat,3*Nat,3*Nat,1e-6,V,3*Nat,FAIL,NOUT) |
1015 |
!DEALLOCATE(V) |
1016 |
ALLOCATE(BBT_inv(NCoord,NCoord)) |
1017 |
! GenInv is in Mat_util.f90 |
1018 |
Call GenInv(NCoord,BBT,BBT_inv,NCoord) |
1019 |
|
1020 |
! Calculation of (B^T)^-1 = (BB^T)^-1B: |
1021 |
!ALLOCATE(BTransInv(3*Nat-6,3*Nat)) ! allocated in Path.f90 |
1022 |
BTransInv = 0.d0 |
1023 |
DO I=1, 3*Nat |
1024 |
DO J=1, NCoord |
1025 |
BTransInv(:,I) = BTransInv(:,I) + BBT_inv(:,J)*BMat_BakerT(I,J) |
1026 |
END DO |
1027 |
END DO |
1028 |
BTransInv_local = BTransInv |
1029 |
UMat_local = UMat |
1030 |
DEALLOCATE(BBT,BBT_inv,UMat_tmp) |
1031 |
DEALLOCATE(x_refGeom,y_refGeom,z_refGeom,x,y,z) |
1032 |
IF (debug) WRITE(*,*) "DBG Calc_baker over." |
1033 |
END SUBROUTINE Calc_baker |