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