-
Notifications
You must be signed in to change notification settings - Fork 0
/
hme2.f90
1457 lines (1363 loc) · 50.3 KB
/
hme2.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
!program hme_wrapper
!implicit none
!integer(kind=4)::natom,ngrid,qtot
!real(kind=8),allocatable::gridcrd(:,:),atomcrd(:,:)
!real(kind=8),allocatable::esp(:),mpinit(:),dpinit(:,:),qpinit(:,:)
!integer(kind=4)::level(3)
!character(len=30)::espfile
!character(len=10)::fitmethod
!integer(kind=4)::istrnt
!real(kind=8)::strength
!integer(kind=4)::nlgn
!integer(kind=4),allocatable::lgnatomidx(:,:)
!real(kind=8),allocatable::lgncons(:)
!integer(kind=4)::nlgnatom_in_this_cnstrnt
!integer(kind=4)::ilgn
!integer(kind=4)::nsymmconstraint,isymmconstraint
!integer(kind=4),allocatable::isymmconstraintpair(:,:)
!integer(kind=4),allocatable::atmidx(:)
!integer(kind=4),allocatable::coef(:)
!real(kind=8)::ssres,crosscorr
!integer(kind=4)::i,j,k,m,n
!espfile='aceticacid.esp'
!open(11,file=espfile)
!read(11,*)natom,ngrid
!if(allocated(gridcrd))deallocate(gridcrd)
!allocate(gridcrd(3,ngrid))
!if(allocated(atomcrd))deallocate(atomcrd)
!allocate(atomcrd(3,natom))
!if(allocated(esp))deallocate(esp)
!allocate(esp(ngrid))
!if(allocated(mpinit))deallocate(mpinit)
!allocate(mpinit(natom))
!if(allocated(dpinit))deallocate(dpinit)
!allocate(dpinit(3,natom))
!if(allocated(qpinit))deallocate(qpinit)
!allocate(qpinit(5,natom))
!if(allocated(atmidx))deallocate(atmidx)
!allocate(atmidx(natom))
!if(allocated(coef))deallocate(coef)
!allocate(coef(natom))
!read(11,'(16X,3E16.7)')((atomcrd(i,j),i=1,3),j=1,natom)
!do j = 1, natom
!! write(*,'(A,I7,16X,3ES16.7)')'Atom ', j, atomcrd(:,j)
!end do
!read(11,'(4E16.7)')((esp(j),(gridcrd(i,j),i=1,3)),j=1,ngrid)
!do j = 1, ngrid
!! write(*,'(A,I7,4ES16.7)')'Grid ', j, esp(j),gridcrd(:,j)
!end do
!close(11)
!qtot=0
!mpinit(1)= -0.602012
!mpinit(2)= 0.224271
!mpinit(3)= 0.224271
!mpinit(4)= 0.224271
!mpinit(5)= 0.762176
!mpinit(6)= -0.701809
!mpinit(7)= 0.459737
!mpinit(8)= -0.533535
!mpinit=0.d0
!dpinit=0.d0
!qpinit=0.d0
!fitmethod='resp'
!istrnt=1
!strength=0.0005d0
!level(1)=1
!level(2)=1
!level(3)=1
!read*,nlgn
!allocate(lgnatomidx(natom,nlgn))
!allocate(lgncons(nlgn))
!do ilgn = 1, nlgn
! ! parse the input for Lagrange constraints
! read(*,*)nlgnatom_in_this_cnstrnt,lgncons(ilgn)
! lgnatomidx(:,ilgn)=0
! coef=0
! read(*,'(8I5)')((atmidx(j),Coef(j)),j=1,nlgnatom_in_this_cnstrnt)
! do j = 1, nlgnatom_in_this_cnstrnt
! lgnatomidx(atmidx(j),ilgn)=Coef(j)
! end do
!end do
!read*,nsymmconstraint
!allocate(isymmconstraintpair(2,nsymmconstraint))
!do isymmconstraint = 1, nsymmconstraint
! read(*,'(2I5)')isymmconstraintpair(:,isymmconstraint)
!end do
!call hme(natom,ngrid,gridcrd,atomcrd,level,esp,qtot,mpinit,dpinit,qpinit,fitmethod,istrnt,strength,nlgn,lgnatomidx,lgncons,nsymmconstraint,isymmconstraintpair,40,ssres,crosscorr)
!deallocate(gridcrd,atomcrd,esp,mpinit,dpinit,qpinit,atmidx,coef,lgnatomidx,lgncons,isymmconstraintpair)
!
!end program hme_wrapper
subroutine hme(natom,ngrid,gridcrd,atomcrd,level,espinit,qtot,mpinit,dpinit,qpinit,fitmethod0, &
istrnt,strength,nlgn0,lgnatomidx0,lgncons0,nsymmconstraint,isymmconstraintpair,outid,ssres,crosscorr)
implicit none
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!! Restrained Electrostatic Potential Fitting of Atomic Charges !!
!! Written by Ye Mei, East China Normal University !!
!! During his stay at Dr. Bernie Brooks' group as Special Volunteer !!
!! Ver. 1.0, Feb. 12, 2014 !!
!! All rights Reserved !!
!! References: !!
!! 1. C. I. Bayly, P. Cieplak, W. Cornell, P. A. Kollman, !!
!! J. Phys. Chem., 1993, 97 (40), pp 10269–10280 !!
!! 2. J. Zeng, L. L., Duan, J. Z. H., Zhang, Y. Mei, !!
!! J. Comput. Chem., 34, 847-853 (2013) !!
!! Input variables !!
!! natom: Number of atoms (multipole sites) !!
!! ngrid: Number of ESP grid points !!
!! gridcrd(3,ngrid): Coordinates of ESP grid points !!
!! atomcrd(3,natoms): Coordinates of atoms (multipole sites) !!
!! level(3): if level(i)=1, multipole L=i-1 will be fitted. if 0, No. =2? No Kidding! !!
!! espinit(ngrid): Standard (QM) ESP on grids !!
!! qtot : total charge !!
!! mpinit(natom): Initial guess of monopole on sites, the fitted monopole on return !!
!! on output, the fitted monopole if level(1)=1 !!
!! dpinit(natom): Initial guess of dipole on sites, the fitted dipole on return !!
!! on output, the fitted dipole if level(2)=1 !!
!! qpinit(natom): Initial guess of quadrupole on sites, the fitted quadrupole on return !!
!! on output, the fitted quadrupole if level(3)=1 !!
!! fitmethod: Fitting method (RESP, DRESP, to be added) !!
!! istrnt: restraint function (0: Harmonic; 1: hyperbolic !!
!! strength: strength of restraint !!
!! nlgn0: number of Lagrange constraints (excluding the constraint of total charge) !!
!! lgnatomidx0: coefficients for the atoms in each Lagrange constraints !!
!! lgncons0: constraint for each Lagrange constraints !!
!! nsymmconstraint: number of symmetry constraints !!
!! isymmconstraintpair: site pair in symmetry !!
!! Output variables !!
!! outid: file pipe for output !!
!! ssres: residual sum of squares !!
!! crosscorr: correlation between fitted and target ESP !!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
interface
subroutine getlgn(natom,nlgn0,lgnatomidx0,lgncons0,nlgn,lgnatomidx,lgncons,qresidual,qbase,outid)
integer(kind=4), intent(in) :: natom
integer(kind=4), intent(in) :: nlgn0
integer(kind=4), intent(in) :: lgnatomidx0(natom,nlgn0)
real(kind=8), intent(in) :: lgncons0(nlgn0)
integer(kind=4), intent(out) :: nlgn
integer(kind=4), allocatable, intent(out) :: lgnatomidx(:,:)
real(kind=8), allocatable, intent(out) :: lgncons(:)
real(kind=8), intent(in) :: qresidual
real(kind=8), intent(in) :: qbase(natom)
integer(kind=4), intent(in) :: outid
end subroutine getlgn
subroutine gendismat(natom,ngrid,atomcrd,gridcrd,r_atm_grid,r_atm_grid_vec,ibigmem,outid)
integer(kind=4), intent(in) :: natom, ngrid
real(kind=8), intent(in) :: atomcrd(3,natom), gridcrd(3,ngrid)
real(kind=8), allocatable, intent(out) :: r_atm_grid(:,:), r_atm_grid_vec(:,:,:)
integer(kind=4), intent(out) :: ibigmem
integer(kind=4), intent(in) :: outid
end subroutine gendismat
end interface
integer(kind=4), intent(in) :: natom, ngrid
real(kind=8), intent(in) :: gridcrd(3,ngrid), atomcrd(3,natom), espinit(ngrid)
integer(kind=4), intent(in) :: level(3)
integer(kind=4), intent(in) :: qtot
real(kind=8), intent(in out) :: mpinit(natom), dpinit(3,natom), qpinit(5,natom)
character(len=*), intent(in) :: fitmethod0
integer(kind=4), intent(in) :: istrnt
real(kind=8), intent(in) :: strength
integer(kind=4), intent(in) :: nlgn0
integer(kind=4), intent(in) :: lgnatomidx0(natom,nlgn0)
real(kind=8), intent(in) :: lgncons0(nlgn0)
integer(kind=4), intent(in) :: nsymmconstraint
integer(kind=4), intent(in) :: isymmconstraintpair(2,nsymmconstraint)
integer(kind=4), intent(in) :: outid
real(kind=8), intent(out) :: ssres, crosscorr
character(len=10)::fitmethod
real(kind=8) :: espfitted(ngrid), espresidual(ngrid)
integer(kind=4) :: nlgn
integer(kind=4), allocatable :: lgnatomidx(:,:)
real(kind=8), allocatable :: lgncons(:)
integer(kind=4) :: isymmconstraintpairlocal(2,nsymmconstraint)
real(kind=8), allocatable :: q(:), qcopy(:)
real(kind=8), allocatable :: awork(:,:), bwork(:)
real(kind=8), allocatable :: acopy(:,:), bcopy(:)
real(kind=8) :: q0(natom), qbase(natom), qresidual
real(kind=8) :: monopoleout(natom), dipoleout(3,natom), quadrupoleout(5,natom)
integer(kind=4) :: neffecsize
real(kind=8), allocatable :: qwt(:)
real(kind=8), allocatable :: r_atm_grid(:,:), r_atm_grid_vec(:,:,:)
real(kind=8) :: sstot, ssreg
integer(kind=4) :: ibigmem
logical :: isconverged
integer(kind=4) :: istep
integer(kind=4) :: icycle
real(kind=8), parameter :: criterion=1.d-4
integer(kind=4),allocatable :: iskilledsite(:) ! some charge sites are not freely variable due to the existence of constraints
integer(kind=4) :: deallocatestatus
integer(kind=4) :: AllocateStatus
integer(kind=4) :: i, j, k, m, n
fitmethod=trim(adjustl(fitmethod0))
do i = 1, len(trim(fitmethod))
if(ichar(fitmethod(i:i))<=ichar('z').and.ichar(fitmethod(i:i))>=ichar('a'))then
fitmethod(i:i)=char(ichar(fitmethod(i:i))-ichar('a')+ichar('A'))
end if
end do
if(fitmethod/='RESP'.and.fitmethod/='DRESP')then
write(outid,'(A)')' Error: Only RESP and DRESP methods are allowed. '
write(outid,'(A)')' Error: Please check your input for "fitmethod" and run me again.'
stop
else
write(outid,'(A)')' Running '//trim(fitmethod)//' multipole fitting method'
end if
write(outid,'(A)')' Fit multipole moments using '//trim(fitmethod)//' fitting method'
write(outid,'(A,I7)')' Number of atoms ', natom
write(outid,'(A,I7)')' Number of grids ', ngrid
write(outid,*)
if(fitmethod/='RESP'.and.fitmethod/='DRESP')then
write(outid,'(A)')' Initial multipole moments'
write(outid,'(A)')' Monopoles'
write(outid,'(8F10.6)')(mpinit(i),i=1,natom)
write(outid,*)
write(outid,'(A)')' Dipoles'
do i = 1, natom
write(outid,'(3F10.6)')dpinit(:,i)
end do
write(outid,*)
write(outid,'(A)')' Quadrupoles'
do i = 1, natom
write(outid,'(5F10.6)')qpinit(:,i)
end do
write(outid,*)
! give a score to the initial multipole
call multipoleESP(ngrid,gridcrd,natom,atomcrd,mpinit,dpinit,qpinit,espfitted)
call qualitymetric(ngrid,espinit,espfitted,sstot,ssreg,ssres,crosscorr)
write(outid,'(A)')' Statistics of the initial multipoles:'
write(outid,'(A,F10.3)')' The initial sum of squares (ssvpot) ', sstot
write(outid,'(A,F10.3)')' The residual sum of squares (chipot) ', ssres
write(outid,'(A,F10.5)')' The std err of estimate (sqrt(chipot/N)) ', sqrt(ssres/ngrid)
write(outid,'(A,F10.5)')' ESP relative RMS (SQRT(chipot/ssvpot)) ', sqrt(ssres/sstot)
write(outid,'(A,F10.5)')' Correlation coefficient ', 1.d0-ssres/sstot
write(outid,'(A,F10.5)')' Cross correlation ', crosscorr
write(outid,*)
end if
! generate distance matrix between atoms and grids if memory allows
call gendismat(natom,ngrid,atomcrd,gridcrd,r_atm_grid,r_atm_grid_vec,ibigmem,outid)
! initialize
if(allocated(qwt))then
deallocate(qwt,stat=deallocatestatus)
if(deallocatestatus/=0)then
write(outid,'(A)')'deallocating qwt failed. deallocatestatus=',deallocatestatus
stop
end if
end if
allocate(qwt(natom),stat=AllocateStatus)
if(AllocateStatus/=0)then
write(outid,'(A)')'Failed when allocating qwt'
stop
end if
call initialize(natom,ngrid,gridcrd,atomcrd,espinit,espresidual,mpinit,qbase,qwt,r_atm_grid,ibigmem,fitmethod,outid)
write(outid,'(A,I8)')' Total charge of this molecule:', qtot
write(outid,'(A,F16.5)')' Sum of initial charge:', sum(mpinit)
qresidual=qtot-sum(qbase)
if(fitmethod=='DRESP')then
write(outid,'(A,F14.5)')' Residueal charge to fit:', qresidual
end if
write(outid,*)
if(allocated(iskilledsite))then
deallocate(iskilledsite,stat=deallocatestatus)
if(deallocatestatus/=0)then
write(outid,'(A)')'deallocating iskilledsite failed. deallocatestatus=',deallocatestatus
stop
end if
end if
allocate(iskilledsite(natom),stat=AllocateStatus)
if(AllocateStatus/=0)then
write(outid,'(A)')'Failed when allocating iskilledsite'
stop
end if
iskilledsite=0
espfitted = 0.d0
istep = 0
if(level(1)<=0)then
if(abs(qresidual)>1.D-4)then
write(outid,'(A)')'Error: Monopole fitting is disabled. However, I noticed that the sum of the initial charge does not equal to qtot.'
write(outid,'(A)')'Error: This is abnormal.'
stop
end if
else
istep = istep + 1
write(outid,'(A,I2,A)')' Step ',istep,': To fit monopoles'
! get info for Lagrange multiplier
call getlgn(natom,nlgn0,lgnatomidx0,lgncons0,nlgn,lgnatomidx,lgncons,qresidual,qbase,outid)
isymmconstraintpairlocal=isymmconstraintpair
write(outid,'(A,I5)')' Total number of symmetry constraints:',nsymmconstraint
do i = 1, nsymmconstraint
write(outid,'(2I3)')isymmconstraintpairlocal(:,i)
end do
write(outid,*)
if(allocated(awork))then
deallocate(awork,stat=deallocatestatus)
if(deallocatestatus/=0)then
write(outid,'(A)')'deallocating awork failed. deallocatestatus=',deallocatestatus
stop
end if
end if
allocate(awork(natom+nlgn,natom+nlgn),stat=AllocateStatus)
if(AllocateStatus/=0)then
write(outid,'(A)')'Failed when allocating awork'
stop
end if
if(allocated(bwork))then
deallocate(bwork,stat=deallocatestatus)
if(deallocatestatus/=0)then
write(outid,'(A)')'deallocating bwork failed. deallocatestatus=',deallocatestatus
stop
end if
end if
allocate(bwork(natom+nlgn),stat=AllocateStatus)
if(AllocateStatus/=0)then
write(outid,'(A)')'Failed when allocating bwork'
stop
end if
if(allocated(q))then
deallocate(q,stat=deallocatestatus)
if(deallocatestatus/=0)then
write(outid,'(A)')'deallocating q failed. deallocatestatus=',deallocatestatus
stop
end if
end if
allocate(q(natom),stat=AllocateStatus)
if(AllocateStatus/=0)then
write(outid,'(A)')'Failed when allocating q'
stop
end if
neffecsize=natom+nlgn-nsymmconstraint
if(allocated(acopy))then
deallocate(acopy,stat=deallocatestatus)
if(deallocatestatus/=0)then
write(outid,'(A)')'deallocating acopy failed. deallocatestatus=',deallocatestatus
stop
end if
end if
allocate(acopy(neffecsize,neffecsize),stat=AllocateStatus)
if(AllocateStatus/=0)then
write(outid,'(A)')'Failed when allocating acopy'
stop
end if
if(allocated(bcopy))then
deallocate(bcopy,stat=deallocatestatus)
if(deallocatestatus/=0)then
write(outid,'(A)')'deallocating bcopy failed. deallocatestatus=',deallocatestatus
stop
end if
end if
allocate(bcopy(neffecsize),stat=AllocateStatus)
if(AllocateStatus/=0)then
write(outid,'(A)')'Failed when allocating bcopy'
stop
end if
if(allocated(qcopy))then
deallocate(qcopy,stat=deallocatestatus)
if(deallocatestatus/=0)then
write(outid,'(A)')'deallocating bcopy failed. deallocatestatus=',deallocatestatus
stop
end if
end if
allocate(qcopy(neffecsize),stat=AllocateStatus)
if(AllocateStatus/=0)then
write(outid,'(A)')'Failed when allocating qcopy'
stop
end if
monopoleout=0.d0
dipoleout=0.d0
quadrupoleout=0.d0
q0=mpinit-qbase
isconverged=.false.
icycle=0
do while(.not.isconverged)
icycle = icycle + 1
write(outid,'(A,I3)')' Cycle:', icycle
! generate least-square matrix
call bldmatmonopole(natom,nlgn,atomcrd,ngrid,gridcrd,espresidual,r_atm_grid,ibigmem,awork,bwork)
! write(outid,'(<natom>F10.5,2X,F10.4)')((awork(i,j),j=1,natom),bwork(i),i=1,natom)
! write(outid,*)
! modify working matrices with charge restraints
call getrestraints(natom,nlgn,awork,bwork,istrnt,strength,q0,qwt)
! write(outid,'(<natom>F10.5,2X,F10.4)')((awork(i,j),j=1,natom),bwork(i),i=1,natom)
! write(outid,*)
! add Lagrange condition to least-square matrix
call appendlgn2mat(natom,nlgn,awork,bwork,lgnatomidx,lgncons)
! write(outid,'(<natom+nlgn>F10.5,2X,F10.4)')((awork(i,j),j=1,natom+nlgn),bwork(i),i=1,natom+nlgn)
! write(outid,*)
if(nsymmconstraint>0)then
call sortsymmconstraint(nsymmconstraint,isymmconstraintpairlocal)
call applysymmetry(natom,nlgn,awork,bwork,nsymmconstraint,isymmconstraintpair,iskilledsite)
end if
! write(outid,'(<natom+nlgn>F10.5,2X,F10.4)')((awork(i,j),j=1,natom+nlgn),bwork(i),i=1,natom+nlgn)
! write(outid,*)
! backups for working matrices. Essential for iterations
m = 0
do i = 1, natom+nlgn
if(i<=natom)then
if(iskilledsite(i)/=0)cycle
end if
m=m+1
n=0
do j = 1, natom+nlgn
if(j<=natom)then
if(iskilledsite(j)/=0)cycle
end if
n=n+1
acopy(m,n)=awork(i,j)
end do
bcopy(m)=bwork(i)
end do
! write(outid,'(<neffecsize>F10.5,2X,F10.4)')((acopy(i,j),j=1,neffecsize),bcopy(i),i=1,neffecsize)
! write(outid,*)
do i = 1, neffecsize
if(abs(acopy(i,i))<1.D-10)acopy(i,i)=1.D-10
end do
! if(allocated(awork))deallocate(awork)
! allocate(awork(neffecsize,neffecsize))
! if(allocated(bwork))deallocate(bwork)
! allocate(bwork(neffecsize))
! awork=acopy
! bwork=bcopy
! SVD fit of charges
call LESVD(neffecsize,neffecsize,acopy,bcopy,qcopy,outid)
! call LELUD(neffecsize,neffecsize,acopy,bcopy,qcopy,outid)
! print*,'killed charge:'
! do i = 1, natom
! if(iskilledsite(i)>0)print*,i
! end do
do i = 1, min(natom,neffecsize)
q(i) = qcopy(i)
end do
call rebuildchargearray(natom,q,nsymmconstraint,isymmconstraintpairlocal,iskilledsite)
forall(i=1:natom)
monopoleout(i)=qbase(i)+q(i)
end forall
write(outid,'(A)')' Fitted monopoles'
write(outid,'(8F10.6)')(q(i),i=1,natom)
write(outid,*)
if(fitmethod=='DRESP')then
write(outid,'(A)')' Fitted + initial monopoles'
write(outid,'(8F10.6)')(monopoleout(i),i=1,natom)
write(outid,*)
end if
! check if converged
if(istrnt==1)then
call checkconvergence(natom,q,q0,criterion,isconverged)
forall(i=1:natom)
q0(i)=q(i)
end forall
else
isconverged=.true.
endif
end do
!
! call multipoleESP(ngrid,gridcrd,natom,atomcrd,monopoleout,dipoleout,quadrupoleout,espfitted)
! call qualitymetric(ngrid,espinit,espfitted,sstot,ssreg,ssres,crosscorr)
! write(outid,'(A)')' Statistics of the fitted multipoles (up to L=0):'
! write(outid,'(A,F10.3)')' The initial sum of squares (ssvpot) ', sstot
! write(outid,'(A,F10.3)')' The residual sum of squares (chipot) ', ssres
! write(outid,'(A,F10.5)')' The std err of estimate (sqrt(chipot/N)) ', sqrt(ssres/ngrid)
! write(outid,'(A,F10.5)')' ESP relative RMS (SQRT(chipot/ssvpot)) ', sqrt(ssres/sstot)
! write(outid,'(A,F10.5)')' Correlation coefficient (1-chipot/ssvpot)', 1.d0-ssres/sstot
! write(outid,'(A,F10.5)')' Cross correlation ', crosscorr
! write(outid,*)
end if
if(level(2)>0)then
istep = istep + 1
write(outid,'(A,I2,A)')' Step ',istep,': To fit dipoles'
if(allocated(awork))deallocate(awork)
allocate(awork(3*natom,3*natom))
if(allocated(bwork))deallocate(bwork)
allocate(bwork(3*natom))
espresidual=espinit-espfitted
call bldmatdipole(natom,atomcrd,ngrid,gridcrd,espresidual,r_atm_grid,r_atm_grid_vec,ibigmem,awork,bwork)
call LESVD(3*natom,3*natom,awork,bwork,dipoleout,outid)
write(outid,'(A)')' Fitted dipoles'
do i = 1, natom
write(outid,'(3F10.6)')dipoleout(:,i)
end do
write(outid,*)
call multipoleESP(ngrid,gridcrd,natom,atomcrd,monopoleout,dipoleout,quadrupoleout,espfitted)
call qualitymetric(ngrid,espinit,espfitted,sstot,ssreg,ssres,crosscorr)
write(outid,'(A)')' Statistics of the fitted multipoles (up to L=1):'
write(outid,'(A,F10.3)')' The initial sum of squares (ssvpot) ', sstot
write(outid,'(A,F10.3)')' The residual sum of squares (chipot) ', ssres
write(outid,'(A,F10.5)')' The std err of estimate (sqrt(chipot/N)) ', sqrt(ssres/ngrid)
write(outid,'(A,F10.5)')' ESP relative RMS (SQRT(chipot/ssvpot)) ', sqrt(ssres/sstot)
write(outid,'(A,F10.5)')' Correlation coefficient (1-chipot/ssvpot)', 1.d0-ssres/sstot
write(outid,'(A,F10.5)')' Cross correlation ', crosscorr
write(outid,*)
end if
if(level(3)>0)then
istep = istep + 1
write(outid,'(A,I2,A)')' Step ',istep,': To fit quadrupoles'
if(allocated(awork))deallocate(awork)
allocate(awork(5*natom,5*natom))
if(allocated(bwork))deallocate(bwork)
allocate(bwork(5*natom))
espresidual=espinit-espfitted
call bldmatquadrupole(natom,atomcrd,ngrid,gridcrd,espresidual,r_atm_grid,r_atm_grid_vec,ibigmem,awork,bwork)
call LESVD(5*natom,5*natom,awork,bwork,quadrupoleout,outid)
write(outid,'(A)')' Fitted quadrupoles'
do i = 1, natom
write(outid,'(5F10.6)')quadrupoleout(:,i)
end do
write(outid,*)
call multipoleESP(ngrid,gridcrd,natom,atomcrd,monopoleout,dipoleout,quadrupoleout,espfitted)
call qualitymetric(ngrid,espinit,espfitted,sstot,ssreg,ssres,crosscorr)
write(outid,'(A)')' Statistics of the fitted multipoles (up to L=2):'
write(outid,'(A,F10.3)')' The initial sum of squares (ssvpot) ', sstot
write(outid,'(A,F10.3)')' The residual sum of squares (chipot) ', ssres
write(outid,'(A,F10.5)')' The std err of estimate (sqrt(chipot/N)) ', sqrt(ssres/ngrid)
write(outid,'(A,F10.5)')' ESP relative RMS (SQRT(chipot/ssvpot)) ', sqrt(ssres/sstot)
write(outid,'(A,F10.5)')' Correlation coefficient (1-chipot/ssvpot)', 1.d0-ssres/sstot
write(outid,'(A,F10.5)')' Cross correlation ', crosscorr
write(outid,*)
end if
write(outid,'(A)')' Summary of multipoles'
write(outid,'(A)')' Fitted monopoles'
write(outid,'(8F10.6)')(monopoleout(i),i=1,natom)
write(outid,*)
write(outid,'(A)')' Fitted dipoles'
do i = 1, natom
write(outid,'(3F10.6)')dipoleout(:,i)
end do
write(outid,*)
write(outid,'(A)')' Fitted quadrupoles'
do i = 1, natom
write(outid,'(5F10.6)')quadrupoleout(:,i)
end do
write(outid,*)
if(level(1)==1)then
mpinit=monopoleout
else if (level(2)==1) then
dpinit=dipoleout
else if (level(3)==1)then
qpinit=quadrupoleout
end if
if(allocated(awork))deallocate(awork)
if(allocated(bwork))deallocate(bwork)
if(allocated(acopy))deallocate(acopy,STAT=deallocatestatus)
if(allocated(bcopy))deallocate(bcopy)
if(allocated(bcopy))deallocate(qcopy)
if(allocated(q))deallocate(q)
if(allocated(qwt))deallocate(qwt)
if(allocated(lgnatomidx))deallocate(lgnatomidx)
if(allocated(lgncons))deallocate(lgncons)
if(allocated(r_atm_grid))deallocate(r_atm_grid)
if(allocated(r_atm_grid_vec))deallocate(r_atm_grid_vec)
end subroutine hme
subroutine bldmatquadrupole(natom,atomcrd,ngrid,gridcrd,esp,rnorm,r,ibigmem,awork,bwork)
implicit none
integer(kind=4), intent(in) :: natom, ngrid
real(kind=8), intent(in) :: atomcrd(3,natom), gridcrd(3,ngrid)
real(kind=8), intent(in) :: esp(ngrid)
real(kind=8), intent(in) :: rnorm(ngrid,natom),r(3,ngrid,natom)
integer(kind=4), intent(in) :: ibigmem
real(kind=8), intent(out) :: awork(5*natom,5*natom), bwork(5*natom)
integer(kind=4) :: idx1, idx2
integer(kind=4) :: i1, j1, k1
integer(kind=4) :: i2, j2, k2
real(kind=8) :: rji1, rji2
real(kind=8) :: rji1_v(3), rji2_v(3)
real(kind=8) :: factor
integer(kind=4) :: i, j, k, m, n
if(ibigmem<4)then
!do i = 1, natom
! do j = 1, ngrid
! do k = 1, 3
! r(k,j,i)=gridcrd(k,j)-atomcrd(k,i)
! end do
! rnorm(j,i)=sqrt(r(1,j,i)**2+r(2,j,i)**2+r(3,j,i)**2)
! end do
!end do
awork=0.d0
bwork=0.d0
do i1 = 1, natom
idx1=5*(i1-1)+1
do j = 1, ngrid
do k = 1, 3
rji1_v(k)=gridcrd(k,j)-atomcrd(k,i1)
end do
rji1=sqrt(rji1_v(1)**2+rji1_v(2)**2+rji1_v(3)**2)
factor=(3.d0*rji1_v(3)**2-rji1**2)/(2.d0*rji1**5)
do i2 = 1, natom
do k = 1, 3
rji2_v(k)=gridcrd(k,j)-atomcrd(k,i2)
end do
rji2=sqrt(rji2_v(1)**2+rji2_v(2)**2+rji2_v(3)**2)
idx2=5*(i2-1)+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*(3.d0*rji2_v(3)**2-rji2**2)/(2.d0*rji2**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(1)*rji2_v(3)/rji2**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(2)*rji2_v(3)/rji2**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*0.5d0*(rji2_v(1)**2-rji2_v(2)**2)/(2.d0*rji2**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(1)*rji2_v(2)/rji2**5
end do
bwork(idx1)=bwork(idx1)+esp(j)*factor
end do
idx1=idx1+1
do j = 1, ngrid
do k = 1, 3
rji1_v(k)=gridcrd(k,j)-atomcrd(k,i1)
end do
rji1=sqrt(rji1_v(1)**2+rji1_v(2)**2+rji1_v(3)**2)
factor=sqrt(3.d0)*rji1_v(1)*rji1_v(3)/rji1**5
do i2 = 1, natom
do k = 1, 3
rji2_v(k)=gridcrd(k,j)-atomcrd(k,i2)
end do
rji2=sqrt(rji2_v(1)**2+rji2_v(2)**2+rji2_v(3)**2)
idx2=5*(i2-1)+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*(3.d0*rji2_v(3)**2-rji2**2)/(2.d0*rji2**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(1)*rji2_v(3)/rji2**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(2)*rji2_v(3)/rji2**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*0.5d0*(rji2_v(1)**2-rji2_v(2)**2)/(2.d0*rji2**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(1)*rji2_v(2)/rji2**5
end do
bwork(idx1)=bwork(idx1)+esp(j)*factor
end do
idx1=idx1+1
do j = 1, ngrid
do k = 1, 3
rji1_v(k)=gridcrd(k,j)-atomcrd(k,i1)
end do
rji1=sqrt(rji1_v(1)**2+rji1_v(2)**2+rji1_v(3)**2)
factor=sqrt(3.d0)*rji1_v(2)*rji1_v(3)/rji1**5
do i2 = 1, natom
do k = 1, 3
rji2_v(k)=gridcrd(k,j)-atomcrd(k,i2)
end do
rji2=sqrt(rji2_v(1)**2+rji2_v(2)**2+rji2_v(3)**2)
idx2=5*(i2-1)+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*(3.d0*rji2_v(3)**2-rji2**2)/(2.d0*rji2**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(1)*rji2_v(3)/rji2**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(2)*rji2_v(3)/rji2**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*0.5d0*(rji2_v(1)**2-rji2_v(2)**2)/(2.d0*rji2**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(1)*rji2_v(2)/rji2**5
end do
bwork(idx1)=bwork(idx1)+esp(j)*factor
end do
idx1=idx1+1
do j = 1, ngrid
do k = 1, 3
rji1_v(k)=gridcrd(k,j)-atomcrd(k,i1)
end do
rji1=sqrt(rji1_v(1)**2+rji1_v(2)**2+rji1_v(3)**2)
factor=sqrt(3.d0)*0.5d0*(rji1_v(1)**2-rji1_v(2)**2)/(2.d0*rji1**5)
do i2 = 1, natom
do k = 1, 3
rji2_v(k)=gridcrd(k,j)-atomcrd(k,i2)
end do
rji2=sqrt(rji2_v(1)**2+rji2_v(2)**2+rji2_v(3)**2)
idx2=5*(i2-1)+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*(3.d0*rji2_v(3)**2-rji2**2)/(2.d0*rji2**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(1)*rji2_v(3)/rji2**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(2)*rji2_v(3)/rji2**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*0.5d0*(rji2_v(1)**2-rji2_v(2)**2)/(2.d0*rji2**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(1)*rji2_v(2)/rji2**5
end do
bwork(idx1)=bwork(idx1)+esp(j)*factor
end do
idx1=idx1+1
do j = 1, ngrid
do k = 1, 3
rji1_v(k)=gridcrd(k,j)-atomcrd(k,i1)
end do
rji1=sqrt(rji1_v(1)**2+rji1_v(2)**2+rji1_v(3)**2)
factor=sqrt(3.d0)*rji1_v(1)*rji1_v(2)/rji1**5
do i2 = 1, natom
do k = 1, 3
rji2_v(k)=gridcrd(k,j)-atomcrd(k,i2)
end do
rji2=sqrt(rji2_v(1)**2+rji2_v(2)**2+rji2_v(3)**2)
idx2=5*(i2-1)+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*(3.d0*rji2_v(3)**2-rji2**2)/(2.d0*rji2**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(1)*rji2_v(3)/rji2**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(2)*rji2_v(3)/rji2**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*0.5d0*(rji2_v(1)**2-rji2_v(2)**2)/(2.d0*rji2**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*rji2_v(1)*rji2_v(2)/rji2**5
end do
bwork(idx1)=bwork(idx1)+esp(j)*factor
end do
end do
else
awork=0.d0
bwork=0.d0
do i1 = 1, natom
idx1=5*(i1-1)+1
do j = 1, ngrid
factor=(3.d0*r(3,j,i1)**2-rnorm(j,i1)**2)/(2.d0*rnorm(j,i1)**5)
do i2 = 1, natom
idx2=5*(i2-1)+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*(3.d0*r(3,j,i2)**2-rnorm(j,i2)**2)/(2.d0*rnorm(j,i2)**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(1,j,i2)*r(3,j,i2)/rnorm(j,i2)**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(2,j,i2)*r(3,j,i2)/rnorm(j,i2)**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*0.5d0*(r(1,j,i2)**2-r(2,j,i2)**2)/(2.d0*rnorm(j,i2)**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(1,j,i2)*r(2,j,i2)/rnorm(j,i2)**5
end do
bwork(idx1)=bwork(idx1)+esp(j)*factor
end do
idx1=idx1+1
do j = 1, ngrid
factor=sqrt(3.d0)*r(1,j,i1)*r(3,j,i1)/rnorm(j,i1)**5
do i2 = 1, natom
idx2=5*(i2-1)+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*(3.d0*r(3,j,i2)**2-rnorm(j,i2)**2)/(2.d0*rnorm(j,i2)**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(1,j,i2)*r(3,j,i2)/rnorm(j,i2)**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(2,j,i2)*r(3,j,i2)/rnorm(j,i2)**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*0.5d0*(r(1,j,i2)**2-r(2,j,i2)**2)/(2.d0*rnorm(j,i2)**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(1,j,i2)*r(2,j,i2)/rnorm(j,i2)**5
end do
bwork(idx1)=bwork(idx1)+esp(j)*factor
end do
idx1=idx1+1
do j = 1, ngrid
factor=sqrt(3.d0)*r(2,j,i1)*r(3,j,i1)/rnorm(j,i1)**5
do i2 = 1, natom
idx2=5*(i2-1)+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*(3.d0*r(3,j,i2)**2-rnorm(j,i2)**2)/(2.d0*rnorm(j,i2)**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(1,j,i2)*r(3,j,i2)/rnorm(j,i2)**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(2,j,i2)*r(3,j,i2)/rnorm(j,i2)**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*0.5d0*(r(1,j,i2)**2-r(2,j,i2)**2)/(2.d0*rnorm(j,i2)**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(1,j,i2)*r(2,j,i2)/rnorm(j,i2)**5
end do
bwork(idx1)=bwork(idx1)+esp(j)*factor
end do
idx1=idx1+1
do j = 1, ngrid
factor=sqrt(3.d0)*0.5d0*(r(1,j,i1)**2-r(2,j,i1)**2)/(2.d0*rnorm(j,i1)**5)
do i2 = 1, natom
idx2=5*(i2-1)+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*(3.d0*r(3,j,i2)**2-rnorm(j,i2)**2)/(2.d0*rnorm(j,i2)**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(1,j,i2)*r(3,j,i2)/rnorm(j,i2)**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(2,j,i2)*r(3,j,i2)/rnorm(j,i2)**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*0.5d0*(r(1,j,i2)**2-r(2,j,i2)**2)/(2.d0*rnorm(j,i2)**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(1,j,i2)*r(2,j,i2)/rnorm(j,i2)**5
end do
bwork(idx1)=bwork(idx1)+esp(j)*factor
end do
idx1=idx1+1
do j = 1, ngrid
factor=sqrt(3.d0)*r(1,j,i1)*r(2,j,i1)/rnorm(j,i1)**5
do i2 = 1, natom
idx2=5*(i2-1)+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*(3.d0*r(3,j,i2)**2-rnorm(j,i2)**2)/(2.d0*rnorm(j,i2)**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(1,j,i2)*r(3,j,i2)/rnorm(j,i2)**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(2,j,i2)*r(3,j,i2)/rnorm(j,i2)**5
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*0.5d0*(r(1,j,i2)**2-r(2,j,i2)**2)/(2.d0*rnorm(j,i2)**5)
idx2=idx2+1
awork(idx1,idx2)=awork(idx1,idx2)+factor*sqrt(3.d0)*r(1,j,i2)*r(2,j,i2)/rnorm(j,i2)**5
end do
bwork(idx1)=bwork(idx1)+esp(j)*factor
end do
end do
end if
end subroutine bldmatquadrupole
subroutine bldmatdipole(natom,atomcrd,ngrid,gridcrd,esp,rnorm,r,ibigmem,awork,bwork)
implicit none
integer(kind=4), intent(in) :: natom, ngrid
real(kind=8), intent(in) :: atomcrd(3,natom), gridcrd(3,ngrid)
real(kind=8), intent(in) :: esp(ngrid)
real(kind=8), intent(in) :: rnorm(ngrid,natom),r(3,ngrid,natom)
integer(kind=4), intent(in) :: ibigmem
real(kind=8), intent(out) :: awork(3*natom,3*natom), bwork(3*natom)
integer(kind=4) :: idx1, idx2
integer(kind=4) :: i1, j1, k1
integer(kind=4) :: i2, j2, k2
real(kind=8) :: rji1, rji2
integer(kind=4) :: i, j, k, m, n
if(ibigmem==4)then
awork=0.d0
bwork=0.d0
idx1=0
do i1 = 1, natom
do k1 = 1, 3
idx1=idx1+1
idx2=0
do i2 = 1, natom
do k2 = 1, 3
idx2=idx2+1
do j = 1, ngrid
awork(idx1,idx2)=awork(idx1,idx2)+r(k1,j,i1)*r(k2,j,i2)/(rnorm(j,i1)*rnorm(j,i2))**3
end do
end do
end do
do j = 1, ngrid
bwork(idx1)=bwork(idx1)+esp(j)*r(k1,j,i1)/rnorm(j,i1)**3
end do
end do
end do
else
awork=0.d0
bwork=0.d0
idx1=0
do i1 = 1, natom
do k1 = 1, 3
idx1=idx1+1
idx2=0
do i2 = 1, natom
do k2 = 1, 3
idx2=idx2+1
do j = 1, ngrid
rji1=sqrt((gridcrd(1,j)-atomcrd(1,i1))**2+(gridcrd(2,j)-atomcrd(2,i1))**2+(gridcrd(3,j)-atomcrd(3,i1))**2)
rji2=sqrt((gridcrd(1,j)-atomcrd(1,i2))**2+(gridcrd(2,j)-atomcrd(2,i2))**2+(gridcrd(3,j)-atomcrd(3,i2))**2)
awork(idx1,idx2)=awork(idx1,idx2)+(gridcrd(k1,j)-atomcrd(k1,i1))*(gridcrd(k2,j)-atomcrd(k2,i2))/(rji1*rji2)**3
end do
end do
end do
do j = 1, ngrid
rji1=sqrt((gridcrd(1,j)-atomcrd(1,i1))**2+(gridcrd(2,j)-atomcrd(2,i1))**2+(gridcrd(3,j)-atomcrd(3,i1))**2)
bwork(idx1)=bwork(idx1)+esp(j)*(gridcrd(k1,j)-atomcrd(k1,i1))/rji1**3
end do
end do
end do
end if
end subroutine bldmatdipole
subroutine qualitymetric(ngrid,esp0,esp1,sstot,ssreg,ssres,crosscorr)
implicit none
integer(kind=4), intent(in) :: ngrid
real(kind=8), intent(in) :: esp0(ngrid), esp1(ngrid)
real(kind=8), intent(out) :: sstot, ssreg, ssres, crosscorr
real(kind=8) :: esp0mean, esp1mean
real(kind=8) :: esp0width, esp1width
integer(kind=4) :: i, j, k, m, n
esp0mean = sum(esp0)/max(1,ngrid)
esp1mean = sum(esp1)/max(1,ngrid)
sstot = 0.d0
ssreg = 0.d0
ssres = 0.d0
do i = 1, ngrid
sstot = sstot + (esp0(i)-esp0mean)**2
ssreg = ssreg + (esp1(i)-esp0mean)**2
ssres = ssres + (esp1(i)-esp0(i))**2
end do
crosscorr=0.d0
esp0width=0.d0
esp1width=0.d0
do i = 1, ngrid
crosscorr=crosscorr+(esp0(i)-esp0mean)*(esp1(i)-esp1mean)
esp0width=esp0width+(esp0(i)-esp0mean)**2
esp1width=esp1width+(esp1(i)-esp1mean)**2
end do
crosscorr = crosscorr/sqrt(esp0width*esp1width+1.D-10)
end subroutine qualitymetric
subroutine multipoleESP(ngrid,gridcrd,natom,atomcrd,monopole,dipole,quadrupole,esp)
implicit none
integer(kind=4), intent(in) :: ngrid, natom
real(kind=8), intent(in) :: gridcrd(3,ngrid), atomcrd(3,natom), monopole(natom), dipole(3,natom), quadrupole(5,natom)
real(kind=8), intent(out) :: esp(ngrid)
real(kind=8) :: thisdipole(3), thisquadrupole(5)
real(kind=8) :: dr(3), rnorm
real(kind=8) :: dipoleesp, quadrupoleesp
integer(kind=4) :: i, j, k, m, n
esp=0.d0
do i =1, ngrid
do j = 1, natom
do k =1, 3
dr(k) = gridcrd(k,i)-atomcrd(k,j)
end do
rnorm=sqrt(dot_product(dr,dr))
thisdipole(:)=dipole(:,j)
thisquadrupole(:)=quadrupole(:,j)
esp(i)=esp(i)+monopole(j)/rnorm
esp(i)=esp(i)+dipoleesp(thisdipole,dr)
esp(i)=esp(i)+quadrupoleesp(thisquadrupole,dr)
end do
end do
end subroutine multipoleESP
function dipoleesp(dipole,r)
implicit none
real(kind=8), intent(in) :: dipole(3),r(3)
real(kind=8) :: dipoleesp
real(kind=8) :: rnorm
rnorm=sqrt(dot_product(r,r))
dipoleesp=dot_product(dipole,r)/rnorm**3
end function dipoleesp
function quadrupoleesp(quadrupole,r)
implicit none
real(kind=8),intent(in) :: quadrupole(5),r(3)
real(kind=8) :: quadrupoleesp
real(kind=8) :: rnorm
rnorm=sqrt(dot_product(r,r))
quadrupoleesp=(quadrupole(1)*(1.5d0*r(3)**2-0.5d0*rnorm**2) &
+quadrupole(2)*sqrt(3.d0)*r(1)*r(3) &
+quadrupole(3)*sqrt(3.d0)*r(2)*r(3) &
+quadrupole(4)*0.5d0*sqrt(3.d0)*(r(1)**2-r(2)**2) &
+quadrupole(5)*sqrt(3.d0)*r(1)*r(2) &
)/rnorm**5
end function quadrupoleesp
subroutine checkconvergence(natom,q,q0,criterion,isconverged)
implicit none
integer(kind=4), intent(in) :: natom