-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathqm2_get_qm_forces.F90
1614 lines (1459 loc) · 75.6 KB
/
qm2_get_qm_forces.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
! <compile=optimized>
#include "copyright.h"
#include "../include/dprec.fh"
subroutine qm2_get_qm_forces(dxyzqm)
!Current code maintained by: Ross Walker (TSRI 2004)
!This routine calculates the derivatives of the energy for QM-QM
!interactions.
! qmmm_struct%qm_coords - QM Atom Coordinates
! dxyzqm - Returned with the forces in for each QM atom.
use constants , only : EV_TO_KCAL
use ElementOrbitalIndex, only: MaxValenceOrbitals
use qmmm_module , only : qmmm_nml,qmmm_struct, qm2_struct, qm2_params, qmmm_mpi
use qm2_pm6_hof_module
use dh_correction_module, only : dh_correction_grad
implicit none
_REAL_, parameter :: change=2.0D-6, halfChange=change/2.0D0, oneChange=1.0D0/change
_REAL_, parameter :: delAdj =1.0D-8, twoOnedelAdj= 0.5D0/delAdj
!Passed in
_REAL_, intent(out) :: dxyzqm(3,qmmm_struct%nquant_nlink)
!Local
_REAL_ e_repul(22) !Used when qmqm_erep_incore = false
_REAL_ pair_force(3)
integer loop_count !Keeps track of number of times through nquant * (nquant-1)/2 loop
! _REAL_ psum(36) !36 = max combinations with heavy and heavy = 4 orbs * 4 orbs (Note, no d orb support)
_REAL_ psum(MaxValenceOrbitals**2*3)
_REAL_ xyz_qmi(3), xyz_qmj(3), vec_qm_qm1, vec_qm_qm2, vec_qm_qm3
integer natqmi, natqmj, qmitype, qmjtype
integer ii, iif, iil, jj, jjf, jjl, ij
integer i,j,k,l
integer n_atomic_orbi, n_atomic_orbj
integer jstart, jend
_REAL_ aa,ee,deriv,angle,refh,heat,sum
_REAL_ corei, corej
_REAL_ betasas, betasap, betapas, betapap
_REAL_ bdd1i, bdd2i, bdd3i, bdd1j, bdd2j, bdd3j
_REAL_ qqi, qqi2, qqj, qqj2, ddi,ddj
_REAL_ htype, fqmii(3)
integer natom
!#define change 1.D-4
!#define halfChange 5.D-5
!!one/change = 10000
!#define onechange 10000
!#define delAdj 1.0D-8
!#define TWOONEdelAdj 50000000
if (qmmm_nml%qmqm_analyt) then !We do analytical derivatives
!RCW: Note: there is a lot of repeated code in the two options below
!but this is necessary to factor this if statement out of the inner loop.
loop_count = 0
#ifdef MPI
do ii = qmmm_mpi%nquant_nlink_istart, qmmm_mpi%nquant_nlink_iend
jstart = qmmm_mpi%nquant_nlink_jrange(1,ii)
jend = qmmm_mpi%nquant_nlink_jrange(2,ii)
#else
do II=2,qmmm_struct%nquant_nlink
jstart = 1
jend = ii-1
#endif
!Loop over all pairs of quantum atoms
! GET FIRST ATOM INFO
iif=qm2_params%orb_loc(1,II)
iil=qm2_params%orb_loc(2,II)
! n_atomic_orbi = iil - iif + 1
n_atomic_orbi = qm2_params%natomic_orbs(ii)
natqmi=qmmm_struct%iqm_atomic_numbers(II)
corei=qm2_params%core_chg(ii)
ddi = qm2_params%multip_2c_elec_params(1,ii)
qqi = qm2_params%multip_2c_elec_params(2,ii)
qqi2 = qqi*qqi
bdd1i = qm2_params%multip_2c_elec_params(3,ii)
bdd2i = qm2_params%multip_2c_elec_params(4,ii)
bdd3i = qm2_params%multip_2c_elec_params(5,ii)
xyz_qmi(1:3)=qmmm_struct%qm_coords(1:3,II)
qmitype = qmmm_struct%qm_atom_type(ii)
fqmii(1:3) = 0.0d0
do JJ=jstart,jend !jj=1,ii-1
! GET SECOND ATOM INFO
jjf=qm2_params%orb_loc(1,JJ)
jjl=qm2_params%orb_loc(2,JJ)
! n_atomic_orbj = jjl - jjf + 1
n_atomic_orbj = qm2_params%natomic_orbs(jj)
natqmj=qmmm_struct%iqm_atomic_numbers(JJ)
corej=qm2_params%core_chg(jj)
ddj = qm2_params%multip_2c_elec_params(1,jj)
qqj = qm2_params%multip_2c_elec_params(2,jj)
qqj2 = qqj*qqj
bdd1j = qm2_params%multip_2c_elec_params(3,jj)
bdd2j = qm2_params%multip_2c_elec_params(4,jj)
bdd3j = qm2_params%multip_2c_elec_params(5,jj)
xyz_qmj(1:3)=qmmm_struct%qm_coords(1:3,JJ)
vec_qm_qm1=xyz_qmi(1) - xyz_qmj(1)
vec_qm_qm2=xyz_qmi(2) - xyz_qmj(2)
vec_qm_qm3=xyz_qmi(3) - xyz_qmj(3)
qmjtype = qmmm_struct%qm_atom_type(jj)
betasas = qm2_params%betasas(qmitype,qmjtype)
betasap = qm2_params%betasap(qmitype,qmjtype)
betapas = qm2_params%betasap(qmjtype,qmitype)
betapap = qm2_params%betapap(qmitype,qmjtype)
IJ=0
! FORM DIATOMIC MATRICES
do I=jjf,jjl
K=qm2_params%pascal_tri1(i)+jjf-1
do J=jjf,I
IJ=IJ+1
K=K+1
psum(IJ)=qm2_struct%den_matrix(K)
end do
end do
! GET SECOND ATOM FIRST ATOM INTERSECTION
do I=iif,iil
L=qm2_params%pascal_tri1(i)
K=L+jjf-1
do J=jjf,jjl
IJ=IJ+1
K=K+1
psum(IJ)=qm2_struct%den_matrix(K)
end do
K=L+iif-1
do L=iif,I
K=K+1
IJ=IJ+1
psum(IJ)=qm2_struct%den_matrix(K)
end do
end do
loop_count=loop_count+1
if (qmmm_nml%qmqm_erep_incore) then
call qm2_deriv_qm_analyt(ii,jj,qm2_struct%qm_qm_e_repul(1:22,loop_count), &
psum,n_atomic_orbj,n_atomic_orbi, &
corei,corej,betasas,betasap,betapas,betapap,vec_qm_qm1,vec_qm_qm2, &
vec_qm_qm3, pair_force, qqi, qqi2, qqj, qqj2, ddi, ddj, &
bdd1i, bdd2i, bdd3i, bdd1j, bdd2j, bdd3j, &
qm2_params%atom_orb_zz_sxs_over_sas(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_ss_eqn_adb(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_zz_sxp_over_sap(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_zz_sxp_over_sap(1:6,1:6,qmjtype,qmitype), &
qm2_params%atom_orb_sp_eqn_xx1(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_sp_eqn_xx2(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_sp_eqn_xx1(1:6,1:6,qmjtype,qmitype), &
qm2_params%atom_orb_sp_eqn_xx2(1:6,1:6,qmjtype,qmitype), &
qm2_params%atom_orb_sp_eqn_xy(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_sp_eqn_xy(1:6,1:6,qmjtype,qmitype), &
qm2_params%atom_orb_zz_pxp_over_pap(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_pp_eqn_xxy1(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_pp_eqn_xxy2(1:6,1:6,qmitype,qmjtype))
else
!Same call as above, just qm2_struct%qm_qm_e_repul(1,loop_count) replaced with local e_repul
call qm2_deriv_qm_analyt(ii,jj,e_repul, &
psum,n_atomic_orbj,n_atomic_orbi, &
corei,corej,betasas,betasap,betapas,betapap,vec_qm_qm1,vec_qm_qm2, &
vec_qm_qm3, pair_force, qqi, qqi2, qqj, qqj2, ddi, ddj, &
bdd1i, bdd2i, bdd3i, bdd1j, bdd2j, bdd3j, &
qm2_params%atom_orb_zz_sxs_over_sas(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_ss_eqn_adb(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_zz_sxp_over_sap(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_zz_sxp_over_sap(1:6,1:6,qmjtype,qmitype), &
qm2_params%atom_orb_sp_eqn_xx1(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_sp_eqn_xx2(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_sp_eqn_xx1(1:6,1:6,qmjtype,qmitype), &
qm2_params%atom_orb_sp_eqn_xx2(1:6,1:6,qmjtype,qmitype), &
qm2_params%atom_orb_sp_eqn_xy(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_sp_eqn_xy(1:6,1:6,qmjtype,qmitype), &
qm2_params%atom_orb_zz_pxp_over_pap(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_pp_eqn_xxy1(1:6,1:6,qmitype,qmjtype), &
qm2_params%atom_orb_pp_eqn_xxy2(1:6,1:6,qmitype,qmjtype))
end if
fqmii(1:3) = fqmii(1:3) + pair_force(1:3)
dxyzqm(1:3,JJ)=dxyzqm(1:3,JJ)+pair_force(1:3)
end do
dxyzqm(1:3,II)=dxyzqm(1:3,II)-fqmii(1:3)
end do
!************** end ANALYTICAL DERIVATIVES **************
else !We will do (pseudo numerical derivatives)
!************** PSEUDO NUMERICAL DERIVATIVES **************
#ifdef MPI
do ii = qmmm_mpi%nquant_nlink_istart, qmmm_mpi%nquant_nlink_iend
jstart = qmmm_mpi%nquant_nlink_jrange(1,ii)
jend = qmmm_mpi%nquant_nlink_jrange(2,ii)
#else
do II=2,qmmm_struct%nquant_nlink
jstart = 1
jend = ii-1
#endif
!Loop over all pairs of quantum atoms
iif=qm2_params%orb_loc(1,II)
iil=qm2_params%orb_loc(2,II)
qmitype = qmmm_struct%qm_atom_type(ii)
natqmi=qmmm_struct%iqm_atomic_numbers(II)
xyz_qmi(1)=qmmm_struct%qm_coords(1,II)
xyz_qmi(2)=qmmm_struct%qm_coords(2,II)
xyz_qmi(3)=qmmm_struct%qm_coords(3,II)
do JJ=jstart,jend !jj=1,ii-1
! FORM DIATOMIC MATRICES
jjf=qm2_params%orb_loc(1,JJ)
jjl=qm2_params%orb_loc(2,JJ)
! GET FIRST ATOM
qmjtype = qmmm_struct%qm_atom_type(jj)
natqmj=qmmm_struct%iqm_atomic_numbers(JJ)
xyz_qmj(1)=qmmm_struct%qm_coords(1,JJ)
xyz_qmj(2)=qmmm_struct%qm_coords(2,JJ)
xyz_qmj(3)=qmmm_struct%qm_coords(3,JJ)
IJ=0
do I=jjf,jjl
K=qm2_params%pascal_tri1(i)+jjf-1
do J=jjf,I
IJ=IJ+1
K=K+1
psum(IJ)=qm2_struct%den_matrix(K)
end do
end do
! GET SECOND ATOM FIRST ATOM INTERSECTION
do I=iif,iil
L=qm2_params%pascal_tri1(i)
K=L+jjf-1
do J=jjf,jjl
IJ=IJ+1
K=K+1
psum(IJ)=qm2_struct%den_matrix(K)
end do
K=L+iif-1
do L=iif,I
K=K+1
IJ=IJ+1
psum(IJ)=qm2_struct%den_matrix(K)
end do
end do
do K=1,3
xyz_qmi(K)=xyz_qmi(K)+halfChange
call qm2_dhc(psum,ii,jj,qmitype,qmjtype,xyz_qmi,xyz_qmj,iif,iil,jjf, &
jjl,AA)
xyz_qmi(K)=xyz_qmi(K)-change
call qm2_dhc(psum,ii,jj,qmitype,qmjtype,xyz_qmi,xyz_qmj,iif,iil,jjf, &
jjl,EE)
xyz_qmi(K)=xyz_qmi(K)+halfChange
DERIV=(EE-AA)*EV_TO_KCAL*onechange
dxyzqm(K,II)=dxyzqm(K,II)-DERIV
dxyzqm(K,JJ)=dxyzqm(K,JJ)+DERIV
end do
end do
end do
!************** end PSEUDO NUMERICAL DERIVATIVES **************
end if
! --------------------------------------------
! PM6: Gradient of PM6 correction to HOF
! (nitrogen non-planarity energy penalty)
! --------------------------------------------
if (qmmm_mpi%commqmmm_master) then
! this is not parallelized - do only on the master
if (qmmm_nml%qmtheory%PM6) then
natom = qmmm_struct%nquant_nlink
call hofCorrectionGradient(natom, dxyzqm)
end if
if (qmmm_nml%qmtheory%DISPERSION .or. qmmm_nml%qmtheory%DISPERSION_HYDROGENPLUS) then
call dh_correction_grad(qmmm_struct%nquant_nlink,qmmm_struct%qm_coords, &
qmmm_struct%iqm_atomic_numbers,qmmm_nml%qmtheory,dxyzqm)
endif
end if
if(qmmm_nml%peptide_corr) then
! NOW ADD IN MOLECULAR-MECHANICS CORRECTION TO THE H-N-C=O TORSION
if (qmmm_nml%qmtheory%PM3 .OR. qmmm_nml%qmtheory%PDDGPM3 .OR. qmmm_nml%qmtheory%PM3CARB1 &
.OR. qmmm_nml%qmtheory%PM3ZNB .OR. qmmm_nml%qmtheory%PDDGPM3_08) then
htype = 7.1853D0
elseif (qmmm_nml%qmtheory%AM1 .OR. qmmm_nml%qmtheory%AM1DHFR .OR. qmmm_nml%qmtheory%RM1) then
htype = 3.3191D0
else !Assume MNDO
htype = 6.1737D0
end if
!Parallel
do I=qmmm_mpi%mytaskid+1,qm2_struct%n_peptide_links,qmmm_mpi%numthreads !1,n_peptide_links
do J=1,4
do K=1,3
qmmm_struct%qm_coords(K,qm2_struct%peptide_links(J,I))= &
qmmm_struct%qm_coords(K,qm2_struct%peptide_links(J,I))-delAdj
call qm2_dihed(qmmm_struct%qm_coords,qm2_struct%peptide_links(1,I), &
qm2_struct%peptide_links(2,I),qm2_struct%peptide_links(3,I), &
qm2_struct%peptide_links(4,I),ANGLE)
REFH=HTYPE*SIN(ANGLE)**2
qmmm_struct%qm_coords(K,qm2_struct%peptide_links(J,I))= &
qmmm_struct%qm_coords(K,qm2_struct%peptide_links(J,I))+delAdj*2.D0
call qm2_dihed(qmmm_struct%qm_coords,qm2_struct%peptide_links(1,I),qm2_struct%peptide_links(2,I), &
qm2_struct%peptide_links(3,I),qm2_struct%peptide_links(4,I),ANGLE)
qmmm_struct%qm_coords(K,qm2_struct%peptide_links(J,I))= &
qmmm_struct%qm_coords(K,qm2_struct%peptide_links(J,I))-delAdj
HEAT=HTYPE*SIN(ANGLE)**2
SUM=(REFH-HEAT)*TWOONEdelAdj
dxyzqm(K,qm2_struct%peptide_links(J,I))=dxyzqm(K,qm2_struct%peptide_links(J,I))-SUM
end do
end do
end do
end if
end subroutine qm2_get_qm_forces
subroutine qm2_deriv_qm_analyt(iqm,jqm,qm_qm_e_repul,PSUM, &
n_atomic_orbj,n_atomic_orbi,corei,corej,betasas,betasap,betapas,betapap, &
vec_qm_qm1,vec_qm_qm2, vec_qm_qm3, pair_force, &
qqi, qqi2, qqj, qqj2, ddi, ddj, bdd1i, bdd2i, bdd3i, bdd1j, bdd2j, &
bdd3j,zz_sxs_over_sas,ss_eqn_adb,zz_sxp_over_sapij, &
zz_sxp_over_sapji,sp_eqn_xx1ij,sp_eqn_xx2ij,sp_eqn_xx1ji,sp_eqn_xx2ji, &
sp_eqn_xyij, sp_eqn_xyji,zz_pxp_over_pap,pp_eqn_xxy1,pp_eqn_xxy2)
!************************************************************************
!* *
!* CALCULATION OF ANALYTICAL DERIVATIVES *
!* *
!* Current routine maintained by Ross Walker (TSRI, 2004) *
!* Inlining and optimisation by Ross Walker (TSRI, 2004) *
!* *
!************************************************************************
!*** Variable Definitions:
!
! qm_qm_e_repul = QM-QM electron repulsion integrals for this QM-QM pair
! psum = Density matrix elements for orbitals centered on this
! QM-QM pair. - For coulomb term
! n_atomic_orbj = number of atomic orbitals on atom j
! n_atomic_orbi = number of atomic orbitals on atom i
! vec_qm_qmx = xyz_qmi(x) - xyz_qmj(x)
! pair_force = Returned with cartesian force on this QM-QM pair
use constants , only: A_TO_BOHRS, A3_TO_BOHRS3, A2_TO_BOHRS2, AU_TO_EV, zero, &
half, one, two, four, fourth, eighth, sixteenth, thirtysecond, &
EV_TO_KCAL
use ElementOrbitalIndex, only: NumberElements, MaxValenceOrbitals, MaxValenceDimension
use qmmm_module , only: qmmm_nml, qm2_params, AXIS_TOL, &
OVERLAP_CUTOFF, EXPONENTIAL_CUTOFF
implicit none
!Passed in
integer, intent(in) :: iqm,jqm
_REAL_, intent(inout) :: qm_qm_e_repul(22) !for qmqm_erep_incore=true it is in, for =false it is out.
_REAL_, intent(in) :: psum(MaxValenceOrbitals**2*3)
integer, intent(in) :: n_atomic_orbj, n_atomic_orbi
_REAL_, intent(in) :: corei,corej,betasas,betasap,betapas,betapap
_REAL_, intent(in) :: vec_qm_qm1, vec_qm_qm2, vec_qm_qm3
_REAL_, intent(in) :: qqi, qqi2, qqj, qqj2, ddi, ddj
_REAL_, intent(in) :: bdd1i, bdd2i, bdd3i, bdd1j, bdd2j, bdd3j
_REAL_, intent(out):: pair_force(3)
_REAL_, intent(in) :: zz_sxs_over_sas(6,6),ss_eqn_adb(6,6)
_REAL_, intent(in) :: zz_sxp_over_sapij(6,6), zz_sxp_over_sapji(6,6)
_REAL_, intent(in) :: sp_eqn_xx1ij(6,6), sp_eqn_xx2ij(6,6)
_REAL_, intent(in) :: sp_eqn_xx1ji(6,6), sp_eqn_xx2ji(6,6)
_REAL_, intent(in) :: sp_eqn_xyij(6,6), sp_eqn_xyji(6,6)
_REAL_, intent(in) :: zz_pxp_over_pap(6,6)
_REAL_, intent(in) :: pp_eqn_xxy1(6,6), pp_eqn_xxy2(6,6)
!Local
_REAL_ FAAX, FAAY, FAAZ,FABX, FABY, FABZ,FNUCX, FNUCY, FNUCZ
_REAL_ dsx,dsy,dsz,DSPX(3), DSPY(3), DSPZ3
_REAL_ DPXPX(3), DPYPY(3), DPZPZ(3), DPCROSS
_REAL_ dgx(22), dgy(22), dgz(22)
_REAL_ drx(MaxValenceDimension**2), dry(MaxValenceDimension**2), drz(MaxValenceDimension**2)
_REAL_ r2, rij, onerij, rr, rr2
_REAL_ bb, aa
integer total_atomic_orb, qm_atomi_orb_start, orb_offset
integer isp, k, l
integer m, n, mn, kk, ll, kl, mk, nk, ml, nl
_REAL_ vec_qm_qm_onerij1, vec_qm_qm_onerij2, vec_qm_qm_onerij3
_REAL_ vec2_qm_qm1, vec2_qm_qm2, vec2_qm_qm3
_REAL_ vec_qm_qm123, vec_qm_qm12, vec_qm_qm23, vec_qm_qm13
integer i, j
_REAL_ ABN, ADBR2, SS
_REAL_ temp_real, temp_real2
_REAL_ SQRTAEE, bdd1ij
!Originally in delri
_REAL_ termx, termy, termz, ee, ade, aqe, dze, qzze, qxxe
_REAL_ aed, aeq, edz, eqzz, eqxx
_REAL_ add, adq, aqd, aqq, dxdx, dzdz,dzqxx,qxxdz,dzqzz,qzzdz
_REAL_ qxxqxx, qxxqyy, qxxqzz, qzzqxx, qzzqzz, dxqxz, qxzdx, qxzqxz
!Originally in delmol
integer mm, nn
_REAL_ ::xTDX(MaxValenceOrbitals-1),xTDY(MaxValenceOrbitals-1),xTDZ(MaxValenceOrbitals-1)
_REAL_ ::yTDX(MaxValenceOrbitals-1),yTDY(MaxValenceOrbitals-1),yTDZ(MaxValenceOrbitals-1)
_REAL_ ::zTDX(MaxValenceOrbitals-1),zTDY(MaxValenceOrbitals-1),zTDZ(MaxValenceOrbitals-1)
_REAL_ ::TX(MaxValenceOrbitals-1),TY(MaxValenceOrbitals-1),TZ(MaxValenceOrbitals-1)
_REAL_ ::TZ3i, TZ3i2
_REAL_ ::xtemp1, ytemp1, ztemp1, xtemp2, ytemp2, ztemp2
!Originally for rotat
_REAL_ RXY2, RYZ2, RZX2, oneRXY
logical LRXY2, LRYZ2, LRZX2
!For calls to core repulsion gradient routine
_REAL_ :: xyzij(3), dxyz(3)
_REAL_ :: gij, dgij(3)
!AWG: Now needed for call to qm2_repp
_REAL_ :: core(10,2)
FAAX=zero; FAAY=zero; FAAZ=zero
FABX=zero; FABY=zero; FABZ=zero
FNUCX = zero
FNUCY = zero
FNUCZ = zero
DGX(:) = ZERO
DGY(:) = ZERO
DGZ(:) = ZERO
if (n_atomic_orbi > 1 .OR. n_atomic_orbj >1) then
vec_qm_qm123 = vec_qm_qm1*vec_qm_qm2*vec_qm_qm3
vec_qm_qm12 = vec_qm_qm1*vec_qm_qm2
vec_qm_qm23 = vec_qm_qm2*vec_qm_qm3
vec_qm_qm13 = vec_qm_qm1*vec_qm_qm3
else
vec_qm_qm123 = ZERO
vec_qm_qm12 = ZERO
vec_qm_qm23 = ZERO
vec_qm_qm13 = ZERO
end if
vec2_qm_qm1 = vec_qm_qm1*vec_qm_qm1 ! (xi - xj)**2
vec2_qm_qm2 = vec_qm_qm2*vec_qm_qm2 ! (yi - yj)**2
vec2_qm_qm3 = vec_qm_qm3*vec_qm_qm3 ! (zi - zj)**2
r2=vec2_qm_qm1+vec2_qm_qm2+vec2_qm_qm3 ! (Ri-Rj)**2
rr2=r2 * A2_TO_BOHRS2
oneRIJ = one/sqrt(r2) !1/sqrt is faster than doing sqrt.
RIJ = r2*oneRIJ !one/oneRIJ
RR=RIJ * A_TO_BOHRS
bdd1ij = bdd1i+bdd1j
bdd1ij = bdd1ij*bdd1ij
SQRTAEE=1.0d0/sqrt(RR2+bdd1ij)
vec_qm_qm_onerij1 = vec_qm_qm1 * oneRIJ ! (xi - xj) / Rij
vec_qm_qm_onerij2 = vec_qm_qm2 * oneRIJ ! (yi - yj) / Rij
vec_qm_qm_onerij3 = vec_qm_qm3 * oneRIJ ! (zi - zj) / Rij
total_atomic_orb=n_atomic_orbi+n_atomic_orbj
qm_atomi_orb_start=n_atomic_orbj+1
!If we don't have the 1-e repul integrals for this QM-QM pair in memory we need
!to calculate them now.
if (.NOT. qmmm_nml%qmqm_erep_incore) then
call qm2_repp(iqm,jqm,rr,rr2,qm_qm_e_repul,core,SQRTAEE)
end if
! THE FIRST DERIVATIVES OF OVERLAP INTEGRALS
!Loop over atomic orbitals for QM atom i
!for PX, PY, PZ - So K=0 = S orbital, K=1 = PX, 2 = PY, 3 = PZ
! CALCULATE OVERLAP DERIVATIVES, STORE RESULTS IN DS
if (RR2 < OVERLAP_CUTOFF) then !10A cutoff on overlaps
!ALL Atom pairs must have at least an S-S interaction - Min 1 orbital per atom.
! (S/S) Terms
dsx=zero ; dsy=zero; dsz = zero
do I=1,6 !1 to NGAUSS
do J=1,6
ADBR2=zz_sxs_over_sas(i,j)*RR2
!Check overlap is non-zero before starting
if(ADBR2 < EXPONENTIAL_CUTOFF) then
SS=ss_eqn_adb(i,j)*EXP(-ADBR2)
DSx=DSx+vec_qm_qm1*SS
DSy=DSy+vec_qm_qm2*SS
DSz=DSz+vec_qm_qm3*SS
end if
end do
end do
! COMBINE TOGETHER THE OVERLAP DERIVATIVE PARTS
orb_offset=1+qm2_params%pascal_tri1(qm_atomi_orb_start)
temp_real=betasas*PSUM(orb_offset)
FABx=FABx+temp_real*DSx
FABy=FABy+temp_real*DSy
FABz=FABz+temp_real*DSz
! end OF S/S terms, if both QM atoms have only 1 orbital then we don't do the loops below.
!NOW DO S with P orbitals if necessary (K=0)
if (n_atomic_orbj>1) then
DSPX=zero ; DSPY=zero !Zeros entire arrays of 3
DSPZ3 = zero
do j=1,6
do i=1,6
ADBR2=zz_sxp_over_sapij(i,j)*RR2
!Check overlap is non-zero - otherwise we can skip this bit.
if (ADBR2<EXPONENTIAL_CUTOFF) then
ADBR2 = EXP(-ADBR2)
! (S/PX-x) TERM
ABN=sp_eqn_xx1ij(i,j) - (sp_eqn_xx2ij(i,j)*vec2_qm_qm1)
DSPX(1)=DSPX(1)+ABN*ADBR2
! (S/PY-y) TERM
ABN=sp_eqn_xx1ij(i,j) - (sp_eqn_xx2ij(i,j)*vec2_qm_qm2)
DSPY(2)=DSPY(2)+ABN*ADBR2
! (S/PZ-z) TERM
ABN=sp_eqn_xx1ij(i,j) - (sp_eqn_xx2ij(i,j)*vec2_qm_qm3)
DSPZ3=DSPZ3+ABN*ADBR2
ABN=ADBR2*sp_eqn_xyij(i,j)
! (PX-y/S) TERM
DSPX(2)=DSPX(2)+ABN*vec_qm_qm12
! (PX-z/S) TERM
DSPX(3)=DSPX(3)+ABN*vec_qm_qm13
! (PY-z/S) TERM
DSPY(3)=DSPY(3)+ABN*vec_qm_qm23
! (PY-x/S) TERM
!DSPY(1)=DSPY(1)-ABN*vec_qm_qm12 DSPY(1)=DSPX(2)
! (PZ-x/S) TERM
!DSPZ(1)=DSPZ(1)-ABN*vec_qm_qm13 DSPZ(1)=DSPX(3)
! (PZ-y/S) TERM
!DSPZ(2)=DSPZ(2)-ABN*vec_qm_qm23 DSPZ(2)=DSPY(3)
end if !(ADBR2<EXPONENTIAL_CUTOFF)
end do
end do
! COMBINE TOGETHER THE OVERLAP DERIVATIVE PARTS
temp_real = betasap*PSUM(orb_offset+1)
FABx=FABx+temp_real*DSPX(1)
FABy=FABy+temp_real*DSPX(2)
FABz=FABz+temp_real*DSPX(3)
temp_real = betasap*PSUM(orb_offset+2)
FABx=FABx+temp_real*DSPX(2) !DSPY(1)=DSPX(2)
FABy=FABy+temp_real*DSPY(2)
FABz=FABz+temp_real*DSPY(3)
temp_real = betasap*PSUM(orb_offset+3)
FABx=FABx+temp_real*DSPX(3) !DSPZ(1)=DSPX(3)
FABy=FABy+temp_real*DSPY(3) !DSPZ(2)=DSPY(3)
FABz=FABz+temp_real*DSPZ3
end if ! if (n_atomic_orbj>1)
!Now do P orbitals of atom 1 with S of atom 2 (K>0 and L=0)
if (n_atomic_orbi>1) then
DSPX=zero ; DSPY=zero !Zeros entire arrays of 3
DSPZ3=zero
do I=1,6
do J=1,6
ADBR2=zz_sxp_over_sapji(j,i)*RR2
if (ADBR2<EXPONENTIAL_CUTOFF) then !Only do the rest if the overlap is not zero
ADBR2 = EXP(-ADBR2)
! (PX-x/S) TERM
ABN=-sp_eqn_xx1ji(j,i) + (sp_eqn_xx2ji(j,i)*vec2_qm_qm1)
DSPX(1)=DSPX(1)+ABN*ADBR2
! (PY-y/S) TERM
ABN=-sp_eqn_xx1ji(j,i) + (sp_eqn_xx2ji(j,i)*vec2_qm_qm2)
DSPY(2)=DSPY(2)+ABN*ADBR2
! (PZ-z/S) TERM
ABN=-sp_eqn_xx1ji(j,i) + (sp_eqn_xx2ji(j,i)*vec2_qm_qm3)
DSPZ3=DSPZ3+ABN*ADBR2
ABN=ADBR2*sp_eqn_xyji(j,i)
! (PX-y/S) TERM
DSPX(2)=DSPX(2)-ABN*vec_qm_qm12
! (PX-z/S) TERM
DSPX(3)=DSPX(3)-ABN*vec_qm_qm13
! (PY-z/S) TERM
DSPY(3)=DSPY(3)-ABN*vec_qm_qm23
! (PY-x/S) TERM
!DSPY(1)=DSPY(1)-ABN*vec_qm_qm12 DSPY(1)=DSPX(2)
! (PZ-x/S) TERM
!DSPZ(1)=DSPZ(1)-ABN*vec_qm_qm13 DSPZ(1)=DSPX(3)
! (PZ-y/S) TERM
!DSPZ(2)=DSPZ(2)-ABN*vec_qm_qm23 DSPZ(2)=DSPY(3)
end if !(ADBR2<EXPONENTIAL_CUTOFF)
end do
end do
temp_real = betapas*PSUM(1+qm2_params%pascal_tri1(qm_atomi_orb_start+1))
FABx=FABx+temp_real*DSPX(1)
FABy=FABy+temp_real*DSPX(2)
FABz=FABz+temp_real*DSPX(3)
temp_real = betapas*PSUM(1+qm2_params%pascal_tri1(qm_atomi_orb_start+2))
FABx=FABx+temp_real*DSPX(2) !DSPY(1)=DSPX(2)
FABy=FABy+temp_real*DSPY(2)
FABz=FABz+temp_real*DSPY(3)
temp_real = betapas*PSUM(1+qm2_params%pascal_tri1(qm_atomi_orb_start+3))
FABx=FABx+temp_real*DSPX(3) !DSPZ(1)=DSPX(3)
FABy=FABy+temp_real*DSPY(3) !DSPZ(2)=DSPY(3)
FABz=FABz+temp_real*DSPZ3
end if !if (n_atomic_orbi>1)
!Ross Walker - PP combinations that are the same have been condensed
!to one variable for speed and to save memory.
!Finally do P's of atom 1 with P's of atom 2 (if necessary)
if (n_atomic_orbi > 1 .and. n_atomic_orbj > 1) then
DPXPX=zero; DPYPY=zero; DPZPZ=zero !Zeros entire arrays of 3
DPCROSS=zero
do I=1,6
do J=1,6
ADBR2=zz_pxp_over_pap(i,j)*RR2
if (ADBR2<EXPONENTIAL_CUTOFF) then !Only do the next bit if the overlap is not zero
ADBR2=EXP(-ADBR2)
! (PX / PX) - x
ABN=((3.0D0*pp_eqn_xxy1(i,j)) + pp_eqn_xxy2(i,j) * vec2_qm_qm1)*ADBR2
DPXPX(1)=DPXPX(1)+ABN*vec_qm_qm1
! (PY / PY) - y
ABN=((3.0D0*pp_eqn_xxy1(i,j)) + pp_eqn_xxy2(i,j) * vec2_qm_qm2)*ADBR2
DPYPY(2)=DPYPY(2)+ABN*vec_qm_qm2
! (PZ / PZ) - z
ABN=((3.0D0*pp_eqn_xxy1(i,j)) + pp_eqn_xxy2(i,j) * vec2_qm_qm3)*ADBR2
DPZPZ(3)=DPZPZ(3)+ABN*vec_qm_qm3
ABN = (pp_eqn_xxy1(i,j) + pp_eqn_xxy2(i,j) * vec2_qm_qm1)*ADBR2
! (PX / PX) - y
DPXPX(2)=DPXPX(2)+ABN*vec_qm_qm2
! (PX / PX) - z
DPXPX(3)=DPXPX(3)+ABN*vec_qm_qm3
! (PX / PY) - x
!DPXPY(1)=DPXPY(1)+ABN*vec_qm_qm2 PXPY(1)=PXPX(2)
! (PX / PZ) - x
!DPXPZ(1)=DPXPZ(1)+ABN*vec_qm_qm3 PXPZ(1)=PXPX(3)
! (PY / PX) - x
!DPYPX(1)=DPYPX(1)+ABN*vec_qm_qm2 PXPY=PYPX
! (PZ / PX) - x
!DPZPX(1)=DPZPX(1)+ABN*vec_qm_qm3 PXPZ=PZPX
ABN=(pp_eqn_xxy1(i,j) + pp_eqn_xxy2(i,j) * vec2_qm_qm2)*ADBR2
! (PY / PY) - x
DPYPY(1)=DPYPY(1)+ABN*vec_qm_qm1
! (PY / PY) - z
DPYPY(3)=DPYPY(3)+ABN*vec_qm_qm3
! (PX / PY) - y
! DPXPY(2)=DPXPY(2)+ABN*vec_qm_qm1 PXPY(2)=PYPX(2)=PYPY(1)
! (PY / PZ) - y
! DPYPZ(2)=DPYPZ(2)+ABN*vec_qm_qm3 PYPZ(2)=PZPY(2)=PYPY(3)
! (PY / PX) - y
!DPYPX(2)=DPYPX(2)+ABN*vec_qm_qm1 PXPY=PYPX
! (PZ / PY) - y
!DPZPY(2)=DPZPY(2)+ABN*vec_qm_qm3 PYPZ=PZPY
ABN=(pp_eqn_xxy1(i,j) + pp_eqn_xxy2(i,j) * vec2_qm_qm3)*ADBR2
! (PZ / PZ) - x
DPZPZ(1)=DPZPZ(1)+ABN*vec_qm_qm1
! (PZ / PZ) - y
DPZPZ(2)=DPZPZ(2)+ABN*vec_qm_qm2
! (PX / PZ) - z
! DPXPZ(3)=DPXPZ(3)+ABN*vec_qm_qm1 PXPZ(3)=PZPX(3)=PZPZ(1)
! (PY / PZ) - z
! DPYPZ(3)=DPYPZ(3)+ABN*vec_qm_qm2 PYPZ(3)=PZPY(3)=PZPZ(2)
! (PZ / PY) - z
!DPZPY(3)=DPZPY(3)+ABN*vec_qm_qm2 PYPZ=PZPY
! (PZ / PX) - z
!DPZPX(3)=DPZPX(3)+ABN*vec_qm_qm1 PXPZ=PZPX
ABN=pp_eqn_xxy2(i,j) * ADBR2*vec_qm_qm123
DPCROSS=DPCROSS+ABN
! (PX / PY) - z
!DPXPY(3)=DPXPY(3)+ABN !PXPY(3)=PYPX(3)=DPCROSS
! (PX / PZ) - y
!DPXPZ(2)=DPXPZ(2)+ABN !PXPZ(2)=PZPX(2)=DPCROSS
! (PY / PZ) - x
!DPYPZ(1)=DPYPZ(1)+ABN !PYPZ(1)=PZPY(1)=DPCROSS
! (PZ / PY) - x
!DPZPY(1)=DPZPY(1)+ABN PYPZ=PZPY
! (PY / PX) - z
!DPYPX(3)=DPYPX(3)+ABN PXPY=PYPX
! (PZ / PX) - y
!DPZPX(2)=DPZPX(2)+ABN PXPZ=PZPX
end if !ADBR2>EXPONENTIAL_CUTOFF
end do
end do
! COMBINE TOGETHER THE OVERLAP DERIVATIVE PARTS
orb_offset=2+qm2_params%pascal_tri1(qm_atomi_orb_start+1)
temp_real = betapap*PSUM(orb_offset)
FABx=FABx+temp_real*DPXPX(1)
FABy=FABy+temp_real*DPXPX(2)
FABz=FABz+temp_real*DPXPX(3)
temp_real = betapap*PSUM(orb_offset+1)
FABx=FABx+temp_real*DPXPX(2) !PYPX(1)=PXPY(1)=PXPX(2)
FABy=FABy+temp_real*DPYPY(1) !PXPY(2)=PYPX(2)=PYPY(1)
FABz=FABz+temp_real*DPCROSS !PXPY(3)=PYPX(3)=DPCROSS
temp_real = betapap*PSUM(orb_offset+2)
FABx=FABx+temp_real*DPXPX(3) !PXPZ(1)=PZPX(1)=PXPX(3)
FABy=FABy+temp_real*DPCROSS !PXPZ(2)=PZPX(2)=DPCROSS
FABz=FABz+temp_real*DPZPZ(1) !PXPZ(3)=PZPX(3)=PZPZ(1)
orb_offset=2+qm2_params%pascal_tri1(qm_atomi_orb_start+2)
temp_real=betapap*PSUM(orb_offset)
FABx=FABx+temp_real*DPXPX(2) !PYPX(1)=PXPY(1)=PXPX(2)
FABy=FABy+temp_real*DPYPY(1) !PXPY(2)=PYPX(2)=PYPY(1)
FABz=FABz+temp_real*DPCROSS !PXPY(3)=PYPX(3)=DPCROSS
temp_real=betapap*PSUM(orb_offset+1)
FABx=FABx+temp_real*DPYPY(1)
FABy=FABy+temp_real*DPYPY(2)
FABz=FABz+temp_real*DPYPY(3)
temp_real=betapap*PSUM(orb_offset+2)
FABx=FABx+temp_real*DPCROSS !PYPZ(1)=PZPY(1)=DPCROSS
FABy=FABy+temp_real*DPYPY(3) !PYPZ(2)=PZPY(2)=PYPY(3)
FABz=FABz+temp_real*DPZPZ(2) !PYPZ(3)=PZPY(3)=PZPZ(2)
orb_offset=2+qm2_params%pascal_tri1(qm_atomi_orb_start+3)
temp_real = betapap*PSUM(orb_offset)
FABx=FABx+temp_real*DPXPX(3) !PXPZ(1)=PZPX(1)=PXPX(3)
FABy=FABy+temp_real*DPCROSS !PXPZ(2)=PZPX(2)=DPCROSS
FABz=FABz+temp_real*DPZPZ(1) !PXPZ(3)=PZPX(3)=PZPZ(1)
temp_real = betapap*PSUM(orb_offset+1)
FABx=FABx+temp_real*DPCROSS !PYPZ(1)=PZPY(1)=DPCROSS
FABy=FABy+temp_real*DPYPY(3) !PYPZ(2)=PZPY(2)=PYPY(3)
FABz=FABz+temp_real*DPZPZ(2) !PYPZ(3)=PZPY(3)=PZPZ(2)
temp_real = betapap*PSUM(orb_offset+2)
FABx=FABx+temp_real*DPZPZ(1)
FABy=FABy+temp_real*DPZPZ(2)
FABz=FABz+temp_real*DPZPZ(3)
end if !n_atomic_orbi >1 .and. n_atomic_orbj > 1
end if !if (RR2 < 100.D0 * A2_TO_BOHRS2) then !10A cutoff on overlaps
!Code is linear from this point on for a given pair.
!NOTE: Ross Walker - In the equations below the only thing that varies for this pair
!during a simulation is the value of RR so we may want to consider pre-computing
!a lot of this. E.g. For a given pair AEE is constant.
! S-orbital of QM-QM pairs - all QM pairs have s-s interactions.
TERMX=vec_qm_qm_onerij1*AU_TO_EV*A_TO_BOHRS
TERMY=vec_qm_qm_onerij2*AU_TO_EV*A_TO_BOHRS
TERMZ=vec_qm_qm_onerij3*AU_TO_EV*A_TO_BOHRS
EE =-RR*(SQRTAEE)**3
DGX(1)=TERMX*EE
DGY(1)=TERMY*EE
DGZ(1)=TERMZ*EE
if(n_atomic_orbi > 2) then
! SP-ATOM-HYDROGEN
ADE=(bdd2i+bdd1j)**2
AQE=(bdd3i+bdd1j)**2
DZE = (RR+ddi)/(SQRT((RR+ddi)**2+ADE))**3 &
-(RR-ddi)/(SQRT((RR-ddi)**2+ADE))**3
temp_real = (two*RR)/(SQRT(RR2+AQE))**3
QZZE =-(RR+two*qqi)/(SQRT((RR+two*qqi)**2+AQE))**3 &
-(RR-two*qqi)/(SQRT((RR-two*qqi)**2+AQE))**3 &
+temp_real
QXXE =-(two*RR)/(SQRT(RR2+four*qqi2+AQE))**3 &
+temp_real
temp_real = -DZE*half
DGX(2)=TERMX*temp_real
DGY(2)=TERMY*temp_real
DGZ(2)=TERMZ*temp_real
temp_real = EE+QZZE*fourth
DGX(3)=TERMX*temp_real
DGY(3)=TERMY*temp_real
DGZ(3)=TERMZ*temp_real
temp_real = EE+QXXE*fourth
DGX(4)=TERMX*temp_real
DGY(4)=TERMY*temp_real
DGZ(4)=TERMZ*temp_real
end if
if (n_atomic_orbj > 2) then
! HYDROGEN-SP-ATOM
AED=(bdd1i+bdd2j)**2
AEQ=(bdd1i+bdd3j)**2
EDZ = (RR-ddj)/(SQRT((RR-ddj)**2+AED))**3 &
-(RR+ddj)/(SQRT((RR+ddj)**2+AED))**3
temp_real = (two*RR)/(SQRT(RR2+AEQ))**3
EQZZ =-(RR-two*qqj)/(SQRT((RR-two*qqj)**2+AEQ))**3 &
-(RR+two*qqj)/(SQRT((RR+two*qqj)**2+AEQ))**3 &
+temp_real
EQXX =-(two*RR)/(SQRT(RR2+four*qqj2+AEQ))**3 &
+temp_real
temp_real = -EDZ*half
DGX(5)=TERMX*temp_real
DGY(5)=TERMY*temp_real
DGZ(5)=TERMZ*temp_real
temp_real = EE+EQZZ*fourth
DGX(11)=TERMX*temp_real
DGY(11)=TERMY*temp_real
DGZ(11)=TERMZ*temp_real
temp_real = EE+EQXX*fourth
DGX(12)=TERMX*temp_real
DGY(12)=TERMY*temp_real
DGZ(12)=TERMZ*temp_real
end if
if (n_atomic_orbi > 2 .and. n_atomic_orbj > 2) then
! SP-ATOM-SP-ATOM
ADD=(bdd2i+bdd2j)**2
ADQ=(bdd2i+bdd3j)**2
AQD=(bdd3i+bdd2j)**2
AQQ=(bdd3i+bdd3j)**2
DXDX =-(two*RR)/(SQRT(RR2+(ddi-ddj)**2+ADD))**3 &
+(two*RR)/(SQRT(RR2+(ddi+ddj)**2+ADD))**3
DZDZ =-(RR+ddi-ddj)/(SQRT((RR+ddi-ddj)**2+ADD))**3 &
-(RR-ddi+ddj)/(SQRT((RR-ddi+ddj)**2+ADD))**3 &
+(RR-ddi-ddj)/(SQRT((RR-ddi-ddj)**2+ADD))**3 &
+(RR+ddi+ddj)/(SQRT((RR+ddi+ddj)**2+ADD))**3
DZQXX = two*(RR+ddi)/(SQRT((RR+ddi)**2+four*qqj2+ADQ))**3 &
-two*(RR-ddi)/(SQRT((RR-ddi)**2+four*qqj2+ADQ))**3 &
-two*(RR+ddi)/(SQRT((RR+ddi)**2+ADQ))**3 &
+two*(RR-ddi)/(SQRT((RR-ddi)**2+ADQ))**3
QXXDZ = two*(RR-ddj)/(SQRT((RR-ddj)**2+four*qqi2+AQD))**3 &
-two*(RR+ddj)/(SQRT((RR+ddj)**2+four*qqi2+AQD))**3 &
-two*(RR-ddj)/(SQRT((RR-ddj)**2+AQD))**3 &
+two*(RR+ddj)/(SQRT((RR+ddj)**2+AQD))**3
DZQZZ = (RR+ddi-two*qqj)/(SQRT((RR+ddi-two*qqj)**2+ADQ))**3 &
-(RR-ddi-two*qqj)/(SQRT((RR-ddi-two*qqj)**2+ADQ))**3 &
+(RR+ddi+two*qqj)/(SQRT((RR+ddi+two*qqj)**2+ADQ))**3 &
-(RR-ddi+two*qqj)/(SQRT((RR-ddi+two*qqj)**2+ADQ))**3 &
+two*(RR-ddi)/(SQRT((RR-ddi)**2+ADQ))**3 &
-two*(RR+ddi)/(SQRT((RR+ddi)**2+ADQ))**3
QZZDZ = (RR+two*qqi-ddj)/(SQRT((RR+two*qqi-ddj)**2+AQD))**3 &
-(RR+two*qqi+ddj)/(SQRT((RR+two*qqi+ddj)**2+AQD))**3 &
+(RR-two*qqi-ddj)/(SQRT((RR-two*qqi-ddj)**2+AQD))**3 &
-(RR-two*qqi+ddj)/(SQRT((RR-two*qqi+ddj)**2+AQD))**3 &
-two*(RR-ddj)/(SQRT((RR-ddj)**2+AQD))**3 &
+two*(RR+ddj)/(SQRT((RR+ddj)**2+AQD))**3
QXXQXX=-(two*RR)/(SQRT(RR2+four*(qqi-qqj)**2+AQQ))**3 &
-(two*RR)/(SQRT(RR2+four*(qqi+qqj)**2+AQQ))**3 &
+(four*RR)/(SQRT(RR2+four*qqi2+AQQ))**3 &
+(four*RR)/(SQRT(RR2+four*qqj2+AQQ))**3 &
-(four*RR)/(SQRT(RR2+AQQ))**3
QXXQYY=-(four*RR)/(SQRT(RR2+four*qqi2+four*qqj2+AQQ))**3 &
+(four*RR)/(SQRT(RR2+four*qqi2+AQQ))**3 &
+(four*RR)/(SQRT(RR2+four*qqj2+AQQ))**3 &
-(four*RR)/(SQRT(RR2+AQQ))**3
QXXQZZ= &
-two*(RR-two*qqj)/(SQRT((RR-two*qqj)**2+four*qqi2+AQQ))**3 &
-two*(RR+two*qqj)/(SQRT((RR+two*qqj)**2+four*qqi2+AQQ))**3 &
+two*(RR-two*qqj)/(SQRT((RR-two*qqj)**2+AQQ))**3 &
+two*(RR+two*qqj)/(SQRT((RR+two*qqj)**2+AQQ))**3 &
+(four*RR)/(SQRT(RR2+four*qqi2+AQQ))**3 &
-(four*RR)/(SQRT(RR2+AQQ))**3
QZZQXX= &
-two*(RR+two*qqi)/(SQRT((RR+two*qqi)**2+four*qqj2+AQQ))**3 &
-two*(RR-two*qqi)/(SQRT((RR-two*qqi)**2+four*qqj2+AQQ))**3 &
+two*(RR+two*qqi)/(SQRT((RR+two*qqi)**2+AQQ))**3 &
+two*(RR-two*qqi)/(SQRT((RR-two*qqi)**2+AQQ))**3 &
+(four*RR)/(SQRT(RR2+four*qqj2+AQQ))**3 &
-(four*RR)/(SQRT(RR2+AQQ))**3
QZZQZZ= &
-(RR+two*qqi-two*qqj)/(SQRT((RR+two*qqi-two*qqj)**2+AQQ))**3 &
-(RR+two*qqi+two*qqj)/(SQRT((RR+two*qqi+two*qqj)**2+AQQ))**3 &
-(RR-two*qqi-two*qqj)/(SQRT((RR-two*qqi-two*qqj)**2+AQQ))**3 &
-(RR-two*qqi+two*qqj)/(SQRT((RR-two*qqi+two*qqj)**2+AQQ))**3 &
+two*(RR-two*qqi)/(SQRT((RR-two*qqi)**2+AQQ))**3 &
+two*(RR+two*qqi)/(SQRT((RR+two*qqi)**2+AQQ))**3 &
+two*(RR-2.D0*qqj)/(SQRT((RR-2.D0*qqj)**2+AQQ))**3 &
+two*(RR+2.D0*qqj)/(SQRT((RR+2.D0*qqj)**2+AQQ))**3 &
-(four*RR)/(SQRT(RR2+AQQ))**3
DXQXZ = two*(RR-qqj)/(SQRT((RR-qqj)**2+(ddi-qqj)**2+ADQ))**3 &
-two*(RR+qqj)/(SQRT((RR+qqj)**2+(ddi-qqj)**2+ADQ))**3 &
-two*(RR-qqj)/(SQRT((RR-qqj)**2+(ddi+qqj)**2+ADQ))**3 &
+two*(RR+qqj)/(SQRT((RR+qqj)**2+(ddi+qqj)**2+ADQ))**3
QXZDX = two*(RR+qqi)/(SQRT((RR+qqi)**2+(qqi-ddj)**2+AQD))**3 &
-two*(RR-qqi)/(SQRT((RR-qqi)**2+(qqi-ddj)**2+AQD))**3 &
-two*(RR+qqi)/(SQRT((RR+qqi)**2+(qqi+ddj)**2+AQD))**3 &
+two*(RR-qqi)/(SQRT((RR-qqi)**2+(qqi+ddj)**2+AQD))**3
QXZQXZ=-two*(RR+qqi-qqj)/(SQRT((RR+qqi-qqj)**2+(qqi-qqj)**2+AQQ))**3 &
+two*(RR+qqi+qqj)/(SQRT((RR+qqi+qqj)**2+(qqi-qqj)**2+AQQ))**3 &
+two*(RR-qqi-qqj)/(SQRT((RR-qqi-qqj)**2+(qqi-qqj)**2+AQQ))**3 &
-two*(RR-qqi+qqj)/(SQRT((RR-qqi+qqj)**2+(qqi-qqj)**2+AQQ))**3 &
+two*(RR+qqi-qqj)/(SQRT((RR+qqi-qqj)**2+(qqi+qqj)**2+AQQ))**3 &
-two*(RR+qqi+qqj)/(SQRT((RR+qqi+qqj)**2+(qqi+qqj)**2+AQQ))**3 &
-two*(RR-qqi-qqj)/(SQRT((RR-qqi-qqj)**2+(qqi+qqj)**2+AQQ))**3 &
+two*(RR-qqi+qqj)/(SQRT((RR-qqi+qqj)**2+(qqi+qqj)**2+AQQ))**3
temp_real = DZDZ*fourth
DGX(6)=TERMX*temp_real
DGY(6)=TERMY*temp_real
DGZ(6)=TERMZ*temp_real
temp_real = DXDX*fourth
DGX(7)=TERMX*temp_real
DGY(7)=TERMY*temp_real
DGZ(7)=TERMZ*temp_real
temp_real = -(EDZ*half+QZZDZ*eighth)
DGX(8)=TERMX*temp_real
DGY(8)=TERMY*temp_real
DGZ(8)=TERMZ*temp_real
temp_real = -(EDZ*half+QXXDZ*eighth)
DGX(9)=TERMX*temp_real
DGY(9)=TERMY*temp_real
DGZ(9)=TERMZ*temp_real
temp_real = -QXZDX*eighth
DGX(10)=TERMX*temp_real
DGY(10)=TERMY*temp_real
DGZ(10)=TERMZ*temp_real
temp_real = -(DZE*half+DZQZZ*eighth)
DGX(13)=TERMX*temp_real
DGY(13)=TERMY*temp_real
DGZ(13)=TERMZ*temp_real
temp_real = -(DZE*half+DZQXX*eighth)
DGX(14)=TERMX*temp_real
DGY(14)=TERMY*temp_real
DGZ(14)=TERMZ*temp_real
temp_real = -DXQXZ*eighth
DGX(15)=TERMX*temp_real
DGY(15)=TERMY*temp_real
DGZ(15)=TERMZ*temp_real
temp_real2 = EE+EQZZ*fourth
temp_real = temp_real2+QZZE*fourth+QZZQZZ*sixteenth
DGX(16)=TERMX*temp_real
DGY(16)=TERMY*temp_real
DGZ(16)=TERMZ*temp_real
temp_real = temp_real2+QXXE*fourth+QXXQZZ*sixteenth
DGX(17)=TERMX*temp_real
DGY(17)=TERMY*temp_real
DGZ(17)=TERMZ*temp_real
temp_real2 = EE+EQXX*fourth
temp_real = temp_real2+QZZE*fourth+QZZQXX*sixteenth
DGX(18)=TERMX*temp_real
DGY(18)=TERMY*temp_real
DGZ(18)=TERMZ*temp_real
temp_real = temp_real2+QXXE*fourth+QXXQXX*sixteenth
DGX(19)=TERMX*temp_real
DGY(19)=TERMY*temp_real
DGZ(19)=TERMZ*temp_real
temp_real = QXZQXZ*sixteenth
DGX(20)=TERMX*temp_real
DGY(20)=TERMY*temp_real
DGZ(20)=TERMZ*temp_real
temp_real = temp_real2+QXXE*fourth+QXXQYY*sixteenth
DGX(21)=TERMX*temp_real
DGY(21)=TERMY*temp_real
DGZ(21)=TERMZ*temp_real
temp_real = (QXXQXX-QXXQYY)*thirtysecond
DGX(22)=TERMX*temp_real
DGY(22)=TERMY*temp_real
DGZ(22)=TERMZ*temp_real
end if
if(n_atomic_orbi>1.OR.n_atomic_orbj>1) then
RXY2 = vec2_qm_qm1 + vec2_qm_qm2
RYZ2 = vec2_qm_qm2 + vec2_qm_qm3
RZX2 = vec2_qm_qm3 + vec2_qm_qm1
LRXY2 = RXY2 < AXIS_TOL
LRYZ2 = RYZ2 < AXIS_TOL
LRZX2 = RZX2 < AXIS_TOL
XTDX=zero; YTDX=zero; ZTDX=zero
XTDY=zero; YTDY=zero; ZTDY=zero
XTDZ=zero; YTDZ=zero; ZTDZ=zero
if (.NOT.(LRXY2 .OR. LRYZ2 .OR. LRZX2)) then
TERMX = vec_qm_qm_onerij1 * oneRIJ !vec_qm_qm1/(RIJ*RIJ)
TERMY = vec_qm_qm_onerij2 * oneRIJ
TERMZ = vec_qm_qm_onerij3 * oneRIJ
oneRXY = one/sqrt(RXY2)
!Ross Walker - rearranged order here slightly for speed.
TZ(3) = oneRIJ/oneRXY
TZ3i = RIJ * oneRXY !Inverse of TZ(3) to avoid other divisions.
TZ3i2 = TZ3i*TZ3i !Square of 1/TZ(3)
TX(1) = vec_qm_qm_onerij1
TX(2) = vec_qm_qm_onerij2
TX(3) = vec_qm_qm_onerij3
TY(1) = -TX(2)*SIGN(one,TX(1)) * TZ3i
TY(2) = ABS(TX(1) * TZ3i)
TY(3) = zero
TZ(1) = -TX(1)*TX(3)*TZ3i
TZ(2) = -TX(2)*TX(3)*TZ3i
!TZ(3) = see above
XTDX(1)=oneRIJ-TX(1)*TERMX
YTDX(1)=-TX(1)*TERMY
ZTDX(1)=-TX(1)*TERMZ
XTDX(2)=-TX(2)*TERMX
YTDX(2)=oneRIJ-TX(2)*TERMY
ZTDX(2)=-TX(2)*TERMZ
XTDX(3)=-TX(3)*TERMX
YTDX(3)=-TX(3)*TERMY
ZTDX(3)=oneRIJ-TX(3)*TERMZ
XTDZ(3)=TX(1)*oneRXY-TZ(3)*TERMX
YTDZ(3)=TX(2)*oneRXY-TZ(3)*TERMY