forked from lgcrego/Dynemol
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathf_intra.f
935 lines (777 loc) · 39.3 KB
/
f_intra.f
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
module F_intra_m
use type_m
use constants_m
use MM_input , only: MM_input_format
use parameters_m , only: PBC , QMMM , driver , n_part
use Allocation_m , only: Allocate_Structures
use setup_m , only: offset
use for_force , only: rcut, vrecut, frecut, pot_INTER, bdpot, angpot, dihpot, Morspot, &
vscut, fscut, KAPPA, LJ_14, LJ_intra, Coul_14, Coul_intra, pot_total, &
Dihedral_Potential_Type, forcefield, rcutsq, ryck_dih, proper_dih, &
harm_dih, imp_dih, harm_bond, morse_bond
use MD_read_m , only: atom , molecule , MM
use MM_types , only: MM_system , MM_molecular , MM_atomic , debug_MM
use gmx2mdflex , only: SpecialPairs , SpecialPairs14 , SpecialMorse
use CSDM_master , only: Ehrenfest_master
use Ehrenfest_Builder , only: EhrenfestForce
use Surface_Hopping , only: SH_Force
use decoherence_m , only: DecoherenceForce
private
public :: ForceINTRA, pot_INTRA, BcastQMArgs, ForceQMMM
! module variables ...
real*8 , dimension (3):: rij , rjk , rkl , rik , rijk , rjkl , rijkl , f1 , f2 , f3 , f4
real*8 :: rijq , rjkq , rklq , rijsq , rjksq , rklsq , fxyz , riju , riku , rijkj , rijkj2 , rjkkl , rjkkl2 , &
rijkl2 , rjksq2 , rijkll , f1x , f1y , f1z , f2x , f2y , f2z , f3x , f3y , f3z , f4x , f4y , f4z , &
sr2 , sr6 , sr12 , fs , phi , cosphi , sinphi , rsinphi , coephi , gamma , KRIJ , expar , eme , dphi , &
term , chrgi , chrgj , freal , sig , eps , pterm , A0 , A1 , A2 , A3 , rtwopi , qterm , qterm0 , rterm, &
sterm , tterm , C0 , C1 , C2 , C3 , C4 , C5 , coephi0 , rterm0 , rikq , riksq , term1 , term2 , term3 , &
term4 , dphi1 , dphi2
real*8 :: pot_INTRA
integer :: i , j , k , l , m , n , ati , atj , atk , atl , loop , ati1 , atj1
logical :: flag1, flag2, flag3, flag4, flag5
logical :: there_are_NB_SpecialPairs = .false.
logical :: there_are_NB_SpecialPairs14 = .false.
integer , allocatable :: species_offset(:)
! imported module variables ...
integer :: PST(2)
real*8 :: t_rate
character(2) :: mode
type(structure) :: system
type(R_eigen) :: QM
type(STO_basis) , allocatable , dimension(:) :: basis(:)
Complex*16 , allocatable , dimension(:,:) :: MO_bra , MO_ket , MO_TDSE_bra , MO_TDSE_ket
interface BcastQMArgs
module procedure AllocateQMArgs
module procedure StoreQMArgs_AO_and_FSSH
module procedure StoreQMArgs_CSDM
end interface
contains
!
!
!
!====================
subroutine ForceINTRA
!====================
implicit none
! local_variables ...
rtwopi = 1.d0/twopi
do j = 1 , MM % N_of_atoms
atom(j) % fnonbd14(:) = D_zero ! Non-bonded 1-4
atom(j) % fnonbd(:) = D_zero ! Non-bonded Intramolecular
atom(j) % fbond(:) = D_zero ! Stretching/Bonding
atom(j) % fang(:) = D_zero ! Bending/Angular
atom(j) % fdihed(:) = D_zero ! Dihedral
atom(j) % fnonch14(:) = D_zero ! Non-bonded coulomb 1-4
atom(j) % fnonch(:) = D_zero ! Non-bonded coulomb Intramolecular
atom(j) % fMorse(:) = D_zero ! Non-bonded Morse
end do
bdpot = D_zero
angpot = D_zero
dihpot = D_zero
Morspot = D_zero
harm_bond = D_zero
morse_bond = D_zero
proper_dih = D_zero
harm_dih = D_zero
ryck_dih = D_zero
imp_dih = D_zero
LJ_14 = D_zero
LJ_intra = D_zero
Coul_14 = D_zero
Coul_intra = D_zero
CALL offset( species_offset )
! Bonding - stretching potential ...
do i = 1 , MM % N_of_molecules
do j = 1 , molecule(i) % Nbonds
ati = molecule(i) % bonds(j,1)
atj = molecule(i) % bonds(j,2)
if( atom(atj) % flex .OR. atom(ati) % flex ) then
rij(:) = atom(atj) % xyz(:) - atom(ati) % xyz(:)
rij(:) = rij(:) - MM % box(:) * DNINT( rij(:) * MM % ibox(:) ) * PBC(:)
rijq = rij(1)*rij(1) + rij(2)*rij(2) + rij(3)*rij(3)
rijsq = SQRT(rijq)
select case ( molecule(i) % bond_type(j) )
case ( "harm" )
! harmonic potential ...
qterm = HALF * molecule(i) % kbond0(j,1) * ( rijsq - molecule(i) % kbond0(j,2) ) * ( rijsq - molecule(i) % kbond0(j,2) )
coephi = molecule(i) % kbond0(j,1)*( rijsq - molecule(i) % kbond0(j,2) )/rijsq
harm_bond = qterm + harm_bond
case ( "Mors" )
! Morse potential ...
qterm0 = exp( -molecule(i) % kbond0(j,3) * ( rijsq - molecule(i) % kbond0(j,2) ) )
qterm = molecule(i) % kbond0(j,1) * ( D_ONE - qterm0 )*( D_ONE - qterm0 )
coephi = TWO * molecule(i) % kbond0(j,1) * molecule(i) % kbond0(j,3) * qterm0 * ( D_ONE - qterm0 ) / rijsq
morse_bond = qterm + morse_bond
end select
atom(atj) % fbond(:) = atom(atj) % fbond(:) - coephi*rij(:)
atom(ati) % fbond(:) = atom(ati) % fbond(:) + coephi*rij(:)
bdpot = qterm + bdpot
end if
end do
end do
!====================================================================
! Angle - bending potential ...
do i = 1 , MM % N_of_molecules
do j = 1 , molecule(i) % Nangs
atj = molecule(i) % angs(j,1)
ati = molecule(i) % angs(j,2)
atk = molecule(i) % angs(j,3)
if( atom(atj) % flex .OR. atom(ati) % flex .OR. atom(atk) % flex ) then
rij(:) = atom(atj) % xyz(:) - atom(ati) % xyz(:)
rij(:) = rij(:) - MM % box(:) * DNINT( rij(:) * MM % ibox(:) ) * PBC(:)
rijq = rij(1)*rij(1) + rij(2)*rij(2) + rij(3)*rij(3)
rijsq = SQRT(rijq)
rjk(:) = atom(atk) % xyz(:) - atom(ati) % xyz(:)
rjk(:) = rjk(:) - MM % box(:)*DNINT( rjk(:) * MM % ibox(:) ) * PBC(:)
rjkq = rjk(1)*rjk(1) + rjk(2)*rjk(2) + rjk(3)*rjk(3)
rjksq = SQRT(rjkq)
phi = ACOS( (rij(1)*rjk(1) + rij(2)*rjk(2) + rij(3)*rjk(3) ) / ( rijsq * rjksq ) )
select case ( molecule(i) % angle_type(j) )
case( "harm" , "urba" )
! Harmonic and Urey-Bradley potentials ...
coephi = 0.d0
if (phi < 1.d-12 .OR. abs(pi - phi) < 1.d-12) then
coephi = 0.d0
rterm = 0.0d0
else
coephi = ( phi - molecule(i) % kang0(j,2) ) / SIN(phi)
rterm = HALF * molecule(i) % kang0(j,1) * ( phi - molecule(i) % kang0(j,2) )**2
end if
angpot = rterm + angpot
do l = 1, 3
if (l == 1) atl = atj
if (l == 2) atl = ati
if (l == 3) atl = atk
do loop = 1, 3 ! X,Y,Z axis (n = 1, 2 or 3)
fxyz = 0.d0
riju = rij(loop)
riku = rjk(loop)
fxyz = ( molecule(i) % kang0(j,1) * coephi ) * &
( (DEL(atl,atj)-DEL(atl,ati))*riku/(rijsq*rjksq) + &
(DEL(atl,atk)-DEL(atl,ati))*riju/(rijsq*rjksq) - &
COS(phi)*( (DEL(atl,atj)-DEL(atl,ati))*riju/(rijsq*rijsq) + &
(DEL(atl,atk)-DEL(atl,ati))*riku/(rjksq*rjksq)) )
atom(atl) % fang(loop) = atom(atl) % fang(loop) + fxyz
end do
end do
! Urey-Bradley bonding term ...
if( molecule(i) % angle_type(j) == "urba" ) then
rik(:) = atom(atk) % xyz(:) - atom(atj) % xyz(:)
rik(:) = rik(:) - MM % box(:) * DNINT( rik(:) * MM % ibox(:) ) * PBC(:)
rikq = rik(1)*rik(1) + rik(2)*rik(2) + rik(3)*rik(3)
riksq = SQRT(rikq)
coephi0 = molecule(i) % kang0(j,3) * ( riksq - molecule(i) % kang0(j,4) )/riksq
rterm0 = HALF * molecule(i) % kang0(j,3) * &
( riksq - molecule(i) % kang0(j,4) ) * ( riksq - molecule(i) % kang0(j,4) )
atom(atk) % fang(:) = atom(atk) % fang(:) - coephi0*rik(:)
atom(atj) % fang(:) = atom(atj) % fang(:) + coephi0*rik(:)
angpot = rterm0 + angpot
end if
end select
end if
end do
end do
!====================================================================
! Dihedral Potential Angle ...
! USING IUPAC first convention for dihedral definitions: trans = 180 deg ...
do i = 1 , MM % N_of_molecules
do j = 1 , molecule(i) % Ndiheds
ati = molecule(i) % diheds(j,1)
atj = molecule(i) % diheds(j,2)
atk = molecule(i) % diheds(j,3)
atl = molecule(i) % diheds(j,4)
if ( atom(atj) % flex .OR. atom(ati) % flex .OR. atom(atk) % flex .OR. atom(atl) % flex ) then
! Definition of vector rij = ri - rj
rij(:) = atom(ati) % xyz(1:3) - atom(atj) % xyz(1:3)
rij(:) = rij(1:3) - MM % box(1:3) * DNINT( rij(1:3) * MM % ibox(1:3) ) * PBC(1:3)
! Definition of vector rjk = rj - rk
rjk(:) = atom(atj) % xyz(1:3) - atom(atk) % xyz(1:3)
rjk(:) = rjk(1:3) - MM % box(1:3) * DNINT( rjk(1:3) * MM % ibox(1:3) ) * PBC(1:3)
rjkq = rjk(1)*rjk(1) + rjk(2)*rjk(2) + rjk(3)*rjk(3)
rjksq = 1.d0 / SQRT(rjkq)
rjksq2 = rjksq * rjksq
! Definition of vector rkl = rk - rl
rkl(:) = atom(atk) % xyz(1:3) - atom(atl) % xyz(1:3)
rkl(:) = rkl(1:3) - MM % box(1:3) * DNINT( rkl(1:3) * MM % ibox(1:3) ) * PBC(1:3)
! Cross Product M = | rij X rjk | :: First dihedral vector ...
rijk(1) = rij(2) * rjk(3) - rij(3) * rjk(2)
rijk(2) = rij(3) * rjk(1) - rij(1) * rjk(3)
rijk(3) = rij(1) * rjk(2) - rij(2) * rjk(1)
rijkj = rijk(1)*rijk(1) + rijk(2)*rijk(2) + rijk(3)*rijk(3)
rijkj = 1.d0 / SQRT(rijkj)
rijkj2 = rijkj * rijkj
! Cross Product N = | rjk X rkl | :: Second dihedral vector ...
rjkl(1) = rjk(2) * rkl(3) - rjk(3) * rkl(2)
rjkl(2) = rjk(3) * rkl(1) - rjk(1) * rkl(3)
rjkl(3) = rjk(1) * rkl(2) - rjk(2) * rkl(1)
rjkkl = rjkl(1)*rjkl(1) + rjkl(2)*rjkl(2) + rjkl(3)*rjkl(3)
rjkkl = 1.d0 / SQRT(rjkkl)
rjkkl2 = rjkkl * rjkkl
! Cross Product O = | rjk X rkl | X | rij X rjk |
rijkl(1) = rjkl(2) * rijk(3) - rjkl(3) * rijk(2)
rijkl(2) = rjkl(3) * rijk(1) - rjkl(1) * rijk(3)
rijkl(3) = rjkl(1) * rijk(2) - rjkl(2) * rijk(1)
rijkl2 = rijkl(1)*rijkl(1)+rijkl(2)*rijkl(2)+rijkl(3)*rijkl(3)
rijkll = SQRT(rijkl2)
! PHI is the dihedral angle defined by PHI = ACOS(B), where B(rij,rjk,rkl) = [ (rij X rjk).(rjk X rkl) / |rij X rjk||rjk X rkl| ]
! and the sign of PHI is positive if the vector O is in the same direction as the bond vector rjk
coephi = sum( rijk(:)*rjkl(:) )
cosphi = coephi * rijkj * rjkkl
If( Abs(cosphi) > 1.0d0 ) cosphi=Sign(1.0d0,cosphi)
sinphi = sum( rjk(:)*rijkl(:) ) * rjksq * rijkj * rjkkl
phi = ATAN2( sinphi,cosphi )
! Avoid singularity in sinphi
sinphi = sign( max(1.d-10 , abs(sinphi) ) , sinphi )
rsinphi = 1.d0 / sinphi
! selection of potential energy function type
if( MM_input_format == "GMX" ) then
CALL gmx
else
CALL not_gmx
end if
! Calculate atomic forces ...
f1x = gamma * ( (-rjkl(2) * rjk(3) + rjkl(3) * rjk(2)) - (-rijk(2) * rjk(3) + rijk(3) * rjk(2)) * coephi * rijkj2)
f1y = gamma * ( ( rjkl(1) * rjk(3) - rjkl(3) * rjk(1)) - ( rijk(1) * rjk(3) - rijk(3) * rjk(1)) * coephi * rijkj2)
f1z = gamma * ( (-rjkl(1) * rjk(2) + rjkl(2) * rjk(1)) - (-rijk(1) * rjk(2) + rijk(2) * rjk(1)) * coephi * rijkj2)
f3x = gamma * ( (-rjkl(2) * rij(3) + rjkl(3) * rij(2)) - (-rijk(2) * rij(3) + rijk(3) * rij(2)) * coephi * rijkj2)
f3y = gamma * ( ( rjkl(1) * rij(3) - rjkl(3) * rij(1)) - ( rijk(1) * rij(3) - rijk(3) * rij(1)) * coephi * rijkj2)
f3z = gamma * ( (-rjkl(1) * rij(2) + rjkl(2) * rij(1)) - (-rijk(1) * rij(2) + rijk(2) * rij(1)) * coephi * rijkj2)
f2x = gamma * ( (-rijk(2) * rkl(3) + rijk(3) * rkl(2)) - (-rjkl(2) * rkl(3) + rjkl(3) * rkl(2)) * coephi * rjkkl2)
f2y = gamma * ( ( rijk(1) * rkl(3) - rijk(3) * rkl(1)) - ( rjkl(1) * rkl(3) - rjkl(3) * rkl(1)) * coephi * rjkkl2)
f2z = gamma * ( (-rijk(1) * rkl(2) + rijk(2) * rkl(1)) - (-rjkl(1) * rkl(2) + rjkl(2) * rkl(1)) * coephi * rjkkl2)
f4x = gamma * ( (-rijk(2) * rjk(3) + rijk(3) * rjk(2)) - (-rjkl(2) * rjk(3) + rjkl(3) * rjk(2)) * coephi * rjkkl2)
f4y = gamma * ( ( rijk(1) * rjk(3) - rijk(3) * rjk(1)) - ( rjkl(1) * rjk(3) - rjkl(3) * rjk(1)) * coephi * rjkkl2)
f4z = gamma * ( (-rijk(1) * rjk(2) + rijk(2) * rjk(1)) - (-rjkl(1) * rjk(2) + rjkl(2) * rjk(1)) * coephi * rjkkl2)
atom(ati) % fdihed(1) = atom(ati) % fdihed(1) + f1x
atom(ati) % fdihed(2) = atom(ati) % fdihed(2) + f1y
atom(ati) % fdihed(3) = atom(ati) % fdihed(3) + f1z
atom(atj) % fdihed(1) = atom(atj) % fdihed(1) - f1x - f3x + f2x
atom(atj) % fdihed(2) = atom(atj) % fdihed(2) - f1y - f3y + f2y
atom(atj) % fdihed(3) = atom(atj) % fdihed(3) - f1z - f3z + f2z
atom(atk) % fdihed(1) = atom(atk) % fdihed(1) - f2x - f4x + f3x
atom(atk) % fdihed(2) = atom(atk) % fdihed(2) - f2y - f4y + f3y
atom(atk) % fdihed(3) = atom(atk) % fdihed(3) - f2z - f4z + f3z
atom(atl) % fdihed(1) = atom(atl) % fdihed(1) + f4x
atom(atl) % fdihed(2) = atom(atl) % fdihed(2) + f4y
atom(atl) % fdihed(3) = atom(atl) % fdihed(3) + f4z
dihpot = dihpot + pterm
end if
end do
end do
!====================================================================
! Non-bonded 1,4 intramolecular interactions ...
If( allocated(SpecialPairs14) ) there_are_NB_SpecialPairs14 = .true.
do i = 1 , MM % N_of_molecules
do j = 1 , molecule(i) % Nbonds14
ati = molecule(i) % bonds14(j,1)
atj = molecule(i) % bonds14(j,2)
if ( atom(atj) % flex .OR. atom(ati) % flex ) then
chrgi = atom(ati) % charge
chrgj = atom(atj) % charge
rij(:) = atom(ati) % xyz(:) - atom(atj) % xyz(:)
rij(:) = rij(:) - MM % box(:) * DNINT( rij(:) * MM % ibox(:) ) * PBC(:)
rklq = rij(1)*rij(1) + rij(2)*rij(2) + rij(3)*rij(3)
! Lennard Jones ...
select case ( MM % CombinationRule )
case (2)
! AMBER FF :: GMX COMB-RULE 2
sr2 = ( ( atom(ati) % sig14 + atom(atj) % sig14 ) * ( atom(ati) % sig14 + atom(atj) % sig14 ) ) / rklq
case (3)
! OPLS FF :: GMX COMB-RULE 3
sr2 = ( ( atom(ati) % sig14 * atom(atj) % sig14 ) * ( atom(ati) % sig14 * atom(atj) % sig14 ) ) / rklq
end select
eps = atom(ati) % eps14 * atom(atj) % eps14
If( there_are_NB_SpecialPairs14 ) then ! <== check whether (I,J) is a SpecialPair ...
read_loop1: do n = 1, size(SpecialPairs14)
flag1 = ( adjustl( SpecialPairs14(n) % MMSymbols(1) ) == adjustl( atom(ati) % MMSymbol ) ) .AND. &
( adjustl( SpecialPairs14(n) % MMSymbols(2) ) == adjustl( atom(atj) % MMSymbol ) )
flag2 = ( adjustl( SpecialPairs14(n) % MMSymbols(2) ) == adjustl( atom(ati) % MMSymbol ) ) .AND. &
( adjustl( SpecialPairs14(n) % MMSymbols(1) ) == adjustl( atom(atj) % MMSymbol ) )
if ( flag1 .OR. flag2 ) then ! <== apply SpecialPair parms ...
sr2 = ( SpecialPairs14(n)%Parms(1) * SpecialPairs14(n)%Parms(1) ) / rklq
eps = SpecialPairs14(n) % Parms(2)
exit read_loop1
end if
end do read_loop1
end if
rklsq = sqrt(rklq)
sr6 = sr2 * sr2 * sr2
sr12 = sr6 * sr6
fs = 24.d0 * eps * ( TWO * sr12 - sr6 )
fs = (fs / rklq) * MM % fudgeLJ
atom(ati) % fnonbd14(1:3) = atom(ati) % fnonbd14(1:3) + fs * rij(1:3)
atom(atj) % fnonbd14(1:3) = atom(atj) % fnonbd14(1:3) - fs * rij(1:3)
! factor used to compensate the factor1 and factor2 factors ...
! factor3 = 1.0d-20
sterm = 4.d0 * eps * factor3 * ( sr12 - sr6 )
! alternative cutoff formula ...
! sterm = sterm - vscut(ati,atj) + fscut(ati,atj) * ( rklsq - rcut )
sterm = sterm * MM % fudgeLJ
! Real part (Numero de cargas igual ao numero de sitios)
sr2 = 1.d0 / rklq
KRIJ = KAPPA * rklsq
expar = EXP( -(KRIJ*KRIJ) )
freal = coulomb * chrgi * chrgj * ( sr2/rklsq )
freal = freal * ( ERFC(KRIJ) + TWO * rsqpi * KAPPA * rklsq * expar ) * MM % fudgeQQ
atom(ati) % fnonch14(1:3) = atom(ati) % fnonch14(1:3) + freal * rij(1:3)
atom(atj) % fnonch14(1:3) = atom(atj) % fnonch14(1:3) - freal * rij(1:3)
! factor used to compensate the factor1 and factor2 factors ...
! factor3 = 1.0d-20
tterm = coulomb * factor3 * chrgi * chrgj * ERFC(KRIJ)/rklsq * MM % fudgeQQ
! alternative cutoff formula ...
! tterm = tterm - vrecut * chrgi * chrgj + frecut * chrgi * chrgj * ( rklsq-rcut ) * factor3
LJ_14 = LJ_14 + sterm
Coul_14 = Coul_14 + tterm
end if
end do
end do
!====================================================================
call LENNARD_JONES()
!====================================================================
! Morse Intra/Inter potential for H transfer ...
If( allocated(SpecialMorse) ) then
do k = 1 , MM % N_of_atoms - 1
do l = k , MM % N_of_atoms
read_loop2: do n = 1, size(SpecialMorse)
flag1 = ( adjustl( SpecialMorse(n) % MMSymbols(1) ) == adjustl( atom(k) % MMSymbol ) ) .AND. &
( adjustl( SpecialMorse(n) % MMSymbols(2) ) == adjustl( atom(l) % MMSymbol ) )
flag2 = ( adjustl( SpecialMorse(n) % MMSymbols(2) ) == adjustl( atom(k) % MMSymbol ) ) .AND. &
( adjustl( SpecialMorse(n) % MMSymbols(1) ) == adjustl( atom(l) % MMSymbol ) )
if ( flag1 .OR. flag2 ) then
atk = atom(k) % my_id
atl = atom(l) % my_id
rkl(:) = atom(atk) % xyz(:) - atom(atl) % xyz(:)
rkl(:) = rkl(:) - MM % box(:) * DNINT( rkl(:) * MM % ibox(:) ) * PBC(:)
rklq = rkl(1)*rkl(1) + rkl(2)*rkl(2) + rkl(3)*rkl(3)
rklsq = SQRT(rklq)
! Morse potential ...
qterm0 = exp( -SpecialMorse(n) % Parms(3) * ( rklsq - SpecialMorse(n) % Parms(2) ) )
qterm = SpecialMorse(n) % Parms(1) * ( D_ONE - qterm0 )*( D_ONE - qterm0 )
coephi = TWO * SpecialMorse(n) % Parms(1) * SpecialMorse(n) % Parms(3) * qterm0 * ( 1.d0 - qterm0 ) / rklsq
atom(atk) % fMorse(:) = atom(atk) % fMorse(:) - coephi*rkl(:)
atom(atl) % fMorse(:) = atom(atl) % fMorse(:) + coephi*rkl(:)
Morspot = qterm + Morspot
end if
end do read_loop2
end do
end do
end If
!
!====================================================================
! factor used to compensate the factor1 and factor2 factors ...
! factor3 = 1.0d-20
pot_INTRA = ( bdpot + angpot + dihpot )*factor3 + LJ_14 + LJ_intra + Coul_14 + Coul_intra
pot_total = pot_INTER + pot_INTRA
pot_total = pot_total * mol * micro / MM % N_of_molecules
if( QMMM ) then
select case (driver)
case( "slice_CSDM" )
CALL Ehrenfest_master( system , size(basis) )
CALL DecoherenceForce( system , MO_bra , MO_ket , QM%erg , PST )
case( "slice_AO")
CALL EhrenfestForce( system , basis , QM , MO_bra , MO_ket )
case("slice_FSSH")
CALL SH_Force( system , basis , MO_bra , MO_ket , QM , t_rate )
end select
endif
! Get total MM force; force units = J/mts = Newtons ...
do i = 1 , MM % N_of_atoms
atom(i)% f_MM(:) = atom(i)% f_MM(:) + (atom(i) % fbond(:) + &
atom(i) % fang(:) + &
atom(i) % fdihed(:) + &
atom(i) % fnonbd14(:) + &
atom(i) % fnonch14(:) + &
atom(i) % fnonbd(:) + &
atom(i) % fMorse(:) + &
atom(i) % fnonch(:) &
) * Angs_2_mts
atom(i)% ftotal(:) = atom(i)% f_MM(:)
enddo
end subroutine ForceINTRA
!
!
!
!==============
subroutine gmx
!==============
implicit none
! local variables ...
real*8 :: psi , cos_Psi , dtheta
select case( adjustl(molecule(i) % Dihedral_Type(j)) )
case ('cos') ! V = k_phi * [ 1 + cos( n * phi - phi_s ) ]
! Eq. 4.60 (GMX 5.0.5 manual)
term = int(molecule(i) % kdihed0(j,3)) * phi - molecule(i) % kdihed0(j,1)
pterm = molecule(i) % kdihed0(j,2) * ( 1.d0 + cos(term) )
proper_dih = proper_dih + pterm
gamma = - molecule(i) % kdihed0(j,2) * int(molecule(i) % kdihed0(j,3)) * sin(term) * rsinphi * rijkj * rjkkl
case ('harm') ! V = 1/2.k( xi - xi_0 )²
! Eq. 4.59 (GMX 5.0.5 manual)
dtheta = ( phi - molecule(i) % kdihed0(j,1) )
dtheta = dtheta - Dnint( dtheta * 1.d0/TWOPI ) * TWOPI
term = molecule(i) % kdihed0(j,2) * dtheta
pterm = 0.5d0 * term * dtheta
harm_dih = harm_dih + pterm
gamma = term * rsinphi * rijkj * rjkkl
case('cos3') ! V = C0 + C1*cos(phi - 180) + C2*cos^2(phi - 180) + C3*cos^3(phi - 180) + C4*cos^4(phi - 180) + C5*cos(phi - 180)
! Eq. 4.61 (GMX 5.0.5 manual)
psi = phi - PI
cos_Psi = cos(psi)
pterm = molecule(i) % kdihed0(j,1) + &
molecule(i) % kdihed0(j,2) * cos_Psi + &
molecule(i) % kdihed0(j,3) * cos_Psi * cos_Psi + &
molecule(i) % kdihed0(j,4) * cos_Psi * cos_Psi * cos_Psi + &
molecule(i) % kdihed0(j,5) * cos_Psi * cos_Psi * cos_Psi * cos_Psi + &
molecule(i) % kdihed0(j,6) * cos_Psi * cos_Psi * cos_Psi * cos_Psi * cos_Psi
ryck_dih = ryck_dih + pterm
gamma = - sin(psi) * ( molecule(i) % kdihed0(j,2) + &
2.0d0 * molecule(i) % kdihed0(j,3) * cos_Psi + &
3.0d0 * molecule(i) % kdihed0(j,4) * cos_Psi * cos_Psi + &
4.0d0 * molecule(i) % kdihed0(j,5) * cos_Psi * cos_Psi * cos_Psi + &
5.0d0 * molecule(i) % kdihed0(j,6) * cos_Psi * cos_Psi * cos_Psi * cos_PSi ) * rsinphi * rijkj * rjkkl
case ('imp') ! V = k_phi * [ 1 + cos( n * phi - phi_s ) ] (improper)
! Eq. 4.60 (GMX 5.0.5 manual)
term = int(molecule(i) % kdihed0(j,3)) * phi - molecule(i) % kdihed0(j,1)
pterm = molecule(i) % kdihed0(j,2) * ( 1.d0 + cos(term) )
imp_dih = imp_dih + pterm
gamma = - molecule(i) % kdihed0(j,2) * int(molecule(i) % kdihed0(j,3)) * sin(term) * rsinphi * rijkj * rjkkl
case ('chrm') ! V = k_phi * [ 1 + cos( n * phi - phi_s ) ] (multiple)
! Eq. 4.60 (GMX 5.0.5 manual)
term = int(molecule(i) % kdihed0(j,3)) * phi - molecule(i) % kdihed0(j,1)
term1 = int(molecule(i) % kdihed0(j,6)) * phi - molecule(i) % kdihed0(j,4)
term2 = int(molecule(i) % kdihed0(j,9)) * phi - molecule(i) % kdihed0(j,7)
term3 = int(molecule(i) % kdihed0(j,12)) * phi - molecule(i) % kdihed0(j,10)
term4 = int(molecule(i) % kdihed0(j,15)) * phi - molecule(i) % kdihed0(j,13)
pterm = molecule(i) % kdihed0(j,2) * ( 1.d0 + cos(term) )
pterm = pterm + molecule(i) % kdihed0(j,5) * ( 1.d0 + cos(term1) )
pterm = pterm + molecule(i) % kdihed0(j,8) * ( 1.d0 + cos(term2) )
pterm = pterm + molecule(i) % kdihed0(j,11) * ( 1.d0 + cos(term3) )
pterm = pterm + molecule(i) % kdihed0(j,14) * ( 1.d0 + cos(term4) )
proper_dih = proper_dih + pterm
gamma = - molecule(i) % kdihed0(j,2) * int(molecule(i) % kdihed0(j,3)) * sin(term) * rsinphi * rijkj * rjkkl
gamma = gamma - molecule(i) % kdihed0(j,5) * int( molecule(i) % kdihed0(j,6)) * sin(term1) * rsinphi * rijkj * rjkkl
gamma = gamma - molecule(i) % kdihed0(j,8) * int( molecule(i) % kdihed0(j,9)) * sin(term2) * rsinphi * rijkj * rjkkl
gamma = gamma - molecule(i) % kdihed0(j,11) * int(molecule(i) % kdihed0(j,12)) * sin(term3) * rsinphi * rijkj * rjkkl
gamma = gamma - molecule(i) % kdihed0(j,14) * int(molecule(i) % kdihed0(j,15)) * sin(term4) * rsinphi * rijkj * rjkkl
end select
end subroutine gmx
!
!
!
!==================
subroutine not_gmx
!==================
implicit none
! local variables ...
select case( adjustl(molecule(i) % Dihedral_Type(j)) )
case ('cos') ! V = k[1 + cos(n.phi - theta)]
term = int(molecule(i) % kdihed0(j,3)) * phi - molecule(i) % kdihed0(j,1)
pterm = molecule(i) % kdihed0(j,2) * ( 1.d0 + cos(term) )
proper_dih = proper_dih + pterm
gamma = - molecule(i) % kdihed0(j,2) * int(molecule(i) % kdihed0(j,3)) * sin(term) * rsinphi * rijkj * rjkkl
case('harm') ! V = 1/2.k(phi - phi0)²
dphi = phi - molecule(i) % kdihed0(j,1)
dphi = dphi - DNINT(dphi * rtwopi) * twopi
dphi1 = phi - molecule(i) % kdihed0(j,3)
dphi1 = dphi1 - DNINT(dphi1 * rtwopi) * twopi
dphi2 = phi - molecule(i) % kdihed0(j,5)
dphi2 = dphi2 - DNINT(dphi2 * rtwopi) * twopi
term = molecule(i) % kdihed0(j,2) * dphi
term = term + molecule(i) % kdihed0(j,4) * dphi1
term = term + molecule(i) % kdihed0(j,6) * dphi2
pterm = HALF * term * dphi
pterm = pterm + HALF * term1 * dphi1
pterm = pterm + HALF * term2 * dphi2
harm_dih = harm_dih + pterm
gamma = term * rsinphi * rijkj * rjkkl
case('hcos') ! V = 1/2.k[cos(phi) - cos(phi0)]²
dphi = cos(phi) - cos( molecule(i) % kdihed0(j,2) )
term = molecule(i) % kdihed0(j,1) * dphi
pterm = 0.5d0 * term * dphi
gamma = term * rijkj * rjkkl
case('cos3') ! V = 1/2.A1[1 + cos(phi)] + 1/2.A2[1 - cos(2.phi)] + 1/2.A3[1 + cos(3.phi)]
pterm = 0.5d0 * ( molecule(i) % kdihed0(j,1) * ( 1.d0 + cos(phi) ) + &
molecule(i) % kdihed0(j,2) * ( 1.d0 - cos(2.d0*phi) ) + &
molecule(i) % kdihed0(j,3) * ( 1.d0 + cos(3.d0*phi) ) )
gamma = 0.5d0 * ( molecule(i) % kdihed0(j,1) * sin(phi) - &
2.0d0 * molecule(i) % kdihed0(j,2) * sin(2.d0*phi) + &
3.0d0 * molecule(i) % kdihed0(j,3) * sin(3.d0*phi) ) * rsinphi * rijkj * rjkkl
case('ryck') ! V = sum_i^5 Ci.[cos(phi)]^i
eme = cos(phi)
C0 = molecule(i) % kdihed0(j,1) ; C1 = molecule(i) % kdihed0(j,2) ; C2 = molecule(i) % kdihed0(j,3)
C3 = molecule(i) % kdihed0(j,4) ; C4 = molecule(i) % kdihed0(j,5) ; C5 = molecule(i) % kdihed0(j,6)
pterm = C0 - C1 * eme + C2 * eme**2 - C3 * eme**3 + C4 * eme**4 - C5 * eme**5
ryck_dih = ryck_dih + pterm
gamma = - ( -C1 + 2.d0 * C2 * eme - 3.d0 * C3 * eme**2 + 4.d0 * C4 * eme**3 - 5.d0 * C5 * eme**4 ) * rijkj * rjkkl
case('opls') ! V = A0 + 1/2{A1[1 + cos(phi)] + A2[1 - cos(2.phi)] + A3[1 + cos(3.phi)]}
dphi = phi - molecule(i) % kdihed0(j,5)
A0 = molecule(i) % kdihed0(j,1) ; A1 = molecule(i) % kdihed0(j,2) ; A2 = molecule(i) % kdihed0(j,3)
A3 = molecule(i) % kdihed0(j,4)
pterm = A0 + 0.5d0 * ( A1 * (1.d0 + cos(dphi)) + A2 * (1.d0 - cos(2.d0*dphi)) + A3 * (1.d0 + cos(3.d0*dphi)) )
gamma = 0.5d0 * ( A1 * sin(dphi) - 2.d0 * A2 * sin(2.d0*dphi) + 3.d0 * A3 * sin(3.d0*dphi) ) * rsinphi * rijkj * rjkkl
case ('chrm') ! V = k_phi * [ 1 + cos( n * phi - phi_s ) ] (multiple)
! Eq. 4.60 (GMX 5.0.5 manual)
term = int(molecule(i) % kdihed0(j,3)) * phi - molecule(i) % kdihed0(j,1)
term1 = int(molecule(i) % kdihed0(j,6)) * phi - molecule(i) % kdihed0(j,4)
term2 = int(molecule(i) % kdihed0(j,9)) * phi - molecule(i) % kdihed0(j,7)
term3 = int(molecule(i) % kdihed0(j,12)) * phi - molecule(i) % kdihed0(j,10)
term4 = int(molecule(i) % kdihed0(j,15)) * phi - molecule(i) % kdihed0(j,13)
pterm = molecule(i) % kdihed0(j,2) * ( 1.d0 + cos(term) )
pterm = pterm + molecule(i) % kdihed0(j,5) * ( 1.d0 + cos(term1) )
pterm = pterm + molecule(i) % kdihed0(j,8) * ( 1.d0 + cos(term2) )
pterm = pterm + molecule(i) % kdihed0(j,11) * ( 1.d0 + cos(term3) )
pterm = pterm + molecule(i) % kdihed0(j,14) * ( 1.d0 + cos(term4) )
proper_dih = proper_dih + pterm
gamma = - molecule(i) % kdihed0(j,2) * int(molecule(i) % kdihed0(j,3)) * sin(term) * rsinphi * rijkj * rjkkl
gamma = gamma - molecule(i) % kdihed0(j,5) * int( molecule(i) % kdihed0(j,6)) * sin(term1) * rsinphi * rijkj * rjkkl
gamma = gamma - molecule(i) % kdihed0(j,8) * int( molecule(i) % kdihed0(j,9)) * sin(term2) * rsinphi * rijkj * rjkkl
gamma = gamma - molecule(i) % kdihed0(j,11) * int(molecule(i) % kdihed0(j,12)) * sin(term3) * rsinphi * rijkj * rjkkl
gamma = gamma - molecule(i) % kdihed0(j,14) * int(molecule(i) % kdihed0(j,15)) * sin(term4) * rsinphi * rijkj * rjkkl
case('none')
pterm = 0.d0
gamma = 0.d0
end select
end subroutine not_gmx
!
!
!
!===================
function ERFC ( X )
!===================
real*8 :: ERFC
real*8 :: A1, A2, A3, A4, A5, P, T, X, XSQ, TP
parameter ( A1 = 0.254829592, A2 = -0.284496736 )
parameter ( A3 = 1.421413741, A4 = -1.453122027 )
parameter ( A5 = 1.061405429, P = 0.3275911 )
T = 1.0 / ( 1.0 + P * X )
XSQ = X * X
TP = T * (A1 + T * (A2 + T * (A3 + T * (A4 + T * A5))))
ERFC = TP * EXP ( -XSQ )
end function ERFC
!
!
!
!====================
function DEL ( X,Y )
!====================
real*8 :: DEL
integer :: X, Y
if (X == Y) then
DEL = 1.
else if (X /= Y) then
DEL = 0.
end if
end function DEL
!
!
!
! Get the flags status based on the current pair number and analysed atoms
!
! @param[in] pair_number
! @param[in] ati
! @param[in] atj
! @param[out] flag1
! @param[out] flag2
subroutine GET_FLAGS(flag1, flag2, pair_number, ati, atj)
implicit none
logical, intent(out) :: flag1, flag2
integer, intent(in) :: pair_number, ati, atj
character(4) :: mm_symbols1, mm_symbols2, symbol_i, symbol_j
symbol_i = atom(ati) % MMSymbol
symbol_j = atom(atj) % MMSymbol
mm_symbols1 = SpecialPairs(pair_number)%MMSymbols(1)
mm_symbols2 = SpecialPairs(pair_number)%MMSymbols(2)
flag1 = adjustl(mm_symbols1) == adjustl(symbol_i)
flag1 = flag1 .AND. (adjustl(mm_symbols2) == adjustl(symbol_j))
flag2 = adjustl(mm_symbols2) == adjustl(symbol_i)
flag2 = flag2 .AND. (adjustl(mm_symbols1) == adjustl(symbol_j))
end subroutine GET_FLAGS
!
!
!
!
!==========================
subroutine LENNARD_JONES()
!==========================
implicit none
real*8 :: local_fnonbd(size(atom), 3)
real*8 :: local_fnonch(size(atom), 3)
If( allocated(SpecialPairs) ) there_are_NB_SpecialPairs = .true.
local_fnonbd = 0
local_fnonch = 0
!$OMP parallel &
!$OMP default(shared) &
!$OMP private(i) &
!$OMP reduction(+ : LJ_intra, Coul_intra, local_fnonbd, local_fnonch)
do i = 1 , MM % N_of_molecules
!$OMP do &
!$OMP schedule(static) &
!$OMP private(j, n, ati, ati1, atj, atj1, chrgi, chrgj, rij, rklq, rklsq, sr2, sr6, sr12, &
!$OMP eps, flag1, flag2, fs, sterm, tterm, expar, freal, KRIJ)
do j = 1 , molecule(i) % NintraLJ
ati = molecule(i) % IntraLJ(j,1)
ati1 = atom(ati) % my_intra_id + species_offset( atom(ati)%my_species )
atj = molecule(i) % IntraLJ(j,2)
atj1 = atom(atj) % my_intra_id + species_offset( atom(atj)%my_species )
if (.NOT. (atom(atj)%flex .OR. atom(ati)%flex)) then
cycle
end if
chrgi = atom(ati)%charge
chrgj = atom(atj)%charge
rij(:) = atom(ati)%xyz(:) - atom(atj)%xyz(:)
rij(:) = rij(:) - MM%box(:) * DNINT(rij(:)*MM % ibox(:)) * PBC(:)
rklq = SUM(rij(:) ** 2)
if (.NOT. (rklq < rcutsq)) then
cycle
end if
! Lennard Jones ...
! AMBER FF :: GMX COMB-RULE 2
if (MM%CombinationRule == 2) then
sr2 = (atom(ati)%sig + atom(atj)%sig) ** 2 / rklq
! OPLS FF :: GMX COMB-RULE 3
else if (MM%CombinationRule == 3) then
sr2 = (atom(ati)%sig ** 2) * (atom(atj)%sig ** 2) / rklq
end if
eps = atom(ati)%eps * atom(atj)%eps
If( there_are_NB_SpecialPairs ) then ! <== check whether (I,J) is a SpecialPair ...
do n = 1, SIZE(SpecialPairs)
call GET_FLAGS(flag1, flag2, n, ati, atj)
if (.NOT. (flag1 .OR. flag2)) cycle
sr2 = ( SpecialPairs(n)%Parms(1) * SpecialPairs(n)%Parms(1) ) / rklq
eps = SpecialPairs(n) % Parms(2)
exit
end do
end if
rklsq = SQRT(rklq)
fs = 24.0 * eps * (2.0 * sr2 ** 6.0 - sr2 ** 3.0)
! with force cut-off ...
fs = (fs / rklq) - fscut(ati1, atj1) / rklsq
! factor3 used to compensate factor1 ...
sterm = 4.0 * eps * factor3 * (sr2 ** 6.0 - sr2 ** 3.0)
! alternative formula with cutoff ...
sterm = sterm - vscut(ati1, atj1) + fscut(ati1, atj1) * (rklsq - rcut)
! Real part (Number of charges equal to the number of sites)
sr2 = 1.0 / rklq
KRIJ = KAPPA * rklsq
expar = EXP(-(KRIJ ** 2.0))
freal = coulomb * chrgi * chrgj * (sr2 / rklsq)
freal = freal * (ERFC(KRIJ) + 2.0 * rsqpi * KAPPA * rklsq * expar)
! with force cut-off ...
freal = freal - frecut / rklsq * chrgi * chrgj
! factor used to compensate factor1 ...
! factor3 = 1.0d-20
tterm = (coulomb * factor3 * chrgi * chrgj * ERFC(KRIJ)) / rklsq
! alternative formula with cutoff ...
tterm = tterm - vrecut * chrgi * chrgj
tterm = tterm + frecut * chrgi * chrgj * (rklsq - rcut) * factor3
! atom(ati)%fnonbd(:) = atom(ati)%fnonbd(:) + fs * rij(:)
! atom(atj)%fnonbd(:) = atom(atj)%fnonbd(:) - fs * rij(:)
! atom(ati)%fnonch(:) = atom(ati)%fnonch(:) + freal * rij(:)
! atom(atj)%fnonch(:) = atom(atj)%fnonch(:) - freal * rij(:)
local_fnonbd(ati,:) = local_fnonbd(ati,:) + fs * rij(:)
local_fnonbd(atj,:) = local_fnonbd(atj,:) - fs * rij(:)
local_fnonch(ati,:) = local_fnonch(ati,:) + freal * rij(:)
local_fnonch(atj,:) = local_fnonch(atj,:) - freal * rij(:)
LJ_intra = LJ_intra + sterm
Coul_intra = Coul_intra + tterm
end do
!$OMP end do
end do
!$OMP end parallel
!$OMP parallel do schedule(static) default(shared)
do n = 1, size(atom)
atom(n)%fnonbd(:) = atom(n)%fnonbd(:) + local_fnonbd(n,:)
atom(n)%fnonch(:) = atom(n)%fnonch(:) + local_fnonch(n,:)
end do
!$OMP end parallel do
end subroutine LENNARD_JONES
!
!
!
!=====================
subroutine ForceQMMM
!=====================
implicit none
! local variables ...
integer :: i
If( driver == "slice_CSDM" ) then
do i = 1 , MM % N_of_atoms
atom(i)% f_QM(:) = atom(i)% Ehrenfest(:) + atom(i)% f_CSDM(:)
atom(i)% ftotal(:) = atom(i)% f_MM(:) + atom(i)% f_QM(:)
end do
else
do i = 1 , MM % N_of_atoms
atom(i)% f_QM(:) = atom(i)% Ehrenfest(:)
atom(i)% ftotal(:) = atom(i)% f_MM(:) + atom(i)% f_QM(:)
end do
endif
end subroutine ForceQMMM
!
!
!
!====================================================================
subroutine StoreQMArgs_CSDM( sys , vec , mtx1 , mtx2 , Eigen , PSE )
!====================================================================
implicit none
type(structure), intent(inout):: sys
type(STO_basis), intent(in) :: vec(:)
complex*16 , intent(in) :: mtx1(:,:)
complex*16 , intent(in) :: mtx2(:,:)
type(R_eigen) , intent(in) :: Eigen
integer , intent(in) :: PSE(2)
! local variables ...
basis = vec
system = sys
PST = PSE
MO_bra = mtx1
MO_ket = mtx2
QM% erg = Eigen% erg
QM% L = Eigen% L
QM% R = Eigen% R
QM% Fermi_state = Eigen% Fermi_state
end subroutine StoreQMArgs_CSDM
!
!
!
!===================================================================================
subroutine StoreQMArgs_AO_and_FSSH( sys , vec , mtx1 , mtx2 , Eigen , delta , txt )
!===================================================================================
implicit none
type(structure) , intent(inout):: sys
type(STO_basis) , intent(in) :: vec(:)
complex*16 , intent(in) :: mtx1(:,:)
complex*16 , intent(in) :: mtx2(:,:)
type(R_eigen) , intent(in) :: Eigen
real*8 , optional , intent(in) :: delta
character(len=*) , optional , intent(in) :: txt
! local variables ...
basis = vec
system = sys
MO_bra = mtx1
MO_ket = mtx2
if(present(delta)) t_rate = delta
if(present(txt )) mode = txt
QM% erg = Eigen% erg
QM% L = Eigen% L
QM% R = Eigen% R
QM% Fermi_state = Eigen% Fermi_state
end subroutine StoreQMArgs_AO_and_FSSH
!
!
!
!===================================================
subroutine AllocateQMArgs( BasisSize , SystemSize )
!===================================================
implicit none
integer , intent(in) :: BasisSize
integer , intent(in) :: SystemSize
allocate(basis (BasisSize) )
allocate(MO_bra (BasisSize,n_part) )
allocate(MO_ket (BasisSize,n_part) )
CALL Allocate_Structures( SystemSize , System )
allocate(QM%erg(BasisSize) )
allocate(QM%L (BasisSize,BasisSize))
allocate(QM%R (BasisSize,BasisSize))
end subroutine AllocateQMArgs
!
!
end module F_intra_m