-
Notifications
You must be signed in to change notification settings - Fork 3
/
ED_FIT_CHI2.f90
115 lines (102 loc) · 4.06 KB
/
ED_FIT_CHI2.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
MODULE ED_FIT_CHI2
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
USE ED_FIT_REPLICA
USE ED_FIT_GENERAL
implicit none
private
interface ed_chi2_fitgf
module procedure chi2_fitgf_generic_normal
#if __GFORTRAN__ && __GNUC__ > 8
!RDMFT_WRAPPER
module procedure chi2_fitgf_lattice_normal
#endif
end interface ed_chi2_fitgf
public :: ed_chi2_fitgf
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
!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
!+----------------------------------------------------------------------+
!PURPOSE : Chi^2 fit of the G0/Delta
!+----------------------------------------------------------------------+
subroutine chi2_fitgf_generic_normal(fg,bath)
complex(8),dimension(:,:,:,:,:,:,:) :: fg ![Nlat][Nlat][Nspin][Nspin][Norb][Norb][Niw]
real(8),dimension(:) :: bath
!
call assert_shape(fg,[Nlat,Nlat,Nspin,Nspin,Norb,Norb,size(fg,7)],"chi2_fitgf_generic_normal","fg")
!
select case(cg_method)
case default
stop "ED Error: cg_method > 1"
case (0)
if(ed_verbose>2)write(LOGfile,"(A,I1,A,A)")"\Chi2 fit with CG-nr and CG-weight: ",cg_weight," on: ",cg_scheme
case (1)
if(ed_verbose>2)write(LOGfile,"(A,I1,A,A)")"\Chi2 fit with CG-minimize and CG-weight: ",cg_weight," on: ",cg_scheme
end select
!
select case(bath_type)
case('replica')
call chi2_fitgf_replica(fg,bath)
case('general')
call chi2_fitgf_general(fg,bath)
end select
!
!set trim_state_list to true after the first fit has been done: this
!marks the ends of the cycle of the 1st DMFT loop.
trim_state_list=.true.
!
end subroutine chi2_fitgf_generic_normal
#if __GFORTRAN__ && __GNUC__ > 8
!+----------------------------------------------------------------------!
! PURPOSE: given a number of independent baths, evaluate N independent
! Delta/G0 functions and fit them to update the effective baths for ED.
!+----------------------------------------------------------------------!
!RDMFT WRAPPER:
subroutine chi2_fitgf_lattice_normal(fg,bath)
real(8),dimension(:,:) :: bath
complex(8),dimension(:,:,:,:,:,:,:,:) :: fg
!MPI auxiliary vars
integer :: isites
integer :: Nsites
character(len=5) :: tmp_suffix
!
! Check dimensions !
Nsites=size(bath,1)
call assert_shape(fg,[Nsites,Nlat,Nlat,Nspin,Nspin,Norb,Norb,size(fg,8)],"chi2_fitgf_generic_normal","fg")
!
!
do isites = 1, Nsites
!
ed_file_suffix=reg(ineq_site_suffix)//str(isites,site_indx_padding)
!
call chi2_fitgf_generic_normal(fg(isites,:,:,:,:,:,:,:),bath(isites,:))
!
end do
!
!
ed_file_suffix=""
end subroutine chi2_fitgf_lattice_normal
#endif
end MODULE ED_FIT_CHI2