-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathED_OBSERVABLES.f90
1128 lines (1052 loc) · 41.5 KB
/
ED_OBSERVABLES.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
MODULE ED_OBSERVABLES
USE SF_CONSTANTS, only:zero,pi,xi
USE SF_IOTOOLS, only:free_unit,reg,txtfy
USE SF_ARRAYS, only: arange
USE SF_INTEGRATE, only: quad
USE SF_LINALG
USE SF_TIMER
USE ED_INPUT_VARS
USE ED_VARS_GLOBAL
USE ED_EIGENSPACE
USE ED_SETUP
USE ED_HAMILTONIAN
USE ED_BATH
USE ED_AUX_FUNX
USE ED_IO
USE ED_BATH_FUNCTIONS
implicit none
private
!
interface add_custom_observable
module procedure :: add_custom_observable_local
module procedure :: add_custom_observable_kdep
end interface add_custom_observable
public :: observables_impurity
public :: local_energy_impurity
public :: density_matrix_impurity
!
public :: init_custom_observables
public :: add_custom_observable
public :: get_custom_observables
public :: clear_custom_observables
logical,save :: iolegend=.true.
real(8),dimension(:,:),allocatable :: dens,dens_up,dens_dw
real(8),dimension(:,:),allocatable :: docc
real(8),dimension(:,:),allocatable :: magz
real(8),dimension(:,:,:,:),allocatable :: sz2,n2
real(8),dimension(:),allocatable :: s2tot
real(8) :: Egs
real(8) :: Ei
real(8) :: integrationR
complex(8),allocatable,dimension(:,:,:) :: sij
complex(8),allocatable,dimension(:,:,:) :: Hk
!
integer :: iorb,jorb,iorb1,jorb1
integer :: ispin,jspin
integer :: ilat,jlat
integer :: ibath,jbath
integer :: r,m,k,k1,k2
integer :: iup,idw
integer :: jup,jdw
integer :: mup,mdw
real(8) :: sgn,sgn1,sgn2,sg1,sg2
real(8) :: gs_weight
!
real(8) :: peso
real(8) :: norm
!
integer :: i,j,ii
integer :: isector,jsector
integer :: idim,idimUP,idimDW
!
complex(8),dimension(:),allocatable :: state_cvec
logical :: Jcondition
!
contains
!+-------------------------------------------------------------------+
!PURPOSE : Evaluate and print out many interesting physical qties
!+-------------------------------------------------------------------+
subroutine observables_impurity()
write(LOGfile,"(A)")" "
write(LOGfile,"(A)")"Get observables:"
call lanc_observables()
end subroutine observables_impurity
subroutine local_energy_impurity()
write(LOGfile,"(A)")" "
write(LOGfile,"(A)")"Get local energy:"
call lanc_local_energy()
end subroutine local_energy_impurity
!+-------------------------------------------------------------------+
!PURPOSE : Lanc method
!+-------------------------------------------------------------------+
subroutine lanc_observables()
integer :: istate,Nud(2,Ns),iud(2),jud(2),is,js
integer,dimension(2*Ns_Ud) :: Indices
integer,dimension(Ns_Ud,Ns_Orb) :: Nups,Ndws ![1,Ns]-[Norb,1+Nbath]
integer,dimension(Ns_Ud) :: iDimUps,iDimDws
integer,dimension(Ns) :: IbUp,IbDw ![Ns]
real(8),dimension(Nlat,Norb) :: nup,ndw,Sz,nt
type(sector_map),dimension(2*Ns_Ud) :: HI
!
!LOCAL OBSERVABLES:
allocate(dens(Nlat,Norb),dens_up(Nlat,Norb),dens_dw(Nlat,Norb))
allocate(docc(Nlat,Norb),s2tot(Nlat))
allocate(magz(Nlat,Norb),sz2(Nlat,Nlat,Norb,Norb),n2(Nlat,Nlat,Norb,Norb))
!
Egs = state_list%emin
dens = 0.d0
dens_up = 0.d0
dens_dw = 0.d0
docc = 0.d0
magz = 0.d0
sz2 = 0.d0
n2 = 0.d0
s2tot = 0.d0
!
do istate=1,state_list%size
isector = es_return_sector(state_list,istate)
Ei = es_return_energy(state_list,istate)
!
#ifdef _MPI
if(MpiStatus)then
call es_return_cvector(MpiComm,state_list,istate,state_cvec)
else
call es_return_cvector(state_list,istate,state_cvec)
endif
#else
call es_return_cvector(state_list,istate,state_cvec)
#endif
!
!
peso = 1.d0 ; if(finiteT)peso=exp(-beta*(Ei-Egs))
peso = peso/zeta_function
!
iDim = getdim(isector)
call get_DimUp(isector,iDimUps)
call get_DimDw(isector,iDimDws)
iDimUp = product(iDimUps)
iDimDw = product(iDimDws)
!
if(MpiMaster)then
call build_sector(isector,HI)
do i = 1,iDim
call state2indices(i,[iDimUps,iDimDws],Indices)
do ii=1,Ns_Ud
mup = HI(ii)%map(Indices(ii))
mdw = HI(ii+Ns_Ud)%map(Indices(ii+Ns_ud))
IbUp = Bdecomp(mup,Ns_Orb) ![Ns_Orb = Ns = Nlat*Norb*(Nbath+1) in CDMFT code]
IbDw = Bdecomp(mdw,Ns_Orb)
enddo
!
gs_weight=peso*abs(state_cvec(i))**2
!
!Get operators:
do ilat=1,Nlat
do iorb=1,Norb
nup(ilat,iorb)= ibup(imp_state_index(ilat,iorb))
ndw(ilat,iorb)= ibdw(imp_state_index(ilat,iorb))
sz(ilat,iorb) = (nup(ilat,iorb) - ndw(ilat,iorb))/2d0
nt(ilat,iorb) = nup(ilat,iorb) + ndw(ilat,iorb)
enddo
enddo
!
!
!Evaluate averages of observables:
!>TODO:
!add non-local averges like spin-spin, density-density etc...
!<TODO
do ilat=1,Nlat
do iorb=1,Norb
dens(ilat,iorb) = dens(ilat,iorb) + nt(ilat,iorb)*gs_weight
dens_up(ilat,iorb) = dens_up(ilat,iorb) + nup(ilat,iorb)*gs_weight
dens_dw(ilat,iorb) = dens_dw(ilat,iorb) + ndw(ilat,iorb)*gs_weight
docc(ilat,iorb) = docc(ilat,iorb) + nup(ilat,iorb)*ndw(ilat,iorb)*gs_weight
magz(ilat,iorb) = magz(ilat,iorb) + (nup(ilat,iorb)-ndw(ilat,iorb))*gs_weight
enddo
s2tot(ilat) = s2tot(ilat) + (sum(sz(ilat,:)))**2*gs_weight
enddo
!
!
do ilat=1,Nlat
do iorb=1,Norb
sz2(ilat,ilat,iorb,iorb) = sz2(ilat,ilat,iorb,iorb) + (sz(ilat,iorb)*sz(ilat,iorb))*gs_weight
n2(ilat,ilat,iorb,iorb) = n2(ilat,ilat,iorb,iorb) + (nt(ilat,iorb)*nt(ilat,iorb))*gs_weight
do jlat=1,Nlat
do jorb=iorb+1,Norb
sz2(ilat,jlat,iorb,jorb) = sz2(ilat,jlat,iorb,jorb) + (sz(ilat,iorb)*sz(jlat,jorb))*gs_weight
sz2(ilat,jlat,jorb,iorb) = sz2(ilat,jlat,jorb,iorb) + (sz(ilat,jorb)*sz(jlat,iorb))*gs_weight
n2(ilat,jlat,iorb,jorb) = n2(ilat,jlat,iorb,jorb) + (nt(ilat,iorb)*nt(jlat,jorb))*gs_weight
n2(ilat,jlat,jorb,iorb) = n2(ilat,jlat,jorb,iorb) + (nt(ilat,jorb)*nt(jlat,iorb))*gs_weight
enddo
enddo
enddo
enddo
enddo
call delete_sector(isector,HI)
endif
!
!
if(allocated(state_cvec))deallocate(state_cvec)
!
enddo
!
!
if(MPIMASTER)then
if(iolegend)call write_legend
call write_observables()
!
do ilat=1,Nlat
write(LOGfile,"(A,10f18.12,f18.12,A)")"dens site "//str(ilat)//" "//reg(ed_file_suffix)//"=",(dens(ilat,iorb),iorb=1,Norb),sum(dens(ilat,:))
enddo
write(LOGfile,"(A,10f18.12,f18.12,A)")"dens avg "//reg(ed_file_suffix)//" =",(sum(dens(:,iorb))/Nlat,iorb=1,Norb),sum(dens)/Nlat
!
write(LOGfile,"(A,10f18.12,A)")"docc "//reg(ed_file_suffix)//" =",(sum(docc(:,iorb))/Nlat,iorb=1,Norb)
if(Nspin==2)write(LOGfile,"(A,10f18.12,A)") "mag "//reg(ed_file_suffix)//"=",(sum(magz(:,iorb))/Nlat,iorb=1,Norb)
!
write(LOGfile,"(A)")" "
!
ed_dens_up=dens_up
ed_dens_dw=dens_dw
ed_dens =dens
ed_docc =docc
endif
#ifdef _MPI
if(MpiStatus)then
call Bcast_MPI(MpiComm,ed_dens_up)
call Bcast_MPI(MpiComm,ed_dens_dw)
call Bcast_MPI(MpiComm,ed_dens)
call Bcast_MPI(MpiComm,ed_docc)
endif
#endif
!
deallocate(dens,docc,dens_up,dens_dw,magz,sz2,n2,s2tot)
end subroutine lanc_observables
!+-------------------------------------------------------------------+
!PURPOSE : Get internal energy from the Impurity problem.
!+-------------------------------------------------------------------+
subroutine lanc_local_energy()
integer :: istate,iud(2),jud(2),is,js
integer,dimension(2*Ns_Ud) :: Indices
integer,dimension(Ns_Ud) :: iDimUps,iDimDws
integer,dimension(Ns_Ud,Ns_Orb) :: Nups,Ndws ![1,Ns]-[Norb,1+Nbath]
real(8),dimension(Ns) :: Nup,Ndw
type(sector_map),dimension(2*Ns_Ud) :: H
!
Egs = state_list%emin
ed_Ehartree= 0.d0
ed_Eknot = 0.d0
ed_Epot = 0.d0
ed_Dust = 0.d0
ed_Dund = 0.d0
ed_Dse = 0.d0
ed_Dph = 0.d0
!
!
do istate=1,state_list%size
isector = es_return_sector(state_list,istate)
Ei = es_return_energy(state_list,istate)
#ifdef _MPI
if(MpiStatus)then
call es_return_cvector(MpiComm,state_list,istate,state_cvec)
else
call es_return_cvector(state_list,istate,state_cvec)
endif
#else
call es_return_cvector(state_list,istate,state_cvec)
#endif
!
iDim = getdim(isector)
call get_DimUp(isector,iDimUps)
call get_DimDw(isector,iDimDws)
iDimUp = product(iDimUps)
iDimDw = product(iDimDws)
!
peso = 1.d0 ; if(finiteT)peso=exp(-beta*(Ei-Egs))
peso = peso/zeta_function
!
!Master:
if(MpiMaster)then
call build_sector(isector,H)
do i=1,iDim
call state2indices(i,[iDimUps,iDimDws],Indices)
do ii=1,Ns_Ud
mup = H(ii)%map(Indices(ii))
mdw = H(ii+Ns_Ud)%map(Indices(ii+Ns_ud))
Nups(ii,:) = Bdecomp(mup,Ns_Orb) ![Ns_Orb = Ns = Nlat*Norb*(Nbath+1) in CDMFT code]
Ndws(ii,:) = Bdecomp(mdw,Ns_Orb)
enddo
Nup = Breorder(Nups) !Actually, they are already reordered in CDMFT code...
Ndw = Breorder(Ndws) !...and breorder() corresponds to an identity: look in ED_SETUP!
!
gs_weight=peso*abs(state_cvec(i))**2
!
!> H_Imp: Diagonal Elements, i.e. local part
do ilat=1,Nlat
do iorb=1,Norb
is = imp_state_index(ilat,iorb)
ed_Eknot = ed_Eknot + impHloc(ilat,ilat,1,1,iorb,iorb)*Nup(is)*gs_weight
ed_Eknot = ed_Eknot + impHloc(ilat,ilat,Nspin,Nspin,iorb,iorb)*Ndw(is)*gs_weight
enddo
enddo
! !> H_imp: Off-diagonal elements, i.e. non-local part.
iup = Indices(1) ; idw = Indices(2)
mup = H(1)%map(iup) ; mdw = H(2)%map(idw)
do ilat=1,Nlat
do jlat=1,Nlat
do iorb=1,Norb
do jorb=1,Norb
is = imp_state_index(ilat,iorb)
js = imp_state_index(jlat,jorb)
!UP
Jcondition = &
(impHloc(ilat,jlat,1,1,iorb,jorb)/=0d0) .AND. &
(Nup(js)==1) .AND. (Nup(is)==0)
if (Jcondition) then
call c(js,mup,k1,sg1)
call cdg(is,k1,k2,sg2)
jup = binary_search(H(1)%map,k2)
j = jup + (idw-1)*iDimUp
ed_Eknot = ed_Eknot + &
impHloc(ilat,jlat,1,1,iorb,jorb)*sg1*sg2*state_cvec(i)*conjg(state_cvec(j))*peso
endif
!
!DW
Jcondition = &
(impHloc(ilat,jlat,Nspin,Nspin,iorb,jorb)/=0d0) .AND. &
(ndw(js)==1) .AND. (ndw(is)==0)
if (Jcondition) then
call c(js,mdw,k1,sg1)
call cdg(is,k1,k2,sg2)
jdw = binary_search(H(2)%map,k2)
j = iup + (jdw-1)*iDimUp
ed_Eknot = ed_Eknot + &
impHloc(ilat,jlat,Nspin,Nspin,iorb,jorb)*sg1*sg2*state_cvec(i)*conjg(state_cvec(j))*peso
endif
enddo
enddo
enddo
enddo
!
!
!DENSITY-DENSITY INTERACTION: SAME ORBITAL, OPPOSITE SPINS
!Euloc=\sum=i U_i*(n_u*n_d)_i
!ed_Epot = ed_Epot + dot_product(uloc,nup*ndw)*gs_weight
do ilat=1,Nlat
do iorb=1,Norb
is = imp_state_index(ilat,iorb)
ed_Epot = ed_Epot + Uloc(iorb)*nup(is)*ndw(is)*gs_weight
enddo
enddo
!
!DENSITY-DENSITY INTERACTION: DIFFERENT ORBITALS, OPPOSITE SPINS
!Eust=\sum_ij Ust*(n_up_i*n_dn_j + n_up_j*n_dn_i)
! "="\sum_ij (Uloc - 2*Jh)*(n_up_i*n_dn_j + n_up_j*n_dn_i)
if(Norb>1)then
do ilat=1,Nlat
do iorb=1,Norb
do jorb=iorb+1,Norb
is = imp_state_index(ilat,iorb)
js = imp_state_index(ilat,jorb)
ed_Epot = ed_Epot + Ust*(nup(is)*ndw(js) + nup(js)*ndw(is))*gs_weight
ed_Dust = ed_Dust + (nup(is)*ndw(js) + nup(js)*ndw(is))*gs_weight
enddo
enddo
enddo
endif
!
!DENSITY-DENSITY INTERACTION: DIFFERENT ORBITALS, PARALLEL SPINS
!Eund = \sum_ij Und*(n_up_i*n_up_j + n_dn_i*n_dn_j)
! "="\sum_ij (Ust-Jh)*(n_up_i*n_up_j + n_dn_i*n_dn_j)
! "="\sum_ij (Uloc-3*Jh)*(n_up_i*n_up_j + n_dn_i*n_dn_j)
if(Norb>1)then
do ilat=1,Nlat
do iorb=1,Norb
do jorb=iorb+1,Norb
is = imp_state_index(ilat,iorb)
js = imp_state_index(ilat,jorb)
ed_Epot = ed_Epot + (Ust-Jh)*(nup(is)*nup(js) + ndw(is)*ndw(js))*gs_weight
ed_Dund = ed_Dund + (nup(is)*nup(js) + ndw(is)*ndw(js))*gs_weight
enddo
enddo
enddo
endif
!
!HARTREE-TERMS CONTRIBUTION:
if(hfmode)then
!ed_Ehartree=ed_Ehartree - 0.5d0*dot_product(uloc,nup+ndw)*gs_weight + 0.25d0*sum(uloc)*gs_weight
do ilat=1,Nlat
do iorb=1,Norb
is =imp_state_index(ilat,iorb)
ed_Ehartree=ed_Ehartree - 0.5d0*uloc(iorb)*(nup(is)+ndw(is))*gs_weight + 0.25d0*uloc(is)*gs_weight
enddo
enddo
if(Norb>1)then
do ilat=1,Nlat
do iorb=1,Norb
do jorb=iorb+1,Norb
is=imp_state_index(ilat,iorb)
js=imp_state_index(ilat,jorb)
ed_Ehartree=ed_Ehartree - 0.5d0*Ust*(nup(is)+ndw(is)+nup(js)+ndw(js))*gs_weight + 0.25d0*Ust*gs_weight
ed_Ehartree=ed_Ehartree - 0.5d0*(Ust-Jh)*(nup(is)+ndw(is)+nup(js)+ndw(js))*gs_weight + 0.25d0*(Ust-Jh)*gs_weight
enddo
enddo
enddo
endif
endif
enddo
call delete_sector(isector,H)
endif
!
if(allocated(state_cvec))deallocate(state_cvec)
!
enddo
!
!
#ifdef _MPI
if(MpiStatus)then
call Bcast_MPI(MpiComm,ed_Epot)
call Bcast_MPI(MpiComm,ed_Eknot)
call Bcast_MPI(MpiComm,ed_Ehartree)
call Bcast_MPI(MpiComm,ed_Dust)
call Bcast_MPI(MpiComm,ed_Dund)
endif
#endif
!
ed_Epot = ed_Epot + ed_Ehartree
!
if(ed_verbose>=3)then
write(LOGfile,"(A,10f18.12)")"<Hint> =",ed_Epot
write(LOGfile,"(A,10f18.12)")"<V> =",ed_Epot-ed_Ehartree
write(LOGfile,"(A,10f18.12)")"<E0> =",ed_Eknot
write(LOGfile,"(A,10f18.12)")"<Ehf> =",ed_Ehartree
write(LOGfile,"(A,10f18.12)")"Dust =",ed_Dust
write(LOGfile,"(A,10f18.12)")"Dund =",ed_Dund
write(LOGfile,"(A)")" "
endif
!
if(MPIMASTER)then
call write_energy_info()
call write_energy()
endif
!
!
end subroutine lanc_local_energy
!+-------------------------------------------------------------------+
!PURPOSE : Build the impurity density matrices
!+-------------------------------------------------------------------+
subroutine density_matrix_impurity()
integer :: istate
integer :: iUP,iDW,jUP,jDW
integer :: IimpUp,IimpDw,JimpUp,JimpDw
integer :: IbathUp,IbathDw,bUP,bDW
integer :: lenBATHup,lenBATHdw
integer,allocatable :: BATHup(:),BATHdw(:)
integer :: Nup, Ndw
integer,dimension(Ns_Ud,Ns_Orb) :: Nups,Ndws ![1,Ns]-[Norb,1+Nbath]
integer,dimension(Ns_Ud) :: iDimUps,iDimDws ![1]-[Norb]
integer,dimension(Ns) :: IbUp,IbDw ![Ns]
type(sector_map),dimension(2*Ns_Ud) :: HI ![2]-[2*Norb]
integer,dimension(2*Ns_Ud) :: Indices,Jndices ![2]-[2*Norb]
integer :: Nud(2,Ns),iud(2),jud(2),is,js,io,jo
!
! Here we build two different density matrices related to the impurity problem
! > The CLUSTER-Reduced density matrix: \rho_IMP = Tr_BATH(\rho)
! > The SINGLE-PARTICLE-Reduced density matrix: \rho_sp = <C^+_a C_b>
!
write(LOGfile,"(A)")"Get cluster density matrix:"
!CLUSTER DENSITY MATRIX [\rho_IMP = Tr_BATH(\rho)]
cluster_density_matrix=zero
!
if(MpiMaster) call start_timer()
do istate=1,state_list%size
isector = es_return_sector(state_list,istate)
Ei = es_return_energy(state_list,istate)
#ifdef _MPI
if(MpiStatus)then
call es_return_cvector(MpiComm,state_list,istate,state_cvec)
else
call es_return_cvector(state_list,istate,state_cvec)
endif
#else
call es_return_cvector(state_list,istate,state_cvec)
#endif
!Finite temperature weighting factor
peso = 1.d0 ; if(finiteT)peso=exp(-beta*(Ei-Egs))
peso = peso/zeta_function
!
!Sector sizes (per spin)
call get_DimUp(isector,iDimUps)
call get_DimDw(isector,iDimDws)
iDimUp = product(iDimUps)
iDimDw = product(iDimDws)
!
if(MpiMaster)then
call build_sector(isector,HI,itrace=.true.)
!
do IimpUp=0,2**Nimp-1
do JimpUp=0,2**Nimp-1
!
!Finding the unique bath states connecting IimpUp and JimpUp -> BATHup(:)
call sp_return_intersection(HI(1)%sp,IimpUp,JimpUp,BATHup,lenBATHup)
if(lenBATHup==0)cycle !there are no bath states intersecting IimpUp,JimpUp
!
do IimpDw=0,2**Nimp-1
do JimpDw=0,2**Nimp-1
!
!Finding the unique bath states connecting IimpDw and JimpDw -> BATHdw(:)
call sp_return_intersection(HI(2)%sp,IimpDw,JimpDw,BATHdw,lenBATHdw)
if(lenBATHdw==0)cycle !there are no bath states intersecting IimpDw,JimpDw
!------------------------------------------------------------------------------------
!Arrived here we have finally determined the {IbathUp,IbathDw} states to trace on.
!We shall use them to build back the Fock space states (IimpSIGMA,IbathSIGMA) and
!through those search the map to retrieve the formal labels i=(iUP,iDW), j=(jUP,jDW)
!------------------------------------------------------------------------------------
!=== >>> TRACE over bath states <<< =================================================
do bUP=1,lenBATHup
IbathUp = BATHup(bUP)
do bDW=1,lenBATHdw
IbathDw = BATHdw(bDW)
!-----------------------------------------------------
!Allowed spin Fock space Istates:
!Iup = IimpUp + 2^Nimp * IbathUp
!Idw = IimpDw + 2^Nimp * IbathDw
!
!Corresponding sector indices per spin:
iUP= binary_search(HI(1)%map,IimpUp + 2**Nimp*IbathUp)
iDW= binary_search(HI(2)%map,IimpDw + 2**Nimp*IbathDw)
!
!Global sector index:
i = iUP + (iDW-1)*iDimUp
!-----------------------------------------------------
!Allowed spin Fock space Jstates:
!Jup = JimpUp + 2^Nimp * IbathUp
!Jdw = JimpDw + 2^Nimp * IbathDw
!
!Corresponding sector jndices per spin:
jUP= binary_search(HI(1)%map,JimpUp + 2**Nimp*IbathUp)
jDW= binary_search(HI(2)%map,JimpDw + 2**Nimp*IbathDw)
!
!Global sector jndex:
j = jUP + (jDW-1)*iDimUp
!-----------------------------------------------------
!Final composition for the impurity:
io = (IimpUp + 2**Nimp*IimpDw) + 1
jo = (JimpUp + 2**Nimp*JimpDw) + 1
!-----------------------------------------------------------------
!(i,j)_th contribution to the (io,jo)_th element of \rho_IMP
cluster_density_matrix(io,jo) = cluster_density_matrix(io,jo) + &
state_cvec(i)*conjg(state_cvec(j))*peso
!-----------------------------------------------------------------
enddo
enddo !=============================================================================
!
enddo
enddo
!
enddo
enddo
!
call delete_sector(isector,HI)
endif
!
if(allocated(state_cvec))deallocate(state_cvec)
!
enddo
if(MpiMaster) call stop_timer()
!
!
!
!
!
!
!
!
!
!
write(LOGfile,"(A)")"Get single-particle density matrix (no print)"
!SINGLE PARTICLE DENSITY MATRIX [<C^+_a C_b>]
do istate=1,state_list%size
isector = es_return_sector(state_list,istate)
Ei = es_return_energy(state_list,istate)
#ifdef _MPI
if(MpiStatus)then
call es_return_cvector(MpiComm,state_list,istate,state_cvec)
else
call es_return_cvector(state_list,istate,state_cvec)
endif
#else
call es_return_cvector(state_list,istate,state_cvec)
#endif
!
peso = 1.d0 ; if(finiteT)peso=exp(-beta*(Ei-Egs))
peso = peso/zeta_function
!
idim = getdim(isector)
call get_DimUp(isector,iDimUps)
call get_DimDw(isector,iDimDws)
iDimUp = product(iDimUps)
iDimDw = product(iDimDws)
!
if(MpiMaster)then
call build_sector(isector,HI)
do i=1,iDim
call state2indices(i,[iDimUps,iDimDws],Indices)
do ii=1,Ns_Ud
mup = HI(ii)%map(Indices(ii))
mdw = HI(ii+Ns_Ud)%map(Indices(ii+Ns_ud))
Nups(ii,:) = Bdecomp(mup,Ns_Orb) ![Ns_Orb = Ns = Nlat*Norb*(Nbath+1) in CDMFT code]
Ndws(ii,:) = Bdecomp(mdw,Ns_Orb)
enddo
Nud(1,:) = Breorder(Nups) !Actually, they are already reordered in CDMFT code...
Nud(2,:) = Breorder(Ndws) !...and breorder() corresponds to an identity: look in ED_SETUP!
!
!Diagonal densities
do ilat=1,Nlat
do ispin=1,Nspin
do iorb=1,Norb
is = imp_state_index(ilat,iorb)
single_particle_density_matrix(ilat,ilat,ispin,ispin,iorb,iorb) = &
single_particle_density_matrix(ilat,ilat,ispin,ispin,iorb,iorb) + &
peso*nud(ispin,is)*conjg(state_cvec(i))*state_cvec(i)
enddo
enddo
enddo
!
!Off-diagonal
do ispin=1,Nspin
do ilat=1,Nlat
do jlat=1,Nlat
do iorb=1,Norb
do jorb=1,Norb
is = imp_state_index(ilat,iorb)
js = imp_state_index(jlat,jorb)
if((Nud(ispin,js)==1).and.(Nud(ispin,is)==0))then
iud(1) = HI(1)%map(Indices(1))
iud(2) = HI(2)%map(Indices(2))
call c(js,iud(ispin),r,sgn1)
call cdg(is,r,k,sgn2)
Jndices = Indices
Jndices(1+(ispin-1)*Ns_Ud) = &
binary_search(HI(1+(ispin-1)*Ns_Ud)%map,k)
call indices2state(Jndices,[iDimUps,iDimDws],j)
!
single_particle_density_matrix(ilat,jlat,ispin,ispin,iorb,jorb) = &
single_particle_density_matrix(ilat,jlat,ispin,ispin,iorb,jorb) + &
peso*sgn1*state_cvec(i)*sgn2*conjg(state_cvec(j))
endif
enddo
enddo
enddo
enddo
enddo
!endif
!
!
enddo
call delete_sector(isector,HI)
endif
!
if(allocated(state_cvec))deallocate(state_cvec)
!
enddo
!
!
!
!
!
!
end subroutine density_matrix_impurity
!####################################################################
! COMPUTATIONAL ROUTINES
!####################################################################
subroutine init_custom_observables(N,Hk)
integer :: N
complex(8),dimension(:,:,:) :: Hk
!
if(MpiMaster)then
custom_o%N_filled=0
custom_o%N_asked=N
allocate(custom_o%Hk(size(Hk,1),size(Hk,2),size(Hk,3)))
custom_o%Hk=Hk
allocate(custom_o%item(N))
custom_o%init=.true.
endif
!
end subroutine init_custom_observables
subroutine add_custom_observable_local(o_name,sij)
integer :: i
complex(8),dimension(:,:) :: sij
character(len=*) :: o_name
!
if(MpiMaster .and. custom_o%init)then
if(custom_o%N_filled .gt. custom_o%N_asked)then
STOP "add_custom_observable: too many observables given"
call clear_custom_observables
endif
!
custom_o%N_filled=custom_o%N_filled+1
custom_o%item(custom_o%N_filled)%o_name=o_name
custom_o%item(custom_o%N_filled)%o_value=0.d0
!
allocate(custom_o%item(custom_o%N_filled)%sij(size(custom_o%Hk,1),size(custom_o%Hk,2),size(custom_o%Hk,3)))
do i=1,size(custom_o%item(custom_o%N_filled)%sij,3)
custom_o%item(custom_o%N_filled)%sij(:,:,i)=sij
enddo
else
STOP "add_custom_observable: custom observables not initialized"
endif
end subroutine add_custom_observable_local
subroutine add_custom_observable_kdep(o_name,sijk)
integer :: i
complex(8),dimension(:,:,:) :: sijk
character(len=*) :: o_name
!
if(MpiMaster .and. custom_o%init)then
if(custom_o%N_filled .gt. custom_o%N_asked)then
STOP "add_custom_observable: too many observables given"
call clear_custom_observables
endif
!
custom_o%N_filled=custom_o%N_filled+1
custom_o%item(custom_o%N_filled)%o_name=o_name
custom_o%item(custom_o%N_filled)%o_value=0.d0
!
allocate(custom_o%item(custom_o%N_filled)%sij(size(custom_o%Hk,1),size(custom_o%Hk,2),size(custom_o%Hk,3)))
custom_o%item(custom_o%N_filled)%sij=sijk
else
STOP "add_custom_observable: custom observables not initialized"
endif
end subroutine add_custom_observable_kdep
subroutine get_custom_observables()
integer :: i
!
if(MpiMaster .and. custom_o%init)then
if(custom_o%N_filled .eq. 0)then
write(Logfile,*)"WARNING! Custom observables initialized but none given."
RETURN
endif
!
write(LOGfile,*)"Calculating custom observables"
!
allocate(sij(size(custom_o%Hk,1),size(custom_o%Hk,2),size(custom_o%Hk,3)))
allocate(Hk(size(custom_o%Hk,1),size(custom_o%Hk,2),size(custom_o%Hk,3)))
sij=zero
Hk=zero
!
Hk=custom_o%Hk
do i=1,custom_o%N_filled
sij=custom_o%item(i)%sij
if(finiteT) then
custom_o%item(i)%o_value=calculate_observable_integral_finite_t()
else
custom_o%item(i)%o_value=calculate_observable_integral_zero_t()
endif
write(LOGfile,"(A,10f18.12,A)")reg(custom_o%item(i)%o_name)//" = ",custom_o%item(i)%o_value
enddo
call write_custom_legend()
call write_custom_observables()
deallocate(sij,Hk)
endif
!
end subroutine get_custom_observables
subroutine clear_custom_observables()
integer :: i
if(MpiMaster .and. custom_o%init)then
do i=1,custom_o%N_filled
deallocate(custom_o%item(i)%sij)
custom_o%item(i)%o_name=""
custom_o%item(i)%o_value=0.d0
enddo
deallocate(custom_o%Hk)
custom_o%N_asked=0
custom_o%N_filled=0
custom_o%init=.false.
endif
end subroutine clear_custom_observables
!+---------------------------------------------------------------------------------+
!PURPOSE : Evaluate and print out custom observable
!+---------------------------------------------------------------------------------+
!T=0
function calculate_observable_integral_zero_t() result(out_2)
integer :: inf
real(8) :: out_2,spin_multiplicity
!
out_2=0.d0
spin_multiplicity=3.d0-Nspin
!
call quad(sum_observable_kmesh,a=0.0d0,inf=1,verbose=(ED_VERBOSE>=3),result=out_2,strict=.false.)
!
out_2=spin_multiplicity*out_2/pi
return
end function calculate_observable_integral_zero_t
!T finite
function calculate_observable_integral_finite_t() result(out_2)
integer :: inf,Nmax,ii
real(8) :: out_2,spin_multiplicity,omegamax,integralpart
!
!1) Find the real omegamax
nmax=int(2*(abs(max_exc)+2.d0*hwband)*beta/pi)
if (mod(nmax,2)==0)then
nmax=nmax/2
else
nmax=(nmax+1)/2
endif
integrationR=2*(nmax+1)*pi/beta
!2) Evaluate discrete sum
!
out_2=0.d0
do ii=0,Nmax
out_2=out_2+dreal(sum_observable_kmesh_complex(xi*(2*ii+1)*pi/beta))
enddo
!
out_2=2.d0*(1/beta)*out_2
!
!3) Evaluate integral part
integralpart=0.d0
call quad(integral_contour,a=-pi,b=pi,verbose=(ED_VERBOSE>=3),key=6,result=integralpart,strict=.false.)
!
!4) Sum all
out_2=out_2+integralpart
!5) Spin trick
spin_multiplicity=3.d0-Nspin
out_2=spin_multiplicity*out_2
return
end function calculate_observable_integral_finite_t
!
!
!+-------------------------------------------------------------------+
!PURPOSE : complex integration
!+-------------------------------------------------------------------+
function integral_contour(zeta) result(f)
real(8) :: zeta,f
complex(8) :: w,fermi
!
w=integrationR*exp(xi*zeta)
if(dreal((w-XMU)*beta)>= 100)then
fermi=0.d0
else
fermi=(1/(exp(beta*(w-XMU))+1))
endif
!
f=dreal((1.d0/pi)*w*fermi*sum_observable_kmesh_complex(w))
end function integral_contour
!
!
!+-------------------------------------------------------------------+
!PURPOSE : sum on k-vectors
!+-------------------------------------------------------------------+
!
function sum_observable_kmesh(omega) result(out_1)
integer :: ii,jj,kk
real(8) :: omega
real(8) :: out_1
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb) :: g,invg0
complex(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: invg_lso,invg0_lso,sigma,Gk_lso
!
out_1=0.d0
!
g=zero
invg0=zero
invg_lso=zero
invg0_lso=zero
Gk_lso=zero
Sigma=zero
!
!Obtain Sigma(iw)
call ed_gf_cluster(dcmplx(0.d0,omega),g)
invg_lso=nnn2lso_reshape(g,Nlat,Nspin,Norb)
call inv(invg_lso)
invg0=invg0_bath(dcmplx(0.d0,omega))
invg0_lso=nnn2lso_reshape(invg0,Nlat,Nspin,Norb)
Sigma=invg0_lso-invg_lso
!
!Do the k-sum
do ii=1,size(Hk,3)
Gk_lso=(dcmplx(0d0,omega)+xmu)*eye(Nlat*Nspin*Norb) - Hk(:,:,ii) - Sigma
call inv(Gk_lso)
out_1=out_1+DREAL(trace(matmul(sij(:,:,ii),Gk_lso)) - trace(sij(:,:,ii))/(dcmplx(-1.1d0,omega)))
enddo
!
out_1=out_1/(size(Hk,3))
!
return
!
end function sum_observable_kmesh
!
function sum_observable_kmesh_complex(omega) result(out_1)
integer :: ii,jj,kk
complex(8) :: omega
complex(8) :: out_1
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb) :: g,invg0
complex(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: invg_lso,invg0_lso,sigma,Gk_lso
!
out_1=0.d0
!
g=zero
invg0=zero
invg_lso=zero
invg0_lso=zero
Gk_lso=zero
Sigma=zero
!
!
call ed_gf_cluster(omega,g)
invg_lso=nnn2lso_reshape(g,Nlat,Nspin,Norb)
call inv(invg_lso)
invg0=invg0_bath(omega)
invg0_lso=nnn2lso_reshape(invg0,Nlat,Nspin,Norb)
Sigma=invg0_lso-invg_lso
!
do ii=1,size(Hk,3)
Gk_lso=(xi*omega+xmu)*eye(Nlat*Nspin*Norb) - Hk(:,:,ii)- Sigma
call inv(Gk_lso)
out_1=out_1+DREAL(trace(matmul(sij(:,:,ii),Gk_lso)))
enddo
out_1=out_1/(size(Hk,3))
!
return
!
end function sum_observable_kmesh_complex
!####################################################################
! PRINTING ROUTINES
!####################################################################
!+-------------------------------------------------------------------+
!PURPOSE : write legend, i.e. info about columns
!+-------------------------------------------------------------------+
subroutine write_legend()
integer :: unit,iorb,jorb,ispin
unit = free_unit()
open(unit,file="observables_info.ed")
write(unit,"(A1,90(A10,6X))")"#",&
(reg(txtfy(iorb))//"dens_"//reg(txtfy(iorb)),iorb=1,Norb),&
(reg(txtfy(Norb+iorb))//"docc_"//reg(txtfy(iorb)),iorb=1,Norb),&
(reg(txtfy(2*Norb+iorb))//"nup_"//reg(txtfy(iorb)),iorb=1,Norb),&
(reg(txtfy(3*Norb+iorb))//"ndw_"//reg(txtfy(iorb)),iorb=1,Norb),&
(reg(txtfy(4*Norb+iorb))//"mag_"//reg(txtfy(iorb)),iorb=1,Norb),&
reg(txtfy(5*Norb+1))//"s2",&
reg(txtfy(5*Norb+2))//"egs",&
((reg(txtfy(5*Norb+2+(iorb-1)*Norb+jorb))//"sz2_"//reg(txtfy(iorb))//reg(txtfy(jorb)),jorb=1,Norb),iorb=1,Norb),&
((reg(txtfy((5+Norb)*Norb+2+(iorb-1)*Norb+jorb))//"n2_"//reg(txtfy(iorb))//reg(txtfy(jorb)),jorb=1,Norb),iorb=1,Norb)
close(unit)
!
unit = free_unit()
open(unit,file="parameters_info.ed")
write(unit,"(A1,90(A14,1X))")"#","1xmu","2beta",&
(reg(txtfy(2+iorb))//"U_"//reg(txtfy(iorb)),iorb=1,Norb),&
reg(txtfy(2+Norb+1))//"U'",reg(txtfy(2+Norb+2))//"Jh"
close(unit)
!
iolegend=.false.
end subroutine write_legend
subroutine write_custom_legend()
integer :: unit,i
unit = free_unit()
open(unit,file="custom_observables_info.ed")
write(unit,"(A1,90(A10,6X))")"#",(reg(txtfy(i))//reg(custom_o%item(i)%o_name),i=1,custom_o%N_filled)
close(unit)