Statistiques
| Révision :

root / src / IntCoord_der.f90

Historique | Voir | Annoter | Télécharger (51,01 ko)

1
!----------------------------------------------------------------------
2
!  Copyright 2003-2014 Ecole Normale Supérieure de Lyon, 
3
!  Centre National de la Recherche Scientifique,
4
!  Université Claude Bernard Lyon 1. All rights reserved.
5
!
6
!  This work is registered with the Agency for the Protection of Programs 
7
!  as IDDN.FR.001.100009.000.S.P.2014.000.30625
8
!
9
!  Authors: P. Fleurat-Lessard, P. Dayal
10
!  Contact: optnpath@gmail.com
11
!
12
! This file is part of "Opt'n Path".
13
!
14
!  "Opt'n Path" is free software: you can redistribute it and/or modify
15
!  it under the terms of the GNU Affero General Public License as
16
!  published by the Free Software Foundation, either version 3 of the License,
17
!  or (at your option) any later version.
18
!
19
!  "Opt'n Path" is distributed in the hope that it will be useful,
20
!  but WITHOUT ANY WARRANTY; without even the implied warranty of
21
!
22
!  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
23
!  GNU Affero General Public License for more details.
24
!
25
!  You should have received a copy of the GNU Affero General Public License
26
!  along with "Opt'n Path". If not, see <http://www.gnu.org/licenses/>.
27
!
28
! Contact The Office of Technology Licensing, valorisation@ens-lyon.fr,
29
! for commercial licensing opportunities.
30
!----------------------------------------------------------------------
31

    
32
      SUBROUTINE VPRD(X,Y,Z)
33
!     ******************************************************************
34
!     **  CROSS PRODUCT (VECTOR PRODUCT OF TWO VECTORS)               **
35
!     ******************************************************************
36
        Use VarTypes
37
      IMPLICIT NONE
38
      REAL(KREAL),INTENT(IN) :: X(3)
39
      REAL(KREAL),INTENT(IN) :: Y(3)
40
      REAL(KREAL),INTENT(OUT):: Z(3)
41
!     ******************************************************************
42
      Z(1)=X(2)*Y(3)-X(3)*Y(2)
43
      Z(2)=X(3)*Y(1)-X(1)*Y(3)
44
      Z(3)=X(1)*Y(2)-X(2)*Y(1)
45
      RETURN
46
      END
47

    
48
!     ..................................................................
49
      SUBROUTINE CONSTRAINTS_BONDLENGTH_DER(NAT,IAT1,IAT2,R0,B)
50
!     ******************************************************************
51
!     **                                                              **
52
!     **  BOND-LENGTH CONSTRAINT : SET UP COORDINATES                 **   
53
!     **                                                              **
54
!     ******************************************************************
55
        Use VarTypes
56
      IMPLICIT NONE
57
      INTEGER(KINT),INTENT(IN) :: NAT
58
      INTEGER(KINT),INTENT(IN) :: IAT1
59
      INTEGER(KINT),INTENT(IN) :: IAT2
60
      REAL(KREAL)   ,INTENT(IN) :: R0(3,NAT)
61
      REAL(KREAL)   ,INTENT(OUT):: B(3,NAT)
62

    
63
      INTEGER(KINT)            :: I
64
      REAL(KREAL)               :: DIS, DI
65
      real(KREAL)               :: dr(3)
66
!     ******************************************************************
67
      DR(:)=R0(:,IAT2)-R0(:,IAT1)
68
!
69
!     ==================================================================
70
!     == CALCULATE VALUE AND FIRST TWO DERIVATIVES OF D=|R_2-R_2|     ==
71
!     ==================================================================
72
      DIS=DSQRT(DOT_PRODUCT(DR,DR))
73
      DO I=1,3
74
         DI=DR(I)/DIS
75
         B(I,IAT2)=DI
76
         B(I,IAT1)=-DI
77
      ENDDO
78
      RETURN
79
      END SUBROUTINE CONSTRAINTS_BONDLENGTH_DER
80

    
81
      SUBROUTINE CONSTRAINTS_GLC_DER(NAT,VEC,B)
82
!     ******************************************************************
83
!     **  CONSTRAINTS_SETGLC                                          **
84
!     **  SET GENERAL LINEAR CONSTRAINT (GLC) CONST+VEC*R=0           **
85
!     **                                                              **
86
!     **  APPROXIMATES CONSTRAINT FUNCTION (WITHOUT VALUE)            **
87
!     **  BY AN QUADRATIC POLYNOM IN A+B*R+0.5*R*C*R                  **
88
!     **  THE CONSTRAINT IS THE POLYNOMIAL MINUS THE VALUE FROM THE   **
89
!     **  REFERENCE STRUCTURE                                         **
90
!     ******************************************************************
91
        Use VarTypes
92
      IMPLICIT NONE
93
      INTEGER(KINT),INTENT(IN) :: NAT
94
      REAL(KREAL)   ,INTENT(IN) :: VEC(3,NAT)
95
      REAL(KREAL)   ,INTENT(OUT):: B(3,NAT)
96
      INTEGER(KINT)            :: IAT,I
97
!     ******************************************************************
98
      DO IAT=1,NAT
99
         DO I=1,3
100
            B(I,IAT)=VEC(I,IAT)
101
         ENDDO
102
      ENDDO
103
      RETURN
104
      END SUBROUTINE CONSTRAINTS_GLC_DER
105

    
106

    
107
      SUBROUTINE CONSTRAINTS_SUBSTITUTION_DER(NAT,IAT1,IAT2,IAT3,R0,B)
108
!     ******************************************************************
109
!     **  CONSTRAINTS_SETsubstitution                                 **
110
!     **  APPROXIMATES CONSTRAINT FUNCTION (WITHOUT VALUE)            **
111
!     **  BY AN QUADRATIC POLYNOMial IN A+B*R+0.5*R*C*R               **
112
!     **  THE CONSTRAINT IS THE POLYNOMIAL MINUS THE VALUE FROM THE   **
113
!     **  REFERENCE STRUCTURE                                         **
114
!     **  CONSTRAINS THE difference OF THE BOND DISTANCES BETWEEN     **
115
!     **  R1-R2 AND R2-R3                                             **
116
!     ******************************************************************
117
      Use Vartypes
118
      implicit none
119
      INTEGER(KINT),INTENT(IN) :: NAT
120
      INTEGER(KINT),INTENT(IN) :: IAT1
121
      INTEGER(KINT),INTENT(IN) :: IAT2
122
      INTEGER(KINT),INTENT(IN) :: IAT3
123
      REAL(KREAL)   ,INTENT(OUT):: B(3,NAT)
124
      REAL(KREAL)   ,INTENT(IN) :: R0(3,NAT)
125
      REAL(KREAL)               :: X(9)
126
      real(KREAL)               :: phi
127
      REAL(KREAL)               :: X1X2,Y1Y2,Z1Z2,X2X3,Y2Y3,Z2Z3,D12,D22,  &
128
                                   D1,D2,D13,D23
129
      integer(KINT)            :: i
130
!     ******************************************************************
131
!
132
! --- TRANSFER R0 TO X ARRAY TO MAKE HANDLING AND PASSING EASIER
133
      DO I=1,3
134
       X(I)  =R0(I,IAT1)
135
       X(I+3)=R0(I,IAT2)
136
       X(I+6)=R0(I,IAT3)
137
      ENDDO
138
!
139
! --- FOR THE FOLLOWING, THE CONSTRAINT ITSELF WILL BE CALLED PHI
140
! --- CALCULATE PHI,D PHI/D X_I, D^2 PHI/D X_I D X_J
141
!
142
      X1X2=X(1)-X(4)
143
      Y1Y2=X(2)-X(5)
144
      Z1Z2=X(3)-X(6)
145
      X2X3=X(4)-X(7)
146
      Y2Y3=X(5)-X(8)
147
      Z2Z3=X(6)-X(9)
148
!
149
      D12=(X1X2)**2 + (Y1Y2)**2 + (Z1Z2)**2
150
      D22=(X2X3)**2 + (Y2Y3)**2 + (Z2Z3)**2
151
      D1= SQRT(D12)
152
      D2= SQRT(D22)
153
      D13= D1**3
154
      D23= D2**3
155
      PHI= D1 - D2
156
!
157
      B(1,IAT1)= (X1X2)/D1
158
      B(1,IAT2)= -((X1X2)/D1) - (X2X3)/D2
159
      B(1,IAT3)= (X2X3)/D2
160
      B(2,IAT1)= (Y1Y2)/D1
161
      B(2,IAT2)= -((Y1Y2)/D1) - (Y2Y3)/D2
162
      B(2,IAT3)= (Y2Y3)/D2
163
      B(3,IAT1)= (Z1Z2)/D1
164
      B(3,IAT2)= -((Z1Z2)/D1) - (Z2Z3)/D2
165
      B(3,IAT3)= (Z2Z3)/D2
166
!
167
      RETURN
168
      END SUBROUTINE CONSTRAINTS_SUBSTITUTION_DER
169

    
170
!     ..................................................................
171
      SUBROUTINE CONSTRAINTS_BONDANGLE_DER(NAT,IAT1,IAT2,IAT3,R0,B)
172
!     ******************************************************************
173
!     **
174
!     **  CONSTRAINS THE ANGLE ENCLOSED BETWEEN ATOMS IAT1--IAT2--IAT3
175
!     **  In the followin, IAT1 is also denoted by A, IAT2 by B and IAT3 by C
176
!     **
177
!     ******************************************************************
178
        Use VarTypes
179
      IMPLICIT none
180
      INTEGER(KINT),INTENT(IN) :: NAT
181
      INTEGER(KINT),INTENT(IN) :: IAT1
182
      INTEGER(KINT),INTENT(IN) :: IAT2
183
      INTEGER(KINT),INTENT(IN) :: IAT3
184
      REAL(KREAL)   ,INTENT(OUT):: B(3,NAT)
185
      REAL(KREAL)   ,INTENT(IN) :: R0(3,NAT)
186
      REAL(KREAL)               :: X(9), DI(9)
187
      real(KREAL)               :: phi
188
      integer(KINT)            :: i
189
      REAL(KREAL)               :: UU    ! U*U
190
      REAL(KREAL)               :: VV    ! V*V
191
      REAL(KREAL)               :: UV    ! U*V
192
      REAL(KREAL)               :: D    ! 1/SQRT(A*B)
193
      REAL(KREAL)               :: DADX(9)
194
      REAL(KREAL)               :: DBDX(9)
195
      REAL(KREAL)               :: DCDX(9)
196
      REAL(KREAL)               :: DDDX(9)
197
      REAL(KREAL)               :: DJDX(9)
198
      REAL(KREAL)               :: DGDX(9)    ! D(A*B)/DX
199
      REAL(KREAL)               :: AB,AB2N,V,PREFAC
200
      REAL(KREAL)               :: COSPHI,SINPHI
201
!
202
! Variable used in case of a linear molecule with three atoms call A,B and C
203
! Phi is the ABC angle
204
! Threshold use to decide that sin(Phi)=0.
205
      REAL(KREAL), PARAMETER    :: ThreshSin=1e-7
206
! Sine and Cosine of the euler Phi angle between ABC and x axis
207
      REAL(KREAL)               :: CosP,SinP
208
! Sine and Cosine of the euler theta angle between ABC and z axis
209
      REAL(KREAL)               :: CosTheta,SinTheta
210
! Vectors CB, CA,AB
211
      REAL(KREAL)               ::vCB(3),vCA(3),vAB(3)
212
! distance AB in the xAy plane
213
      REAL(KREAL)               :: dABxy
214
! cosine of the BAC angle
215
      REAL(KREAL)               :: CosBAC
216
! cosine of the BCA angle
217
      REAL(KREAL)               :: CosBCA
218

    
219
!
220
! For debugginh purposes:
221
      LOGICAL, PARAMETER :: Debug=.False.
222

    
223

    
224
      IF (Debug) WRITE(*,*) "================ Entering  CONSTRAINTS_BONDANGLE_DER =========="
225
!     ******************************************************************
226
!
227
!     ==================================================================
228
!     ==  SET UP CONSTRAINT IN COORDINATE SPACE                       ==
229
!     ==================================================================
230
!
231
!     ==================================================================
232
!     == CALCULATE VALUE AND FIRST TWO DERIVATIVES OF                 ==
233
!     == PHI=<X1-X2|X3-X2>/SQRT(<X1-X2|X1-X2><X3-X2|X3-X2>)           ==
234
!     ==================================================================
235
! --- TRANSFER R0 TO X ARRAY TO MAKE HANDLING AND PASSING EASIER
236
      DO I=1,3
237
       X(I)  =R0(I,IAT1)
238
       X(I+3)=R0(I,IAT2)
239
       X(I+6)=R0(I,IAT3)
240
      ENDDO
241

    
242
      if (debug) THEN
243
         WRITE(*,*) "X=",X
244
         WRITE(*,*) "A=",X(1)-X(4),X(2)-X(5),X(3)-X(6)
245
         WRITE(*,*) "B=",X(7)-X(4),X(8)-X(5),X(9)-X(6)
246
      END IF
247

    
248
!     ----------------------------------------------------------------
249
!     --- R1=X(1:3)   R2=X(4:6)  R3=X(7:9)
250
!     --- U=R1-R2 V=R3-R2
251
!     --- A=U*U=d(At1-At2)**2
252
!     --- B=V*V=d(At2-At3)**2
253
!     --- C=U*V=vec(At2-At1)*vec(At2-At3)
254
!     ----------------------------------------------------------------
255
      UU=(X(1)-X(4))**2 + (X(2)-X(5))**2 + (X(3)-X(6))**2
256
      VV=(X(7)-X(4))**2 + (X(8)-X(5))**2 + (X(9)-X(6))**2
257
      UV=(X(1)-X(4))*(X(7)-X(4))  &
258
           +(X(2)-X(5))*(X(8)-X(5))   &
259
           +(X(3)-X(6))*(X(9)-X(6))
260
      AB=UU*VV
261
      D=SQRT(AB)
262
      COSPHI=UV/D       !=COS(PHI)
263
      PHI=DACOS(COSPHI)
264
      SINPHI=DSQRT(1.D0-COSPHI**2)
265

    
266
      IF (SINPHI.GE.ThreshSin) THEN
267
! The following evaluation of the derivative involves
268
! 1/sinphi... so we test that it is not 0 :-)
269
!     ---  CALCULATE DADX
270
      DADX(1)= 2.D0*(X(1)-X(4))
271
      DADX(2)= 2.D0*(X(2)-X(5))
272
      DADX(3)= 2.D0*(X(3)-X(6))
273
      DADX(4)=-2.D0*(X(1)-X(4))
274
      DADX(5)=-2.D0*(X(2)-X(5))
275
      DADX(6)=-2.D0*(X(3)-X(6))
276
      DADX(7)= 0.D0
277
      DADX(8)= 0.D0
278
      DADX(9)= 0.D0
279
!     ---  CALCULATE DBDX
280
      DBDX(1)= 0.D0
281
      DBDX(2)= 0.D0
282
      DBDX(3)= 0.D0
283
      DBDX(4)=-2.D0*(X(7)-X(4))
284
      DBDX(5)=-2.D0*(X(8)-X(5))
285
      DBDX(6)=-2.D0*(X(9)-X(6))
286
      DBDX(7)= 2.D0*(X(7)-X(4))
287
      DBDX(8)= 2.D0*(X(8)-X(5))
288
      DBDX(9)= 2.D0*(X(9)-X(6))
289
! --- CALCULATE DCDX
290
      DCDX(1)= X(7)-X(4)
291
      DCDX(2)= X(8)-X(5)
292
      DCDX(3)= X(9)-X(6)
293
      DCDX(4)=-X(7)+2.D0*X(4)-X(1)
294
      DCDX(5)=-X(8)+2.D0*X(5)-X(2)
295
      DCDX(6)=-X(9)+2.D0*X(6)-X(3)
296
      DCDX(7)= X(1)-X(4)
297
      DCDX(8)= X(2)-X(5)
298
      DCDX(9)= X(3)-X(6)
299

    
300
!       WRITE(*,*) 'PFL DER ANGLE'
301
!       WRITE(*,'(15(1X,F10.6))') DADX
302
!       WRITE(*,'(15(1X,F10.6))') DBDX
303
!       WRITE(*,'(15(1X,F10.6))') DCDX
304

    
305
!     ==================================================================
306
!     ==  CALCULATE PHI ETC  ===========================================
307
!     ==================================================================
308
!     --- CALCULATE DERIVATIVE OF ARCCOS
309
      V=-1.D0/SINPHI                 !=-SIN(PHI)=D(ACOS(PHI))/D(COSPHI)
310
! --- CALCULATE FIRST DERIVATIVES OF G=A*B
311
      DO I=1,9
312
        DGDX(I)=DADX(I)*VV+UU*DBDX(I)
313
      ENDDO
314
! --- CALCULATE DDDX
315
      AB2N=0.5D0/D
316
      DO I=1,9
317
        DDDX(I)=AB2N*DGDX(I)
318
      ENDDO
319
! --- CALCULATE FIRST DERIVATIVES OF J
320
      PREFAC=1.D0/AB
321
      DO I=1,9
322
        DJDX(I)=(DCDX(I)*D-UV*DDDX(I))*PREFAC
323
      ENDDO
324
! --- CALCULATE FIRST DERIVATIVE
325
      DO I=1,9
326
        DI(I)=DJDX(I)*V
327
      ENDDO
328
!
329
! Store it back to B
330
      DO I=1,3
331
        B(I,IAT1)=DI(I)
332
        B(I,IAT2)=DI(I+3)
333
        B(I,IAT3)=DI(I+6)
334
      ENDDO
335

    
336
   ELSE
337
! First calculate cosP, sinP, CosTheta, SinTheta
338
         dABxy=sqrt((X(4)-X(1))**2+(X(5)-X(2))**2)
339
         cosP=abs((X(7)-X(1)))/sqrt((X(7)-X(1))**2+(X(8)-X(2))**2)
340
!         SinP=(X(8)-X(2))/sqrt((X(7)-X(1))**2+(X(8)-X(2))**2)
341
         sinP=sqrt(1.d0-cosP**2)
342
         cosTheta=(X(6)-X(3))/(sqrt(UU))
343
         SinTheta=dABxy/(sqrt(UU))
344
         IF (SinTheta.le.ThreshSin) THEN
345
! Theta=0, so that the euler Phi angle is ill-defined. We put it to 0:
346
            CosP=1.d0
347
            SinP=0.d0
348
         END IF
349
         IF (Debug)   WRITE(*,*) "CosP,SinP ",CosP,SinP
350
         if (debug)   WRITE(*,*) "CosTheta,SinTheta",CosTheta,SinTheta
351
! We calculate cos(BAC) and cos(BCA) because they are needed to calculate
352
! the derivatives for the centre atom (B)
353
         DO I=1,3
354
            vCB(I)=X(3+I)-X(6+I)
355
            vCA(I)=X(I)-X(6+I)
356
            vAB(I)=x(3+I)-x(I)
357
         END DO
358
         CosBCA=DOT_PRODUCT(vCB,vCA)/sqrt(DOT_PRODUCT(vCB,vCB)*DOT_PRODUCT(vCA,vCA))
359
         CosBAC=-DOT_PRODUCT(vAB,vCA)/sqrt(DOT_PRODUCT(vAB,vAB)*DOT_PRODUCT(vCA,vCA))
360
         if (debug) WRITE(*,*) "CosBAC,CosBCA:",CosBAC,CosBCA
361
! Sign of the derivatives are -1 for Phi=Pi and 1 for Phi=0
362
! ie is cos(Phi).
363
         DI(1)=COSPHI/sqrt(UU)
364
         DI(2)=DI(1)
365
         DI(3)=0.d0
366
         DI(7)=COSPHI/sqrt(VV)
367
         DI(8)=DI(7)
368
         DI(9)=0.d0
369
         DI(4)=CosBCA*DI(1)+cosBAC*DI(7)
370
         DI(5)=DI(4)
371
         DI(6)=0.D0
372
         if (debug) WRITE(*,'(A,9(1X,F12.6))') 'DI=',DI(1:9)
373
         if (debug) THEN
374
            WRITE(*,'(3(1X,F12.6))') cosTheta*CosP,sinP,-cosP*sinTheta
375
            WRITE(*,'(3(1X,F12.6))') -cosTheta*sinP,cosP,sinP*sinTheta
376
            WRITE(*,'(3(1X,F12.6))') sinTheta,0.,cosTheta
377
         END IF
378
! We now take care of the fact that the three fragment are
379
! not along the Z-axis
380
         B(1,Iat1)=cosTheta*CosP*DI(1)+sinP*DI(2)-cosP*sinTheta*DI(3)
381
         B(1,Iat2)=cosTheta*CosP*DI(4)+sinP*DI(5)-cosP*sinTheta*DI(6)
382
         B(1,Iat3)=cosTheta*CosP*DI(7)+sinP*DI(8)-cosP*sinTheta*DI(9)
383
         B(2,Iat1)=-cosTheta*SinP*DI(1)+cosP*DI(2)+sinP*sinTheta*DI(3)
384
         B(2,Iat2)=-cosTheta*SinP*DI(4)+cosP*DI(5)+sinP*sinTheta*DI(6)
385
         B(2,Iat3)=-cosTheta*SinP*DI(7)+cosP*DI(8)+sinP*sinTheta*DI(9)
386
         B(3,Iat1)=sinTheta*DI(1)+cosTheta*DI(3)
387
         B(3,Iat2)=sinTheta*DI(4)+cosTheta*DI(6)
388
         B(3,Iat3)=sinTheta*DI(7)+cosTheta*DI(9)
389
      END IF
390

    
391
      IF (Debug) WRITE(*,*) "================ Exiting  CONSTRAINTS_BONDANGLE_DER =========="
392
      END SUBROUTINE CONSTRAINTS_BONDANGLE_DER
393

    
394

    
395
!     ..................................................................
396
      SUBROUTINE CONSTRAINTS_TORSION_DER(NAT,IAT1,IAT2,IAT3,IAT4,R0,B)
397
!     ******************************************************************
398
!     **
399
!     **  CONSTRAINS THE TORSION OF IAT1-IAT2---IAT3-IAT4
400
!     **
401
!     ******************************************************************
402
        Use VarTypes
403
      IMPLICIT none
404
      INTEGER(KINT),INTENT(IN) :: NAT
405
      INTEGER(KINT),INTENT(IN) :: IAT1
406
      INTEGER(KINT),INTENT(IN) :: IAT2
407
      INTEGER(KINT),INTENT(IN) :: IAT3
408
      INTEGER(KINT),INTENT(IN) :: IAT4
409
      REAL(KREAL)   ,INTENT(IN) :: R0(3,NAT)
410
      REAL(KREAL)   ,INTENT(OUT):: B(3,NAT)
411
      REAL(KREAL)               :: X(12), DI(12)
412
      REAL(KREAL)               :: P(3),Q(3),T(3),U(3),R(3),W(3)
413
      integer(KINT)            :: i, j
414

    
415
      REAL(KREAL)   :: TORS, PP,QQ, PQ
416
      REAL(KREAL) :: DADX(12),DBDX(12),DCDX(12)
417
      REAL(KREAL) :: DDDX(12), DJDX(12), DGDX(12), DPDX(3, 12), DQDX(3, 12)
418
      real(KREAL) :: d, coverd, abyb2n, abyb, prefac, v
419

    
420
!     ******************************************************************
421
!
422
! --- TRANSFER R0 TO X ARRAY TO MAKE HANDLING AND PASSING EASIER
423
      DO I=1,3
424
        X(I)=R0(I,IAT1)
425
        X(I+3)=R0(I,IAT2)
426
        X(I+6)=R0(I,IAT3)
427
        X(I+9)=R0(I,IAT4)
428
      ENDDO
429
!
430
! --- CALCULATE PHI,D PHI/D X_I, D^2 PHI/D X_I D X_J
431
      T(:)=X(1:3)  -X(4:6)
432
      U(:)=X(7:9)  -X(4:6)
433
      R(:)=X(10:12)-X(7:9)
434
      W(:)=X(4:6)  -X(7:9)
435
      CALL VPRD(T,U,P)
436
      CALL VPRD(W,R,Q)
437
      PQ=DOT_PRODUCT(P,Q)
438
      PP=DOT_PRODUCT(P,P)
439
      QQ=DOT_PRODUCT(Q,Q)
440
!
441
      ABYB=PP*QQ
442
      D=ABYB**0.5D0
443
      COVERD=PQ/D
444
      ABYB2N=0.5D0/D
445
!
446
      TORS=DACOS(COVERD)
447
      WRITE(*,*) "Torsion derivs  Phi=", Tors*180./3.1415926,Iat1,Iat2,Iat3,Iat4
448
!
449
! ----------------------------------------------------------------
450
! --- CALCULATE FIRST DERIVATIVES
451
! ----------------------------------------------------------------
452
! --- CALCULATE DERIVATIVE OF ARCCOS
453
      IF (dabs(DABS(COVERD)-1.d0) < 1.D-9) THEN
454
        WRITE(*,*)   'WARNING: DERIVATIVE OF TORSION AT ZERO AND PI IS UNDEFINED!'
455
        COVERD=sign(0.999d0,COVERD)
456
      TORS=DACOS(COVERD)
457
      !WRITE(*,*) "Torsion derivs, Phi=", Tors
458

    
459
!        V=0.D0
460
!      ELSE
461
     END IF
462
        V=-1.D0/DSQRT(1.D0-(COVERD)**2.D0)
463
!      ENDIF
464
! --- CALCULATE DPDX
465
! --- ZERO OUT SPARSE MATRIX
466
      DO J=1,3
467
       DO I=1,12
468
        DPDX(J,I)=0.D0
469
        DQDX(J,I)=0.D0
470
       ENDDO
471
      ENDDO
472
      DPDX(1,2)=X(9)-X(6)
473
      DPDX(1,3)=X(5)-X(8)
474
      DPDX(1,5)=X(3)-X(9)
475
      DPDX(1,6)=X(8)-X(2)
476
      DPDX(1,8)=X(6)-X(3)
477
      DPDX(1,9)=X(2)-X(5)
478
      DPDX(2,1)=X(6)-X(9)
479
      DPDX(2,3)=X(7)-X(4)
480
      DPDX(2,4)=X(9)-X(3)
481
      DPDX(2,6)=X(1)-X(7)
482
      DPDX(2,7)=X(3)-X(6)
483
      DPDX(2,9)=X(4)-X(1)
484
      DPDX(3,1)=X(8)-X(5)
485
      DPDX(3,2)=X(4)-X(7)
486
      DPDX(3,4)=X(2)-X(8)
487
      DPDX(3,5)=X(7)-X(1)
488
      DPDX(3,7)=X(5)-X(2)
489
      DPDX(3,8)=X(1)-X(4)
490
! --- CALCULATE DQDX
491
      DQDX(1,5)=X(12)-X(9)
492
      DQDX(1,6)=X(8)-X(11)
493
      DQDX(1,8)=X(6)-X(12)
494
      DQDX(1,9)=X(11)-X(5)
495
      DQDX(1,11)=X(9)-X(6)
496
      DQDX(1,12)=X(5)-X(8)
497
      DQDX(2,4)=X(9)-X(12)
498
      DQDX(2,6)=X(10)-X(7)
499
      DQDX(2,7)=X(12)-X(6)
500
      DQDX(2,9)=X(4)-X(10)
501
      DQDX(2,10)=X(6)-X(9)
502
      DQDX(2,12)=X(7)-X(4)
503
      DQDX(3,4)=X(11)-X(8)
504
      DQDX(3,5)=X(7)-X(10)
505
      DQDX(3,7)=X(5)-X(11)
506
      DQDX(3,8)=X(10)-X(4)
507
      DQDX(3,10)=X(8)-X(5)
508
      DQDX(3,11)=X(4)-X(7)
509

    
510
!       WRITE(*,*) 'DP/DX(1) '
511
!       WRITE(*,'(1X,12(F9.5,1X))') DPDX(1,:)
512
!       WRITE(*,*) 'DP/DX(2) '
513
!       WRITE(*,'(1X,12(F9.5,1X))') DPDX(2,:)
514
!       WRITE(*,*) 'DP/DX(3) '
515
!       WRITE(*,'(1X,12(F9.5,1X))') DPDX(3,:)
516
!       WRITE(*,*) 'DQ/DX(1) '
517
!       WRITE(*,'(1X,12(F9.5,1X))') DQDX(1,:)
518
!       WRITE(*,*) 'DQ/DX(2) '
519
!       WRITE(*,'(1X,12(F9.5,1X))') DQDX(2,:)
520
!       WRITE(*,*) 'DQ/DX(3) '
521
!       WRITE(*,'(1X,12(F9.5,1X))') DQDX(3,:)
522

    
523
! ---  CALCULATE DADX
524
      DO I=1,12
525
       DADX(I)=0.D0
526
       DBDX(I)=0.D0
527
       DCDX(I)=0.D0
528
       DO J=1,3
529
        DADX(I)=DADX(I)+2.D0*P(J)*DPDX(J,I)
530
        DBDX(I)=DBDX(I)+2.D0*Q(J)*DQDX(J,I)
531
        DCDX(I)=DCDX(I)+DPDX(J,I)*Q(J)+P(J)*DQDX(J,I)
532
       ENDDO
533
      ENDDO
534

    
535
! --- CALCULATE FIRST DERIVATIVES OF G
536
      DO I=1,12
537
       DGDX(I)=DADX(I)*QQ+PP*DBDX(I)
538
      ENDDO
539
! --- CALCULATE DDDX
540
      DO I=1,12
541
       DDDX(I)=ABYB2N*DGDX(I)
542
      ENDDO
543
! --- CALCULATE FIRST DERIVATIVES OF J
544
      PREFAC=1.D0/ABYB
545
      DO I=1,12
546
       DJDX(I)=(DCDX(I)*D-PQ*DDDX(I))*PREFAC
547
      ENDDO
548
! --- CALCULATE FIRST DERIVATIVE
549
      DO I=1,12
550
       DI(I)=DJDX(I)*V
551
      ENDDO
552

    
553
!        WRITE(*,*) 'PFL 1st DER TORSION'
554
!        WRITE(*,'(15(1X,F10.6))') DADX(:)
555
!        WRITE(*,'(15(1X,F10.6))') DBDX(:)
556
!        WRITE(*,'(15(1X,F10.6))') DCDX(:)
557
!        WRITE(*,'(15(1X,F10.6))') DGDX(:)
558
!        WRITE(*,'(15(1X,F10.6))') DDDX(:)
559
!        WRITE(*,'(15(1X,F10.6))') DJDX(:)
560
!        WRITE(*,'(15(1X,F10.6))') DI(:)
561

    
562

    
563
! -------------------------------------------------------------
564
! --- CALCULATE B
565
! -------------------------------------------------------------
566

    
567
      WRITE(*,*) 'DER TOR 1'
568
      WRITE(*,*) DI(1:3)
569
      WRITE(*,*) DI(4:6)
570
      WRITE(*,*) DI(7:9)
571
      WRITE(*,*) DI(10:12)
572
 
573
      DO I=1,3
574
         B(I,IAT1)=DI(I)
575
         B(I,IAT2)=DI(I+3)
576
         B(I,IAT3)=DI(I+6)
577
         B(I,IAT4)=DI(I+9)
578
      END DO
579
      END SUBROUTINE CONSTRAINTS_TORSION_DER
580

    
581

    
582
!     ..................................................................
583
      SUBROUTINE CONSTRAINTS_TORSION_DER3(NAT,IAT1,IAT2,IAT3,IAT4,R0,B)
584
!     ******************************************************************
585
!     **
586
!     **  CONSTRAINS THE TORSION OF IAT1-IAT2---IAT3-IAT4
587
!     **
588
!     ******************************************************************
589
      use VarTypes
590
      IMPLICIT none
591
      INTEGER(KINT),INTENT(IN) :: NAT
592
      INTEGER(KINT),INTENT(IN) :: IAT1
593
      INTEGER(KINT),INTENT(IN) :: IAT2
594
      INTEGER(KINT),INTENT(IN) :: IAT3
595
      INTEGER(KINT),INTENT(IN) :: IAT4
596
      REAL(KREAL)   ,INTENT(IN) :: R0(3,NAT)
597
      REAL(KREAL)   ,INTENT(OUT):: B(3,NAT)
598
      REAL(KREAL)               :: X(12), DI(12)
599
      REAL(KREAL)               :: T(3)
600
      REAL(KREAL)  :: Rij(3), rkj(3), rkl(3),rmj(3),rnk(3)
601
      REAL(KREAL)  :: Nrkj, Nrmj, Nrnk
602
      REAL(KREAL) :: SignPhi
603
      integer(KINT)            :: i
604

    
605
      REAL(KREAL)   :: TORS
606
      real(KREAL) :: coverd
607

    
608
      LOGICAL :: Debug=.TRUE.
609

    
610
!     ******************************************************************
611
!
612
      if (debug)    WRITE(*,*) "============= Entering CONSTRAINTS_TORSION_DER3 ==============="
613

    
614
! --- TRANSFER R0 TO X ARRAY TO MAKE HANDLING AND PASSING EASIER
615
    
616
      DO I=1,3
617
        X(I)=R0(I,IAT1)
618
        X(I+3)=R0(I,IAT2)
619
        X(I+6)=R0(I,IAT3)
620
        X(I+9)=R0(I,IAT4)
621
      ENDDO
622

    
623
      if (debug) THEN
624
         WRITE(*,'(1X,I5,3(1X,F10.6))') IAt1,R0(:,IAT1)
625
         WRITE(*,'(1X,I5,3(1X,F10.6))') IAt2,R0(:,IAT2)
626
         WRITE(*,'(1X,I5,3(1X,F10.6))') IAt3,R0(:,IAT3)
627
         WRITE(*,'(1X,I5,3(1X,F10.6))') IAt4,R0(:,IAT4)
628
      END IF
629
!
630
      Rij=X(4:6)-X(1:3)
631
      Rkj=X(4:6)-X(7:9)
632
      Rkl=X(10:12)-X(7:9)
633
      
634
      Call VPRD(Rij,Rkj,Rmj)
635
      Call Vprd(Rkj,Rkl,Rnk)
636
      NRmj=DOT_PRODUCT(Rmj,Rmj)
637
      NRnk=DOT_PRODUCT(Rnk,Rnk)
638
      NRkj=sqrt(DOT_PRODUCT(Rkj,Rkj))
639

    
640
      Coverd=DOT_PRODUCT(Rmj,Rnk)/sqrt(NRmj*NRnk)
641
      CALL VPRD(Rmj,Rnk,T)
642
      SignPhi=Sign(1.d0,DOT_PRODUCT(Rkj,T))
643
!      SignPhi=Sign(1.d0,DOT_PRODUCT(Rkj,T))
644
!      WRITE(*,*) "DBG Tor 2 signphi, dacos(coverd), sign(1.d0,signphi)", &
645
!           signphi, dacos(coverd), sign(1.d0,signphi)
646
!
647
      TORS=-1.d0*DACOS(COVERD)*signPhi
648
      
649
      if (debug) WRITE(*,*) "Torsion derivs 3, Phi=", Tors*180./3.1415926,Iat1,Iat2,Iat3,Iat4
650
!
651
! ----------------------------------------------------------------
652
! --- CALCULATE FIRST DERIVATIVES
653
! ----------------------------------------------------------------
654
      DI(1:3)=NRkj/NRmj*Rmj
655
      DI(10:12)=-NRkj/NRnk*Rnk
656
      DI(4:6)=(DOT_PRODUCT(Rij,Rkj)/Nrkj**2-1.d0)*DI(1:3)-(DOT_PRODUCT(Rkl,Rkj)/NRkj**2)*DI(10:12)
657
      DI(7:9)=(DOT_PRODUCT(Rkl,Rkj)/Nrkj**2-1.d0)*DI(10:12)-(DOT_PRODUCT(Rij,Rkj)/NRkj**2)*DI(1:3)
658

    
659

    
660
      DI=DI*SignPhi
661

    
662
      if (debug) THEN
663
         WRITE(*,*) 'DER TOR 3'
664
         WRITE(*,*) DI(1:3)
665
         WRITE(*,*) DI(4:6)
666
         WRITE(*,*) DI(7:9)
667
         WRITE(*,*) DI(10:12)
668
      END IF
669
 
670
      DO I=1,3
671
         B(I,IAT1)=DI(I)
672
         B(I,IAT2)=DI(I+3)
673
         B(I,IAT3)=DI(I+6)
674
         B(I,IAT4)=DI(I+9)
675
      END DO
676

    
677
      if (debug)    WRITE(*,*) "============= CONSTRAINTS_TORSION_DER3  over ==============="
678

    
679
    END SUBROUTINE CONSTRAINTS_TORSION_DER3
680

    
681

    
682
!     ..................................................................
683
      SUBROUTINE CONSTRAINTS_TORSION_DER2(NAT,IAT1,IAT2,IAT3,IAT4,R0,DZDC)
684
!     ******************************************************************
685
!     **
686
!     **  CONSTRAINS THE TORSION OF IAT1-IAT2---IAT3-IAT4
687
!
688
!          corresponding to P - Q - R - S atoms iun the following
689
! This comes from ADF, and in turns comes from formulas from: "Vibrational States", S.Califano (Wiley, 1976)
690
    !       derivatives w.r.t P: (a*b)B/abs(a*b)**2)
691
    !       w.r.t.  S : (b*c)B/abs(b*c)**2
692
    !       w.r.t.  Q : 1/B**2  x ( C*dS + (a-b)*dP ) * b
693
    !       where * stands for vector product and capitals denote lengths
694
    !       -------------------------------------------------------------
695
    !       a=P-Q
696
    !       b=R-Q
697
    !       theta=angle(a,b)
698
    !       phi (negative of the dihedral angle)
699
    !       c=S-R
700

    
701
!     **
702
!     ******************************************************************
703
      use VarTypes
704
      IMPLICIT none
705
      INTEGER(KINT),INTENT(IN) :: NAT
706
      INTEGER(KINT),INTENT(IN) :: IAT1
707
      INTEGER(KINT),INTENT(IN) :: IAT2
708
      INTEGER(KINT),INTENT(IN) :: IAT3
709
      INTEGER(KINT),INTENT(IN) :: IAT4
710
      REAL(KREAL)   ,INTENT(IN) :: R0(3,NAT)
711
      REAL(KREAL)   ,INTENT(OUT):: DZDC(3,NAT)
712
      REAL(KREAL)               :: A(3),B(3),C(3),Vec(3)
713
      REAL(KREAL)               :: U(3),V(3)
714

    
715
      REAL(KREAL)   :: norm1,norm2,norm3
716
      real(KREAL) :: s,bleng
717

    
718
      REAL(KREAL), PARAMETER :: eps=1e-10
719

    
720
      LOGICAL :: Debug,Failed,valid
721
      EXTERNAL :: Valid
722

    
723
!     ******************************************************************
724
!
725
      debug=valid('TORSION_DER') 
726
      if (debug)    WRITE(*,*) "============= Entering CONSTRAINTS_TORSION_DER2 ==============="
727

    
728
      if (debug) THEN
729
         WRITE(*,'(1X,I5,3(1X,F10.6))') IAt1,R0(:,IAT1)
730
         WRITE(*,'(1X,I5,3(1X,F10.6))') IAt2,R0(:,IAT2)
731
         WRITE(*,'(1X,I5,3(1X,F10.6))') IAt3,R0(:,IAT3)
732
         WRITE(*,'(1X,I5,3(1X,F10.6))') IAt4,R0(:,IAT4)
733
      END IF
734
!
735
      A=R0(:,IAT1)-R0(:,IAT2)
736
      B=R0(:,IAT3)-R0(:,IAT2)
737
      C=R0(:,IAT4)-R0(:,IAT3)
738
      DZDC=0.
739
   
740
      
741
! ----------------------------------------------------------------
742
! --- We check that PQR and QRS are not aligned.
743
! If they are, we return 0. as derivatives instead of just crashing.
744
! ----------------------------------------------------------------
745
    call produit_vect (b(1),b(2),b(3),c(1),c(2),c(3),vec(1),vec(2),vec(3),Norm3)
746
     norm1=sqrt(dot_product(b,b))
747
     norm2=sqrt(dot_product(c,c))
748
    Failed=.FALSE.
749
     
750
     if (Norm3 <= eps) then
751
           WRITE(*,*) '*** WARNING: ILL-DEFINED DIHEDRAL ANGLE(S)'
752
           write (*,*) 'QRS aligned',iat2,iat3,iat4
753
           Failed=.TRUE.
754
     end if
755
     call produit_vect (b(1),b(2),b(3),a(1),a(2),a(3),vec(1),vec(2),vec(3),Norm3)
756
     if (norm3 <= eps) then
757
           WRITE(*,*) '*** WARNING: ILL-DEFINED DIHEDRAL ANGLE(S)'
758
           write (*,*) 'PQR aligned',iat1,iat2,iat3
759
           Failed=.TRUE.
760
     end if
761

    
762
! ----------------------------------------------------------------
763
! --- CALCULATE FIRST DERIVATIVES
764
! ----------------------------------------------------------------
765

    
766
     bleng = sqrt (b(1)**2 + b(2)**2 + b(3)**2)
767
     call produit_vect (a(1),a(2),a(3),b(1),b(2),b(3),vec(1),vec(2),vec(3),Norm3)
768
     s = bleng / (vec(1)**2 + vec(2)**2 + vec(3)**2)
769
!     dzdc(:,ip,3,ip) = s * vec
770
      dzdc(:,iat1) = s * vec
771

    
772
     call produit_vect (b(1),b(2),b(3),c(1),c(2),c(3),vec(1),vec(2),vec(3),Norm3)
773
     s = bleng / (vec(1)**2 + vec(2)**2 + vec(3)**2)
774
!    dzdc(:,is,3,ip) = s * vec
775
     dzdc(:,iat4) = s * vec
776

    
777

    
778
     call produit_vect (c(1),c(2),c(3),   &
779
          dzdc(1,iat4),  dzdc(2,iat4),  dzdc(3,iat4), &
780
          u(1),u(2),u(3),Norm3)
781
     vec = a - b
782
     call produit_vect (vec(1),vec(2),vec(3),    &
783
          dzdc(1,iat1),  dzdc(2,iat1),  dzdc(3,iat1),  &
784
          v(1),v(2),v(3),Norm3)
785
     vec = u + v
786
     call produit_vect (vec(1),vec(2),vec(3),    &
787
          b(1), b(2), b(3), u(1),u(2),u(3),Norm3)
788

    
789
!     dzdc(:,iq,3,ip) = u / bleng**2
790
!     dzdc(:,ir,3,ip) = - (dzdc(:,ip,3,ip) + dzdc(:,iq,3,ip) + dzdc(:,is,3,ip))
791
     dzdc(:,iat2) = u / bleng**2
792
     dzdc(:,iat3) = - (dzdc(:,iat1) + dzdc(:,iat2) + dzdc(:,iat4))
793

    
794
      if (debug) THEN
795
         WRITE(*,*) 'DER TOR 2'
796
         WRITE(*,*) DZDC(:,IAT1)
797
         WRITE(*,*) DZDC(:,IAT2)
798
         WRITE(*,*) DZDC(:,IAT3)
799
         WRITE(*,*) DZDC(:,IAT4)
800
      END IF
801
 
802
      if (Failed) DZDC=0.
803

    
804
      if (debug)    WRITE(*,*) "============= CONSTRAINTS_TORSION_DER2  over ==============="
805

    
806
    END SUBROUTINE CONSTRAINTS_TORSION_DER2
807

    
808
! Not yet implemented
809

    
810

    
811
! !     ..................................................................
812
!       SUBROUTINE CONSTRAINTS_COGSEP_DER(NAT,GROUP1,GROUP2,R0,RMass,B)
813
! !     ******************************************************************
814
! !     **  SET UP CONSTRAINT IN COORDINATE SPACE                       **
815
! !     **  CONSTRAINT IS THE DISTANCE BETWEEN THE CENTER OF GRAVITIES  **
816
! !     **  OF TWO NAMED GROUPS                                         **
817
! !     **                                                              **
818
! !     **  THE CONSTRAINT IS EXPRESSED AS A+B*R+0.5*R*C*R=0            **
819
! !     **                                                              **
820
! !     ******************************************************************
821
!        use Geometry_module
822

    
823
!        Use VarTypes
824
!       IMPLICIT NONE
825
!       INTEGER(KINT)  ,INTENT(IN) :: NAT
826
!       CHARACTER(*),INTENT(IN) :: GROUP1  ! NAME OF FIRST GROUP
827
!       CHARACTER(*),INTENT(IN) :: GROUP2  ! NAME OF SECOND GROUP
828
!       REAL(KREAL)     ,INTENT(OUT):: B(3,NAT)
829
!       REAL(KREAL)     ,INTENT(IN) :: R0(3,NAT)
830
!       REAL(KREAL)                 :: RMASS(NAT),A
831
!       LOGICAL(4)              :: TMEMBER1(NAT)
832
!       LOGICAL(4)              :: TMEMBER2(NAT)
833
!       REAL(KREAL)                 :: SUMMASS1,SUMMASS2
834
!       REAL(KREAL)                 :: COG1(3)   ! COG OF FIRST GROUP
835
!       REAL(KREAL)                 :: COG2(3)   ! COG OF SECOND GROUP
836
!       REAL(KREAL)                 :: X(3)      ! COG1-COG2
837
!       REAL(KREAL)                 :: COGSEP    
838
!       REAL(KREAL)                 :: DADX(3)   
839
!       INTEGER(KINT)              :: I,J,IAT,IAT1,IAT2
840
! !     ******************************************************************
841

    
842
!       CALL GROUPLIST_MEMBERS(GROUP1,TMEMBER1)
843
!       CALL GROUPLIST_MEMBERS(GROUP2,TMEMBER2)
844
! !
845
! !     ==================================================================
846
! !     ==  CALCULATE CENTERS OF GRAVITY                                ==
847
! !     ==================================================================
848
!       SUMMASS1=0.D0
849
!       SUMMASS2=0.D0
850
!       COG1(:)=0.D0
851
!       COG2(:)=0.D0
852
!       DO IAT=1,NAT
853
!         IF(TMEMBER1(IAT)) THEN
854
!           SUMMASS1=SUMMASS1+RMASS(IAT)
855
!           COG1(:)=COG1(:)+RMASS(IAT)*R0(:,IAT)
856
!         ENDIF
857
!         IF(TMEMBER2(IAT)) THEN
858
!           SUMMASS2=SUMMASS2+RMASS(IAT)
859
!           COG2(:)=COG2(:)+RMASS(IAT)*R0(:,IAT)
860
!         ENDIF
861
!       ENDDO
862
!       COG1(:)=COG1(:)/SUMMASS1
863
!       COG2(:)=COG2(:)/SUMMASS2
864
! !
865
! !     =============================================================
866
! !     ==  SEPARATION BETWEEN THE CENTERS OF GRAVITY              ==
867
! !     =============================================================
868
!       X(:)=COG1(:)-COG2(:)
869
! !
870
! !     =============================================================
871
! !     ==  CALCULATE A(X) AND ITS DERIVATIVES                     ==
872
! !     =============================================================
873
!       A=SQRT(DOT_PRODUCT(X(:),X(:)))
874
!       DADX(:)=X(:)/A
875

    
876
!       DO IAT1=1,NAT
877
!          IF(TMEMBER1(IAT1)) THEN
878
!             B(:,IAT1)=B(:,IAT1)+DADX(:)*RMASS(IAT1)/SUMMASS1
879
!          END IF
880
!          IF(TMEMBER2(IAT1)) THEN
881
!             B(:,IAT1)=B(:,IAT1)-DADX(:)*RMASS(IAT1)/SUMMASS2
882
!          END IF
883
!       END DO
884

    
885
!       END SUBROUTINE CONSTRAINTS_COGSEP_DER
886

    
887

    
888
! !     ..................................................................
889
!       SUBROUTINE CONSTRAINTS_COGANG_DER(NAT,GROUP1,GROUP2,GROUP3,R0,RMass,B)
890
! !     ******************************************************************
891
! !     **  SET UP CONSTRAINT IN COORDINATE SPACE                       **
892
! !     **  CONSTRAINT IS THE ANGLE    BETWEEN THE CENTER OF GRAVITIES  **
893
! !     **  OF THREE NAMED GROUPS                                       **
894
! !     **                                                              **
895
! !     **  THE CONSTRAINT IS EXPRESSED AS A+B*R+0.5*R*C*R=0            **
896
! ! This makes heavy use of BONDANGLEDERIVS
897
! !     **                                                              **
898
! !     ******************************************************************
899
!       IMPLICIT NONE
900
!       INTEGER(KINT)  ,INTENT(IN) :: NAT
901
!       CHARACTER(*),INTENT(IN) :: GROUP1  ! NAME OF FIRST GROUP
902
!       CHARACTER(*),INTENT(IN) :: GROUP2  ! NAME OF SECOND GROUP
903
!       CHARACTER(*),INTENT(IN) :: GROUP3  ! NAME OF THIRD GROUP
904
!       REAL(KREAL)     ,INTENT(OUT):: B(3,NAT)
905
!       REAL(KREAL)     ,INTENT(IN) :: R0(3,NAT)
906
!       REAL(KREAL)                 :: RMASS(NAT)
907
!       LOGICAL(4)              :: TMEMBER1(NAT)
908
!       LOGICAL(4)              :: TMEMBER2(NAT)
909
!       LOGICAL(4)              :: TMEMBER3(NAT)
910
!       REAL(KREAL)                 :: SUMMASS1,SUMMASS2, SUMMASS3
911
!       REAL(KREAL)                 :: COG1(3)   ! COG OF FIRST GROUP
912
!       REAL(KREAL)                 :: COG2(3)   ! COG OF SECOND GROUP
913
!       REAL(KREAL)                 :: COG3(3)   ! COG OF THIRD GROUP
914
!       REAL(KREAL)            :: DCOG1(3,NAT), DCOG2(3,NAT), DCOG3(3,NAT)
915
!       REAL(KREAL)                 :: vecU(3), vecV(3)
916
!       REAL(KREAL)                 :: UU, VV, UV
917
!       INTEGER(KINT)              :: I,J,IAT,IAT1,IAT2,JAT
918
!       REAL(KREAL)           :: Pi
919
! ! From BONDANGLEDERIVS
920
!       REAL(KREAL)               :: D    ! 1/SQRT(A*B)
921
!       REAL(KREAL)               :: DUUDX(3,NAT)
922
!       REAL(KREAL)               :: DVVDX(3,NAT)
923
!       REAL(KREAL)               :: DUVDX(3,NAT)
924
!       REAL(KREAL)               :: DDDX(3,NAT)
925
!       REAL(KREAL)               :: DJDX(3,NAT)
926
!       REAL(KREAL)               :: DVDX(3,NAT)
927
!       REAL(KREAL)               :: DGDX(3,NAT)    ! D(A*B)/DX
928
!       REAL(KREAL)               :: AB,AB2N,V,PREFAC
929
!       REAL(KREAL)               :: PHI,COSPHI,SINPHI
930
!       REAL(KREAL)               :: DIVDX,AI
931

    
932

    
933
! !     ******************************************************************
934
!       CALL GROUPLIST_MEMBERS(GROUP1,TMEMBER1)
935
!       CALL GROUPLIST_MEMBERS(GROUP2,TMEMBER2)
936
!       CALL GROUPLIST_MEMBERS(GROUP3,TMEMBER3)
937
! !
938
! !     ==================================================================
939
! !     ==  CALCULATE CENTERS OF GRAVITY                                ==
940
! !     ==================================================================
941
!       SUMMASS1=0.D0
942
!       SUMMASS2=0.D0
943
!       SUMMASS3=0.D0
944
!       COG1(:)=0.D0
945
!       COG2(:)=0.D0
946
!       COG3(:)=0.D0
947
!       DO IAT=1,NAT
948
!         IF(TMEMBER1(IAT)) THEN
949
!           SUMMASS1=SUMMASS1+RMASS(IAT)
950
!           COG1(:)=COG1(:)+RMASS(IAT)*R0(:,IAT)
951
!         ENDIF
952
!         IF(TMEMBER2(IAT)) THEN
953
!           SUMMASS2=SUMMASS2+RMASS(IAT)
954
!           COG2(:)=COG2(:)+RMASS(IAT)*R0(:,IAT)
955
!         ENDIF
956
!         IF(TMEMBER3(IAT)) THEN
957
!           SUMMASS3=SUMMASS3+RMASS(IAT)
958
!           COG3(:)=COG3(:)+RMASS(IAT)*R0(:,IAT)
959
!         ENDIF
960

    
961
!       ENDDO
962
!       COG1(:)=COG1(:)/SUMMASS1
963
!       COG2(:)=COG2(:)/SUMMASS2
964
!       COG3(:)=COG3(:)/SUMMASS3
965
! !
966
! !
967
!       vecU(:)=COG1(:)-COG2(:)
968
!       vecV(:)=COG3(:)-COG2(:)
969

    
970
!       UV=DOT_PRODUCT(vecU,vecV)
971
!       UU=DOT_PRODUCT(vecU,vecU)
972
!       VV=DOT_PRODUCT(vecV,vecV)
973
! ! We now calculate the first derivatives of the COG vs the atomic
974
! ! Coords... 
975
!       DCOG1=0.
976
!       DCOG2=0.
977
!       DCOG3=0.
978
!       DO IAT=1,NAT
979
!         IF (TMEMBER1(IAT))  DCOG1(:,IAT)=RMASS(IAT)/SUMMASS1
980
!         IF (TMEMBER2(IAT))  DCOG2(:,IAT)=RMASS(IAT)/SUMMASS2
981
!         IF (TMEMBER3(IAT))  DCOG3(:,IAT)=RMASS(IAT)/SUMMASS3
982
!       END DO
983
! ! We now calculate the first derivatives of UU, VV, UV
984
!       DUUDX=0.
985
!       DVVDX=0.
986
!       DUVDX=0.
987
!       DO IAT=1,NAT
988
!         DUUDX(:,IAT)=2.D0*(COG1(:)-COG2(:))*(DCOG1(:,IAT)-DCOG2(:,IAT))
989
!         DVVDX(:,IAT)=2.D0*(COG3(:)-COG2(:))*(DCOG3(:,IAT)-DCOG2(:,IAT))
990
!         DUVDX(:,IAT)=(COG3(:)-COG2(:))*(DCOG1(:,IAT)-DCOG2(:,IAT)) +    &
991
!         (DCOG3(:,IAT)-DCOG2(:,IAT))*(COG1(:)-COG2(:))
992
!       END DO
993

    
994
! !       WRITE(*,*) 'PFL DER COGANGLE'
995
! !       WRITE(*,'(15(1X,F10.6))') DUUDX(:,:)
996
! !       WRITE(*,'(15(1X,F10.6))') DVVDX(:,:)
997
! !       WRITE(*,'(15(1X,F10.6))') DUVDX(:,:)
998

    
999

    
1000
! ! We now go for the hard part: calculate Phi, dPhi/dx, d2Phi/dxidxj
1001
!       Pi=DACOS(-1.D0)
1002
!       AB=UU*VV
1003
!       D=SQRT(AB)
1004
!       COSPHI=UV/D       !=COS(PHI)
1005
!       PHI=DACOS(COSPHI)
1006
! !      WRITE(*,*) 'COGANGLE Phi=',Phi,Phi*180/Pi
1007
! !     --- CALCULATE DERIVATIVE OF ARCCOS
1008
!       SINPHI=DSQRT(1.D0-COSPHI**2)
1009
!       V=-1.D0/SINPHI                 !=-SIN(PHI)=D(ACOS(PHI))/D(COSPHI)
1010
! ! --- CALCULATE FIRST DERIVATIVES OF G=A*B
1011
!       DGDX(:,:)=DUUDX(:,:)*VV+DVVDX(:,:)*UU
1012
! ! --- CALCULATE DDDX
1013
!       AB2N=0.5D0/D
1014
!       DDDX(:,:)=AB2N*DGDX(:,:)
1015
! ! --- CALCULATE FIRST DERIVATIVES OF J = COS(PHI)
1016
!       PREFAC=1.D0/AB
1017
!       DJDX(:,:)=(DUVDX(:,:)*D-UV*DDDX(:,:))*PREFAC
1018
! ! --- CALCULATE FIRST DERIVATIVE
1019
!       B(:,:)=DJDX(:,:)*V
1020

    
1021
!       END SUBROUTINE CONSTRAINTS_COGANG_DER
1022

    
1023
! !     ..................................................................
1024
!       SUBROUTINE CONSTRAINTS_COGTORSION_DER(NAT,GROUP1,GROUP2, GROUP3,GROUP4, R0,RMAss,B)
1025
! !     ******************************************************************
1026
! !     **  SET UP CONSTRAINT IN COORDINATE SPACE                       **
1027
! !     **  CONSTRAINT IS THE ANGLE    BETWEEN THE CENTER OF GRAVITIES  **
1028
! !     **  OF THREE NAMED GROUPS                                       **
1029
! !     **                                                              **
1030
! !     **  THE CONSTRAINT IS EXPRESSED AS A+B*R+0.5*R*C*R=0            **
1031
! ! This makes heavy use of BONDANGLEDERIVS
1032
! !     **                                                              **
1033
! !     ******************************************************************
1034
!       IMPLICIT NONE
1035
!       INTEGER(KINT)  ,INTENT(IN) :: NAT
1036
!       CHARACTER(*),INTENT(IN) :: GROUP1  ! NAME OF FIRST GROUP
1037
!       CHARACTER(*),INTENT(IN) :: GROUP2  ! NAME OF SECOND GROUP
1038
!       CHARACTER(*),INTENT(IN) :: GROUP3  ! NAME OF THIRD GROUP
1039
!       CHARACTER(*),INTENT(IN) :: GROUP4  ! NAME OF FOURTH GROUP
1040
!       REAL(KREAL)     ,INTENT(OUT):: B(3,NAT)
1041
!       REAL(KREAL)     ,INTENT(IN) :: R0(3,NAT)
1042
!       REAL(KREAL)                 :: RMASS(NAT)
1043
!       LOGICAL(4)              :: TMEMBER1(NAT)
1044
!       LOGICAL(4)              :: TMEMBER2(NAT)
1045
!       LOGICAL(4)              :: TMEMBER3(NAT)
1046
!       LOGICAL(4)              :: TMEMBER4(NAT)
1047
!       REAL(KREAL)                 :: SUMMASS1,SUMMASS2, SUMMASS3, SUMMASS4
1048
!       REAL(KREAL)                 :: COG1(3)   ! COG OF FIRST GROUP
1049
!       REAL(KREAL)                 :: COG2(3)   ! COG OF SECOND GROUP
1050
!       REAL(KREAL)                 :: COG3(3)   ! COG OF THIRD GROUP
1051
!       REAL(KREAL)                 :: COG4(3)   ! COG OF FOURTH GROUP
1052
! ! Even thought D COG(i)/DXj =0 (ie COG_x depends only on x-coordinates
1053
! ! of the atoms, not on y or z), we store a 3,3*NAT matrix to make 
1054
! ! handling easier below
1055
!       REAL(KREAL)                 :: DCOG1(3,3,NAT), DCOG2(3,3,NAT),              &
1056
!                                 DCOG3(3,3,NAT), DCOG4(3,3,NAT)
1057
!       REAL(KREAL)                 :: T(3), U(3), R(3), W(3), P(3), Q(3)
1058
!       REAL(KREAL)                 :: PP, QQ, PQ
1059
!       INTEGER(KINT)              :: I,J,IAT,IAT1,IAT2,JAT
1060
!       REAL(KREAL)                 :: Pi, Phi
1061
! ! From BONDANGLEDERIVS
1062
!       REAL(KREAL)               :: D    ! SQRT(PP*QQ)
1063
!       REAL(KREAL)               :: DPDX(3,3,NAT)
1064
!       REAL(KREAL)               :: DQDX(3,3,NAT)
1065
!       REAL(KREAL)               :: DPPDX(3,NAT)
1066
!       REAL(KREAL)               :: DQQDX(3,NAT)
1067
!       REAL(KREAL)               :: DPQDX(3,NAT)
1068
!       REAL(KREAL)               :: DDDX(3,NAT)
1069
!       REAL(KREAL)               :: DJDX(3,NAT) ! D(C/D)/DX*D2
1070
!       REAL(KREAL)               :: DVDX(3,NAT)
1071
!       REAL(KREAL)               :: DGDX(3,NAT)    ! D(A*B)/DX
1072
!       REAL(KREAL)               :: AB,AB2N,V,PREFAC,ABYB,ABYB2N
1073
!       REAL(KREAL)               :: COSPHI,SINPHI
1074

    
1075

    
1076
! !     ******************************************************************
1077

    
1078
! !      WRITE(*,*) 'ENTERING SET_COGTORSION'
1079
! !      WRITE(*,*) 'G1:',TRIM(GROUP1),' G2:',TRIM(GROUP2),' G3:',TRIM(GROUP3),' G4:',TRIM(GROUP4)
1080

    
1081
! !     ==================================================================
1082
! !     == DETERMINE GROUPS AND ATOM-MASSES
1083
! !     ==================================================================
1084

    
1085
!       CALL GROUPLIST_MEMBERS(GROUP1,TMEMBER1)
1086
!       CALL GROUPLIST_MEMBERS(GROUP2,TMEMBER2)
1087
!       CALL GROUPLIST_MEMBERS(GROUP3,TMEMBER3)
1088
!       CALL GROUPLIST_MEMBERS(GROUP4,TMEMBER4)
1089

    
1090
! !
1091
! !     ==================================================================
1092
! !     == CALCULATE CENTER OF GRAVITY FOR GROUPS
1093
! !     ==================================================================
1094
!       SUMMASS1=0.D0
1095
!       SUMMASS2=0.D0
1096
!       SUMMASS3=0.D0
1097
!       SUMMASS4=0.D0
1098

    
1099
!       COG1=0.D0
1100
!       COG2=0.D0
1101
!       COG3=0.D0
1102
!       COG4=0.D0
1103

    
1104
!       DO IAT=1,NAT
1105
! ! One atom can belong to many groups as we are looking at
1106
! ! center of gravity of the groups
1107
!         IF (TMEMBER1(IAT)) THEN
1108
!           SUMMASS1=SUMMASS1+RMASS(IAT)
1109
!           COG1(:)=COG1(:)+RMASS(IAT)*R0(:,IAT)
1110
!         END IF
1111
!         IF (TMEMBER2(IAT)) THEN
1112
!           SUMMASS2=SUMMASS2+RMASS(IAT)
1113
!           COG2(:)=COG2(:)+RMASS(IAT)*R0(:,IAT)
1114
!         END IF
1115
!         IF (TMEMBER3(IAT)) THEN
1116
!           SUMMASS3=SUMMASS3+RMASS(IAT)
1117
!           COG3(:)=COG3(:)+RMASS(IAT)*R0(:,IAT)
1118
!         END IF
1119
!         IF (TMEMBER4(IAT)) THEN
1120
!           SUMMASS4=SUMMASS4+RMASS(IAT)
1121
!           COG4(:)=COG4(:)+RMASS(IAT)*R0(:,IAT)
1122
!         END IF
1123
!       ENDDO
1124

    
1125
!       COG1(:)=COG1(:)/SUMMASS1
1126
!       COG2(:)=COG2(:)/SUMMASS2
1127
!       COG3(:)=COG3(:)/SUMMASS3
1128
!       COG4(:)=COG4(:)/SUMMASS4
1129
!       T(:)=COG1(:)-COG2(:)
1130
!       U(:)=COG3(:)-COG2(:)
1131
!       R(:)=COG4(:)-COG3(:)
1132
!       W(:)=COG2(:)-COG3(:)
1133

    
1134
!       CALL VPRD(T,U,P)                     ! P=CROSS-PRODUCT(T,U)
1135
!       CALL VPRD(W,R,Q)
1136
!       PQ=DOT_PRODUCT(P,Q)
1137
!       PP=DOT_PRODUCT(P,P)
1138
!       QQ=DOT_PRODUCT(Q,Q)
1139
!       ABYB=PP*QQ
1140
!       D=SQRT(ABYB)
1141

    
1142
!       ABYB2N=0.5D0/D
1143
!       COSPHI=PQ/D
1144
!       Phi=DACOS(COSPHI)
1145
!       Pi=DACOS(-1.0D0)
1146
!       SINPHI=DSQRT(1.D0-(COSPHI)**2.D0)
1147

    
1148
! ! --- CALCULATE DERIVATIVE OF ARCCOS
1149
!       IF (dabs(DABS(COSPHI)-1.d0) < 1.D-9) THEN
1150
!         WRITE(*,*)   'WARNING: DERIVATIVE OF TORSION AT ZERO AND PI IS 
1151
!      &UNDEFINED!'
1152
!         V=0.D0
1153
!       ELSE
1154
!         V=-1.D0/SINPHI
1155
!       ENDIF
1156

    
1157
! !      WRITE(*,*)  'Phi et al',Phi,Phi*180./Pi,cosphi,sinphi,V
1158

    
1159
! ! We now calculate the first derivatives of the COG vs the atomic
1160
! ! Coords... 
1161
!       DCOG1=0.
1162
!       DCOG2=0.
1163
!       DCOG3=0.
1164
!       DCOG4=0.
1165
!       DO IAT=1,NAT
1166
!         IF (TMEMBER1(IAT))  THEN
1167
!            DO J=1,3
1168
!               DCOG1(J,J,IAT)=RMASS(IAT)/SUMMASS1
1169
!            END DO
1170
!         END IF
1171
!         IF (TMEMBER2(IAT))  THEN
1172
!            DO J=1,3
1173
!               DCOG2(J,J,IAT)=RMASS(IAT)/SUMMASS2
1174
!            END DO
1175
!         END IF
1176
!         IF (TMEMBER3(IAT))  THEN
1177
!            DO J=1,3
1178
!               DCOG3(J,J,IAT)=RMASS(IAT)/SUMMASS3
1179
!            END DO
1180
!         END IF
1181
!         IF (TMEMBER4(IAT))  THEN
1182
!            DO J=1,3
1183
!               DCOG4(J,J,IAT)=RMASS(IAT)/SUMMASS4
1184
!            END DO
1185
!         END IF
1186
!       END DO
1187

    
1188
! ! We now calculate the first derivatives of P, Q
1189
!       DPDX=0.
1190
!       DQDX=0.
1191
!       DO I=1,NAT
1192
!          DO J=1,3
1193
!             DPDX(1,J,I)=(DCOG1(2,J,I)-DCOG2(2,J,I))*(COG3(3)-COG2(3))     &
1194
!      &                 +(COG1(2)-COG2(2))*(DCOG3(3,J,I)-DCOG2(3,J,I))     &
1195
!      &                 -(DCOG1(3,J,I)-DCOG2(3,J,I))*(COG3(2)-COG2(2))     &
1196
!      &                 -(COG1(3)-COG2(3))*(DCOG3(2,J,I)-DCOG2(2,J,I)) 
1197
!             DPDX(2,J,I)=(DCOG1(3,J,I)-DCOG2(3,J,I))*(COG3(1)-COG2(1))     &
1198
!      &                 +(COG1(3)-COG2(3))*(DCOG3(1,J,I)-DCOG2(1,J,I))     &
1199
!      &                 -(DCOG1(1,J,I)-DCOG2(1,J,I))*(COG3(3)-COG2(3))     &
1200
!      &                 -( COG1(1)-COG2(1))*(DCOG3(3,J,I)-DCOG2(3,J,I)) 
1201
!             DPDX(3,J,I)=(DCOG1(1,J,I)-DCOG2(1,J,I))*(COG3(2)-COG2(2))     &
1202
!      &                 +(COG1(1)-COG2(1))*(DCOG3(2,J,I)-DCOG2(2,J,I))     &
1203
!      &                 -(DCOG1(2,J,I)-DCOG2(2,J,I))*(COG3(1)-COG2(1))     &
1204
!      &                 -(COG1(2)-COG2(2))*(DCOG3(1,J,I)-DCOG2(1,J,I)) 
1205
!             DQDX(1,J,I)=(DCOG2(2,J,I)-DCOG3(2,J,I))*(COG4(3)-COG3(3))     &
1206
!      &                 +(COG2(2)-COG3(2))*(DCOG4(3,J,I)-DCOG3(3,J,I))     &
1207
!      &                 -(DCOG2(3,J,I)-DCOG3(3,J,I))*(COG4(2)-COG3(2))     &
1208
!      &                 -(COG2(3)-COG3(3))*(DCOG4(2,J,I)-DCOG3(2,J,I)) 
1209
!             DQDX(2,J,I)=(DCOG2(3,J,I)-DCOG3(3,J,I))*(COG4(1)-COG3(1))     &
1210
!      &                 +(COG2(3)-COG3(3))*(DCOG4(1,J,I)-DCOG3(1,J,I))     &
1211
!      &                 -(DCOG2(1,J,I)-DCOG3(1,J,I))*(COG4(3)-COG3(3))     &
1212
!      &                 -( COG2(1)-COG3(1))*(DCOG4(3,J,I)-DCOG3(3,J,I)) 
1213
!             DQDX(3,J,I)=(DCOG2(1,J,I)-DCOG3(1,J,I))*(COG4(2)-COG3(2))     &
1214
!      &                 +(COG2(1)-COG3(1))*(DCOG4(2,J,I)-DCOG3(2,J,I))     &
1215
!      &                 -(DCOG2(2,J,I)-DCOG3(2,J,I))*(COG4(1)-COG3(1))     &
1216
!      &                 -(COG2(2)-COG3(2))*(DCOG4(1,J,I)-DCOG3(1,J,I)) 
1217

    
1218
!          END DO
1219
!       END DO
1220

    
1221

    
1222
! !       WRITE(*,*) 'DP/DX(1) '
1223
! !       WRITE(*,'(1X,12(F9.5,1X))') DPDX(1,1:3,5:8)
1224
! !       WRITE(*,*) 'DP/DX(2) '
1225
! !       WRITE(*,'(1X,12(F9.5,1X))') DPDX(2,1:3,5:8)
1226
! !       WRITE(*,*) 'DP/DX(3) '
1227
! !       WRITE(*,'(1X,12(F9.5,1X))') DPDX(3,1:3,5:8)
1228
! !       WRITE(*,*) 'DQ/DX(1) '
1229
! !       WRITE(*,'(1X,12(F9.5,1X))') DQDX(1,1:3,5:8)
1230
! !       WRITE(*,*) 'DQ/DX(2) '
1231
! !       WRITE(*,'(1X,12(F9.5,1X))') DQDX(2,1:3,5:8)
1232
! !       WRITE(*,*) 'DQ/DX(3) '
1233
! !       WRITE(*,'(1X,12(F9.5,1X))') DQDX(3,1:3,5:8)
1234

    
1235
! ! We now calculate the first derivatives of PP, QQ and PQ and DD/DX*2D
1236
      
1237
!       DPPDX=0.
1238
!       DQQDX=0.
1239
!       DPQDX=0.
1240
!       PREFAC=1.d0/ABYB
1241
!       DO IAT=1,NAT
1242
!          DO J=1,3
1243
!             DPPDX(J,IAT)=2.D0*(P(1)*DPDX(1,J,IAT)+P(2)*DPDX(2,J,IAT)+
1244
!      &                P(3)*DPDX(3,J,IAT))
1245
!             DQQDX(J,IAT)=2.D0*(Q(1)*DQDX(1,J,IAT)+Q(2)*DQDX(2,J,IAT)+ 
1246
!      &                Q(3)*DQDX(3,J,IAT))
1247
!             DPQDX(J,IAT)=P(1)*DQDX(1,J,IAT)+DPDX(1,J,IAT)*Q(1)  
1248
!      &               +P(2)*DQDX(2,J,IAT)+DPDX(2,J,IAT)*Q(2)
1249
!      &               +P(3)*DQDX(3,J,IAT)+DPDX(3,J,IAT)*Q(3) 
1250
!             DGDX(J,IAT)=DPPDX(J,IAT)*QQ+PP*DQQDX(J,IAT)
1251
!             DDDX(J,IAT)=DGDX(J,IAT)*ABYB2N
1252
!             DJDX(J,IAT)=(DPQDX(J,IAT)*D-PQ*DDDX(J,IAT))*PREFAC
1253
!             B(J,IAT)=DJDX(J,IAT)*V
1254
!          END DO
1255
!       END DO
1256
      
1257
!       END SUBROUTINE CONSTRAINTS_COGTORSION_DER
1258

    
1259

    
1260
! !     ..................................................................
1261
!       SUBROUTINE CONSTRAINTS_COGDIFF_DER(NAT,GROUP1,GROUP2,   
1262
!      &    GROUP3,R0,RMAss,B)
1263
! !     ******************************************************************
1264
! !     **  SET UP CONSTRAINT IN COORDINATE SPACE                       **
1265
! !     **  CONSTRAINT IS THE DIFFERENCE BETWEEN THE DISTANCES BETWEEN  **
1266
! !     **  THE CENTER OF GRAVITIES  OF THREE NAMED GROUPS              **
1267
! !     **                                                              **
1268
! !     **  THE CONSTRAINT IS EXPRESSED AS A+B*R+0.5*R*C*R=0            **
1269
! !     **                                                              **
1270
! !     ******************************************************************
1271
!       IMPLICIT NONE
1272
!       INTEGER(KINT)  ,INTENT(IN) :: NAT
1273
!       CHARACTER(*),INTENT(IN) :: GROUP1  ! NAME OF FIRST GROUP
1274
!       CHARACTER(*),INTENT(IN) :: GROUP2  ! NAME OF SECOND GROUP
1275
!       CHARACTER(*),INTENT(IN) :: GROUP3  ! NAME OF THIRD GROUP
1276
!       REAL(KREAL)     ,INTENT(OUT):: B(3,NAT)
1277
!       REAL(KREAL)     ,INTENT(IN) :: R0(3,NAT)
1278
!       REAL(KREAL)     ,INTENT(IN) :: RMASS(NAT)
1279
!       LOGICAL(4)              :: TMEMBER1(NAT)
1280
!       LOGICAL(4)              :: TMEMBER2(NAT)
1281
!       LOGICAL(4)              :: TMEMBER3(NAT)
1282
!       REAL(KREAL)                 :: SUMMASS1,SUMMASS2, SUMMASS3
1283
!       REAL(KREAL)                 :: COG1(3)   ! COG OF FIRST GROUP
1284
!       REAL(KREAL)                 :: COG2(3)   ! COG OF SECOND GROUP
1285
!       REAL(KREAL)                 :: COG3(3)   ! COG OF THIRD GROUP
1286
!       REAL(KREAL)            :: DCOG1(3,NAT), DCOG2(3,NAT), DCOG3(3,NAT)
1287
!       REAL(KREAL)                 :: vecU(3), vecV(3)
1288
!       REAL(KREAL)                 :: UU, VV, UV
1289
!       INTEGER(KINT)              :: I,J,IAT,IAT1,IAT2,JAT
1290
!       REAL(KREAL)           :: Pi
1291
! ! From BONDANGLEDERIVS
1292
!       REAL(KREAL)               :: D1,D2 
1293
!       REAL(KREAL)               :: DUUDX(3,NAT)
1294
!       REAL(KREAL)               :: DVVDX(3,NAT)
1295
!       REAL(KREAL)               :: DUVDX(3,NAT)
1296
!       REAL(KREAL)               :: DDDX(3,NAT)
1297
!       REAL(KREAL)               :: DJDX(3,NAT)
1298
!       REAL(KREAL)               :: DVDX(3,NAT)
1299
!       REAL(KREAL)               :: DGDX(3,NAT)    ! D(A*B)/DX
1300
!       REAL(KREAL)               :: DIDJ(3,NAT,3,NAT), DI(3,NAT)
1301
!       REAL(KREAL)               :: AB,AB2N,V,PREFAC
1302
!       REAL(KREAL)               :: PHI,COSPHI,SINPHI
1303
!       REAL(KREAL)               :: DIVDX,AI
1304

    
1305

    
1306
! !     ******************************************************************
1307

    
1308
!       CALL GROUPLIST_MEMBERS(GROUP1,TMEMBER1)
1309
!       CALL GROUPLIST_MEMBERS(GROUP2,TMEMBER2)
1310
!       CALL GROUPLIST_MEMBERS(GROUP3,TMEMBER3)
1311
! !
1312
! !     ==================================================================
1313
! !     ==  CALCULATE CENTERS OF GRAVITY                                ==
1314
! !     ==================================================================
1315
!       SUMMASS1=0.D0
1316
!       SUMMASS2=0.D0
1317
!       SUMMASS3=0.D0
1318
!       COG1(:)=0.D0
1319
!       COG2(:)=0.D0
1320
!       COG3(:)=0.D0
1321
!       DO IAT=1,NAT
1322
!         IF(TMEMBER1(IAT)) THEN
1323
!           SUMMASS1=SUMMASS1+RMASS(IAT)
1324
!           COG1(:)=COG1(:)+RMASS(IAT)*R0(:,IAT)
1325
!         ENDIF
1326
!         IF(TMEMBER2(IAT)) THEN
1327
!           SUMMASS2=SUMMASS2+RMASS(IAT)
1328
!           COG2(:)=COG2(:)+RMASS(IAT)*R0(:,IAT)
1329
!         ENDIF
1330
!         IF(TMEMBER3(IAT)) THEN
1331
!           SUMMASS3=SUMMASS3+RMASS(IAT)
1332
!           COG3(:)=COG3(:)+RMASS(IAT)*R0(:,IAT)
1333
!         ENDIF
1334

    
1335
!       ENDDO
1336
!       COG1(:)=COG1(:)/SUMMASS1
1337
!       COG2(:)=COG2(:)/SUMMASS2
1338
!       COG3(:)=COG3(:)/SUMMASS3
1339
! !
1340
! !
1341
!       vecU(:)=COG1(:)-COG2(:)
1342
!       vecV(:)=COG3(:)-COG2(:)
1343

    
1344
!       UU=DOT_PRODUCT(vecU,vecU)
1345
!       VV=DOT_PRODUCT(vecV,vecV)
1346
! ! We now calculate the first derivatives of the COG vs the atomic
1347
! ! Coords... 
1348
!       DCOG1=0.
1349
!       DCOG2=0.
1350
!       DCOG3=0.
1351
!       DO IAT=1,NAT
1352
!         IF (TMEMBER1(IAT))  DCOG1(:,IAT)=RMASS(IAT)/SUMMASS1
1353
!         IF (TMEMBER2(IAT))  DCOG2(:,IAT)=RMASS(IAT)/SUMMASS2
1354
!         IF (TMEMBER3(IAT))  DCOG3(:,IAT)=RMASS(IAT)/SUMMASS3
1355
!       END DO
1356
! ! We now calculate the first derivatives of UU, VV
1357
!       DUUDX=0.
1358
!       DVVDX=0.
1359
!       DO IAT=1,NAT
1360
!         DUUDX(:,IAT)=2.D0*(COG1(:)-COG2(:))*(DCOG1(:,IAT)-DCOG2(:,IAT))
1361
!         DVVDX(:,IAT)=2.D0*(COG3(:)-COG2(:))*(DCOG3(:,IAT)-DCOG2(:,IAT))
1362
!       END DO
1363
! !       WRITE(*,*) 'PFL DER COGDIFF'
1364
! !       WRITE(*,'(15(1X,F10.6))') DUUDX(:,:)
1365
! !       WRITE(*,'(15(1X,F10.6))') DVVDX(:,:)
1366

    
1367
! ! We now use the fact that Sigma=sqrt(UU) - sqrt(VV) to calculate the 
1368
! ! first  derivatives of Sigma along R
1369

    
1370
! ! Calculate the first derivatives and second derivatives.
1371
!       D1=dsqrt(UU)
1372
!       D2=dsqrt(VV)
1373
!       B(:,:)=0.5d0*( DUUDX(:,:)/D1 - DVVDX(:,:)/D2 )
1374

    
1375
!       END SUBROUTINE CONSTRAINTS_COGDIFF_DER
1376