-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathED_HAMILTONIAN_SPARSE_HxV.f90
327 lines (299 loc) · 9.47 KB
/
ED_HAMILTONIAN_SPARSE_HxV.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
! > BUILD SPARSE HAMILTONIAN of the SECTOR
MODULE ED_HAMILTONIAN_SPARSE_HxV
USE ED_HAMILTONIAN_COMMON
USE ED_SETUP
USE ED_AUX_FUNX
implicit none
private
!>Sparse Matric constructors
public :: ed_buildh_main
!>Sparse Mat-Vec product using stored sparse matrix
public :: spMatVec_main
#ifdef _MPI
public :: spMatVec_MPI_main
#endif
integer :: i,iup,idw
integer :: j,jup,jdw
integer :: m,mup,mdw
integer :: ms,iud
integer :: impi
integer :: ilat,jlat,iorb,jorb,ispin,jspin,ibath,is,js
integer :: k1,k2,k3,k4
integer :: ialfa,ibeta,indx
real(8) :: sg1,sg2,sg3,sg4
complex(8) :: htmp,htmpup,htmpdw
logical :: Jcondition
contains
!####################################################################
! BUILD SPARSE HAMILTONIAN of the SECTOR
!####################################################################
subroutine ed_buildh_main(isector,Hmat)
integer :: isector
complex(8),dimension(:,:),optional :: Hmat
complex(8),dimension(:,:),allocatable :: Htmp_up,Htmp_dw,Hrdx
integer,dimension(Ns) :: ibup,ibdw
integer,dimension(Nlat,Norb) :: Nup,Ndw
real(8),dimension(Nlat,Nspin,Norb,Nbath) :: diag_hybr
real(8),dimension(Nlat,Nspin,Norb,Nbath) :: bath_diag
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb,Nbath) :: Hbath_reconstructed
!
nup=zero
ndw=zero
#ifdef _MPI
if(Mpistatus .AND. MpiComm == MPI_COMM_NULL)return
#endif
!
if(.not.Hstatus)stop "ed_buildh_main ERROR: Hsector NOT set"
isector=Hsector
!
if(present(Hmat))&
call assert_shape(Hmat,[getdim(isector), getdim(isector)],"ed_buildh_main","Hmat")
!
diag_hybr=0d0
bath_diag=0d0
do ibath=1,Nbath
Hbath_reconstructed(:,:,:,:,:,:,ibath)=Hbath_build(dmft_bath%item(ibath)%lambda)
do ilat=1,Nlat
do ispin=1,Nspin
do iorb=1,Norb
diag_hybr(ilat,ispin,iorb,ibath)=dmft_bath%item(ibath)%v(index_stride_lso(ilat,ispin,iorb))
bath_diag(ilat,ispin,iorb,ibath)=DREAL(Hbath_Reconstructed(ilat,ilat,ispin,ispin,iorb,iorb,ibath))
enddo
enddo
enddo
enddo
!
#ifdef _MPI
if(MpiStatus)then
call sp_set_mpi_matrix(MpiComm,spH0d,mpiIstart,mpiIend,mpiIshift)
call sp_init_matrix(MpiComm,spH0d,Dim)
if(Jhflag)then
call sp_set_mpi_matrix(MpiComm,spH0nd,mpiIstart,mpiIend,mpiIshift)
call sp_init_matrix(MpiComm,spH0nd,Dim)
endif
else
call sp_init_matrix(spH0d,Dim)
if(Jhflag)call sp_init_matrix(spH0nd,Dim)
endif
#else
call sp_init_matrix(spH0d,Dim)
if(Jhflag)call sp_init_matrix(spH0nd,Dim)
#endif
call sp_init_matrix(spH0dws(1),DimDw)
call sp_init_matrix(spH0ups(1),DimUp)
!
!-----------------------------------------------!
!LOCAL HAMILTONIAN TERMS
include "ED_HAMILTONIAN/sparse/H_local.f90"
!
!NON-LOCAL HAMILTONIAN TERMS
if(jhflag)then
include "ED_HAMILTONIAN/sparse/H_non_local.f90"
endif
!
!UP TERMS
include "ED_HAMILTONIAN/sparse/H_up.f90"
!
!DW TERMS
include "ED_HAMILTONIAN/sparse/H_dw.f90"
!-----------------------------------------------!
!
if(present(Hmat))then
Hmat = zero
allocate(Htmp_up(DimUp,DimUp));Htmp_up=zero
allocate(Htmp_dw(DimDw,DimDw));Htmp_dw=zero
!
#ifdef _MPI
if(MpiStatus)then
call sp_dump_matrix(MpiComm,spH0d,Hmat)
else
call sp_dump_matrix(spH0d,Hmat)
endif
#else
call sp_dump_matrix(spH0d,Hmat)
#endif
!
if(Jhflag)then
allocate(Hrdx(Dim,Dim));Hrdx=zero
#ifdef _MPI
if(MpiStatus)then
call sp_dump_matrix(MpiComm,spH0nd,Hrdx)
else
call sp_dump_matrix(spH0nd,Hrdx)
endif
#else
call sp_dump_matrix(spH0nd,Hrdx)
#endif
Hmat = Hmat + Hrdx
deallocate(Hrdx)
endif
!
call sp_dump_matrix(spH0ups(1),Htmp_up)
call sp_dump_matrix(spH0dws(1),Htmp_dw)
Hmat = Hmat + kronecker_product(Htmp_dw,one*eye(DimUp))
Hmat = Hmat + kronecker_product(one*eye(DimDw),Htmp_up)
!
deallocate(Htmp_up,Htmp_dw)
endif
!
return
!
end subroutine ed_buildh_main
!####################################################################
! SPARSE MAT-VEC PRODUCT USING STORED SPARSE MATRIX
!####################################################################
!+------------------------------------------------------------------+
!PURPOSE: Perform the matrix-vector product H*v used in the
! - serial
! - MPI
!+------------------------------------------------------------------+
subroutine spMatVec_main(Nloc,v,Hv)
integer :: Nloc
complex(8),dimension(Nloc) :: v
complex(8),dimension(Nloc) :: Hv
complex(8) :: val
integer :: i,iup,idw,j,jup,jdw,jj
!
!
Hv=zero
!
!Local:
do i = 1,Nloc
do j=1,spH0d%row(i)%Size
Hv(i) = Hv(i) + spH0d%row(i)%vals(j)*v(spH0d%row(i)%cols(j))
enddo
enddo
!
!DW:
do iup=1,DimUp
!
do idw=1,DimDw
i = iup + (idw-1)*DimUp
do jj=1,spH0dws(1)%row(idw)%Size
jup = iup
jdw = spH0dws(1)%row(idw)%cols(jj)
val = spH0dws(1)%row(idw)%vals(jj)
j = jup + (jdw-1)*DimUp
Hv(i) = Hv(i) + val*V(j)
enddo
enddo
!
enddo
!
!UP:
do idw=1,DimDw
!
do iup=1,DimUp
i = iup + (idw-1)*DimUp
do jj=1,spH0ups(1)%row(iup)%Size
jup = spH0ups(1)%row(iup)%cols(jj)
jdw = idw
val = spH0ups(1)%row(iup)%vals(jj)
j = jup + (jdw-1)*DimUp
Hv(i) = Hv(i) + val*V(j)
enddo
enddo
!
enddo
!
!Non-Local:
if(jhflag)then
do i = 1,Nloc
do j=1,spH0nd%row(i)%Size
val = spH0nd%row(i)%vals(j)
jj = spH0nd%row(i)%cols(j)
Hv(i) = Hv(i) + val*V(jj)
enddo
enddo
endif
!
end subroutine spMatVec_main
#ifdef _MPI
subroutine spMatVec_mpi_main(Nloc,v,Hv)
integer :: Nloc
complex(8),dimension(Nloc) :: v
complex(8),dimension(Nloc) :: Hv
!
integer :: N
complex(8),dimension(:),allocatable :: vt,Hvt
complex(8),dimension(:),allocatable :: vin
complex(8) :: val
integer :: i,iup,idw,j,jup,jdw,jj
!local MPI
integer :: irank,MpiIerr
integer,allocatable,dimension(:) :: Counts
integer,allocatable,dimension(:) :: Offset
!
! if(MpiComm==Mpi_Comm_Null)return
! if(MpiComm==MPI_UNDEFINED)stop "spMatVec_mpi_cc ERROR: MpiComm = MPI_UNDEFINED"
if(.not.MpiStatus)stop "spMatVec_mpi_cc ERROR: MpiStatus = F"
!
!Evaluate the local contribution: Hv_loc = Hloc*v
Hv=zero
do i=1,Nloc !==spH0%Nrow
do j=1,spH0d%row(i)%Size
Hv(i) = Hv(i) + spH0d%row(i)%vals(j)*v(i)
end do
end do
!
!Non-local terms.
!UP part: contiguous in memory.
do idw=1,MpiQdw
do iup=1,DimUp
i = iup + (idw-1)*DimUp
hxv_up: do jj=1,spH0ups(1)%row(iup)%Size
jup = spH0ups(1)%row(iup)%cols(jj)
jdw = idw
val = spH0ups(1)%row(iup)%vals(jj)
j = jup + (idw-1)*DimUp
Hv(i) = Hv(i) + val*v(j)
end do hxv_up
enddo
end do
!
!DW part: non-contiguous in memory -> MPI transposition
!Transpose the input vector as a whole:
mpiQup=DimUp/MpiSize
if(MpiRank<mod(DimUp,MpiSize))MpiQup=MpiQup+1
!
allocate(vt(mpiQup*DimDw)) ;vt=zero
allocate(Hvt(mpiQup*DimDw));Hvt=zero
call vector_transpose_MPI(DimUp,MpiQdw,v,DimDw,MpiQup,vt)
Hvt=0d0
do idw=1,MpiQup !<= Transposed order: column-wise DW <--> UP
do iup=1,DimDw !<= Transposed order: column-wise DW <--> UP
i = iup + (idw-1)*DimDw
hxv_dw: do jj=1,spH0dws(1)%row(iup)%Size
jup = spH0dws(1)%row(iup)%cols(jj)
jdw = idw
j = jup + (jdw-1)*DimDw
val = spH0dws(1)%row(iup)%vals(jj)
Hvt(i) = Hvt(i) + val*vt(j)
end do hxv_dw
enddo
end do
deallocate(vt) ; allocate(vt(DimUp*mpiQdw)) ; vt=zero
call vector_transpose_MPI(DimDw,mpiQup,Hvt,DimUp,mpiQdw,vt)
Hv = Hv + Vt
deallocate(vt)
!
!
!Non-Local:
if(jhflag)then
N = 0
call AllReduce_MPI(MpiComm,Nloc,N)
!
allocate(vt(N)) ; vt = zero
call allgather_vector_MPI(MpiComm,v,vt)
!
do i=1,Nloc
matmul: do j=1,spH0nd%row(i)%Size
Hv(i) = Hv(i) + spH0nd%row(i)%vals(j)*Vt(spH0nd%row(i)%cols(j))
enddo matmul
enddo
deallocate(Vt)
endif
!
end subroutine spMatVec_mpi_main
#endif
end MODULE ED_HAMILTONIAN_SPARSE_HXV