forked from TinkerTools/tinker
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathechgtrn3.f
922 lines (922 loc) · 27 KB
/
echgtrn3.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
c
c
c ############################################################
c ## COPYRIGHT (C) 2018 by Joshua Rackers & Jay W. Ponder ##
c ## All Rights Reserved ##
c ############################################################
c
c ##################################################################
c ## ##
c ## subroutine echgtrn3 -- charge transfer energy & analysis ##
c ## ##
c ##################################################################
c
c
c "echgtrn3" calculates the charge transfer energy; also partitions
c the energy among the atoms
c
c
subroutine echgtrn3
use limits
implicit none
c
c
c choose method for summing over charge transfer interactions
c
if (use_mlist) then
call echgtrn3c
else if (use_lights) then
call echgtrn3b
else
call echgtrn3a
end if
return
end
c
c
c #############################################################
c ## ##
c ## subroutine echgtrn3a -- double loop chgtrn analysis ##
c ## ##
c #############################################################
c
c
c "echgtrn3a" calculates the charge transfer interaction energy
c and also partitions the energy among the atoms using a pairwise
c double loop
c
c
subroutine echgtrn3a
use action
use analyz
use atomid
use atoms
use bound
use chgpot
use chgtrn
use cell
use couple
use ctrpot
use energi
use group
use inform
use inter
use iounit
use molcul
use mplpot
use mpole
use shunt
use usage
implicit none
integer i,j,k
integer ii,kk
integer jcell
real*8 e,f,fgrp
real*8 r,r2,r3
real*8 r4,r5
real*8 xi,yi,zi
real*8 xr,yr,zr
real*8 chgi,chgk
real*8 chgik
real*8 alphai,alphak
real*8 expi,expk
real*8 expik
real*8 taper
real*8, allocatable :: mscale(:)
logical proceed,usei
logical header,huge
character*6 mode
c
c
c zero out the charge transfer energy and partitioning terms
c
nect = 0
ect = 0.0d0
do i = 1, n
aect(i) = 0.0d0
end do
if (nct .eq. 0) return
c
c perform dynamic allocation of some local arrays
c
allocate (mscale(n))
c
c initialize connected atom exclusion coefficients
c
do i = 1, n
mscale(i) = 1.0d0
end do
c
c set conversion factor, cutoff and switching coefficients
c
f = electric / dielec
mode = 'CHGTRN'
call switch (mode)
c
c print header information if debug output was requested
c
header = .true.
if (debug .and. nct.ne.0) then
header = .false.
write (iout,10)
10 format (/,' Individual Charge Transfer Interactions :',
& //,' Type',14x,'Atom Names',15x,'Distance',
& 8x,'Energy',/)
end if
c
c calculate the charge transfer energy term
c
do ii = 1, npole-1
i = ipole(ii)
xi = x(i)
yi = y(i)
zi = z(i)
chgi = chgct(ii)
alphai = dmpct(ii)
if (alphai .eq. 0.0d0) alphai = 1000.0d0
usei = use(i)
c
c set exclusion coefficients for connected atoms
c
do j = 1, n12(i)
mscale(i12(j,i)) = m2scale
end do
do j = 1, n13(i)
mscale(i13(j,i)) = m3scale
end do
do j = 1, n14(i)
mscale(i14(j,i)) = m4scale
end do
do j = 1, n15(i)
mscale(i15(j,i)) = m5scale
end do
c
c evaluate all sites within the cutoff distance
c
do kk = ii+1, npole
k = ipole(kk)
proceed = .true.
if (use_group) call groups (proceed,fgrp,i,k,0,0,0,0)
if (.not. use_intra) proceed = .true.
if (proceed) proceed = (usei .or. use(k))
if (proceed) then
xr = x(k) - xi
yr = y(k) - yi
zr = z(k) - zi
if (use_bounds) call image (xr,yr,zr)
r2 = xr*xr + yr* yr + zr*zr
if (r2 .le. off2) then
r = sqrt(r2)
chgk = chgct(kk)
alphak = dmpct(kk)
if (alphak .eq. 0.0d0) alphak = 1000.0d0
if (ctrntyp .eq. 'SEPARATE') then
expi = exp(-alphai*r)
expk = exp(-alphak*r)
e = -chgi*expk - chgk*expi
else
chgik = sqrt(abs(chgi*chgk))
expik = exp(-0.5d0*(alphai+alphak)*r)
e = -chgik * expik
end if
e = f * e * mscale(k)
c
c use energy switching if near the cutoff distance
c
if (r2 .gt. cut2) then
r3 = r2 * r
r4 = r2 * r2
r5 = r2 * r3
taper = c5*r5 + c4*r4 + c3*r3
& + c2*r2 + c1*r + c0
e = e * taper
end if
c
c scale the interaction based on its group membership
c
if (use_group) e = e * fgrp
c
c increment the overall charge transfer energy components
c
if (e .ne. 0.0d0) then
nect = nect + 1
ect = ect + e
aect(i) = aect(i) + 0.5d0*e
aect(k) = aect(k) + 0.5d0*e
end if
c
c increment the total intermolecular energy
c
if (molcule(i) .ne. molcule(k)) then
einter = einter + e
end if
c
c print a message if the energy of this interaction is large
c
huge = (e .gt. 10.0d0)
if ((debug.and.e.ne.0.0d0)
& .or. (verbose.and.huge)) then
if (header) then
header = .false.
write (iout,20)
20 format (/,' Individual Charge Transfer',
& ' Interactions :',
& //,' Type',14x,'Atom Names',
& 15x,'Distance',8x,'Energy',/)
end if
write (iout,30) i,name(i),k,name(k),r,e
30 format (' ChgTrn',4x,2(i7,'-',a3),
& 9x,f10.4,2x,f12.4)
end if
end if
end if
end do
c
c reset exclusion coefficients for connected atoms
c
do j = 1, n12(i)
mscale(i12(j,i)) = 1.0d0
end do
do j = 1, n13(i)
mscale(i13(j,i)) = 1.0d0
end do
do j = 1, n14(i)
mscale(i14(j,i)) = 1.0d0
end do
do j = 1, n15(i)
mscale(i15(j,i)) = 1.0d0
end do
end do
c
c for periodic boundary conditions with large cutoffs
c neighbors must be found by the replicates method
c
if (use_replica) then
c
c calculate interaction energy with other unit cells
c
do ii = 1, npole
i = ipole(ii)
xi = x(i)
yi = y(i)
zi = z(i)
chgi = chgct(ii)
alphai = dmpct(ii)
if (alphai .eq. 0.0d0) alphai = 1000.0d0
usei = use(i)
c
c set exclusion coefficients for connected atoms
c
do j = 1, n12(i)
mscale(i12(j,i)) = m2scale
end do
do j = 1, n13(i)
mscale(i13(j,i)) = m3scale
end do
do j = 1, n14(i)
mscale(i14(j,i)) = m4scale
end do
do j = 1, n15(i)
mscale(i15(j,i)) = m5scale
end do
c
c evaluate all sites within the cutoff distance
c
do kk = i, npole
k = ipole(kk)
proceed = .true.
if (use_group) call groups (proceed,fgrp,i,k,0,0,0,0)
if (.not. use_intra) proceed = .true.
if (proceed) proceed = (usei .or. use(k))
if (proceed) then
do jcell = 2, ncell
xr = x(k) - xi
yr = y(k) - yi
zr = z(k) - zi
call imager (xr,yr,zr,jcell)
r2 = xr*xr + yr* yr + zr*zr
if (r2 .le. off2) then
r = sqrt(r2)
chgk = chgct(kk)
alphak = dmpct(kk)
if (alphak .eq. 0.0d0) alphak = 1000.0d0
if (ctrntyp .eq. 'SEPARATE') then
expi = exp(-alphai*r)
expk = exp(-alphak*r)
e = -chgi*expk - chgk*expi
else
chgik = sqrt(abs(chgi*chgk))
expik = exp(-0.5d0*(alphai+alphak)*r)
e = -chgik * expik
end if
e = f * e * mscale(k)
c
c use energy switching if near the cutoff distance
c
if (r2 .gt. cut2) then
r3 = r2 * r
r4 = r2 * r2
r5 = r2 * r3
taper = c5*r5 + c4*r4 + c3*r3
& + c2*r2 + c1*r + c0
e = e * taper
end if
c
c scale the interaction based on its group membership
c
if (use_group) e = e * fgrp
c
c increment the overall charge transfer energy components
c
if (e .ne. 0.0d0) then
nect = nect + 1
if (i .eq. k) then
ect = ect + 0.5d0*e
aect(i) = aect(i) + 0.5d0*e
else
ect = ect + e
aect(i) = aect(i) + 0.5d0*e
aect(k) = aect(k) + 0.5d0*e
end if
end if
c
c increment the total intermolecular energy
c
if (molcule(i) .ne. molcule(k)) then
einter = einter + e
end if
c
c print a message if the energy of this interaction is large
c
huge = (e .gt. 10.0d0)
if ((debug.and.e.ne.0.0d0)
& .or. (verbose.and.huge)) then
if (header) then
header = .false.
write (iout,40)
40 format (/,' Individual Charge Transfer',
& ' Interactions :',
& //,' Type',14x,'Atom Names',
& 15x,'Distance',8x,'Energy',/)
end if
write (iout,50) i,name(i),k,name(k),r,e
50 format (' ChgTrn',4x,2(i7,'-',a3),2x,
& '(XTAL)',1x,f10.4,2x,f12.4)
end if
end if
end do
end if
end do
c
c reset exclusion coefficients for connected atoms
c
do j = 1, n12(i)
mscale(i12(j,i)) = 1.0d0
end do
do j = 1, n13(i)
mscale(i13(j,i)) = 1.0d0
end do
do j = 1, n14(i)
mscale(i14(j,i)) = 1.0d0
end do
do j = 1, n15(i)
mscale(i15(j,i)) = 1.0d0
end do
end do
end if
c
c perform deallocation of some local arrays
c
deallocate (mscale)
return
end
c
c
c ###############################################################
c ## ##
c ## subroutine echgtrn3b -- method lights chgtrn analysis ##
c ## ##
c ###############################################################
c
c
c "echgtrn3b" calculates the charge transfer interaction energy
c and also partitions the energy among the atoms using the method
c of lights
c
c
subroutine echgtrn3b
use action
use analyz
use atomid
use atoms
use bound
use boxes
use chgpot
use chgtrn
use cell
use couple
use ctrpot
use energi
use group
use inform
use inter
use iounit
use light
use molcul
use mplpot
use mpole
use shunt
use usage
implicit none
integer i,j,k
integer ii,kk
integer kgy,kgz
integer start,stop
real*8 e,f,fgrp
real*8 r,r2,r3
real*8 r4,r5
real*8 xi,yi,zi
real*8 xr,yr,zr
real*8 chgi,chgk
real*8 chgik
real*8 alphai,alphak
real*8 expi,expk
real*8 expik
real*8 taper
real*8, allocatable :: mscale(:)
real*8, allocatable :: xsort(:)
real*8, allocatable :: ysort(:)
real*8, allocatable :: zsort(:)
logical proceed,usei
logical unique,repeat
logical header,huge
character*6 mode
c
c
c zero out the charge transfer energy
c
nect = 0
ect = 0.0d0
do i = 1, n
aect(i) = 0.0d0
end do
if (nct .eq. 0) return
c
c perform dynamic allocation of some local arrays
c
allocate (mscale(n))
allocate (xsort(8*n))
allocate (ysort(8*n))
allocate (zsort(8*n))
c
c initialize connected atom exclusion coefficients
c
do i = 1, n
mscale(i) = 1.0d0
end do
c
c set conversion factor, cutoff and switching coefficients
c
f = electric / dielec
mode = 'CHGTRN'
call switch (mode)
c
c print header information if debug output was requested
c
header = .true.
if (debug .and. nct.ne.0) then
header = .false.
write (iout,10)
10 format (/,' Individual Charge Transfer Interactions :',
& //,' Type',14x,'Atom Names',15x,'Distance',
& 8x,'Energy',/)
end if
c
c transfer the interaction site coordinates to sorting arrays
c
do ii = 1, npole
i = ipole(ii)
xsort(i) = x(i)
ysort(i) = y(i)
zsort(i) = z(i)
end do
c
c use the method of lights to generate neighbors
c
unique = .true.
call lights (off,nct,xsort,ysort,zsort,unique)
c
c calculate the charge transfer energy term
c
do ii = 1, npole
i = ipole(ii)
xi = x(i)
yi = y(i)
zi = z(i)
chgi = chgct(ii)
alphai = dmpct(ii)
if (alphai .eq. 0.0d0) alphai = 1000.0d0
usei = use(i)
c
c set exclusion coefficients for connected atoms
c
do j = 1, n12(i)
mscale(i12(j,i)) = m2scale
end do
do j = 1, n13(i)
mscale(i13(j,i)) = m3scale
end do
do j = 1, n14(i)
mscale(i14(j,i)) = m4scale
end do
do j = 1, n15(i)
mscale(i15(j,i)) = m5scale
end do
c
c loop over method of lights neighbors of current atom
c
if (kbx(ii) .le. kex(ii)) then
repeat = .false.
start = kbx(ii) + 1
stop = kex(ii)
else
repeat = .true.
start = 1
stop = kex(ii)
end if
20 continue
do j = start, stop
kk = locx(j)
kgy = rgy(kk)
if (kby(ii) .le. key(ii)) then
if (kgy.lt.kby(ii) .or. kgy.gt.key(ii)) goto 50
else
if (kgy.lt.kby(ii) .and. kgy.gt.key(ii)) goto 50
end if
kgz = rgz(kk)
if (kbz(ii) .le. kez(ii)) then
if (kgz.lt.kbz(ii) .or. kgz.gt.kez(ii)) goto 50
else
if (kgz.lt.kbz(ii) .and. kgz.gt.kez(ii)) goto 50
end if
k = ipole(kk-((kk-1)/npole)*npole)
c
c decide whether to compute the current interaction
c
proceed = .true.
if (use_group) call groups (proceed,fgrp,i,k,0,0,0,0)
if (.not. use_intra) proceed = .true.
if (proceed) proceed = (usei .or. use(k))
if (proceed) then
xr = x(k) - xi
yr = y(k) - yi
zr = z(k) - zi
if (use_bounds) then
if (abs(xr) .gt. xcell2) xr = xr - sign(xcell,xr)
if (abs(yr) .gt. ycell2) yr = yr - sign(ycell,yr)
if (abs(zr) .gt. zcell2) zr = zr - sign(zcell,zr)
if (monoclinic) then
xr = xr + zr*beta_cos
zr = zr * beta_sin
else if (triclinic) then
xr = xr + yr*gamma_cos + zr*beta_cos
yr = yr*gamma_sin + zr*beta_term
zr = zr * gamma_term
end if
end if
r2 = xr*xr + yr* yr + zr*zr
if (r2 .le. off2) then
r = sqrt(r2)
chgk = chgct(kk)
alphak = dmpct(kk)
if (alphak .eq. 0.0d0) alphak = 1000.0d0
if (ctrntyp .eq. 'SEPARATE') then
expi = exp(-alphai*r)
expk = exp(-alphak*r)
e = -chgi*expk - chgk*expi
else
chgik = sqrt(abs(chgi*chgk))
expik = exp(-0.5d0*(alphai+alphak)*r)
e = -chgik * expik
end if
e = f * e * mscale(k)
c
c use energy switching if near the cutoff distance
c
if (r2 .gt. cut2) then
r3 = r2 * r
r4 = r2 * r2
r5 = r2 * r3
taper = c5*r5 + c4*r4 + c3*r3
& + c2*r2 + c1*r + c0
e = e * taper
end if
c
c scale the interaction based on its group membership
c
if (use_group) e = e * fgrp
c
c increment the overall charge transfer energy components
c
if (e .ne. 0.0d0) then
nect = nect + 1
ect = ect + e
aect(i) = aect(i) + 0.5d0*e
aect(k) = aect(k) + 0.5d0*e
end if
c
c increment the total intermolecular energy
c
if (molcule(i) .ne. molcule(k)) then
einter = einter + e
end if
c
c print a message if the energy of this interaction is large
c
huge = (e .gt. 10.0d0)
if ((debug.and.e.ne.0.0d0)
& .or. (verbose.and.huge)) then
if (header) then
header = .false.
write (iout,30)
30 format (/,' Individual Charge Transfer',
& ' Interactions :',
& //,' Type',14x,'Atom Names',
& 15x,'Distance',8x,'Energy',/)
end if
write (iout,40) i,name(i),k,name(k),r,e
40 format (' ChgTrn',4x,2(i7,'-',a3),
& 9x,f10.4,2x,f12.4)
end if
end if
end if
50 continue
end do
if (repeat) then
repeat = .false.
start = kbx(ii) + 1
stop = nlight
goto 20
end if
c
c reset exclusion coefficients for connected atoms
c
do j = 1, n12(i)
mscale(i12(j,i)) = 1.0d0
end do
do j = 1, n13(i)
mscale(i13(j,i)) = 1.0d0
end do
do j = 1, n14(i)
mscale(i14(j,i)) = 1.0d0
end do
do j = 1, n15(i)
mscale(i15(j,i)) = 1.0d0
end do
end do
c
c perform deallocation of some local arrays
c
deallocate (mscale)
deallocate (xsort)
deallocate (ysort)
deallocate (zsort)
return
end
c
c
c ###############################################################
c ## ##
c ## subroutine echgtrn3c -- neighbor list chgtrn analysis ##
c ## ##
c ###############################################################
c
c
c "echgtrn3c" calculates the charge transfer interaction energy
c and also partitions the energy among the atoms using a pairwise
c neighbor list
c
c
subroutine echgtrn3c
use action
use analyz
use atomid
use atoms
use bound
use chgpot
use chgtrn
use cell
use couple
use ctrpot
use energi
use group
use inform
use inter
use iounit
use molcul
use mplpot
use mpole
use neigh
use shunt
use usage
implicit none
integer i,j,k
integer ii,kk,kkk
real*8 e,f,fgrp
real*8 r,r2,r3
real*8 r4,r5
real*8 xi,yi,zi
real*8 xr,yr,zr
real*8 chgi,chgk
real*8 chgik
real*8 alphai,alphak
real*8 expi,expk
real*8 expik
real*8 taper
real*8, allocatable :: mscale(:)
logical proceed,usei
logical header,huge
character*6 mode
c
c zero out the charge transfer energy and partitioning terms
c
nect = 0
ect = 0.0d0
do i = 1, n
aect(i) = 0.0d0
end do
if (nct .eq. 0) return
c
c perform dynamic allocation of some local arrays
c
allocate (mscale(n))
c
c initialize connected atom exclusion coefficients
c
do i = 1, n
mscale(i) = 1.0d0
end do
c
c set conversion factor, cutoff and switching coefficients
c
f = electric / dielec
mode = 'CHGTRN'
call switch (mode)
c
c print header information if debug output was requested
c
header = .true.
if (debug .and. nct.ne.0) then
header = .false.
write (iout,10)
10 format (/,' Individual Charge Transfer Interactions :',
& //,' Type',14x,'Atom Names',15x,'Distance',
& 8x,'Energy',/)
end if
c
c OpenMP directives for the major loop structure
c
!$OMP PARALLEL default(private)
!$OMP& shared(npole,ipole,x,y,z,chgct,dmpct,n12,i12,n13,i13,
!$OMP& n14,i14,n15,i15,m2scale,m3scale,m4scale,m5scale,nelst,
!$OMP& elst,use,use_group,use_intra,use_bounds,ctrntyp,f,off2,
!$OMP& cut2,molcule,c0,c1,c2,c3,c4,c5,name,verbose,debug,
!$OMP& header,iout)
!$OMP& firstprivate(mscale) shared(ect,nect,aect,einter)
!$OMP DO reduction(+:ect,nect,aect,einter) schedule(guided)
c
c compute the charge transfer energy
c
do ii = 1, npole
i = ipole(ii)
xi = x(i)
yi = y(i)
zi = z(i)
chgi = chgct(ii)
alphai = dmpct(ii)
if (alphai .eq. 0.0d0) alphai = 1000.0d0
usei = use(i)
c
c set exclusion coefficients for connected atoms
c
do j = 1, n12(i)
mscale(i12(j,i)) = m2scale
end do
do j = 1, n13(i)
mscale(i13(j,i)) = m3scale
end do
do j = 1, n14(i)
mscale(i14(j,i)) = m4scale
end do
do j = 1, n15(i)
mscale(i15(j,i)) = m5scale
end do
c
c evaluate all sites within the cutoff distance
c
do kkk = 1, nelst(ii)
kk = elst(kkk,ii)
k = ipole(kk)
proceed = .true.
if (use_group) call groups (proceed,fgrp,i,k,0,0,0,0)
if (.not. use_intra) proceed = .true.
if (proceed) proceed = (usei .or. use(k))
if (proceed) then
xr = x(k) - xi
yr = y(k) - yi
zr = z(k) - zi
if (use_bounds) call image (xr,yr,zr)
r2 = xr*xr + yr* yr + zr*zr
if (r2 .le. off2) then
r = sqrt(r2)
chgk = chgct(kk)
alphak = dmpct(kk)
if (alphak .eq. 0.0d0) alphak = 1000.0d0
if (ctrntyp .eq. 'SEPARATE') then
expi = exp(-alphai*r)
expk = exp(-alphak*r)
e = -chgi*expk - chgk*expi
else
chgik = sqrt(abs(chgi*chgk))
expik = exp(-0.5d0*(alphai+alphak)*r)
e = -chgik * expik
end if
e = f * e * mscale(k)
c
c use energy switching if near the cutoff distance
c
if (r2 .gt. cut2) then
r3 = r2 * r
r4 = r2 * r2
r5 = r2 * r3
taper = c5*r5 + c4*r4 + c3*r3
& + c2*r2 + c1*r + c0
e = e * taper
end if
c
c scale the interaction based on its group membership
c
if (use_group) e = e * fgrp
c
c increment the overall charge transfer energy components
c
if (e .ne. 0.0d0) then
nect = nect + 1
ect = ect + e
aect(i) = aect(i) + 0.5d0*e
aect(k) = aect(k) + 0.5d0*e
end if
c
c increment the total intermolecular energy
c
if (molcule(i) .ne. molcule(k)) then
einter = einter + e
end if
c
c print a message if the energy of this interaction is large
c
huge = (e .gt. 10.0d0)
if ((debug.and.e.ne.0.0d0)
& .or. (verbose.and.huge)) then
if (header) then
header = .false.
write (iout,20)
20 format (/,' Individual Charge Transfer',
& ' Interactions :',
& //,' Type',14x,'Atom Names',
& 15x,'Distance',8x,'Energy',/)
end if
write (iout,30) i,name(i),k,name(k),r,e
30 format (' ChgTrn',4x,2(i7,'-',a3),
& 9x,f10.4,2x,f12.4)
end if
end if
end if
end do
c
c reset exclusion coefficients for connected atoms
c
do j = 1, n12(i)
mscale(i12(j,i)) = 1.0d0
end do
do j = 1, n13(i)
mscale(i13(j,i)) = 1.0d0
end do
do j = 1, n14(i)
mscale(i14(j,i)) = 1.0d0
end do
do j = 1, n15(i)
mscale(i15(j,i)) = 1.0d0
end do
end do
c
c OpenMP directives for the major loop structure
c
!$OMP END DO
!$OMP END PARALLEL
c
c perform deallocation of some local arrays
c
deallocate (mscale)
return
end