-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathRUQT.f90
1091 lines (875 loc) · 36.9 KB
/
RUQT.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 RUQT
Use InterfaceMod
Use FunctionMod
!Use InterfaceMod2
Use TypeMod
implicit none
real(8),allocatable,dimension(:,:) :: H_one,H_Two,OneInts,H_Two_le,H_Two_re,H_Two_cen,H_Two_le_trans,H_Two_re_trans
real(8), allocatable, dimension(:) :: TwoIntsCompact,transm,current,energy_list,mo_ener,B0_coeff
real(8), allocatable, dimension(:,:) :: Smat_le,Smat_re,Smat_cen,Smat,coupling_mat,mo_coeff,mo_coeff2
complex(8), allocatable, dimension(:,:) :: gfc_r,gfc_a,current_temp,Sigma_l,Sigma_r,Gamma_L,Gamma_R
complex(8),allocatable,dimension(:) :: voltage
real(8) :: coupling_r,coupling_l,energy,energy_start,energy_end,delta_en,volt_start,volt_end,delta_volt,corr_ener
integer :: size_l,size_r,size_c,size_lc,size_lcr,numatomic,num_threads,numfcore,numfvirt
character(len=100) :: inputfile,outfile,option
logical :: libint,gamess,rdm_flag,invert,doubles,currentflag,cisd_flag,qchem,hf_flag
logical :: dft_flag,pyscf,maple,use_b0,molcas
integer :: i,j,k,counter,counter2,current_values,norb,numact,energy_val,volt_val,ioerror,numocc,numvirt
character(len=40) :: ElectrodeType,CalcType,functional,inputcode,b0_type
real(8) :: KT,current_con,Fermi_enl,Fermi_enR,localden_fermi,localden_fermi_l,localden_fermi_r,temp
complex(8), allocatable, dimension(:,:) :: test
type(B1) :: B1data,l1data
type(B2) :: B2data,l2data
type(energr) :: G_S
real(8) :: time_start,time_end
call cpu_time(time_start)
invert=.true.
currentflag=.false.
write(*,*) "Just getting started"
Call Get_Command_Argument(1,inputfile)
Call ReadInput(inputfile,norb,numfcore,numfvirt,numocc,numvirt,size_l,size_r,size_c,energy_start,energy_end,delta_en,volt_start,volt_end,delta_volt,inputcode,KT,ElectrodeType,Fermi_enl,Fermi_enR,CalcType,localden_fermi_l,localden_fermi_r,doubles,numatomic,functional,num_threads,use_b0,b0_type)
write(*,*) "Input File Read"
write(*,*) "Using the following parameters for the transport calculation"
write(*,*) "Number of OpenMP Threads:",num_threads
write(*,*) "Number of Molecular and Atomic Orbitals:",norb,numatomic
write(*,*) "Number of Active Orbitals:",norb-numfcore-numfvirt
write(*,*) "Number of Occupied Orbitals:",numocc
write(*,*) "Number of Virtual Orbitals:",numvirt
write(*,*) "Orbitals in left electrode:",size_l
write(*,*) "Orbitals in device region:",size_c
write(*,*) "Orbitals in right electrode:",size_r
write(*,*) "Transmission Energy Window(eV) and dE:",energy_start,energy_end,delta_en
write(*,*) "Voltage Window and dV:",volt_start,volt_end,delta_volt
write(*,*) "Fermi Density:",localden_fermi_l,localden_fermi_r
write(*,*) "KT:",KT
size_lc = size_l + size_c
size_lcr = size_l + size_c + size_r
libint=.false.
Call Flag_set(inputcode,functional,cisd_flag,rdm_flag,hf_flag,dft_flag,qchem,gamess,pyscf,maple,molcas)
if(qchem.eqv..true.) then
write(*,*) 'Using Qchem data for this run'
Call Get_HF_Qchem(inputfile,norb,H_Two,Smat)
elseif(molcas.eqv..true.) then
Call Get_HF_Molcas(inputfile,norb,H_Two,Smat)
elseif(libint.eqv..true.) then
Call Get_HF_libint(inputfile,norb,numact,H_one,Smat,mo_coeff,OneInts,TwoIntsCompact)
elseif(gamess.eqv..true.) then
Call Get_HF_GAMESS(inputfile,numatomic,H_Two,Smat,norb)
write(*,*) 'Using GAMESS data for this run'
elseif(pyscf.eqv..true.) then
write(*,*) 'Using PySCF data for this run'
Call Get_HF_PySCF(inputfile,numatomic,H_Two,Smat,norb)
!Call Get_HF_Psi4()
elseif(maple.eqv..true.) then
write(*,*) "Using Maple+QuantumChemistry data for this run"
Call Get_HF_PySCF(inputfile,numatomic,H_Two,Smat,norb)
endif
write(*,*) 'This run using:'
if(rdm_flag.eqv..true.) then
write(*,*) 'The Lehmann representation of a'
if(doubles.eqv..true.) then
write(*,*) 'p2-RDM Greens function'
elseif(doubles.eqv..false.) then
write(*,*) 'HF Greens function'
end if
elseif(qchem.eqv..true.) then
write(*,*) 'QCHEM HF/DFT Greens function'
elseif(molcas.eqv..true.) then
write(*,*) 'Using Molcas Fock Matrix: FOCK_AO'
elseif(pyscf.eqv..true.) then
write(*,*) 'PYSCF HF/DFT Greens Function'
end if
!CALL OMP_SET_NUM_THREADS(num_threads)
!Here we want to parition the Smatrix and H matrix into electrodes and
!the device
if(trim(ElectrodeType).eq."Metal_WBL") then
write(*,*) 'Starting Metal WBL calculation'
allocate(Sigma_l(1:size_c,1:size_c))
allocate(Sigma_r(1:size_c,1:size_c))
Call PartitionHS_MetalWBL(Smat,H_Two,size_l,size_r,size_c,size_lc,size_lcr,Smat_le,Smat_re,Smat_cen,H_Two_le,H_Two_re,H_Two_cen,H_Two_le_trans,H_Two_re_trans)
write(*,*) 'Getting Electrodes'
Call Electrodes_MetalWBL(Sigma_l,Sigma_r,Smat_re,Smat_le,H_Two_le,H_Two_re,localden_fermi_l,localden_fermi_r,size_c,size_lc,size_lcr,size_l,size_r,H_Two_le_trans,H_Two_re_trans)
allocate(Gamma_L(1:size_c,1:size_c))
allocate(Gamma_R(1:size_c,1:size_c))
Gamma_L=-(DIMAG(Sigma_L)-DIMAG(adjoint(Sigma_L,size_c)))!,1E-12,8)
Gamma_R=-(DIMAG(Sigma_R)-DIMAG(adjoint(Sigma_R,size_c)))!,1E-12,8)
if(trim(CalcType).eq."current".or.trim(CalcType).eq."Current") then
write(*,*) "Starting Current Calculation"
volt_val = abs(int((volt_end - volt_start)/delta_volt)+1)
energy_val = abs(int((energy_start-energy_end)/delta_en)+1)
current_values = volt_val
allocate(gfc_r(1:size_c,1:size_c))
allocate(gfc_a(1:size_c,1:size_c))
allocate(current(1:current_values))
allocate(voltage(1:current_values))
allocate(transm(1:energy_val))
allocate(current_temp(1:size_c,1:size_c))
energy = energy_start
current_con = 2*1.6021766E-19*(4.135667E-15)**(-1)
counter = 1
do j=1,volt_val
temp = (j-1)*delta_volt + volt_start
voltage(j) = temp
current(j) = 0
do k=1,energy_val
transm(k) = 0
energy = energy_start + (k-1)*delta_en
current_temp = 0
if((hf_flag.eqv..true.).or.(dft_flag.eqv..true.)) then
gfc_r = 0
gfc_a = 0
gfc_r = energy*Smat_cen-H_Two_cen
gfc_r = gfc_r - Sigma_l - Sigma_r
gfc_r = inv(gfc_r)
gfc_a = adjoint(gfc_r,size_c)
else if((cisd_flag.eqv..true.).or.(rdm_flag.eqv..true.)) then
gfc_r = 0
gfc_a = 0
Call Build_G_SD_Invert(gfc_r,Sigma_l,Sigma_r,energy,size_l,size_c,size_lc,size_lcr,norb,inputfile,numocc,numvirt,counter,B1data,B2data,mo_ener,mo_coeff,mo_coeff2,doubles,currentflag,energy_val,k,G_S,corr_ener,numatomic,B0_coeff,use_b0,gamess,maple,numfcore,numfvirt,b0_type)
gfc_a = adjoint(gfc_r,size_c)
counter=2
end if
current_temp = matmul_zgemm(Gamma_R,gfc_a)
current_temp = matmul_zgemm(gfc_r,current_temp)
current_temp = matmul_zgemm(Gamma_L,current_temp)
do i=1,size_c
transm(k) = transm(k) + real(current_temp(i,i))
end do
transm(k) = transm(k)*(fermi_function(energy-temp*0.5,fermi_enL,KT)-fermi_function(energy+temp*0.5,fermi_enR,KT))
end do
currentflag=.true.
do k=1,energy_val
current(j) = current(j) + delta_en*current_con*transm(k)
end do
write(*,*) 'Done with current at voltage:',real(voltage(j))
end do
do j=1,volt_val
write(*,*) 'IV curve',real(voltage(j)),real(current(j))
end do
outfile = trim(inputfile) // ".dat"
open(unit=7,file=outfile,action='write',iostat=ioerror)
do j=1,volt_val
write(7,*) real(voltage(j)),real(current(j))
end do
close(7)
elseif(trim(CalcType).eq."Transmission".or.trim(CalcType).eq."transmission") then
write(*,*) "Starting Transmission Calculation at Energy"
energy_val = abs(int((energy_start-energy_end)/delta_en))
allocate(gfc_r(1:size_c,1:size_c))
allocate(gfc_a(1:size_c,1:size_c))
allocate(transm(1:energy_val))
allocate(energy_list(1:energy_val))
allocate(current_temp(1:size_c,1:size_c))
energy = energy_start
counter=1
do k=1,energy_val
transm(k) = 0
current_temp = 0
energy = energy_start + (k-1)*delta_en
energy_list(k) = energy
if((hf_flag.eqv..true.).or.(dft_flag.eqv..true.)) then
gfc_r = 0
gfc_a = 0
gfc_r = energy*Smat_cen-H_Two_cen
gfc_r = gfc_r - Sigma_l - Sigma_r
gfc_r = inv(gfc_r)
gfc_a = adjoint(gfc_r,size_c)
else if((cisd_flag.eqv..true.).or.(rdm_flag.eqv..true.)) then
gfc_r = 0
gfc_a = 0
Call Build_G_SD_Invert(gfc_r,Sigma_l,Sigma_r,energy,size_l,size_c,size_lc,size_lcr,norb,inputfile,numocc,numvirt,counter,B1data,B2data,mo_ener,mo_coeff,mo_coeff2,doubles,currentflag,energy_val,k,G_S,corr_ener,numatomic,B0_coeff,use_b0,gamess,maple,numfcore,numfvirt,b0_type)
gfc_a = adjoint(gfc_r,size_c)
counter=2
end if
current_temp = matmul_zgemm(Gamma_L,gfc_r)
current_temp = matmul_zgemm(current_temp,Gamma_R)
current_temp = matmul_zgemm(current_temp,gfc_a)
do i=1,size_c
transm(k) = transm(k) + real(current_temp(i,i))
end do
write(*,*) 'Done with transmission function at energy:',transm(k)
end do
do j=1,energy_val
write(*,*) 'Transm vs Energy curve',real(energy_list(j)),real(transm(j))
end do
outfile = trim(inputfile) // ".dat"
open(unit=7,file=outfile,action='write',iostat=ioerror)
do j=1,energy_val
write(7,*) real(energy_list(j)),real(transm(j))
end do
close(7)
end if
else if(trim(ElectrodeType).eq."Molecule_WBL") then
write(*,*) 'Using Molecule WBL Electrodes'
write(*,*) "***This option only supports DFT/HF calcs and is outdated/buggy***"
write(*,*) "***Use at own risk***"
allocate(Sigma_l(1:norb,1:norb))
allocate(Sigma_r(1:norb,1:norb))
Sigma_l = 0
Sigma_r = 0
write(*,*) 'Calculating Coupling'
Call Calculate_Coupling_MoleculeWBL(Coupling_R,Coupling_L,localden_fermi)
write(*,*) 'Calculating Electrodes'
Call Electrodes_MoleculeWBL(Sigma_l,Sigma_r,Smat,Coupling_R,Coupling_L,size_c,size_lc,size_lcr)
allocate(Gamma_L(1:size_c,1:size_c))
allocate(Gamma_R(1:size_c,1:size_c))
Gamma_L=DIMAG(Sigma_L)-DIMAG(adjoint(Sigma_L,norb))
Gamma_R=DIMAG(Sigma_R)-DIMAG(adjoint(Sigma_R,norb))
if(trim(CalcType).eq."current".or.trim(CalcType).eq."Current") then
write(*,*) "Starting Current Calculation"
volt_val = abs(int((volt_end - volt_start)/delta_volt))
energy_val = abs(int((energy_start-energy_end)/delta_en))
current_values = volt_val
allocate(gfc_r(1:norb,1:norb))
allocate(gfc_a(1:norb,1:norb))
allocate(current(1:current_values))
allocate(voltage(1:current_values))
allocate(transm(1:energy_val))
allocate(current_temp(1:norb,1:norb))
energy = energy_start
current_con = 1.6021766E-19*4.135667E-15**(-1)
counter = 1
do j=1,volt_val
temp = (j-1)*delta_volt + volt_start
voltage(j) = temp
current(j) = 0
do k=1,energy_val
transm(k) = 0
energy = energy_start + (k-1)*delta_en
current_temp = 0
gfc_r = 0
gfc_a = 0
gfc_r = energy*Smat-H_Two
gfc_r = gfc_r - Sigma_l - Sigma_r
gfc_r = inv(gfc_r)
gfc_a = adjoint(gfc_r,norb)
current_temp = matmul_zgemm(Gamma_R,gfc_a)
current_temp = matmul_zgemm(gfc_r,current_temp)
current_temp = matmul_zgemm(Gamma_L,current_temp)
do i=1,norb
transm(k) = transm(k) + real(current_temp(i,i))
end do
transm(k) = transm(k)*(fermi_function(energy+temp*0.5,fermi_enL,KT)-fermi_function(energy-temp*0.5,fermi_enR,KT))
end do
do k=1,energy_val
current(j) = current(j) + delta_en*current_con*transm(k)
end do
end do
do j=1,volt_val
write(*,*) 'IV curve',real(voltage(j)),real(current(j))
end do
outfile = trim(inputfile) // ".dat"
open(unit=7,file=outfile,action='write',iostat=ioerror)
do j=1,volt_val
write(7,*) real(voltage(j)),real(current(j))
end do
close(7)
elseif(trim(CalcType).eq."Transmission".or.trim(CalcType).eq."transmission") then
write(*,*) "Starting Transmission Calculation at Energy"
energy_val = abs(int((energy_start-energy_end)/delta_en))
allocate(gfc_r(1:norb,1:norb))
allocate(gfc_a(1:norb,1:norb))
allocate(transm(1:energy_val))
allocate(energy_list(1:energy_val))
allocate(current_temp(1:norb,1:norb))
energy = energy_start
temp = volt_start
do k=1,energy_val
transm(k) = 0
energy = energy_start + (k-1)*delta_en
energy_list(k) = energy
gfc_r = 0
gfc_a = 0
gfc_r = energy*Smat-H_Two
gfc_r = gfc_r - Sigma_l - Sigma_r
gfc_r = inv(gfc_r)
gfc_a = adjoint(gfc_r,norb)
current_temp = matmul_zgemm(Gamma_R,gfc_a)
current_temp = matmul_zgemm(gfc_r,current_temp)
current_temp = matmul_zgemm(Gamma_L,current_temp)
do i=1,norb
transm(k) = transm(k) + real(current_temp(i,i))
end do
end do
do j=1,energy_val
write(*,*) 'Transm vs Energy curve',real(energy_list(j)),real(transm(j))
end do
outfile = trim(inputfile) // ".dat"
open(unit=7,file=outfile,action='write',iostat=ioerror)
do j=1,energy_val
write(7,*) real(energy_list(j)),real(transm(j))
end do
close(7)
end if
end if
call cpu_time(time_end)
write(*,*) 'CPU Time:',time_end-time_start
end program
function adjoint(A,norb)
implicit none
complex(8),allocatable,dimension(:,:) :: A,adjoint
integer :: i,j,norb
allocate(adjoint(1:norb,1:norb))
do i=1,norb
do j=1,norb
adjoint(j,i) = DCONJG(A(i,j))
end do
end do
end function
subroutine Get_HF_QChem(inputfile,norb,H_two,Smat)
!Use InterfaceMod
implicit none
character(len=100) :: inputfile,datafile
real(8),allocatable,dimension(:,:) :: H_Two,Smat
integer :: norb,ioerror,i,j
20 format(A)
!write(*,*) inputfile
allocate(H_two(1:norb,1:norb))
allocate(Smat(1:norb,1:norb))
H_two=0
Smat=0
datafile = trim(inputfile) // "_Smat"
open(unit=2,file=datafile,action='READ', iostat = ioerror)
do i=1,norb
do j=1,norb
read(2,*) Smat(i,j)
end do
end do
close(2)
write(*,*) 'Start Htwo'
datafile = trim(inputfile) // "_Htwo"
open(unit=3,file=datafile,action='READ', iostat = ioerror)
do i=1,norb
do j=1,norb
read(3,*) H_Two(i,j)
end do
end do
close(3)
write(*,*) 'Done Getting HF Values'
end subroutine
subroutine Get_HF_Molcas(inputfile,norb,H_two,Smat)
!Use InterfaceMod
implicit none
character(len=100) :: inputfile,datafile
real(8),allocatable,dimension(:,:) :: H_Two,Smat
integer :: norb,ioerror,i,j,x,y
real(8) :: readtemp
20 format(A)
!write(*,*) inputfile
allocate(H_two(1:norb,1:norb))
allocate(Smat(1:norb,1:norb))
H_two=0
Smat=0
datafile = "Overlap"
open(unit=2,file=datafile,action='READ', iostat = ioerror)
do i=1,norb
do j=i,norb
read(2,*) x,y,readtemp
Smat(x,y)=readtemp
Smat(y,x)=readtemp
end do
end do
close(2)
write(*,*) 'Start Fock Matrix'
datafile = "FOCK_AO"
open(unit=3,file=datafile,action='READ', iostat = ioerror)
do i=1,norb
do j=1,norb
read(3,*) x,y,readtemp
H_Two(x,y)=readtemp
end do
end do
close(3)
write(*,*) 'Done Getting HF Values'
end subroutine
Subroutine Get_HF_GAMESS(inputfile,numatomic,H_two,Smat,norb)
use FunctionMod
implicit none
character(len=100) :: inputfile,datafile
real(8),allocatable,dimension(:,:) :: H_Two,Smat,mo_data,mo_data2,ECP_m,ECP_a,ECP_temp
integer :: numatomic,ioerror,i,j,norb
20 format(A)
allocate(H_two(1:numatomic,1:numatomic))
allocate(Smat(1:numatomic,1:numatomic))
H_two=0
Smat=0
datafile = trim(inputfile) // "_Smat"
open(unit=2,file=datafile,action='READ', iostat = ioerror)
do i=1,numatomic
do j=1,numatomic
read(2,*) Smat(i,j)
end do
end do
close(2)
write(*,*) 'Start Htwo'
datafile = trim(inputfile) // "_Htwo"
open(unit=3,file=datafile,action='READ', iostat = ioerror)
do i=1,numatomic
do j=1,numatomic
read(3,*) H_Two(i,j)
end do
end do
close(3)
write(*,*) 'Get ECP'
datafile = trim(inputfile) // "_ecp"
open(unit=4,file=datafile,status='OLD',action='READ', iostat = ioerror)
if(ioerror.eq.0) then
write(*,*) 'ECP found'
allocate(ECP_m(1:norb,1:norb))
allocate(ECP_temp(1:numatomic,1:norb))
allocate(ECP_a(1:numatomic,1:numatomic))
allocate(mo_data(1:numatomic,1:norb))
allocate(mo_data2(1:norb,1:numatomic))
do i=1,norb
do j=1,norb
read(4,*) ECP_m(i,j)
end do
end do
close(5)
datafile = trim(inputfile) // ".mo_dat"
open(unit=6,file=datafile,action='READ', iostat = ioerror)
do i=1,numatomic
do j=1,norb
read(6,*) mo_data(i,j)
end do
end do
close(6)
mo_data2 = transpose(mo_data)
ECP_temp=matmul_dgemm(mo_data,ECP_m)
ECP_a=matmul_dgemm(ECP_temp,mo_data2)
H_two=H_two-ECP_a
else
write(*,*) "No removal of ECP"
end if
write(*,*) 'H_Two',size(H_Two)
write(*,*) 'Smat',size(Smat)
!stop
write(*,*) 'Done Getting HF Values from GAMESS'
end subroutine
Subroutine Get_HF_PySCF(inputfile,numatomic,H_two,Smat,norb)
!Use InterfaceMod
use FunctionMod
implicit none
character(len=100) :: inputfile,datafile,readtemp
real(8),allocatable,dimension(:,:) :: H_Two,Smat,mo_data,mo_data2,ECP_m,ECP_a,ECP_temp
integer :: numatomic,ioerror,i,j,norb,tempx,tempy
real(8) :: tempval
20 format(A)
allocate(H_two(1:numatomic,1:numatomic))
allocate(Smat(1:numatomic,1:numatomic))
H_two=0
Smat=0
datafile = trim(inputfile) // ".scf_dat"
open(unit=2,file=datafile,action='READ', iostat = ioerror)
do while(trim(readtemp)/="Overlap Matrix")
read(2,'(A)') readtemp
end do
do i=1,numatomic
do j=1,numatomic
read(2,*) tempx,tempy,tempval
Smat(tempx,tempy)=tempval
end do
end do
do while(trim(readtemp)/="Fock Matrix")
read(2,'(A)') readtemp
end do
do i=1,numatomic
do j=1,numatomic
read(2,*) tempx,tempy,tempval
H_Two(tempx,tempy)=tempval
end do
end do
close(2)
end subroutine
subroutine Get_HF_libint(inputfile,norb,numact,H_one,Smat,mo_coeff,OneInts,TwoIntsCompact)
use TypeMod
use FunctionMod
implicit none
character(len=100) :: inputfile,datafile
real(8),allocatable,dimension(:,:) :: H_one,H_Two,Smat,OneInts,mo_coeff
real(8),allocatable,dimension(:) :: TwoIntsCompact
integer :: norb,numact,ioerror,i,j,k,l
integer(8) :: index1,index2,compindex1
integer, dimension(1:4) :: orbind
real(8) :: integral
allocate(H_one(1:norb,1:norb))
allocate(Smat(1:norb,1:norb))
allocate(mo_coeff(1:norb,1:norb))
allocate(OneInts(1:numact,1:numact))
allocate(TwoIntsCompact(1:numact*(numact+1)/2*(numact*(numact+1)/2+1)/2))
20 format(A)
datafile = inputfile // ".Hone"
open(unit=1,file=datafile,action='READ', iostat = ioerror)
do i=1,norb
do j=1,norb
read(1,20) H_one(i,j)
end do
end do
close(1)
datafile = inputfile // ".Smat"
open(unit=2,file=datafile,action='READ', iostat = ioerror)
do i=1,norb
do j=1,norb
read(2,20) Smat(i,j)
end do
end do
close(2)
open(unit=6,file=datafile,action='READ',iostat= ioerror)
do i=1,norb
do j=1,norb
read(6,20) mo_coeff(i,j)
end do
end do
datafile = inputfile // ".OneInts"
open(unit=4,file=datafile,action='READ', iostat = ioerror)
do i=1,norb
do j=1,norb
read(4,20) OneInts(i,j)
end do
end do
close(4)
datafile = inputfile // ".TwoInts"
open(unit=5,file=datafile,action='READ', iostat = ioerror)
do while(orbind(1).ne.0)
read(5,*) orbind(1),orbind(2),orbind(3),orbind(4),integral
i = orbind(1)
k = orbind(2)
j = orbind(3)
l = orbind(4)
index1 = FirstIndex(i,k)
index2 = FirstIndex(j,l)
compindex1 = CompositeIndex(index1,index2)
TwoIntsCompact(compindex1) = integral
end do
end subroutine
subroutine IntTransform(TwoIntsCompact,mo_coeff)
use TypeMod
implicit none
real(8), allocatable,dimension(:,:) :: mo_coeff
real(8), allocatable,dimension(:) :: TwoIntsCompact
!There is nothing here yet!
end subroutine
subroutine Calculate_Coupling_MoleculeWBL(Coupling_R,Coupling_L,localden_fermi)
implicit none
real(8) :: Coupling_R,Coupling_L,localden_fermi
Coupling_R = 0.2
Coupling_L = 0.2
end subroutine
subroutine Electrodes_MoleculeWBL(Sigma_l,Sigma_r,Smat,Coupling_R,Coupling_L,size_c,size_lc,size_lcr)
implicit none
integer :: size_lc,size_c,size_lcr,i,j
real(8) :: Coupling_R,Coupling_L,temp
complex(8), allocatable, dimension(:,:) :: Sigma_L,Sigma_R
real(8), allocatable, dimension(:,:) :: Smat_LE,Smat_RE,Smat
do i = size_c+1,size_lc
do j = 1,size_c
Sigma_L(i,j) = CMPLX(0,-0.5*Coupling_L*Smat(i,j),8)
Sigma_L(j,i) = CMPLX(0,-0.5*Coupling_L*Smat(j,i),8)
end do
end do
do i=size_lc+1,size_lcr
do j= 1,size_c
Sigma_R(i,j) = CMPLX(0,-0.5*Coupling_R*Smat(i,j),8)
Sigma_R(j,i) = CMPLX(0,-0.5*Coupling_R*Smat(j,i),8)
end do
end do
end subroutine
subroutine Electrodes_MetalWBL(Sigma_l,Sigma_r,Smat_re,Smat_le,H_Two_le,H_Two_re,localden_fermi_l,localden_fermi_r,size_c,size_lc,size_lcr,size_l,size_r,H_Two_le_trans,H_Two_re_trans)
use FunctionMod
implicit none
integer :: size_lc,size_c,size_lcr,i,j,size_l,size_r
real(8) :: Coupling_R,Coupling_L,temp,localden_fermi_l,localden_fermi_r,pi
complex(8), allocatable, dimension(:,:) :: Sigma_L,Sigma_R
real(8), allocatable, dimension(:,:) :: Smat_LE,Smat_RE,Smat,H_Two_re,H_Two_le,Sigma_L_temp,Sigma_R_temp,Sigma_R_temp2,Sigma_L_temp2,sigma_r_temp3,sigma_l_temp3,H_Two_le_trans,H_Two_re_trans
pi = 3.14159265359
allocate(sigma_l_temp(1:size_l,1:size_l))
allocate(sigma_l_temp2(1:size_c,1:size_l))
allocate(sigma_l_temp3(1:size_c,1:size_c))
write(*,*) 'alloc l done'
Sigma_L_temp = -pi*localden_fermi_l*Smat_le
call matmul_dgemm2(H_Two_le_trans,Sigma_L_temp,Sigma_L_temp2)
call matmul_dgemm2(Sigma_L_temp2,H_Two_le,Sigma_L_temp3)
Sigma_L = CMPLX(0,sigma_l_temp3,8)
write(*,*) 'l done'
deallocate(sigma_l_temp)
deallocate(sigma_l_temp2)
deallocate(sigma_l_temp3)
allocate(sigma_r_temp(1:size_r,1:size_r))
allocate(sigma_r_temp2(1:size_c,1:size_r))
allocate(sigma_r_temp3(1:size_c,1:size_c))
Sigma_R_temp = -pi*localden_fermi_r*Smat_re
call matmul_dgemm2(H_Two_re_trans,Sigma_R_temp,Sigma_R_temp2)
call matmul_dgemm2(Sigma_R_temp2,H_Two_re,Sigma_R_temp3)
Sigma_R = CMPLX(0,Sigma_R_temp3,8)
deallocate(Sigma_R_temp)
deallocate(Sigma_R_temp2)
deallocate(Sigma_R_temp3)
write(*,*) 'Sigma Calculated'
end subroutine
subroutine PartitionHS_MetalWBL(Smat,H_Two,size_l,size_r,size_c,size_lc,size_lcr,Smat_le,Smat_re,Smat_cen,H_Two_le,H_Two_re,H_Two_cen,H_Two_le_trans,H_Two_re_trans)
implicit none
real(8), allocatable, dimension(:,:) :: Smat,Smat_le,Smat_re,Smat_cen,H_Two,H_Two_cen,H_Two_re,H_Two_le,temp,H_Two_le_trans,H_Two_re_trans
real(8), allocatable, dimension(:) :: eigen,work
integer :: size_l,size_c,size_r,size_lc,size_lcr,info,lwork
allocate(Smat_le(1:size_l,1:size_l))
allocate(Smat_re(1:size_r,1:size_r))
allocate(Smat_cen(1:size_c,1:size_c))
Smat_le(1:size_l,1:size_l)=Smat(1:size_l,1:size_l)
Smat_re(1:size_r,1:size_r)=Smat(size_lc+1:size_lcr,size_lc+1:size_lcr)
Smat_cen(1:size_c,1:size_c)=Smat(size_l+1:size_lc,size_l+1:size_lc)
write(*,*) 'Done partitioning S Matrix'
deallocate(Smat)
allocate(H_Two_le(1:size_l,1:size_c))
allocate(H_Two_re(1:size_r,1:size_c))
allocate(H_Two_cen(1:size_c,1:size_c))
allocate(H_Two_le_trans(1:size_c,1:size_l))
allocate(H_Two_re_trans(1:size_c,1:size_r))
H_Two_le(1:size_l,1:size_c)=27.2114*H_Two(1:size_l,size_l+1:size_lc)
H_Two_re(1:size_r,1:size_c)=27.2114*H_Two(size_lc+1:size_lcr,size_l+1:size_lc)
H_Two_cen(1:size_c,1:size_c)=27.2114*H_Two(size_l+1:size_lc,size_l+1:size_lc)
H_Two_le_trans(1:size_c,1:size_l) = 27.2114*H_Two(size_l+1:size_lc,1:size_l) !transpose(H_Two_le)
H_Two_re_trans(1:size_c,1:size_r) = 27.2114*H_Two(size_l+1:size_lc,size_lc+1:size_lcr)!transpose(H_Two_re)
write(*,*) 'Done partitioning Fock Matrix'
deallocate(H_Two)
end subroutine
function fermi_function(energy,fermi_energy,KT)
implicit none
real(8) :: energydiff,KT,energy,fermi_energy,fermi_function
energydiff = energy - fermi_energy
fermi_function = 1.00/(exp(energydiff/KT)+1)
end function
function inv(A) result(Ainv)
implicit none
complex(8),allocatable, dimension(:,:), intent(in) :: A
complex(8), allocatable, dimension(:,:) :: Ainv
complex(8),allocatable, dimension(:) :: work ! work array for LAPACK
integer,allocatable, dimension(:) :: ipiv ! pivot indices
integer :: n, info
allocate(Ainv(1:size(A,1),1:size(A,2)))
allocate(work(1:size(A,1)))
allocate(ipiv(1:size(A,1)))
Ainv = A
n = size(A,1)
call ZGETRF(n, n, Ainv, n, ipiv, info)
if (info /= 0) then
write(*,*) info!,A
stop 'Matrix is numerically singular!'
end if
call ZGETRI(n, Ainv, n, ipiv, work, n, info)
if (info /= 0) then
write(*,*) info
stop 'Matrix inversion failed!'
end if
deallocate(work)
deallocate(ipiv)
end function inv
function inv_real(A) result(Ainv)
implicit none
real(8),allocatable, dimension(:,:), intent(in) :: A
real(8), allocatable, dimension(:,:) :: Ainv
real(8),allocatable, dimension(:) :: work ! work array for LAPACK
integer,allocatable, dimension(:) :: ipiv ! pivot indices
integer :: n, info
allocate(Ainv(1:size(A,1),1:size(A,2)))
allocate(work(1:size(A,1)))
allocate(ipiv(1:size(A,1)))
Ainv = A
n = size(A,1)
call DGETRF(n, n, Ainv, n, ipiv, info)
if (info /= 0) then
write(*,*) info!,A
stop 'Matrix is numerically singular!'
end if
call DGETRI(n, Ainv, n, ipiv, work, n, info)
if (info /= 0) then
write(*,*) info
stop 'Matrix inversion failed!'
end if
deallocate(work)
deallocate(ipiv)
end function inv_real
function FirstIndex(i,k)
Implicit None
integer :: i,k
integer(8) :: FirstIndex
if(i.lt.k) then
FirstIndex = (k-1)*k/2 + i
else
FirstIndex = (i-1)*i/2 + k
end if
end function FirstIndex
function matmul_zgemm(leftmatrix,rightmatrix)
Implicit NONE
complex(8), allocatable, dimension(:,:) :: leftmatrix,rightmatrix
complex(8), allocatable, dimension(:,:) :: matmul_zgemm
integer :: lmr, lmc
integer :: rmr, rmc
lmr = size(leftmatrix,1) !left matrix row size
lmc = size(leftmatrix,2) !left matrix col size
rmr = size(rightmatrix,1) !right matrix row size
rmc = size(rightmatrix,2) !right matrix col size
allocate(matmul_zgemm(1:lmr,1:rmc))
matmul_zgemm(1:lmr,1:rmc) = 0.0
if(lmr.ne.0.and.lmc.ne.0.and.rmr.ne.0.and.rmc.ne.0) then
Call ZGEMM('N','N',lmr,rmc,lmc,1.d0,leftmatrix,lmr,rightmatrix,rmr,0.d0,matmul_zgemm,lmr)
end if
end function matmul_zgemm
function matmul_dgemm(leftmatrix,rightmatrix)
Implicit NONE
real(8), allocatable, dimension(:,:) :: leftmatrix,rightmatrix
real(8), allocatable, dimension(:,:) :: matmul_dgemm
integer :: lmr, lmc
integer :: rmr, rmc
lmr = size(leftmatrix,1) !left matrix row size
lmc = size(leftmatrix,2) !left matrix col size
rmr = size(rightmatrix,1) !right matrix row size
rmc = size(rightmatrix,2) !right matrix col size
matmul_dgemm = 0.0
if(lmr.ne.0.and.lmc.ne.0.and.rmr.ne.0.and.rmc.ne.0) then
Call DGEMM('N','N',lmr,rmc,lmc,1.d0,leftmatrix,lmr,rightmatrix,rmr,0.d0,matmul_dgemm,lmr)
end if
end function matmul_dgemm
subroutine matmul_dgemm2(leftmatrix,rightmatrix,outmat)
Implicit NONE
real(8), allocatable, dimension(:,:) :: leftmatrix,rightmatrix
real(8), allocatable, dimension(:,:) :: outmat
integer :: lmr, lmc,cmc
integer :: rmr, rmc,cmr
lmr = size(leftmatrix,1) !left matrix row size
lmc = size(leftmatrix,2) !left matrix col size
rmr = size(rightmatrix,1) !right matrix row size
rmc = size(rightmatrix,2) !right matrix col size
cmr= size(outmat,1)
cmc=size(outmat,2)
outmat = 0.0
if(lmr.ne.0.and.lmc.ne.0.and.rmr.ne.0.and.rmc.ne.0) then
Call DGEMM("N","N",lmr,rmc,lmc,1.d0,leftmatrix,lmr,rightmatrix,rmr,0.d0,outmat,cmr)
end if
end subroutine matmul_dgemm2
function CompositeIndex(ik,jl)
Implicit NONE
integer(8) :: ik, jl
integer(8) :: CompositeIndex
CompositeIndex = 0
if(ik.lt.jl) then
CompositeIndex = (jl-1)*jl/2 +ik
else
CompositeIndex = (ik-1)*ik/2 + jl
end if
end function CompositeIndex
subroutine ReadInput(inputfile,norb,numfcore,numfvirt,numocc,numvirt,size_l,size_r,size_c,energy_start,energy_end,delta_en,volt_start,volt_end,delta_volt,inputcode,KT,Electrode_Type,Fermi_enl,Fermi_enr,CalcType,localden_fermi_l,localden_fermi_r,doubles,numatomic,functional,num_threads,use_b0,b0_type)
implicit none
character(len=100) :: inputfile
character(len=40) :: inputcode,filename,Electrode_Type,CalcType,functional,b0_type
integer :: norb, size_c,size_r,size_l,numfcore,numfvirt,numocc,numvirt,numatomic,num_threads
real(8) :: energy_start,energy_end,delta_en,volt_start,volt_end,delta_volt,KT
logical :: libint,doubles,use_b0
real(8) :: Fermi_enl,Fermi_enr,localden_fermi_l,localden_fermi_r
filename = trim(inputfile)
open(unit=1,file=filename,action="read")
20 format(A)
read(1,20) CalcType
read(1,20) Electrode_Type
read(1,*) Fermi_enl
read(1,*) Fermi_enr
read(1,*) localden_fermi_l
read(1,*) localden_fermi_r
read(1,*) norb
read(1,*) numatomic
read(1,*) numfcore
read(1,*) numfvirt
read(1,*) numocc
read(1,*) numvirt
read(1,*) size_c
read(1,*) size_l
read(1,*) size_r
read(1,*) energy_start
read(1,*) energy_end
read(1,*) delta_en
read(1,*) volt_start
read(1,*) volt_end
read(1,*) delta_volt
read(1,*) KT
read(1,*) inputcode
read(1,*) doubles
read(1,*) functional
read(1,*) use_b0