-
Notifications
You must be signed in to change notification settings - Fork 1
/
VCA_VARS_GLOBAL.f90
341 lines (274 loc) · 13.2 KB
/
VCA_VARS_GLOBAL.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
MODULE VCA_VARS_GLOBAL
USE VCA_INPUT_VARS
USE VCA_SPARSE_MATRIX
!
USE SF_CONSTANTS
USE SF_ARRAYS, only: arange,linspace
USE SF_IOTOOLS, only: str
#ifdef _MPI
USE MPI
USE SF_MPI
#endif
implicit none
!-------------------- CUSTOM OBSERVABLE STRUCTURE ----------------------!
type observable
complex(8),dimension(:,:,:,:,:,:,:),allocatable :: sij
character(len=32) :: o_name
real(8) :: o_value
end type observable
type custom_observables
type(observable),dimension(:),allocatable :: item
integer :: N_asked
integer :: N_filled
logical :: init=.false.
end type custom_observables
!-------------------- EFFECTIVE BATH STRUCTURE ----------------------!
type effective_bath
complex(8),dimension(:,:,:,:,:,:),allocatable :: h !bath hamiltonian [Nlat_bath,Nlat_bath][Nspin][Nspin][Norb_bath][Norb_bath]
complex(8),dimension(:,:,:,:,:,:),allocatable :: v !coupling part [Nlat][Nlat_bath][Nspin][Nspin][Norb][Norb_bath]
logical :: status=.false.
end type effective_bath
!------------------ FULL HAMILTONIAN STRUCTURE ---------------------!
type full_espace
real(8),dimension(:),pointer :: e
complex(8),dimension(:,:),pointer :: M
end type full_espace
!---------------- SECTOR-TO-FOCK SPACE STRUCTURE -------------------!
type sector_map
integer,dimension(:),allocatable :: map
logical :: status=.false.
end type sector_map
interface map_allocate
module procedure :: map_allocate_scalar
module procedure :: map_allocate_vector
end interface map_allocate
interface map_deallocate
module procedure :: map_deallocate_scalar
module procedure :: map_deallocate_vector
end interface map_deallocate
!-------------- GMATRIX FOR FAST EVALUATION OF GF ------------------!
!note that we use a single Qmatrix here which must be intended as
!component corresponding to the poles.
type GFspectrum
complex(8),dimension(:),allocatable :: weight
complex(8),dimension(:),allocatable :: poles
end type GFspectrum
type GFchannel
type(GFspectrum),dimension(:),allocatable :: channel !N_channel = 2 (c,cdag), 4 (c,cdag,c pm cdag)
end type GFchannel
type GFmatrix
type(GFchannel),dimension(:),allocatable :: state !state_list%size = # of state in the spectrum
end type GFmatrix
interface GFmatrix_allocate
module procedure :: allocate_GFmatrix_Nstate
module procedure :: allocate_GFmatrix_Nchan
module procedure :: allocate_GFmatrix_Nexc
end interface GFmatrix_allocate
!------------------ ABTRACT INTERFACES PROCEDURES ------------------!
!SPARSE MATRIX-VECTOR PRODUCTS USED IN ED_MATVEC
!dbleMat*dbleVec
abstract interface
subroutine cc_sparse_HxV(Nloc,v,Hv)
integer :: Nloc
complex(8),dimension(Nloc) :: v
complex(8),dimension(Nloc) :: Hv
end subroutine cc_sparse_HxV
end interface
!-------------------------- ED VARIABLES --------------------------!
!SIZE OF THE PROBLEM
!=========================================================
integer,save :: Ns !Number of levels per spin
integer,save :: Nsectors !Number of sectors
integer,save :: Ns_orb
integer,save :: Ns_ud
!non-interacting cluster Hamiltonian and full system hopping matrix
!=========================================================
complex(8),dimension(:,:,:,:,:,:),allocatable :: impHloc ![Nlat][Nlat][Nspin][Nspin][Norb][Norb]
complex(8),dimension(:,:,:,:,:,:,:),allocatable :: impHk ![Nlat][Nlat][Nspin][Nspin][Norb][Norb][Nktot] (horizontal rectangular)
!Some maps between sectors and full Hilbert space (pointers) CHECK!!!!!!!!!!!!!!!!!!!!!
!=========================================================
integer,allocatable,dimension(:,:) :: getsector
integer,allocatable,dimension(:,:,:) :: getCsector
integer,allocatable,dimension(:,:,:) :: getCDGsector
integer,allocatable,dimension(:) :: getDim,getDimUp,getDimDw
integer,allocatable,dimension(:) :: getNup,getNdw
integer,allocatable,dimension(:,:) :: getBathStride
integer,allocatable,dimension(:,:) :: impIndex
logical,allocatable,dimension(:) :: twin_mask
logical,allocatable,dimension(:) :: sectors_mask
!Effective Bath used in the VCA code (this is opaque to user)
!=========================================================
type(effective_bath) :: vca_bath
!Eigenvalues,Eigenvectors FULL DIAGONALIZATION
!=========================================================
type(full_espace),dimension(:),allocatable :: espace
!Sparse matrix for Lanczos diagonalization.
!=========================================================
type(sparse_matrix_csr) :: spH0d !diagonal part
type(sparse_matrix_csr) :: spH0nd !non-diagonal part
type(sparse_matrix_csr),dimension(:),allocatable :: spH0ups,spH0dws !reduced UP and DW parts
procedure(cc_sparse_HxV),pointer :: spHtimesV_p=>null()
!Variables for DIAGONALIZATION
!=========================================================
integer,allocatable,dimension(:) :: neigen_sector
!--------------- LATTICE WRAP VARIABLES -----------------!
integer,allocatable,dimension(:,:) :: neigen_sectorii
integer,allocatable,dimension(:) :: neigen_totalii
logical :: trim_state_list=.false.
!Partition function, Omega potential, SFT potential,original lattice bandwidth
!=========================================================
real(8) :: zeta_function
real(8) :: omega_potential
real(8) :: sft_potential
!Cluster Green's functions
!(Nlat,Nlat,Nspin,Nspin,Norb,Norb,:)
!=========================================================
complex(8),allocatable,dimension(:,:,:,:,:,:,:) :: impGmats ![Nlat][Nlat][Nspin][Nspin][Norb][Norb][L]
complex(8),allocatable,dimension(:,:,:,:,:,:,:) :: impGreal ![Nlat][Nlat][Nspin][Nspin][Norb][Norb][L]
!
complex(8),allocatable,dimension(:,:,:,:,:,:,:) :: impG0mats ![Nlat][Nlat][Nspin][Nspin][Norb][Norb][L]
complex(8),allocatable,dimension(:,:,:,:,:,:,:) :: impG0real ![Nlat][Nlat][Nspin][Nspin][Norb][Norb][L]
!
complex(8),allocatable,dimension(:,:,:,:,:,:,:) :: impSmats ![Nlat][Nlat][Nspin][Nspin][Norb][Norb][L]
complex(8),allocatable,dimension(:,:,:,:,:,:,:) :: impSreal ![Nlat][Nlat][Nspin][Nspin][Norb][Norb][L]
!
type(GFmatrix),allocatable,dimension(:,:,:,:,:,:) :: impGmatrix
!Spin Susceptibilities
!=========================================================
real(8),allocatable,dimension(:,:) :: spinChi_tau
complex(8),allocatable,dimension(:,:) :: spinChi_w
complex(8),allocatable,dimension(:,:) :: spinChi_iv
!Cluster local observables:
!=========================================================
real(8),dimension(:,:),allocatable :: imp_dens ![Nlat][Norb]
real(8),dimension(:,:),allocatable :: imp_dens_up ![Nlat][Norb]
real(8),dimension(:,:),allocatable :: imp_dens_dw ![Nlat][Norb]
real(8),dimension(:,:),allocatable :: imp_docc ![Nlat][Norb]
!Custom observables:
!=========================================================
type(custom_observables) :: custom_o
!Suffix string attached to the output files.
!=========================================================
character(len=64) :: file_suffix=""
logical :: offdiag_gf_flag=.false.
!This is the internal Mpi Communicator and variables.
!=========================================================
#ifdef _MPI
integer :: MpiComm_Global=MPI_COMM_NULL
integer :: MpiComm=MPI_COMM_NULL
#endif
integer :: MpiGroup_Global=MPI_GROUP_NULL
integer :: MpiGroup=MPI_GROUP_NULL
logical :: MpiStatus=.false.
logical :: MpiMaster=.true.
integer :: MpiRank=0
integer :: MpiSize=1
integer,allocatable,dimension(:) :: MpiMembers
integer :: mpiQup=0
integer :: mpiRup=0
integer :: mpiQdw=0
integer :: mpiRdw=0
integer :: mpiQ=0
integer :: mpiR=0
integer :: mpiIstart
integer :: mpiIend
integer :: mpiIshift
logical :: mpiAllThreads=.true.
!Frequency and time arrays:
!=========================================================
!real(8),dimension(:),allocatable :: wm,tau,wr,vm
!flag for finite temperature calculation
!=========================================================
logical :: finiteT
contains
!> Get stride position in the one-particle many-body space
function index_stride_lso(ilat,ispin,iorb) result(indx)
integer :: ilat
integer :: iorb
integer :: ispin
integer :: indx
indx = iorb + (ilat-1)*Norb + (ispin-1)*Norb*Nlat
end function index_stride_lso
!>Allocate and Deallocate Hilbert space maps (sector<-->Fock)
subroutine map_allocate_scalar(H,N)
type(sector_map) :: H
integer :: N
if(H%status) call map_deallocate_scalar(H)
allocate(H%map(N))
H%status=.true.
end subroutine map_allocate_scalar
!
subroutine map_allocate_vector(H,N)
type(sector_map),dimension(:) :: H
integer,dimension(size(H)) :: N
integer :: i
do i=1,size(H)
call map_allocate_scalar(H(i),N(i))
enddo
end subroutine map_allocate_vector
subroutine map_deallocate_scalar(H)
type(sector_map) :: H
if(.not.H%status)then
write(*,*) "WARNING map_deallocate_scalar: H is not allocated"
return
endif
if(allocated(H%map))deallocate(H%map)
H%status=.false.
end subroutine map_deallocate_scalar
!
subroutine map_deallocate_vector(H)
type(sector_map),dimension(:) :: H
integer :: i
do i=1,size(H)
call map_deallocate_scalar(H(i))
enddo
end subroutine map_deallocate_vector
!=========================================================
subroutine vca_set_MpiComm(comm)
#ifdef _MPI
integer :: comm,ierr
MpiComm_Global = comm
MpiComm = MpiComm_Global
MpiStatus = .true.
MpiSize = get_Size_MPI(MpiComm_Global)
MpiRank = get_Rank_MPI(MpiComm_Global)
MpiMaster = get_Master_MPI(MpiComm_Global)
call Mpi_Comm_group(MpiComm_Global,MpiGroup_Global,ierr)
#else
integer,optional :: comm
#endif
end subroutine vca_set_MpiComm
subroutine vca_del_MpiComm()
#ifdef _MPI
MpiComm_Global = MPI_UNDEFINED
MpiComm = MPI_UNDEFINED
MpiStatus = .false.
MpiSize = 1
MpiRank = 0
MpiMaster = .true.
#endif
end subroutine vca_del_MpiComm
!Allocate the channels in GFmatrix structure
subroutine allocate_gfmatrix_Nstate(self,Nstate)
type(GFmatrix) :: self
integer :: Nstate
if(allocated(self%state))deallocate(self%state)
allocate(self%state(Nstate))
end subroutine allocate_gfmatrix_Nstate
subroutine allocate_gfmatrix_Nchan(self,istate,Nchan)
type(GFmatrix) :: self
integer :: istate,Nchan
if(allocated(self%state(istate)%channel))deallocate(self%state(istate)%channel)
allocate(self%state(istate)%channel(Nchan))
end subroutine allocate_gfmatrix_Nchan
!Allocate the Excitations spectrum at a given channel
subroutine allocate_gfmatrix_Nexc(self,istate,ichan,Nexc)
type(GFmatrix) :: self
integer :: istate,ichan
integer :: Nexc
if(allocated(self%state(istate)%channel(ichan)%weight))deallocate(self%state(istate)%channel(ichan)%weight)
if(allocated(self%state(istate)%channel(ichan)%poles))deallocate(self%state(istate)%channel(ichan)%poles)
allocate(self%state(istate)%channel(ichan)%weight(Nexc))
allocate(self%state(istate)%channel(ichan)%poles(Nexc))
end subroutine allocate_gfmatrix_Nexc
END MODULE VCA_VARS_GLOBAL