-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathVCA_HAMILTONIAN_SPARSE_HxV.f90
289 lines (261 loc) · 8.97 KB
/
VCA_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
! > BUILD SPARSE HAMILTONIAN of the SECTOR
MODULE VCA_HAMILTONIAN_SPARSE_HxV
USE VCA_HAMILTONIAN_COMMON
implicit none
private
!>build sparse hamiltonian of the sector
!public :: buildH_c
public :: vca_buildh_main
!public :: vca_buildh_orbs
!
!
!>Sparse Mat-Vec product using stored sparse matrix
public :: spMatVec_main
!public :: spMatVec_orbs
#ifdef _MPI
public :: spMatVec_MPI_main
!public :: spMatVec_MPI_orbs
#endif
!
!
!> Related auxiliary routines:
!public :: setup_Hv_sector
!public :: delete_Hv_sector
integer :: i,iup,idw
integer :: j,jup,jdw
integer :: m,mup,mdw
integer :: ms,iud
integer :: impi
integer :: ilat,jlat,iorb,jorb,ispin,jspin,is,js,ilat_bath,iorb_bath
integer :: kp,k1,k2,k3,k4
integer :: ialfa,ibeta,indx
real(8) :: sg1,sg2,sg3,sg4
complex(8) :: htmp,htmpup,htmpdw
logical :: Jcondition
integer :: Nfoo,Nfoo2
complex(8),dimension(:,:,:,:,:),allocatable :: diag_hybr ![Nlat,Nlat_bath,Nspin,Norb,Norb_bath]
complex(8),dimension(:,:,:),allocatable :: bath_diag ![Nlat_bath,Nspin,Norb_bath]
contains
!####################################################################
! BUILD SPARSE HAMILTONIAN of the SECTOR
!####################################################################
subroutine vca_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(2*Ns_Ud) :: Indices ![2-2*Norb]
integer,dimension(Ns_Ud,Ns_Orb) :: Nups,Ndws ![1,Ns]-[Norb,1+Nbath]
integer,dimension(Nlat,Norb) :: Nup,Ndw !
!
nup=zero
ndw=zero
#ifdef _MPI
if(Mpistatus .AND. MpiComm == MPI_COMM_NULL)return
#endif
!
if(.not.Hstatus)stop "buildH_c ERROR: Hsector NOT set"
isector=Hsector
!
!
if(present(Hmat))&
call assert_shape(Hmat,[getdim(isector), getdim(isector)],"vca_buildh_main","Hmat")
!
if(spH0d%status)call sp_delete_matrix(spH0d)
!
!Get diagonal hybridization, bath energy
if(Nlat_bath>0 .and. Norb_bath>0)then
include "VCA_HAMILTONIAN/diag_hybr_bath.f90"
endif
!
#ifdef _MPI
if(MpiStatus)then
call sp_set_mpi_matrix(MpiComm,spH0d,mpiIstart,mpiIend,mpiIshift)
call sp_init_matrix(MpiComm,spH0d,Dim)
else
call sp_init_matrix(spH0d,Dim)
endif
#else
call sp_init_matrix(spH0d,Dim)
#endif
call sp_init_matrix(spH0dws(1),DimDw)
call sp_init_matrix(spH0ups(1),DimUp)
!
!
!-----------------------------------------------!
!LOCAL CLUSTER HAMILTONIAN TERMS
include "VCA_HAMILTONIAN/stored/H_cluster_local.f90"
!
!UP CLUSTER TERMS
include "VCA_HAMILTONIAN/stored/H_cluster_up.f90"
!
!DW CLUSTER TERMS
include "VCA_HAMILTONIAN/stored/H_cluster_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
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)) ! kron(Htmp_dw,eye(DimUp))
Hmat = Hmat + kronecker_product(one*eye(DimDw),Htmp_up) ! kron(eye(DimDw),Htmp_up)
!
deallocate(Htmp_up,Htmp_dw)
endif
!
if(allocated(diag_hybr))deallocate(diag_hybr)
if(allocated(bath_diag))deallocate(bath_diag)
return
!if(present(Hmat))then
!call assert_shape(Hmat,[dim,dim],"buildH_c","Hmat")
!Hmat=zero
!call sp_dump_matrix(spH0,Hmat)
!endif
end subroutine vca_buildh_main
!####################################################################
! TODO: ORBS
!####################################################################
!####################################################################
! SPARSE MAT-VEC PRODUCT USING STORED SPARSE MATRIX
!####################################################################
!+------------------------------------------------------------------+
!PURPOSE: Perform the matrix-vector product H*v used in the
!+------------------------------------------------------------------+
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
!
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
!
end subroutine spMatVec_main
!####################################################################
! TODO: ORBS
!####################################################################
!####################################################################
! MPI MATVEC
!####################################################################
#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) :: val
integer :: i,iup,idw,j,jup,jdw,jj
!local MPI
integer :: irank
!
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=zero
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)
end subroutine spMatVec_mpi_main
!####################################################################
! TODO: ORBS
!####################################################################
#endif
END MODULE VCA_HAMILTONIAN_SPARSE_HxV