-
Notifications
You must be signed in to change notification settings - Fork 6
/
gcp.f90
4832 lines (4465 loc) · 154 KB
/
gcp.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
! This file is part of mctc-gcp.
! SPDX-Identifier: LGPL-3.0-or-later
!
! mctc-gcp is free software: you can redistribute it and/or modify it under
! the terms of the Lesser GNU General Public License as published by
! the Free Software Foundation, either version 3 of the License, or
! (at your option) any later version.
!
! mctc-gcp is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! Lesser GNU General Public License for more details.
!
! You should have received a copy of the Lesser GNU General Public License
! along with mctc-gcp. If not, see <https://www.gnu.org/licenses/>.
module gcp
use gcp_version, only : get_gcp_version
contains
!subroutine for interfacing with orca/turbomole/crystal
subroutine gcp_call(n,xyz,lat,iz,gcp_e,gcp_g,gcp_glat,dograd,dohess,pbc,method,echo,parfile)
implicit none
integer n !number of atoms
real*8 xyz(3,n) !xyzcoordinates
real*8 lat(3,3) !lattice matrix
integer iz(n) !element number
real*8 gcp_g(3,n) !xyz gradient
real*8 gcp_glat(3,3)
logical dograd,dohess,echo !calculate gradient, calculate hessian, print informations
logical pbc !periodic system
character*20 method !method string
integer max_elem,max_para
integer nn,nb,i,zz,j
! max elements possible
parameter (max_elem=94)
! actually parameterized
parameter (max_para=36)
real*8 autoang
parameter (autoang =0.5291772083d0)
! restore cardinal numbers of elements
integer iz2(n)
! frozen coordinate
integer ifrez(n)
! #electrons of elements
integer nel(max_elem)
! # basis functions
integer nbas(max_elem)
real*8 emiss(max_elem)
real*8 tmp (max_elem)
! parameters
real*8 p(4)
real*8 xva(n),xvb(n)
real*8 xx(10),gcp_e,r
real*8 t0,t1
logical ex,warn,lib,parm,test,verb
character*20 ctmp
character*80 arg(10)
character*80 atmp,ftmp,dtmp
character*80 basname
!v2.0 stuff
logical tmen,base,damp,stress,parfile,srb,newgcp
character*80 args(90), comment, line
real*8 er,el,ebas,gbas(3,n)
integer nargs, ios
real*8 dmp_scal,dmp_exp,rscal,qscal
!default no srvb
srb=.false.
!damping params
dmp_scal=4.0d0 ! 4
dmp_exp=6.0d0 ! 6
!default
base=.false.
damp=.false.
!************************
!* modify input string *
!* 1. use lower case *
!* 2. delete hyphends *
!************************
do while (index(method,'-').ne.0)
i=scan(method,'-')
ctmp=trim(method(:(i-1)))//trim(method((i+1):))
method=trim(ctmp)
enddo
select case(method)
case('hf3c')
!JGB special method combination 1: HF-3c
if(echo) write(*,*) ' Perform special gCP-SRB correction for HF-3c '
base=.true.
case('pbeh3c','pbeh3cmod')
!JGB special method combination 2: PBEh-3c
if(echo) write(*,*) ' Perform special gCP(damped) correction for PBEh-3c'
damp=.true.
case('hse3c')
!JGB special method combination 3: HSE-3c
if(echo) write(*,*) ' Perform special gCP(damped) correction for HSE-3c'
damp=.true.
!case('b3pbe3c')
!JGB special method combination 3: B3PBE-3c
!if(echo) write(*,*) ' Perform special gCP(damped) correction for B3PBE-3c'
!damp=.true.
case('b973c')
!JGB special method combination 4: B97-3c
!JGB no gCP needed, switch to SRB correction
if(echo) write(*,*) ' Perform special SRB correction for B97-3c '
srb=.true.
rscal=10.00d0
qscal=0.08d0
case('r2scan3c','def2mtzvpp','mtzvpp')
damp=.true.
case default
end select
!skip for B97-3c
if(.not.srb) then
!get parameters
!****************************************************************************
!* read parameter file *
!* format: *
!* basis sigma eta alpha beta *
!* To read in nbas and emiss, generate param file with *
!* with '-parfile', edit copy it to .gcppar: *
!* A '#' and the beginning will tell the program to read a different format *
!****************************************************************************
if(method.eq.''.or.method.eq.'file')then
!* Fortran 2003 *
call get_environment_variable('HOSTNAME',ftmp)
call get_environment_variable('HOME',atmp)
!* uncomment the following line *
!* if your compiler wont support above features *
!
! call system('hostname > .tmpx')
! open(unit=43,file='.tmpx')
! read(43,'(a)')ftmp
! close(43,status='delete')
! call system('echo $HOME > .tmpx')
! open(unit=43,file='.tmpx')
! read(43,'(a)')atmp
! close(43,status='delete')
write(dtmp,'(a,''/.gcppar.'',a)') trim(atmp), trim(ftmp)
inquire(file=trim(dtmp),exist=ex)
if(ex)then
if(echo)write(*,*) 'reading ',trim(dtmp)
open(unit=42,file=dtmp)
read(42,'(a)',end=9)ftmp
call charXsplit(ftmp,basname,1)
call lower_case(basname)
call readl(ftmp,xx,nn)
p(1:4)=xx(1:4)
if(adjustl(trim(basname)).eq.'#') then
if(echo) then
write(*,*) ' found extended param file ! '
endif
parm=.true.
read(42,'(x,F8.4)',end=9) p(1:4)
do i=1,36
read(42,'(x,I3,x,F9.6)',end=9) nbas(i), emiss(i)
enddo
!additional read for new param file
elseif(adjustl(trim(basname)).eq.'##') then
if(echo) then
write(*,*) ' found extended param file v2! '
endif
call readl(ftmp,xx,nn)
parm=.true.
damp=.true.
dmp_scal=xx(5)
dmp_exp=xx(6)
do i=1,36
read(42,'(x,I3,x,F9.6)',end=9) nbas(i), emiss(i)
enddo
endif
9 close(42)
if(.not.parm) then
if(echo) write(*,*) 'loading ',trim(basname),' params'
call setparam(emiss,nbas,p,basname)
endif
else
write(*,'(a)') ''
write(*,'(2a)') 'found not param file in ',trim(dtmp)
error stop '* ERROR *'
endif
method=trim(basname)
!*******************
!* built-in method *
!*******************
else
call setparam(emiss,nbas,p,method)
endif
!open(unit=42,file='~/.gcppar')
!read(42,*) p(1:4),dmp_scal,dmp_exp
!do i=1,36
!read(42,*) emiss(i),tmp(i)
!enddo
!close(42)
!**********************************************************
!* check for missing parameters for current elements *
!* elements >36 will be treated as their lower homologues *
!* eg: Ag -> Cu, In -> Ga *
!* prepare virtual bf array xva,xvb *
!**********************************************************
! backup original
iz2=iz
! re-arrange elements
do i=1,n
zz=iz(i)
select case (zz)
case(37:54)
iz(i)=iz(i)-18
write(*,'(4a)') '** NOTE ** -> element ',esym(zz),' will be treated as ', esym(iz(i))
warn=.true.
case(55:57)
iz(i)=iz(i)-18*2
write(*,'(4a)') '** NOTE ** -> element ',esym(zz),' will be treated as ', esym(iz(i))
warn=.true.
case(58:71,90:94)
iz(i)=21
write(*,'(4a)') '** NOTE ** -> element ',esym(zz),' will be treated as ', esym(iz(i))
warn=.true.
case(72:89)
iz(i)=iz(i)-2*18-14
write(*,'(4a)') '** NOTE ** -> element ',esym(zz),' will be treated as ', esym(iz(i))
warn=.true.
end select
enddo
! set electrons
do i=1,maxval(iz)
nel(i)=i
enddo
!* handle ECP basis set
if(index(method,'vmb').ne.0) then
write(*,'(a)') 'VMB basis has ECPs. Adjusting electrons'
do i=5,10
nel(i)=nel(i)-2
enddo
do i=11,18
nel(i)=nel(i)-10
enddo
endif
! count bf
nb=0
do nn=1,n
nb=nb+nbas(iz(nn))
enddo
if(nb.lt.1) error stop 'Nbf setup gone wrong'
! set virtuals and look for missing parameters
xva=0
do i=1,n
if(method.eq.'def2mtzvpp'.or.method.eq.'mtzvpp'.or.method.eq.'r2scan3c')then !SG
xva(i)=1.0d0 !SG
if(iz(i).eq.6) xva(i)=3.0d0 !SG
if(iz(i).eq.7) xva(i)=0.5d0 !SG
if(iz(i).eq.8) xva(i)=0.5d0 !SG
else
xva(i)=(nbas(iz(i))-0.5d0*nel(iz(i)))
endif
if(emiss(iz(i)).gt.01e-6.and.xva(i).lt.0.0d0) then
write(*,'(a)') 'element/emiss/nvirt/nel/nbas'
write(*,'(a4,f8.5,f8.5,i3,i3)') esym(iz(i)),emiss(iz(i)),xva(i),nel(iz(i)),nbas(iz(i))
error stop 'negative number of virtual orbitals. Something is wrong in the parameters !'
endif
if( emiss(iz(i)).lt.0.1e-6) then
! print*, '** WARNING ** -> element ',esym(iz(i)),' has no parameters (no contribution)!'
! xva(i)=0
cycle
warn=.true.
endif
if(xva(i).lt.0.5d0) then
write(*,'(3a)') '** WARNING ** -> element ',esym(iz(i)),' has no virtual orbitals (no contribution)!'
warn=.true.
endif
enddo
xvb=xva
endif !B97-3c
!*********************************
!* Parameter handling/setup done *
!*********************************
if(echo)then
write(*,*)''
write(*,'(2x,''level '',3x,a12,2x,a)') adjustr(trim(method)),adjustr(trim(getlevel(method)))
write(*,'(2x,''Nbf '',3x,I12)') nb
write(*,'(2x,''Atoms '',3x,I12)') n
write(*,*)''
if(.not.srb) then
write(*,'(2x,a)') 'Parameters: '
write(*,'(2x,a,2x,f10.4)') 'sigma ',p(1)
write(*,'(2x,a,2x,f10.4)') 'eta ',p(2)
write(*,'(2x,a,2x,f10.4)') 'alpha ',p(3)
write(*,'(2x,a,2x,f10.4)') 'beta ',p(4)
if(damp) then
write(*,'(2x,a,2x,f10.4)') 'dmp_scal',dmp_scal
write(*,'(2x,a,2x,f10.4)') 'dmp_exp ',dmp_exp
endif
else
write(*,'(/2x,a22)')'Parameters for SRB: '
write(*,'(2x,a,2x,f10.4)') 'rscal',rscal
write(*,'(2x,a,2x,f10.4)') 'qscal',qscal
endif
write(*,*)''
endif
! **************************
! * print parameter table *
! **************************
if(.not.srb) then
if(echo) then
write(*,*) ' '
write(*,*) 'element parameters ',trim(method)
write(*,'(2x,3(a4,2x,a6,3x,a4,4x))') 'elem','emiss','nbas','elem','emiss','nbas','elem','emiss','nbas'
do i=1,36,3
write(*,'(2x,3(A4,2x,F8.5,2x,I3,4x))') esym(i), emiss(i),nbas(i), &
esym(i+1), emiss(i+1),nbas(i+1),&
esym(i+2), emiss(i+2),nbas(i+2)
enddo
write(*,*) ' '
endif
endif
! **************************
! * write parameter file *
! **************************
if(parfile) then
write(*,'(a)') 'printing extended parameter file ...'
write(*,'(a)') 'modfiy & copy gcp.param to $HOME/.gcppar.$HOSTNAME'
open(unit=22,file='gcp.param')
write(22,*) ' # ',trim(method),' (comment line)'
write(22,'(x,F8.5)') p(1:4)
do i=1,36
write(22,'(x,I3,x,F8.5)') nbas(i), emiss(i)
enddo
close(22)
endif
!**************************
!* calc energy & gradient *
!**************************
gcp_g(1:3,1:n)=0.0d0
gcp_e=0.0d0
if(srb) then
call srb_egrad2(xyz,iz,lat,n,gcp_e,gcp_g,dograd,rscal,qscal,echo,pbc)
else
call gcp_egrad(n,max_elem,emiss,xyz,iz,p,gcp_e,gcp_g,dograd,echo,xva,xvb,pbc,lat,damp,base,dmp_scal,dmp_exp)
endif
!*************************
!* calc num. cell stress *
!*************************
if(dograd.and.pbc) then
call sgcp(n,max_elem,emiss,xyz,iz,p,gcp_e,gcp_g,dograd,echo,xva,xvb,lat,gcp_glat,base,damp,dmp_scal,dmp_exp)
! calc hessian
if(dohess) then
! call hessian(n,max_elem,emiss,xyz,iz,p,gcp_e,xva,xvb,dohess,damp,base,dmp_scal,dmp_exp)
call hessian(n,max_elem,emiss,xyz,iz,p,xva,xvb,pbc,lat,damp,base,dmp_scal,dmp_exp)
endif
! restore original cardinal numbers
iz=iz2
! set cart. frozen coords to zero
if(maxval(ifrez,n).eq.1.and.dograd) then
if(echo) write(*,'(a)') ' setting gradient of frozen cartesians to zero'
do i=1,n
if(ifrez(i).eq.1) gcp_g(1:3,i)=0.0d0
enddo
endif
!output
if(echo) then
write(*,*) ' '
write(*,*) '|S|=',abs(sum(gcp_glat(1:3,1:3)))
if(.not.verb) then
write(*,'(a)') 'Write cell gradient to file gcp_cellgradient'
open(unit=1,file='gcp_cellgradient')
write(1,*) gcp_glat
close(1)
end if
write(*,*) ' '
end if
end if
!for scripting
open(unit=1,file='.CPC')
write(1,*) gcp_e
close(1)
open(unit=1,file='.CP')
write(1,'(F22.16)') gcp_e
close(1)
end subroutine gcp_call
!***************************************************
!* This subroutine calculates the gcp correction *
!***************************************************
subroutine gcp_egrad(n,max_elem,emiss,xyz,iz,p,gcp,g,grad,echo,xva,xvb,pbc,lat,damp,base,dmp_scal,dmp_exp)
implicit none
integer n,iz(n),max_elem,np
integer iat,jat
real*8 xyz(3,n), xyzjat(3)
real*8 g (3,n)
real*8 gcp,tg,t0,t1
real*8 emiss(max_elem),dum22
real*8 p(4),thrR,thrE
real*8 xva(*),xvb(*)
real*8 va,vb
real*8 r,rscal,rscalexp,rscalexpm1,r0abij
real*8 tmp,ea,ecp,dum,tmpa,tmpb,dum2,tmpc,tmpd
real*8 sab,p1,p2,p3,p4,ene_old_num,ene_old_den
real*8 vec(3),gs(3),gab,ebas,gbas(3,n)
real*8 za(36),zb(36),lat(3,3)
logical echo,grad,damp,base,pbc
integer tau_max(3),a,b,c
real*8 rad(94)
data rad /&
0.32D0,0.37D0,1.30D0,0.99D0,0.84D0,0.75D0,0.71D0,0.64D0,0.60D0,&
0.62D0,1.60D0,1.40D0,1.24D0,1.14D0,1.09D0,1.04D0,1.00D0,1.01D0,&
2.00D0,1.74D0,1.59D0,1.48D0,1.44D0,1.30D0,1.29D0,1.24D0,1.18D0,&
1.17D0,1.22D0,1.20D0,1.23D0,1.20D0,1.20D0,1.18D0,1.17D0,1.16D0,&
2.15D0,1.90D0,1.76D0,1.64D0,1.56D0,1.46D0,1.38D0,1.36D0,1.34D0,&
1.30D0,1.36D0,1.40D0,1.42D0,1.40D0,1.40D0,1.37D0,1.36D0,1.36D0,&
2.38D0,2.06D0,1.94D0,1.84D0,1.90D0,1.88D0,1.86D0,1.85D0,1.83D0,&
1.82D0,1.81D0,1.80D0,1.79D0,1.77D0,1.77D0,1.78D0,1.74D0,1.64D0,&
1.58D0,1.50D0,1.41D0,1.36D0,1.32D0,1.30D0,1.30D0,1.32D0,1.44D0,&
1.45D0,1.50D0,1.42D0,1.48D0,1.46D0,2.42D0,2.11D0,2.01D0,1.90D0,&
1.84D0,1.83D0,1.80D0,1.80D0 /
! cut-off radii for all element pairs
real*8 r0ab(max_elem,max_elem),autoang
!real*8 autoang
!for damping
real*8 dampval, grdfirst,grdsecond,ene_old,ene_dmp,grd_dmp,dmp_scal,dmp_exp
parameter (autoang =0.5291772083d0)
! Two threshold. thrR: distance cutoff thrE: numerical noise cutoff
thrR=60 ! 60 bohr
thrE=epsilon(1.d0) ! cut below machine precision rounding
if(echo) then
write(*,*) ' '
write(*,'('' cutoffs: '',F5.1)')
write(*,'('' distance [bohr]'',F5.1)') thrR
write(*,'('' noise [a.u.]'',Es8.1)') thrE
write(*,'('' SR-damping '',L2)') damp
write(*,*) ' '
endif
g=0.0d0
ecp=0.0d0
gcp=0.0d0
dum=0.0d0
tg=0d0
p1=abs(p(1))
p2=abs(p(2))
p3=abs(p(3))
p4=abs(p(4))
!get r0
call setr0ab(max_elem,autoang,r0ab)
if(abs(p(1)-1.0d0).lt.1.d-6.and.abs(p(2)-1.315d0).lt.1.d-6)then !SG dirty way to set the special r2SCAN-3c mode
call setzet(p2,1.15d0,za,zb)
else
call setzet(p2,1.00d0,za,zb)
endif
gs=0
if(echo) write(*,'(2x,a5,2x,a5,2x,a5,2x,a7,2x,a14,4x,a15)') &
'#','ON','sites','Nvirt','Emiss','BSSE [kcal/mol]'
if(pbc) then
!Determine supercell
call criteria(thrR, lat, tau_max)
! Loop over all i atoms
do iat=1,n
va=xva(iat)
ea=0.0d0
np=0
! the BSSE due to atom jat, Loop over all j atoms in supercell
do jat=1,n
! # of bf that are available from jat
vb=xvb(jat)
if(vb.lt.0.5) cycle
do a=-tau_max(1),tau_max(1)
do b=-tau_max(2),tau_max(2)
do c=-tau_max(3),tau_max(3)
!remove selfinteraction
if((iat.eq.jat).and.(abs(a)+abs(b)+abs(c).eq.0)) cycle
xyzjat(1)=xyz(1,jat)+a*lat(1,1)+b*lat(2,1)+c*lat(3,1)
xyzjat(2)=xyz(2,jat)+a*lat(1,2)+b*lat(2,2)+c*lat(3,2)
xyzjat(3)=xyz(3,jat)+a*lat(1,3)+b*lat(2,3)+c*lat(3,3)
! write(*,*) 'a1 a2 a3', lat(1,1),lat(1,2),lat(1,3)
! write(*,*) 'b1 b2 b3', lat(2,1),lat(2,2),lat(2,3)
! write(*,*) 'c1 c2 c3', lat(3,1),lat(3,2),lat(3,3)
! stop
! xyzjat(1)=xyz(1,jat)+a*lat(1,1)+b*lat(1,2)+c*lat(1,3)
! xyzjat(2)=xyz(2,jat)+a*lat(2,1)+b*lat(2,2)+c*lat(2,3)
! xyzjat(3)=xyz(3,jat)+a*lat(3,1)+b*lat(3,2)+c*lat(3,3)
vec(1)=(xyz(1,iat)-xyzjat(1))
vec(2)=(xyz(2,iat)-xyzjat(2))
vec(3)=(xyz(3,iat)-xyzjat(3))
r=sqrt(vec(1)*vec(1)+vec(2)*vec(2)+vec(3)*vec(3))
! distance cutoff
if(r.gt.thrR) cycle
! calulate slater overlap sab
call ssovl(r,iat,jat,iz,za(iz(iat)),zb(iz(jat)),sab)
! noise cutoff(sqrt(sab))
if(sqrt(abs(sab)).lt.thrE) cycle
! evaluate gcp central expression
ene_old_num=exp(-p3*r**p4)
ene_old_den=sqrt(vb*Sab)
ene_old=ene_old_num/ene_old_den
! noise cutoff(damp)
if(abs(ene_old).lt.thrE) cycle
if(damp) then
!D3 r0ab radii
r0abij=r0ab(iz(iat),iz(jat))
rscal=r/r0abij
!covalent radii
! rscal=r/(rad(iz(iat))+rad(iz(jat)))
rscalexp=rscal**dmp_exp
dampval=(1d0-1d0/(1d0+dmp_scal*rscalexp))
ene_dmp=ene_old*dampval
ea=ea+emiss(iz(iat))*ene_dmp
else
ea=ea+emiss(iz(iat))*ene_old
endif
! sites counter (i.e. # atoms contributing to the 'atomic bsse')
np=np+1
! gradient for i,j pair
if(grad)then
call cpu_time(t0)
call gsovl(r,iat,jat,iz,za(iz(iat)),zb(iz(jat)),gab)
gs(1)=gab*vec(1)
gs(2)=gab*vec(2)
gs(3)=gab*vec(3)
dum=exp(-p3*r**p4)*(-1d0/2d0)
dum2=2d0*p3*p4*r**p4*sab/r
dum22=r*sab**(3d0/2d0)
tmpb=dum22*sqrt(vb)
if(damp) then
rscalexpm1=rscal**(dmp_exp-1)
grd_dmp=dmp_scal*dmp_exp*rscalexpm1/r0abij
grd_dmp=grd_dmp/((dmp_scal*rscalexp+1.0d0)**2)
endif
tmpa=dum2*vec(1)+gs(1)
tmp=dum*tmpa/tmpb
if(damp) then
tmp=tmp*dampval+ene_old*grd_dmp*(vec(1)/r)
endif
g(1,iat)=g(1,iat)+tmp*emiss(iz(iat))
tmpa=dum2*vec(2)+gs(2)
tmp=dum*tmpa/tmpb
if(damp) then
tmp=tmp*dampval+ene_old*grd_dmp*(vec(2)/r)
endif
g(2,iat)=g(2,iat)+tmp*emiss(iz(iat))
tmpa=dum2*vec(3)+gs(3)
tmp=dum*tmpa/tmpb
if(damp) then
tmp=tmp*dampval+ene_old*grd_dmp*(vec(3)/r)
endif
g(3,iat)=g(3,iat)+tmp*emiss(iz(iat))
if(va.lt.0.5) cycle
if(damp) then
ene_old_den=sqrt(va*Sab)
ene_old=ene_old_num/ene_old_den
endif
tmpb=dum22*sqrt(va)
tmpa=dum2*(-vec(1))-gs(1)
tmp=dum*tmpa/tmpb
if(damp) then
tmp=tmp*dampval+ene_old*grd_dmp*(-vec(1)/r)
endif
g(1,iat)=g(1,iat)-tmp*emiss(iz(jat))
tmpa=dum2*(-vec(2))-gs(2)
tmp=dum*tmpa/tmpb
if(damp) then
tmp=tmp*dampval+ene_old*grd_dmp*(-vec(2)/r)
endif
g(2,iat)=g(2,iat)-tmp*emiss(iz(jat))
tmpa=dum2*(-vec(3))-gs(3)
tmp=dum*tmpa/tmpb
if(damp) then
tmp=tmp*dampval+ene_old*grd_dmp*(-vec(3)/r)
endif
g(3,iat)=g(3,iat)-tmp*emiss(iz(jat))
call cpu_time(t1)
tg=tg+t1-t0
endif
! end of j-loop
enddo
!end supercell loop
enddo
enddo
enddo
if(echo) &
write(*,'(2x,3(I5,2x),2x,F5.1,2x,F14.4,2x,F14.4,2x,a)') &
iat,iz(iat),np,va,emiss(iz(iat)), ea*627.5099*p(1)
ecp=ecp+ea
! end of i-loop
enddo
else
! Loop over all i atoms
do iat=1,n
va=xva(iat)
ea=0.0d0
np=0
! the BSSE due to atom jat, Loop over all j atoms
do jat=1,n
if(iat.eq.jat) cycle
vec(1)=(xyz(1,iat)-xyz(1,jat))
vec(2)=(xyz(2,iat)-xyz(2,jat))
vec(3)=(xyz(3,iat)-xyz(3,jat))
r=sqrt(vec(1)*vec(1)+vec(2)*vec(2)+vec(3)*vec(3))
! # of bf that are available from jat
vb=xvb(jat)
if(vb.lt.0.5) cycle
! distance cutoff
if(r.gt.thrR) cycle
! calulate slater overlap sab
call ssovl(r,iat,jat,iz,za(iz(iat)),zb(iz(jat)),sab)
! noise cutoff(sqrt(sab))
if(sqrt(abs(sab)).lt.thrE) cycle
! evaluate gcp central expression
ene_old_num=exp(-p3*r**p4)
ene_old_den=sqrt(vb*Sab)
ene_old=ene_old_num/ene_old_den
! noise cutoff(damp)
if(abs(ene_old).lt.thrE) cycle
if(damp) then
!D3 r0ab radii
r0abij=r0ab(iz(iat),iz(jat))
rscal=r/r0abij
!covalent radii
! rscal=r/(rad(iz(iat))+rad(iz(jat)))
rscalexp=rscal**dmp_exp
dampval=(1d0-1d0/(1d0+dmp_scal*rscalexp))
ene_dmp=ene_old*dampval
ea=ea+emiss(iz(iat))*ene_dmp
else
ea=ea+emiss(iz(iat))*ene_old
endif
! sites counter (i.e. # atoms contributing to the 'atomic bsse')
np=np+1
! gradient for i,j pair
if(grad)then
call cpu_time(t0)
call gsovl(r,iat,jat,iz,za(iz(iat)),zb(iz(jat)),gab)
gs(1)=gab*vec(1)
gs(2)=gab*vec(2)
gs(3)=gab*vec(3)
dum=exp(-p3*r**p4)*(-1d0/2d0)
dum2=2d0*p3*p4*r**p4*sab/r
dum22=r*sab**(3d0/2d0)
tmpb=dum22*sqrt(vb)
if(damp) then
rscalexpm1=rscal**(dmp_exp-1)
grd_dmp=dmp_scal*dmp_exp*rscalexpm1/r0abij
grd_dmp=grd_dmp/((dmp_scal*rscalexp+1.0d0)**2)
endif
tmpa=dum2*vec(1)+gs(1)
tmp=dum*tmpa/tmpb
if(damp) then
tmp=tmp*dampval+ene_old*grd_dmp*(vec(1)/r)
endif
g(1,iat)=g(1,iat)+tmp*emiss(iz(iat))
tmpa=dum2*vec(2)+gs(2)
tmp=dum*tmpa/tmpb
if(damp) then
tmp=tmp*dampval+ene_old*grd_dmp*(vec(2)/r)
endif
g(2,iat)=g(2,iat)+tmp*emiss(iz(iat))
tmpa=dum2*vec(3)+gs(3)
tmp=dum*tmpa/tmpb
if(damp) then
tmp=tmp*dampval+ene_old*grd_dmp*(vec(3)/r)
endif
g(3,iat)=g(3,iat)+tmp*emiss(iz(iat))
if(va.lt.0.5) cycle
if(damp) then
ene_old_den=sqrt(va*Sab)
ene_old=ene_old_num/ene_old_den
endif
tmpb=dum22*sqrt(va)
tmpa=dum2*(-vec(1))-gs(1)
tmp=dum*tmpa/tmpb
if(damp) then
tmp=tmp*dampval+ene_old*grd_dmp*(-vec(1)/r)
endif
g(1,iat)=g(1,iat)-tmp*emiss(iz(jat))
tmpa=dum2*(-vec(2))-gs(2)
tmp=dum*tmpa/tmpb
if(damp) then
tmp=tmp*dampval+ene_old*grd_dmp*(-vec(2)/r)
endif
g(2,iat)=g(2,iat)-tmp*emiss(iz(jat))
tmpa=dum2*(-vec(3))-gs(3)
tmp=dum*tmpa/tmpb
if(damp) then
tmp=tmp*dampval+ene_old*grd_dmp*(-vec(3)/r)
endif
g(3,iat)=g(3,iat)-tmp*emiss(iz(jat))
call cpu_time(t1)
tg=tg+t1-t0
endif
! end of j-loop
enddo
if(echo) &
write(*,'(2x,3(I5,2x),2x,F5.1,2x,F14.4,2x,F14.4,2x,a)') &
iat,iz(iat),np,va,emiss(iz(iat)), ea*627.5099*p(1)
ecp=ecp+ea
! end of i-loop
enddo
endif
gcp=ecp*p(1)
g=g*p(1)
!Special HF-3c correction
if(base) then
call basegrad(n,max_elem,iz,xyz,lat,pbc,0.7d0,0.03d0,ebas,gbas,echo)
gcp=gcp+ebas
if(grad) g=g+gbas
endif
if(grad.and.echo) call prtim(6,tg,'t','gradient')
end subroutine gcp_egrad
!*********************************
!* subroutine calc stress tensor *
!*********************************
subroutine sgcp(n,max_elem,emiss,xyz,iz,p,gcp,g,grad,echo,xva,xvb,lat,sigma,base,damp,dmp_scal,dmp_exp)
implicit none
integer n,iz(n),max_elem,np
integer iat,jat
real*8 xyz(3,n), abc(3,n), xyzdum(3,n)
real*8 g(3,n)
real*8 gcp,tg,t0,t1
real*8 emiss(max_elem),dum22
real*8 p(4),thrR,thrE
real*8 xva(*),xvb(*)
real*8 va,vb
real*8 r,vec(3)
real*8 tmp,ea,ecp,dum,tmpa,tmpb,dum2
real*8 sab,p1,p2,p3,p4,dmp_scal,dmp_exp
real*8 gs(3),gab
real*8 za(36),zb(36), detlat
real*8 lat(3,3), xyzjat(3,n), sigma(3,3), delta, gcpl, gcpr, vol, latdum(3,3), gdum(3,n), edum
logical echo,grad,damp,base
integer a, b, c, tau_max(3), abcsum, selftest, i, j,k
delta=10E-6
!cell volume
latdum=0
detlat=0
detlat =detlat+lat(1,1)*lat(2,2)*lat(3,3)+lat(1,2)*lat(2,3)*lat(3,1)+lat(1,3)*lat(2,1)*lat(3,2)
detlat=detlat-lat(3,1)*lat(2,2)*lat(1,3)-lat(3,2)*lat(2,3)*lat(1,1)-lat(3,3)*lat(2,1)*lat(1,1)
vol=abs(detlat)
! Transform to direct system
! simultaniously shilft all atoms
latdum=lat
call xyz2abc(n,xyz,abc,lat)
do i=1,3
do j=1,3
latdum(i,j)=lat(i,j)-delta
call abc2xyz(n,xyzdum,abc,latdum)
call gcp_egrad(n,max_elem,emiss,xyzdum,iz,p,gcpl,gdum,.false.,.false.,xva,xvb,.true.,latdum,damp,base,dmp_scal,dmp_exp)
! call e(n,max_elem,emiss,xyzdum,iz,p,gcpl,gdum,.false.,.false.,xva,xvb,.true.,latdum,damp,dmp_scal,dmp_exp)
latdum(i,j)=lat(i,j)+delta
call abc2xyz(n,xyzdum,abc,latdum)
! call e(n,max_elem,emiss,xyzdum,iz,p,gcpr,gdum,.false.,.false.,xva,xvb,.true.,latdum,damp,dmp_scal,dmp_exp)
call gcp_egrad(n,max_elem,emiss,xyzdum,iz,p,gcpr,gdum,.false.,.false.,xva,xvb,.true.,latdum,damp,base,dmp_scal,dmp_exp)
! symmetric numerical derivative
sigma(i,j)=(gcpr-gcpl)/(2*delta)
!sigma(i,j)=sigma(i,j)/vol
end do
end do
end subroutine sgcp
!*****************
!* print timings *
!*****************
subroutine prtim(io,tt,is,string)
integer io
real*8 tt,t,sec
integer day,hour,min
character*(*) string
character*1 is
t=tt
day=idint(t/86400)
t=t-day*86400
hour=idint(t/3600)
t=t-hour*3600
min=idint(t/60)
t=t-60*min
sec=t
if(day.ne.0)then
if(is.eq.'w')write(io,2)trim(string),day,hour,min,sec
if(is.eq.'t')write(io,22)trim(string),day,hour,min,sec
return
endif
if(hour.ne.0)then
if(is.eq.'w')write(io,3)trim(string),hour,min,sec
if(is.eq.'t')write(io,33)trim(string),hour,min,sec
return
endif
if(min .ne.0)then
if(is.eq.'w')write(io,4)trim(string),min,sec
if(is.eq.'t')write(io,44)trim(string),min,sec
return
endif
if(is.eq.'w')write(io,5)trim(string),sec
if(is.eq.'t')write(io,55)trim(string),sec
return
1 format('======================================')
2 format('wall-time for ',a,2x,i3,' d',i3,' h',i3,' m',f5.1,' s')
3 format('wall-time for ',a,2x,i3,' h',i3,' m',f5.1,' s')
4 format('wall-time for ',a,2x,i3,' m',f5.1,' s')
5 format('wall-time for ',a,2x,f5.1,' s')
22 format('cpu-time for ',a,2x,i3,' d',i3,' h',i3,' m',f5.1,' s')
33 format('cpu-time for ',a,2x,i3,' h',i3,' m',f5.1,' s')
44 format('cpu-time for ',a,2x,i3,' m',f5.1,' s')
55 format('cpu-time for ',a,2x,f5.1,' s')
return
end
!*********************************
!* Load all necessary parameters *
!*********************************
subroutine setparam(emiss,nbas,p,method)
implicit none
integer, parameter :: mpar=86 ! maximal
integer, parameter :: apar=36 ! actual
integer nbas(mpar)
character*(*) method
real(8) emiss(mpar),p(*)
real(8) HFsv(apar),HFminis(apar),HF631gd(apar),HFsvp(apar),HFtz(apar),&
HFvmb(apar),HFminisd(apar),oldHFsvp(apar),HFpobtz(apar),HFpobdzvp(apar),&
HF2gcore(apar),HF2g(apar), HFdef1tzvp(apar),HFccdz(apar),HFaccdz(apar),&
HFdzp(apar),HFhsv(apar),HFdz(apar),HFmsvp(apar),HFdef2mtzvp(apar),HFdef2mtzvpp(apar) !SG
integer BASsv(apar),BASminis(apar),BAS631gd(apar),BAStz(apar),&
BASsvp(apar),BASvmb(apar),BASminisd(apar),oldBASsvp(apar),BASpobtz(apar),BASpobdzvp(apar),&
BAS2gcore(apar),BAS2g(apar),BASdef1tzvp(apar),BASccdz(apar),BASaccdz(apar),&
BASdzp(apar),BAShsv(apar),BASdz(apar),BASmsvp(apar),BASdef2mtzvp(apar),BASdef2mtzvpp(apar)
real(8) HFlanl2(10)
integer BASlanl2(10)
! ***************
! * Emiss data *
! ***************
data HFsv/ & !H-Kr (no 3d)
0.009037,0.008843,& ! He,He
0.204189,0.107747,0.049530,0.055482,0.072823,0.100847,0.134029,0.174222,& ! Li-Ne
0.315616,0.261123,0.168568,0.152287,0.146909,0.168248,0.187882,0.211160,& !Na -Ar
0.374252,0.460972,& ! K-Ca
0.444886,0.404993,0.378406,0.373439,0.361245,0.360014,0.362928,0.243801,0.405299,0.396510,& ! 3d-TM
0.362671,0.360457,0.363355,0.384170,0.399698,0.417307/ !Ga-Kr
data HFminis/ &
0.042400,0.028324,&
0.252661,0.197201,0.224237,0.279950,0.357906,0.479012,0.638518,0.832349, &
1.232920,1.343390,1.448280,1.613360,1.768140,1.992010,2.233110,2.493230, &
3.029640,3.389980,& ! H-Ca
10*0,6*0/
data HF631GD/ &! H-Ca + Br (no 3d)
0.010083,0.008147,&
0.069260,0.030540,0.032736,0.021407,0.024248,0.036746,0.052733,0.075120,&
0.104255,0.070691,0.100260,0.072534,0.054099,0.056408,0.056025,0.057578,&
0.079198,0.161462,&
10*0.0, &
0.000000,0.000000,0.000000,0.000000,0.381049,0.000000/
data oldHFsvp / & ! Li,Na,Mg,K had wrong parameters
0.008107,0.008045,&
0.142957,0.028371,0.049369,0.055376,0.072785,0.100310,0.133273,0.173600,&
0.191109,0.222839,0.167188,0.149843,0.145396,0.164308,0.182990,0.205668,&
0.221189,0.299661,&
0.325995,0.305488,0.291723,0.293801,0.29179,0.296729,0.304603,0.242041,0.354186,0.350715,&
0.350021,0.345779,0.349532,0.367305,0.382008,0.399709/
data HFsvp / & ! H-Kr
0.008107,0.008045,&
0.113583,0.028371,0.049369,0.055376,0.072785,0.100310,0.133273,0.173600,&
0.181140,0.125558,0.167188,0.149843,0.145396,0.164308,0.182990,0.205668,&
0.200956,0.299661, &
0.325995,0.305488,0.291723,0.293801,0.29179,0.296729,0.304603,0.242041,0.354186,0.350715,&
0.350021,0.345779,0.349532,0.367305,0.382008,0.399709/
data HFtz /& ! H-Kr !def2-TZVP
0.007577,0.003312,&
0.086763,0.009962,0.013964,0.005997,0.004731,0.005687,0.006367,0.007511,&
0.077721,0.050003,0.068317,0.041830,0.025796,0.025512,0.023345,0.022734,&
0.097241,0.099167,&
0.219194,0.189098,0.164378,0.147238,0.137298,0.12751,0.118589,0.0318653,0.120985,0.0568313, &
0.090996,0.071820,0.063562,0.064241,0.061848,0.061021/
data HFdef2mtzvp /& ! m def2-TZVP, no f for B-Ne
0.007930,0.003310,&
0.086760,0.009960,0.013960,0.006000,0.003760,0.004430,0.005380,0.006750,&
0.077720,0.050000,0.068320,0.041830,0.025800,0.025510,0.023340,0.022730,&
0.097240,0.099170,&
0.219190,0.189100,0.164380,0.147240,0.137300,0.127510,0.118590,0.031870,0.120990,0.056830,&
0.091000,0.071820,0.063560,0.064240,0.061850,0.061020/
data HFvmb/&
0.042400,0.028324,&
0.252661,0.197201,0.156009,0.164586,0.169273,0.214704,0.729138,0.336072,&
0.262329,0.231722,0.251169,0.287795,0.213590,0.250524,0.728579,0.260658, &
2*0,&
16*0/
data HFminisd/& !Al-Ar MINIS + Ahlrichs "P" funktions
0.0,0.0,&
8*0.0,&
2*0.0,1.446950,1.610980,1.766610,1.988230,2.228450,2.487960,&
2*0.0,&
16*0.0/
data HFlanl2/ & ! LANL2TZ+ vs LANL2DZ (ORCA), only Sc-Zn
0.102545,0.0719529,0.0491798,0.0362976,0.0266369,0.0235484,0.0171578,0.0438906,0.0100259,0.016208/
data HFpobtz / & ! H-Kr, no RG
0.010077,0.000000,&
0.173239,0.101973,0.131181,0.032234,0.011630,0.008475,0.011673,0.000000,&
0.240653,0.661819,0.522306,0.14163,0.052456,0.184547,0.040837,0.000000,&
0.318136,0.564721,&
0.523194,0.767449,0.620122,0.390227,0.237571,0.247940,0.249589,0.117864,0.325725,0.197183,&
0.264489,0.180375,0.111262,0.147239,0.081747,0.000000/
data HF2gcore / & ! only HCNOF yet
0.000539,0.000000,&
0.000000,0.000000,0.000000,0.173663,0.269952,0.364341,0.384923,0.000000,&
0.000000,0.000000,0.000000,0.000000,0.000000,0.000000,0.000000,0.000000,&
0.000000,0.000000,&
0.000000,0.000000,0.000000,0.000000,0.000000,0.000000,0.000000,0.000000,0.000000,0.000000,&