-
Notifications
You must be signed in to change notification settings - Fork 3
/
ED_FIT_REPLICA.f90
969 lines (926 loc) · 41.8 KB
/
ED_FIT_REPLICA.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
MODULE ED_FIT_REPLICA
USE SF_CONSTANTS
USE SF_OPTIMIZE, only:fmin_cg,fmin_cgminimize
USE SF_LINALG, only:eye,zeye,inv,inv_her,trace,operator(.x.) !BLAS xgemm operator overloading
USE SF_IOTOOLS, only:reg,free_unit,str
USE SF_ARRAYS, only:arange
USE SF_MISC, only:assert_shape
USE ED_INPUT_VARS
USE ED_VARS_GLOBAL
USE ED_AUX_FUNX
USE ED_BATH
USE ED_BATH_FUNCTIONS
implicit none
private
public :: chi2_fitgf_replica
integer :: Ldelta
complex(8),dimension(:,:,:,:,:,:,:),allocatable :: FGmatrix
logical(8),dimension(:,:,:,:,:,:),allocatable :: Hmask
complex(8),dimension(:,:),allocatable :: Fdelta
real(8),dimension(:),allocatable :: Xdelta,Wdelta
integer :: totNorb,totNspin
integer,dimension(:),allocatable :: getIorb,getJorb,getIspin,getJspin,getIlat,getJlat
integer :: Orb_indx,Spin_indx,Spin_mask
integer,dimension(:),allocatable :: Nlambdas
!location of the maximum of the chisquare over Nlso.
integer :: maxchi_loc
!
type nsymm_vector
real(8),dimension(:),allocatable :: element
end type nsymm_vector
!
contains
!##################################################################
! THE CALCULATION OF THE \chi^2 FUNCTIONS USE PROCEDURES FURTHER
! BELOW TO EVALUATE INDEPENDENTLY THE ANDERSON MODEL:
! - DELTA,
! -\GRAD DELTA
! - G0
! THE LATTER ARE ADAPTED FROM THE PROCEDURES:
! DELTA_BATH_MATS
! GRAD_DELTA_BATH_MATS
! G0 BATH_MATS
! FOR, YOU NEED TO DECOMPOSE THE a INPUT ARRAY INTO ELEMENTS.
!##################################################################
!+-------------------------------------------------------------+
!PURPOSE : Chi^2 interface for REPLICA BATH
!+-------------------------------------------------------------+
subroutine chi2_fitgf_replica(fg,bath_)
complex(8),dimension(:,:,:,:,:,:,:) :: fg ![Nlat][Nlat][Nspin][Nspin][Norb][Norb][Lmats]
real(8),dimension(:),intent(inout) :: bath_
real(8),dimension(:),allocatable :: array_bath
integer :: i,j,ilat,jlat,iorb,jorb,ispin,jspin,ibath,io,jo
integer :: iter,stride,counter,Asize
real(8) :: chi
logical :: check
character(len=256) :: suffix
integer :: unit
!
if(size(fg,1)/=Nlat) stop "chi2_fitgf_replica error: size[fg,1]!=Nlat"
if(size(fg,2)/=Nlat) stop "chi2_fitgf_replica error: size[fg,2]!=Nlat"
if(size(fg,3)/=Nspin)stop "chi2_fitgf_replica error: size[fg,3]!=Nspin"
if(size(fg,4)/=Nspin)stop "chi2_fitgf_replica error: size[fg,4]!=Nspin"
if(size(fg,5)/=Norb) stop "chi2_fitgf_replica error: size[fg,5]!=Norb"
if(size(fg,6)/=Norb) stop "chi2_fitgf_replica error: size[fg,6]!=Norb"
!
allocate(Nlambdas(Nbath))
!
check= check_bath_dimension(bath_)
if(.not.check)stop "chi2_fitgf_replica error: wrong bath dimensions"
!
if(cg_pow/=2.AND.cg_norm=="frobenius")then
print *, "WARNING: CG_POW must be 2 for a meaningful definition of the Frobenius norm."
print *, " we'll still let you go ahead with the desired input, but please be "
print *, " be aware that CG_POW is not doing what you would expect for a chi^q"
endif
!
call allocate_dmft_bath()
call set_dmft_bath(bath_)
allocate(array_bath(size(bath_)-Nbath))
counter=0
do ibath=1,Nbath
counter=counter+1
Nlambdas(ibath)=NINT(bath_(counter))
enddo
array_bath=bath_(Nbath+1:size(bath_))
!
Ldelta = Lfit ; if(Ldelta>size(fg,7))Ldelta=size(fg,7)
!
!
allocate(FGmatrix(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta))
allocate(Hmask(Nlat,Nlat,Nspin,Nspin,Norb,Norb))
allocate(Xdelta(Ldelta))
allocate(Wdelta(Ldelta))
!
Xdelta = pi/beta*(2*arange(1,Ldelta)-1)
!
select case(cg_weight)
case default
Wdelta=1d0
case(2)
Wdelta=1d0*arange(1,Ldelta)
case(3)
Wdelta=Xdelta
end select
!
!
!----------------------------------------------------------------------------------------
!POSSIBLY TO BE USED AT SOME POINT (but beware that is meaningless for `FROBENIUS` norm)
Hmask=.true.
! Aren't we sure about FG_ij = FG_ji, coming from Hbath being hermitian?
!
! -> Hmask=Hbath_mask(wdiag=.false.,uplo=.true.) NO, this does more than
! that, putting .FALSE. where Hbath is zero, which may lead to tricky
! wrong fits: I have data showing this for a trimer, where at ~ 6.30d-10
! dmft error, hence a very well converged point, we have a nonnegligible
! Im(Weiss) component for ilat=1, jlat=3, despite Hmask(1,3) = Hmask(3,1)
! would result .FALSE. if computed with the Hbath_mask call above.
!
! For now I'd say that we'd better dump everything inside the \chi^2, hence Hmask=.true.,
! but we might want to consider exploiting hermiticity at least (probably not looking at
! Hbath though: who guarantees zeros therein imply zeros here?). Discussion needed...
!----------------------------------------------------------------------------------------
!
!
FGmatrix = fg
!
!
select case(cg_method) !0=NR-CG[default]; 1=CG-MINIMIZE
case default
if(cg_grad==0)then
#if __GNUC__ >= 8 || __INTEL_COMPILER >= 1500
write(LOGfile,*)" Using analytic gradient"
select case (cg_scheme)
case ("weiss")
call fmin_cg(array_bath,&
chi2_weiss_replica,&
grad_chi2_weiss_replica,&
iter,&
chi, &
itmax=cg_niter,&
ftol=cg_Ftol, &
istop=cg_stop, &
iverbose=(ed_verbose>3))
case ("delta")
call fmin_cg(array_bath,&
chi2_delta_replica,&
grad_chi2_delta_replica,&
iter,&
chi, &
itmax=cg_niter,&
ftol=cg_Ftol, &
istop=cg_stop, &
iverbose=(ed_verbose>3))
case default
stop "chi2_fitgf_replica error: cg_scheme != [weiss,delta]"
end select
#else
STOP "analytic gradient not supported for gfortran < 8"
#endif
else
write(LOGfile,*)" Using numerical gradient"
select case (cg_scheme)
case ("weiss")
call fmin_cg(array_bath,chi2_weiss_replica,&
iter,&
chi, &
itmax=cg_niter,&
ftol=cg_Ftol, &
istop=cg_stop, &
iverbose=(ed_verbose>3))
case ("delta")
call fmin_cg(array_bath,chi2_delta_replica,&
iter,&
chi, &
itmax=cg_niter,&
ftol=cg_Ftol, &
istop=cg_stop, &
iverbose=(ed_verbose>3))
case default
stop "chi2_fitgf_replica error: cg_scheme != [weiss,delta]"
end select
endif
!
case (1)
if(cg_grad==0)then
write(*,*) " "
write(*,*) "WARNING: analytic gradient not available with cg-method=1 (minimize f77 routine)"
write(*,*) " > we will force cg_grad=1 (so let the routine estimate the gradient) "
write(*,*) " "
endif
select case (cg_scheme)
case ("weiss")
call fmin_cgminimize(array_bath,chi2_weiss_replica,&
iter,&
chi, &
itmax=cg_niter,&
ftol=cg_Ftol, &
new_version=cg_minimize_ver,&
hh_par=cg_minimize_hh, &
iverbose=(ed_verbose>3))
case ("delta")
call fmin_cgminimize(array_bath,chi2_delta_replica,&
iter,&
chi, &
itmax=cg_niter,&
ftol=cg_Ftol, &
new_version=cg_minimize_ver,&
hh_par=cg_minimize_hh, &
iverbose=(ed_verbose>3))
case default
stop "chi2_fitgf_replica error: cg_scheme != [weiss,delta]"
end select
!
end select
!
write(LOGfile,"(A,ES18.9,A,I5,A)")"chi^2|iter"//reg(ed_file_suffix)//'= ',chi," | ",iter," <-- All Orbs, All Spins"
!
suffix="_ALLorb_ALLspins"//reg(ed_file_suffix)
unit=free_unit()
open(unit,file="chi2fit_results"//reg(suffix)//".ed",position="append")
write(unit,"(ES18.9,1x,I5)") chi,iter
close(unit)
!
bath_(Nbath+1:size(bath_))=array_bath
call set_dmft_bath(bath_) ! *** bath_ --> dmft_bath *** (per write fit result)
call write_dmft_bath(LOGfile)
call save_dmft_bath()
!
call write_fit_result()
!
call get_dmft_bath(bath_) ! *** dmft_bath --> bath_ *** (bath in output)
call deallocate_dmft_bath()
deallocate(FGmatrix,Hmask,Xdelta,Wdelta)
deallocate(array_bath)
deallocate(Nlambdas)
!
contains
!
subroutine write_fit_result()
complex(8) :: fgand(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta)
integer :: i,j,s,l,ilat,jlat,iorb,jorb,ispin,jspin
real(8) :: w
if(cg_scheme=='weiss')then
fgand(:,:,:,:,:,:,:) = g0and_bath(xi*Xdelta(:))
else
fgand(:,:,:,:,:,:,:) = delta_bath(xi*Xdelta(:))
endif
!
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
suffix="_i"//reg(str(ilat))//&
"_j"//reg(str(jlat))//&
"_l"//reg(str(iorb))//&
"_m"//reg(str(jorb))//&
"_s"//reg(str(ispin))//&
"_r"//reg(str(jspin))//reg(ed_file_suffix)
unit=free_unit()
if(cg_scheme=='weiss')then
open(unit,file="fit_weiss"//reg(suffix)//".ed")
else
open(unit,file="fit_delta"//reg(suffix)//".ed")
endif
do i=1,Ldelta
w = Xdelta(i)
write(unit,"(5F24.15)")Xdelta(i),&
dimag(fg(ilat,jlat,ispin,jspin,iorb,jorb,i)),dimag(fgand(ilat,jlat,ispin,jspin,iorb,jorb,i)),&
dreal(fg(ilat,jlat,ispin,jspin,iorb,jorb,i)),dreal(fgand(ilat,jlat,ispin,jspin,iorb,jorb,i))
enddo
close(unit)
enddo
enddo
enddo
enddo
enddo
enddo
end subroutine write_fit_result
!
end subroutine chi2_fitgf_replica
!##################################################################
! THESE PROCEDURES EVALUATE THE \chi^2 FUNCTIONS TO MINIMIZE.
!##################################################################
!
!
!+-----------------------------------------------------------------+
!PURPOSE: Evaluate the \chi^2 distance of \Delta_Anderson function.
!+-----------------------------------------------------------------+
function chi2_delta_replica(a) result(chi2)
real(8),dimension(:) :: a
real(8) :: chi2
!
select case(cg_norm)
case ("elemental")
chi2 = chi2_delta_replica_elemental(a)
case ("frobenius")
chi2 = chi2_delta_replica_frobenius(a)
case default
stop "chi2_fitgf_replica error: cg_norm != [elemental,frobenius]"
end select
!
end function chi2_delta_replica
!
!
!> ELEMENTAL NORM: weighted sum over i\omega for each matrix element, then weighted sum over elements
function chi2_delta_replica_elemental(a) result(chi2)
real(8),dimension(:) :: a
real(8) :: chi2
real(8),dimension(Ldelta) :: chi2_freq
real(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb) :: chi2_mtrx,Wmat
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta) :: Delta
integer :: ilat,jlat,iorb,jorb,ispin,jspin
!
#ifdef _DEBUG
if(ed_verbose>5)write(Logfile,"(A,"//str(size(a))//"ES10.2)")"DEBUG chi2_delta_replica_elemental. a:",a
#endif
!
Delta = delta_replica(a)
!
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
chi2_freq = abs(Delta(ilat,jlat,ispin,jspin,iorb,jorb,1:Ldelta) - FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,1:Ldelta))
chi2_mtrx(ilat,jlat,ispin,jspin,iorb,jorb) = sum(chi2_freq**cg_pow/Wdelta) !Weighted sum over matsubara frqs
select case(cg_matrix)
case(0) !FLAT (all matrix elements weighted equal)
Wmat(ilat,jlat,ispin,jspin,iorb,jorb) = 1d0!0.25d0 !needs to depend on the hopping I think (\Delta=(D/2)^2*Gloc…)
case(1) !SPECTRAL (normalization through A(iw), element by element)
Wmat(ilat,jlat,ispin,jspin,iorb,jorb) = abs(sum(FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,1:Lmats)))/beta
#ifdef _DEBUG
if(ed_verbose>4)then
print*, "ilat: "//str(ilat)//" ispin: "//str(ispin)//" iorb: "//str(iorb)
print*, "jlat: "//str(jlat)//" jspin: "//str(jspin)//" jorb: "//str(jorb)
print*, "> Wmat: ", Wmat(ilat,jlat,ispin,jspin,iorb,jorb)/0.25d0
endif
#endif
end select
enddo
enddo
enddo
enddo
enddo
enddo
!
chi2 = sum(chi2_mtrx / Wmat, Hmask) !Weighted sum over matrix elements
chi2 = chi2 / Ldelta / count(Hmask) !Normalization over {iw} and Hmask
!
#ifdef _DEBUG
if(ed_verbose>3)write(Logfile,"(A,ES10.2)")"DEBUG chi2_delta_replica_elemental. chi2:",chi2
#endif
!
end function chi2_delta_replica_elemental
!
!
!> FROBENIUS NORM: global \chi^2 for all components, only i\omega are weighted
function chi2_delta_replica_frobenius(a) result(chi2)
real(8),dimension(:) :: a
real(8) :: chi2
real(8),dimension(Ldelta) :: chi2_freq
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta) :: Delta
complex(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: Delta_lso
integer :: l
!
#ifdef _DEBUG
if(ed_verbose>5)write(Logfile,"(A,"//str(size(a))//"ES10.2)")"DEBUG chi2_delta_replica_frobenius. a:",a
#endif
!
Delta = delta_replica(a)
!
do l=1,Ldelta
Delta_lso = nnn2lso_reshape(delta(:,:,:,:,:,:,l) - FGmatrix(:,:,:,:,:,:,l),Nlat,Nspin,Norb)
chi2_freq(l) = sqrt(trace(matmul(Delta_lso,conjg(transpose(Delta_lso)))))
enddo
!
chi2 = sum(chi2_freq**cg_pow/Wdelta) !Weighted sum over matsubara frqs
chi2 = chi2/Ldelta/(Nlat*Nspin*Norb) !Normalization over {iw} and Nlso
!
#ifdef _DEBUG
if(ed_verbose>3)write(Logfile,"(A,ES10.2)")"DEBUG chi2_delta_replica_frobenius. chi2:",chi2
#endif
!
end function chi2_delta_replica_frobenius
!
!
!+-------------------------------------------------------------+
!PURPOSE: Evaluate the \chi^2 distance of G0_Anderson function.
!+-------------------------------------------------------------+
function chi2_weiss_replica(a) result(chi2)
real(8),dimension(:) :: a
real(8) :: chi2
!
select case(cg_norm)
case ("elemental")
chi2 = chi2_weiss_replica_elemental(a)
case ("frobenius")
chi2 = chi2_weiss_replica_frobenius(a)
case default
stop "chi2_fitgf_replica error: cg_norm != [elemental,frobenius]"
end select
!
end function chi2_weiss_replica
!
!
!> ELEMENTAL NORM: weighted sum over i\omega for each matrix element, then weighted sum over elements
function chi2_weiss_replica_elemental(a) result(chi2)
real(8),dimension(:) :: a
real(8) :: chi2
real(8),dimension(Ldelta) :: chi2_freq
real(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb) :: chi2_mtrx,Wmat
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta) :: g0and
integer :: ilat,jlat,iorb,jorb,ispin,jspin
!
#ifdef _DEBUG
if(ed_verbose>5)write(Logfile,"(A,"//str(size(a))//"ES10.2)")"DEBUG chi2_weiss_replica_elemental. a:",a
#endif
!
g0and = g0and_replica(a)
!
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
chi2_freq = abs(g0and(ilat,jlat,ispin,jspin,iorb,jorb,1:Ldelta) - FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,1:Ldelta))
chi2_mtrx(ilat,jlat,ispin,jspin,iorb,jorb) = sum(chi2_freq**cg_pow/Wdelta) !Weighted sum over matsubara frqs
select case(cg_matrix)
case(0) !FLAT (all matrix elements weighted equal)
Wmat(ilat,jlat,ispin,jspin,iorb,jorb) = 1d0
case(1) !SPECTRAL (normalization through A(iw), element by element)
Wmat(ilat,jlat,ispin,jspin,iorb,jorb) = abs(sum(FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,1:Lmats)))/beta
#ifdef _DEBUG
if(ed_verbose>4)then
print*, "ilat: "//str(ilat)//"ispin: "//str(ispin)//"iorb: "//str(iorb)
print*, "jlat: "//str(jlat)//"jspin: "//str(jspin)//"jorb: "//str(jorb)
print*, "> Wmat = ", Wmat(ilat,jlat,ispin,jspin,iorb,jorb)
endif
#endif
end select
enddo
enddo
enddo
enddo
enddo
enddo
!
chi2 = sum(chi2_mtrx / Wmat, Hmask) !Weighted sum over matrix elements
chi2 = chi2 / Ldelta / count(Hmask) !Normalization over {iw} and Hmask
!
#ifdef _DEBUG
if(ed_verbose>3)write(Logfile,"(A,ES10.2)")"DEBUG chi2_weiss_replica_elemental. chi2:",chi2
#endif
!
end function chi2_weiss_replica_elemental
!
!
!> FROBENIUS NORM: global \chi^2 for all components, only i\omega are weighted
function chi2_weiss_replica_frobenius(a) result(chi2)
real(8),dimension(:) :: a
real(8) :: chi2
real(8),dimension(Ldelta) :: chi2_freq
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta) :: g0and
complex(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: Delta_lso
integer :: l
!
#ifdef _DEBUG
if(ed_verbose>5)write(Logfile,"(A,"//str(size(a))//"ES10.2)")"DEBUG chi2_weiss_replica_frobenius. a:",a
#endif
!
g0and = g0and_replica(a)
!
do l=1,Ldelta
Delta_lso = nnn2lso_reshape(g0and(:,:,:,:,:,:,l) - FGmatrix(:,:,:,:,:,:,l),Nlat,Nspin,Norb)
chi2_freq(l) = sqrt(trace(matmul(Delta_lso,conjg(transpose(Delta_lso)))))
enddo
!
chi2 = sum(chi2_freq**cg_pow/Wdelta) !Weighted sum over matsubara frqs
chi2 = chi2/Ldelta/(Nlat*Nspin*Norb) !Normalization over {iw} and Nlso
!
#ifdef _DEBUG
if(ed_verbose>3)write(Logfile,"(A,ES10.2)")"DEBUG chi2_weiss_replica_frobenius. chi2:",chi2
#endif
!
end function chi2_weiss_replica_frobenius
!######################################################################
! THESE PROCEDURES EVALUATE THE >GRADIENTS< OF THE \chi^2 TO MINIMIZE.
!######################################################################
!
!
#if __GNUC__ >= 8 || __INTEL_COMPILER >= 1500
!+----------------------------------------------------------------------+
!PURPOSE: Evaluate the gradient \Grad\chi^2 of \Delta_Anderson function.
!+----------------------------------------------------------------------+
function grad_chi2_delta_replica(a) result(dchi2)
real(8),dimension(:) :: a
real(8),dimension(size(a)) :: dchi2
!
select case(cg_norm)
case ("elemental")
dchi2 = grad_chi2_delta_replica_elemental(a)
case ("frobenius")
dchi2 = grad_chi2_delta_replica_frobenius(a)
case default
stop "chi2_fitgf_replica error: cg_norm != [elemental,frobenius]"
end select
!
end function grad_chi2_delta_replica
!
!
!> ELEMENTAL NORM: weighted sum over i\omega for each matrix element, then weighted sum over elements
function grad_chi2_delta_replica_elemental(a) result(dchi2)
real(8),dimension(:) :: a
real(8),dimension(size(a)) :: dchi2
real(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,size(a)) :: df
real(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb) :: Wmat
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta) :: Delta
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta,size(a)) :: dDelta
complex(8),dimension(Ldelta) :: Ftmp
real(8),dimension(Ldelta) :: Ctmp
integer :: ia,ilat,jlat,iorb,jorb,ispin,jspin
!
#ifdef _DEBUG
if(ed_verbose>5)write(Logfile,"(A,"//str(size(a))//"ES10.2)")"DEBUG grad_chi2_delta_replica_elemental. a:",a
#endif
!
Delta = delta_replica(a)
dDelta = grad_delta_replica(a)
!
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
Ftmp = delta(ilat,jlat,ispin,jspin,iorb,jorb,1:Ldelta) - FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,1:Ldelta)
Ctmp = abs(Ftmp)**(cg_pow-2)
do ia = 1,size(a)
df(ilat,jlat,ispin,jspin,iorb,jorb,ia) = & !Weighted sum over matsubara frqs
sum( dreal(Ftmp) * dreal(dDelta(ilat,jlat,ispin,jspin,iorb,jorb,:,ia)) * Ctmp/Wdelta ) + &
sum( dimag(Ftmp) * dimag(dDelta(ilat,jlat,ispin,jspin,iorb,jorb,:,ia)) * Ctmp/Wdelta )
enddo
select case(cg_matrix)
case(0) !FLAT (all matrix elements weighted equal)
Wmat(ilat,jlat,ispin,jspin,iorb,jorb) = 1d0!0.25d0 !needs to depend on the hopping I think (\Delta=(D/2)^2*Gloc…)
case(1) !SPECTRAL (normalization through A(iw), element by element)
Wmat(ilat,jlat,ispin,jspin,iorb,jorb) = abs(sum(FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,1:Lmats)))/beta
end select
enddo
enddo
enddo
enddo
enddo
enddo
!
do ia=1,size(a)
dchi2(ia) = + cg_pow * sum( df(:,:,:,:,:,:,ia) / Wmat, Hmask) !Weighted sum over matrix elements
dchi2(ia) = dchi2(ia) / Ldelta / count(Hmask) !Normalization over {iw} and Hmask
enddo
!
#ifdef _DEBUG
if(ed_verbose>4)write(Logfile,"(A,"//str(size(a))//"ES10.2)")"DEBUG grad_chi2_delta_replica_elemental. dChi2:",dchi2
#endif
!
end function grad_chi2_delta_replica_elemental
!
!
!> FROBENIUS NORM: global \chi^2 for all components, only i\omega are weighted
function grad_chi2_delta_replica_frobenius(a) result(dchi2)
real(8),dimension(:) :: a
real(8),dimension(size(a)) :: dchi2
real(8),dimension(Ldelta,size(a)) :: df
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta) :: Delta
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta,size(a)) :: dDelta
complex(8),dimension(Ldelta) :: Ftmp
real(8),dimension(Ldelta,size(a)) :: dChi_freq
integer :: i,j,idelta,ilat,jlat,iorb,jorb,ispin,jspin
!
Delta = delta_replica(a)
dDelta = grad_delta_replica(a)
Ftmp=zero
df=zero
!
do idelta=1,Ldelta
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
!
Ftmp(idelta) = Ftmp(idelta) + abs(Delta(ilat,jlat,ispin,jspin,iorb,jorb,idelta)-FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,idelta))**2
do j=1,size(a)
df(idelta,j) = df(idelta,j) + &
real(Delta(ilat,jlat,ispin,jspin,iorb,jorb,idelta) - FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,idelta)) * &
real(dDelta(ilat,jlat,ispin,jspin,iorb,jorb,idelta,j)) + &
imag(Delta(ilat,jlat,ispin,jspin,iorb,jorb,idelta) - FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,idelta)) * &
imag(dDelta(ilat,jlat,ispin,jspin,iorb,jorb,idelta,j))
enddo
enddo
enddo
enddo
enddo
enddo
enddo
Ftmp(idelta) = cg_pow * (sqrt(Ftmp(idelta))**(cg_pow-2)) / Wdelta(idelta)
dchi_freq(idelta,:) = Ftmp(idelta) * df(idelta,:)
enddo
!
dchi2 = sum(dchi_freq,1)/Ldelta/(Nlat*Nspin*Norb)
!
#ifdef _DEBUG
if(ed_verbose>4)write(Logfile,"(A,"//str(size(a))//"ES10.2)")"DEBUG grad_chi2_delta_replica_frobenius. dChi2:",dchi2
#endif
!
end function grad_chi2_delta_replica_frobenius
!
!
!+------------------------------------------------------------------+
!PURPOSE: Evaluate the gradient \Grad\chi^2 of G0_Anderson function.
!+------------------------------------------------------------------+
function grad_chi2_weiss_replica(a) result(dchi2)
real(8),dimension(:) :: a
real(8),dimension(size(a)) :: dchi2
!
select case(cg_norm)
case ("elemental")
dchi2 = grad_chi2_weiss_replica_elemental(a)
case ("frobenius")
dchi2 = grad_chi2_weiss_replica_frobenius(a)
case default
stop "chi2_fitgf_replica error: cg_norm != [elemental,frobenius]"
end select
!
end function grad_chi2_weiss_replica
!
!
!> ELEMENTAL NORM: weighted sum over i\omega for each matrix element, then weighted sum over elements
function grad_chi2_weiss_replica_elemental(a) result(dchi2)
real(8),dimension(:) :: a
real(8),dimension(size(a)) :: dchi2
real(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,size(a)) :: df
real(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb) :: Wmat
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta) :: g0and
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta,size(a)) :: dg0and
complex(8),dimension(Ldelta) :: Ftmp
real(8),dimension(Ldelta) :: Ctmp
integer :: ia,ilat,jlat,iorb,jorb,ispin,jspin
!
#ifdef _DEBUG
if(ed_verbose>5)write(Logfile,"(A,"//str(size(a))//"ES10.2)")"DEBUG grad_chi2_weiss_replica_elemental. a:",a
#endif
!
g0and = g0and_replica(a)
dg0and = grad_g0and_replica(a)
!
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
Ftmp = g0and(ilat,jlat,ispin,jspin,iorb,jorb,1:Ldelta) - FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,1:Ldelta)
Ctmp = abs(Ftmp)**(cg_pow-2)
do ia = 1,size(a)
df(ilat,jlat,ispin,jspin,iorb,jorb,ia) = & !Weighted sum over matsubara frqs
sum( dreal(Ftmp) * dreal(dg0and(ilat,jlat,ispin,jspin,iorb,jorb,:,ia)) * Ctmp/Wdelta ) + &
sum( dimag(Ftmp) * dimag(dg0and(ilat,jlat,ispin,jspin,iorb,jorb,:,ia)) * Ctmp/Wdelta )
enddo
select case(cg_matrix)
case(0) !FLAT (all matrix elements weighted equal)
Wmat(ilat,jlat,ispin,jspin,iorb,jorb) = 1d0
case(1) !SPECTRAL (normalization through A(iw), element by element)
Wmat(ilat,jlat,ispin,jspin,iorb,jorb) = abs(sum(FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,1:Lmats)))/beta
end select
enddo
enddo
enddo
enddo
enddo
enddo
!
do ia=1,size(a)
dchi2(ia) = + cg_pow * sum( df(:,:,:,:,:,:,ia) / Wmat, Hmask) !Weighted sum over matrix elements
dchi2(ia) = dchi2(ia) / Ldelta / count(Hmask) !Normalization over {iw} and Hmask
enddo
!
#ifdef _DEBUG
if(ed_verbose>4)write(Logfile,"(A,"//str(size(a))//"ES10.2)")"DEBUG grad_chi2_weiss_replica_elemental. dChi2:",dchi2
#endif
!
end function grad_chi2_weiss_replica_elemental
!
!
!> FROBENIUS NORM: global \chi^2 for all components, only i\omega are weighted
function grad_chi2_weiss_replica_frobenius(a) result(dchi2)
real(8),dimension(:) :: a
real(8),dimension(size(a)) :: dchi2
real(8),dimension(Ldelta,size(a)) :: df
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta) :: g0and
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta,size(a)) :: dg0and
complex(8),dimension(Ldelta) :: Ftmp
real(8),dimension(Ldelta,size(a)) :: dChi_freq
integer :: i,j,idelta,ilat,jlat,iorb,jorb,ispin,jspin
!
!
print*, " "
print*, "WARNING: the analytic gradient of the Weiss field Frobenius distance "
print*, " has been found giving /QUALITATIVELY WRONG/ fitted spectra! "
print*, " > the issue has emerged in some Nlat=Nspin=Norb=1 test runs "
print*, " > while the bug is investigated please switch cg_scheme to "
print*, " 'delta' or cg_norm to 'elemental' (or give up analytic cg)"
print*, " "
!
!
g0and = g0and_replica(a)
dg0and = grad_g0and_replica(a)
Ftmp=zero
df=zero
!
do idelta=1,Ldelta
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
!
Ftmp(idelta) = Ftmp(idelta) + abs(g0and(ilat,jlat,ispin,jspin,iorb,jorb,idelta)-FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,idelta))**2
do j=1,size(a)
df(idelta,j) = df(idelta,j) + &
real(g0and(ilat,jlat,ispin,jspin,iorb,jorb,idelta) - FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,idelta)) * &
real(dg0and(ilat,jlat,ispin,jspin,iorb,jorb,idelta,j)) + &
imag(g0and(ilat,jlat,ispin,jspin,iorb,jorb,idelta) - FGmatrix(ilat,jlat,ispin,jspin,iorb,jorb,idelta)) * &
imag(dg0and(ilat,jlat,ispin,jspin,iorb,jorb,idelta,j))
enddo
enddo
enddo
enddo
enddo
enddo
enddo
Ftmp(idelta) = cg_pow * (sqrt(Ftmp(idelta))**(cg_pow-2)) / Wdelta(idelta)
dchi_freq(idelta,:) = Ftmp(idelta) * df(idelta,:)
enddo
!
dchi2 = sum(dchi_freq,1)/Ldelta/(Nlat*Nspin*Norb)
!
#ifdef _DEBUG
if(ed_verbose>4)write(Logfile,"(A,"//str(size(a))//"ES10.2)")"DEBUG grad_chi2_weiss_replica_frobenius. dChi2:",dchi2
#endif
!
end function grad_chi2_weiss_replica_frobenius
#endif
!##################################################################
! THESE PROCEDURES EVALUATE THE ANDERSON FUNCTIONS:
! - \Delta (hybridization)
! - g0 (weiss field)
!##################################################################
function delta_replica(a) result(Delta)
real(8),dimension(:) :: a
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta) :: Delta
integer :: ilat,jlat,ispin,jspin,iorb,jorb,ibath,isym
integer :: i,io,jo,ndx,stride
complex(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: V_k
complex(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: Haux
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb) :: invH_knnn
real(8),dimension(Nbath) :: dummy_Vbath
type(nsymm_vector),dimension(Nbath) :: dummy_lambda
complex(8) :: iw
!
!ACHTUNG! here the bath was a temporary one, since we removed the possibility to act on other baths we need to replicate the
!function behaviour. Rather ugly...
!Get Hs
stride = 0
do ibath=1,Nbath
allocate(dummy_lambda(ibath)%element(Nlambdas(ibath)))
!Get Vs
stride = stride + 1
dummy_vbath(ibath) = a(stride)
!get Lambdas
dummy_lambda(ibath)%element=a(stride+1:stride+Nlambdas(ibath))
stride=stride+Nlambdas(ibath)
enddo
!
!
Delta=zero
do i=1,Ldelta
iw = xi*Xdelta(i)
do ibath=1,Nbath
invH_knnn = Hbath_build(dummy_lambda(ibath)%element)
Haux = zeye(Nlat*Nspin*Norb)*iw - nnn2lso_reshape(invH_knnn,Nlat,Nspin,Norb)
call inv(Haux)
invH_knnn = lso2nnn_reshape(Haux,Nlat,Nspin,Norb)
Delta(:,:,:,:,:,:,i)=Delta(:,:,:,:,:,:,i)+ dummy_Vbath(ibath)*invH_knnn(:,:,:,:,:,:)*dummy_Vbath(ibath)
enddo
enddo
!
do ibath=1,Nbath
deallocate(dummy_lambda(ibath)%element)
enddo
!
end function delta_replica
!
function g0and_replica(a) result(G0and)
real(8),dimension(:) :: a
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta) :: G0and,Delta
complex(8),dimension(Nlat*Norb*Nspin,Nlat*Norb*Nspin) :: zeta,fgorb
integer :: i,Nlso
integer :: ilat,jlat,iorb,jorb,ispin,jspin,io,jo
!
Nlso = Nlat*Norb*Nspin
!
Delta = delta_replica(a)
!
do i=1,Ldelta
zeta = (xi*Xdelta(i)+xmu)*zeye(Nlso)
FGorb = zeta - nnn2lso_reshape(impHloc + Delta(:,:,:,:,:,:,i), Nlat,Nspin,Norb)
call inv(FGorb)
G0and(:,:,:,:,:,:,i) = lso2nnn_reshape(FGorb,Nlat,Nspin,Norb)
enddo
!
end function g0and_replica
!
!
!
!
!
!
#if __GNUC__ >= 8 || __INTEL_COMPILER >= 1500
!##################################################################
! THESE PROCEDURES EVALUATE THE GRADIENT OF ANDERSON FUNCTIONS:
! - \Delta (hybridization)
! - g0 (weiss field)
!##################################################################
function grad_delta_replica(a) result(dDelta)
real(8),dimension(:) :: a
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta,size(a)) :: dDelta
integer :: ilat,jlat,ispin,iorb,jorb,ibath
integer :: k,l,io,counter
complex(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: H_reconstructed, Htmp,Hbasis_lso
complex(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb,Ldelta) :: Haux
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta) :: invH_knn
real(8),dimension(1,Nbath) :: dummy_Vbath !TODO) to extend to Vup != Vdw: 1->NSPIN
type(nsymm_vector),dimension(Nbath) :: dummy_lambda
!
!
!Get Hs
counter = 0
do ibath=1,Nbath
if(allocated(dummy_lambda(ibath)%element))deallocate(dummy_lambda(ibath)%element)
allocate(dummy_lambda(ibath)%element(Nlambdas(ibath)))
!
counter = counter + 1
do ispin=1,Nspin
dummy_vbath(1,ibath) = a(counter)
! ^ TODO) to extend to Vup != Vdw: 1->NSPIN
enddo
!
dummy_lambda(ibath)%element=a(counter+1:counter+Nlambdas(ibath))
counter=counter+Nlambdas(ibath)
enddo
!
dDelta=zero
counter=0
!
do ibath=1,Nbath
H_reconstructed = nnn2lso_reshape(Hbath_build(dummy_lambda(ibath)%element),Nlat,Nspin,Norb)
do l=1,Ldelta
Haux(:,:,l) = zeye(Nlat*Nspin*Norb)*(xi*Xdelta(l)) - H_reconstructed
call inv(Haux(:,:,l))
invH_knn(:,:,:,:,:,:,l) = lso2nnn_reshape(Haux(:,:,l),Nlat,Nspin,Norb)
enddo
!Derivate_Vp
counter=counter+1
do ispin=1,Nspin
dDelta(:,:,ispin,ispin,:,:,:,counter)=2d0*dummy_Vbath(1,ibath)*invH_knn(:,:,ispin,ispin,:,:,:)
! ^ TODO) to extend to Vup != Vdw: 1->ispin
enddo
!Derivate_lambda_p
do k=1,Nlambdas(ibath)
counter = counter + 1
Hbasis_lso=nnn2lso_reshape(Hbath_basis(k)%O,Nlat,Nspin,Norb)
do l=1,Ldelta
Htmp = ((Haux(:,:,l) .x. Hbasis_lso)) .x. Haux(:,:,l)
dDelta(:,:,:,:,:,:,l,counter)=lso2nnn_reshape((dummy_Vbath(1,ibath)**2)*Htmp,Nlat,Nspin,Norb)
enddo
enddo
enddo
!
do ibath=1,Nbath
if(allocated(dummy_lambda(ibath)%element))deallocate(dummy_lambda(ibath)%element)
enddo
!
end function grad_delta_replica
!
!
function grad_g0and_replica(a) result(dG0and)
real(8),dimension(:) :: a
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta,size(a)) :: dG0and,dDelta
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Ldelta) :: G0and
complex(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: dDelta_lso,dG0and_lso,G0and_lso
integer :: ilat,jlat,ispin,iorb,jorb
integer :: ik,l
!
G0and = g0and_replica(a)
dDelta = grad_delta_replica(a)
!
dG0and = zero
!
do l=1,Ldelta
G0and_lso=nnn2lso_reshape(g0and(:,:,:,:,:,:,l),Nlat,Nspin,Norb)
do ik=1,size(a)
dDelta_lso=nnn2lso_reshape(dDelta(:,:,:,:,:,:,l,ik),Nlat,Nspin,Norb)
dG0and_lso = (G0and_lso .x. dDelta_lso) .x. G0and_lso
dG0and(:,:,:,:,:,:,l,ik)=lso2nnn_reshape(dG0and_lso,Nlat,Nspin,Norb)
enddo
enddo
!
end function grad_g0and_replica
#endif
end MODULE ED_FIT_REPLICA