Line data Source code
1 : !--------------------------------------------------------------------------------
2 : ! Copyright (c) 2016 Peter Grünberg Institut, Forschungszentrum Jülich, Germany
3 : ! This file is part of FLEUR and available as free software under the conditions
4 : ! of the MIT license as expressed in the LICENSE file in more detail.
5 : !--------------------------------------------------------------------------------
6 :
7 : MODULE m_wann_read_umatrix
8 : use m_juDFT
9 : c************************************************************
10 : c Read in the transformation matrix between Bloch
11 : c and Wannier functions which is calculated by
12 : c wannier90.
13 : c Frank Freimuth, October 2006
14 : c************************************************************
15 : CONTAINS
16 0 : SUBROUTINE wann_read_umatrix2(
17 : > fullnkpts,num_wann,num_bands,
18 : > um_format,jspin,wan90version,
19 : < have_disentangled,
20 0 : < lwindow,ndimwin,
21 0 : < u_matrix_opt,u_matrix_tmp,m_matrix)
22 :
23 : USE m_constants
24 :
25 : implicit none
26 :
27 : integer,intent(in) :: fullnkpts,jspin
28 : integer,intent(in) :: num_wann,num_bands
29 : logical,intent(in) :: um_format
30 : integer,intent(in) :: wan90version
31 : logical,intent(out) :: have_disentangled
32 : logical,intent(out) :: lwindow(num_bands,fullnkpts)
33 : integer,intent(out) :: ndimwin(fullnkpts)
34 : complex,intent(out) :: u_matrix_opt(num_bands,num_wann,fullnkpts)
35 : complex,intent(out) :: u_matrix_tmp(num_wann,num_wann,fullnkpts)
36 : complex,intent(out),optional :: m_matrix(:,:,:,:)
37 :
38 0 : integer,allocatable :: exclude_bands(:)
39 : integer :: num_kpts,chk_unit,mp_grid(3),l
40 : integer :: i,j,k,nkp,ntmp,num_nnmax,nntot
41 : logical :: l_chk,l_umdat
42 : character(len=3) :: spin12(2)
43 : character(len=33) :: header
44 : real :: tmp_latt(3,3)
45 0 : real,allocatable :: tmp_kpt_latt(:,:)
46 : character(len=20) :: checkpoint
47 : real :: tmp_omi
48 : real :: omega_invariant
49 : character(len=2) :: spinspin12(2)
50 : data spin12/'WF1' , 'WF2'/
51 : data spinspin12/'.1','.2'/
52 :
53 : c**************************************************************
54 : c read in chk
55 : c**************************************************************
56 0 : call timestart("wann_read_umatrix2")
57 :
58 0 : write(oUnit,*)"read in chk"
59 0 : l_chk=.false.
60 0 : inquire (file=spin12(jspin)//'.chk',exist=l_chk)
61 0 : IF(.NOT.l_chk) then
62 0 : write(*,*)spin12(jspin)//'.chk'
63 : CALL juDFT_error("file chk not found",calledby
64 0 : + ="wann_read_umatrix")
65 : ENDIF
66 :
67 0 : chk_unit=152
68 : open(chk_unit,file=spin12(jspin)//'.chk',status='old',
69 0 : & form='unformatted')
70 :
71 : !write(*,*)'wan90version',wan90version
72 0 : if(wan90version.eq.1)then !wannier90 version 1.1
73 : ! Read comment line
74 0 : read(chk_unit) header
75 0 : write(oUnit,*)header
76 : ! Real and reciprocal lattice (units: Angstroem)
77 0 : read(chk_unit) ((tmp_latt(i,j),i=1,3),j=1,3) ! Real lattice
78 0 : write(oUnit,*)tmp_latt
79 0 : read(chk_unit) ((tmp_latt(i,j),i=1,3),j=1,3) ! Reciprocal lattice
80 0 : write(oUnit,*)tmp_latt
81 :
82 0 : read(chk_unit) num_kpts ! K-points
83 0 : write(oUnit,*)"num_kpts=",num_kpts
84 0 : if (num_kpts.ne.fullnkpts) CALL
85 : & juDFT_error("num_kpts.ne.fullnkpts",calledby
86 0 : + ="wann_read_umatrix")
87 : elseif((wan90version.eq.2).or.(wan90version.eq.3)
88 0 : & .or.(wan90version.eq.30).or.(wan90version.eq.31) ) then
89 : !wannier90 version 1.2 or wannier90 version 2.0
90 :
91 : ! Read comment line
92 0 : read(chk_unit) header
93 0 : write(oUnit,*)header
94 :
95 0 : read(chk_unit) ntmp
96 :
97 0 : read(chk_unit) ntmp
98 0 : allocate(exclude_bands(ntmp))
99 0 : read(chk_unit) (exclude_bands(i),i=1,ntmp) ! Excluded bands
100 :
101 : ! Real and reciprocal lattice (units: Angstroem)
102 0 : read(chk_unit) ((tmp_latt(i,j),i=1,3),j=1,3) ! Real lattice
103 0 : write(oUnit,*)tmp_latt
104 0 : read(chk_unit) ((tmp_latt(i,j),i=1,3),j=1,3) ! Reciprocal lattice
105 0 : write(oUnit,*)tmp_latt
106 :
107 0 : read(chk_unit) num_kpts ! K-points
108 0 : write(oUnit,*)"num_kpts=",num_kpts
109 0 : if (num_kpts.ne.fullnkpts) stop "num_kpts.ne.fullnkpts"
110 0 : read(chk_unit) (mp_grid(i),i=1,3) ! M-P grid
111 : else
112 : CALL juDFT_error("unknown wan90version",calledby
113 0 : + ="wann_read_umatrix")
114 : endif
115 :
116 0 : allocate(tmp_kpt_latt(3,fullnkpts))
117 0 : read(chk_unit) ((tmp_kpt_latt(i,nkp),i=1,3),nkp=1,fullnkpts)
118 0 : deallocate(tmp_kpt_latt)
119 :
120 0 : read(chk_unit) nntot ! nntot
121 0 : read(chk_unit) ntmp ! num_wann
122 0 : IF (ntmp/=num_wann) CALL juDFT_error("mismatch in num_wann"
123 0 : + ,calledby ="wann_read_umatrix")
124 :
125 0 : read(chk_unit) checkpoint ! checkpoint
126 0 : checkpoint=adjustl(trim(checkpoint))
127 0 : write(oUnit,*)checkpoint
128 :
129 :
130 :
131 :
132 :
133 0 : read(chk_unit) have_disentangled
134 : ! whether a disentanglement has been performed
135 :
136 0 : if (have_disentangled) then
137 0 : write(oUnit,*)"You used disentangling"
138 0 : write(oUnit,*)"Reading in disentangling information"
139 0 : read(chk_unit) omega_invariant ! omega invariant
140 0 : write(oUnit,*)omega_invariant
141 :
142 :
143 : ! U matrix opt
144 0 : read(chk_unit) ((lwindow(i,nkp),i=1,num_bands),nkp=1,num_kpts)
145 0 : read(chk_unit) (ndimwin(nkp),nkp=1,num_kpts)
146 0 : read(chk_unit)(((u_matrix_opt(i,j,nkp),i=1,num_bands)
147 0 : & ,j=1,num_wann),nkp=1,num_kpts)
148 :
149 : endif
150 :
151 :
152 :
153 0 : if(.not.wan90version.eq.0)then
154 :
155 0 : read(chk_unit) (((u_matrix_tmp(i,j,k),i=1,num_wann),
156 0 : & j=1,num_wann),k=1,num_kpts)
157 :
158 : endif
159 :
160 0 : if(present(m_matrix))then
161 0 : write(*,*)"nntot=",nntot
162 0 : write(*,*)"num_kpts=",num_kpts
163 0 : write(*,*)"num_wann=",num_wann
164 0 : read(chk_unit) ((((m_matrix(i,j,k,l),i=1,num_wann),
165 0 : & j=1,num_wann),k=1,nntot),l=1,num_kpts)
166 : endif
167 :
168 0 : close(chk_unit)
169 :
170 0 : if(wan90version.eq.0)then
171 : c****************************************************************
172 : c read in _um.dat (old version of wannier90)
173 : c****************************************************************
174 0 : l_umdat=.false.
175 0 : inquire (file=spin12(jspin)//'_um.dat',exist=l_umdat)
176 0 : IF(.NOT.l_umdat) CALL juDFT_error("where is your um_dat?",
177 0 : + calledby ="wann_read_umatrix")
178 :
179 0 : open(111,file=spin12(jspin)//'_um.dat',form='unformatted')
180 0 : write(oUnit,*)"read in um_data"
181 0 : read(111)header
182 0 : read(111)tmp_omi
183 0 : read(111) ntmp,num_kpts,num_nnmax
184 0 : IF(ntmp/=num_wann) CALL judft_error("mismatch in num_wann",
185 0 : + calledby ="wann_read_umatrix")
186 :
187 0 : read(111)(((u_matrix_tmp(i,j,k),i=1,num_wann),j=1,num_wann),
188 0 : & k=1,num_kpts)
189 0 : close(111)
190 : endif !wannier90 version
191 :
192 : c$$$ if (um_format) then
193 : c$$$ write(oUnit,*)"write um_data to formatted file WF1.umn/WF2.umn"
194 : c$$$ open(222,
195 : c$$$ & file=spin12(jspin)//'.umn',form='formatted')
196 : c$$$c write(222,*)header
197 : c$$$c write(222,*)tmp_omi
198 : c$$$c write(222,*)num_wann,num_kpts,num_nnmax
199 : c$$$ write(222,*)"transformation between Bloch and Wannier"
200 : c$$$ write(222,*)num
201 : c$$$ do k=1,num_kpts
202 : c$$$ do i=1,num_wann
203 : c$$$ do j=1,num_wann
204 : c$$$ write(222,'(i3,3x,i3,3x,i3,3x,f20.16,f20.16)')
205 : c$$$ & i,j,k,u_matrix_tmp(i,j,k)
206 : c$$$ enddo
207 : c$$$ enddo
208 : c$$$ enddo
209 : c$$$ close(222)
210 : c$$$ endif !um_format
211 :
212 0 : call timestop("wann_read_umatrix2")
213 0 : END SUBROUTINE wann_read_umatrix2
214 :
215 0 : SUBROUTINE wann_read_umatrix(
216 : > fullnkpts,num_wann,num_bands,
217 : > um_format,jspin,wan90version,
218 : < have_disentangled,
219 0 : < lwindow,ndimwin,u_matrix)
220 :
221 : USE m_constants
222 :
223 : implicit none
224 :
225 : integer,intent(in) :: fullnkpts,jspin
226 : integer,intent(in) :: num_wann,num_bands
227 : logical,intent(in) :: um_format
228 : integer,intent(in) :: wan90version
229 : logical,intent(out) :: have_disentangled
230 : logical,intent(out) :: lwindow(num_bands,fullnkpts)
231 : integer,intent(out) :: ndimwin(fullnkpts)
232 : complex,intent(out) :: u_matrix(num_bands,num_wann,fullnkpts)
233 :
234 0 : complex :: u_matrix_opt(num_bands,num_wann,fullnkpts)
235 0 : complex :: u_matrix_tmp(num_wann,num_wann,fullnkpts)
236 : integer :: i,j,k
237 : character(len=3) :: spin12(2)
238 : data spin12/'WF1' , 'WF2'/
239 :
240 0 : call timestart("wann_read_umatrix")
241 : call wann_read_umatrix2(
242 : > fullnkpts,num_wann,num_bands,
243 : > um_format,jspin,wan90version,
244 : < have_disentangled,
245 : < lwindow,ndimwin,u_matrix_opt,
246 0 : < u_matrix_tmp)
247 :
248 0 : u_matrix(:,:,:)=0.0
249 0 : if (have_disentangled) then
250 : c**************************************************************
251 : c calculate u_matrix for the entangled case
252 : c*************************************************************
253 0 : do k=1,fullnkpts
254 0 : do i=1,num_wann
255 0 : do j=1,num_wann
256 : u_matrix(:,j,k)= u_matrix(:,j,k)+
257 0 : & u_matrix_opt(:,i,k)*u_matrix_tmp(i,j,k)
258 : enddo
259 : enddo
260 : enddo
261 : else
262 0 : IF(num_bands/=num_wann) CALL juDFT_error("num_bands/=num_wann"
263 0 : + ,calledby ="wann_read_umatrix")
264 0 : u_matrix(:,:,:)=u_matrix_tmp(:,:,:)
265 : endif
266 :
267 0 : if (um_format) then
268 : write(oUnit,*)
269 0 : + "write um_data to formatted file WF1.umn/WF2.umn"
270 : open(222,
271 0 : & file=spin12(jspin)//'.umn',form='formatted')
272 : c write(222,*)header
273 : c write(222,*)tmp_omi
274 : c write(222,*)num_wann,num_kpts,num_nnmax
275 0 : write(222,*)"transformation between Bloch and Wannier"
276 0 : write(222,*)num_bands,fullnkpts,num_wann
277 0 : do k=1,fullnkpts
278 0 : do i=1,num_wann
279 0 : do j=1,num_bands
280 : write(222,'(i5,1x,i5,1x,i5,1x,f18.12,1x,f18.12)')
281 0 : & j,i,k,u_matrix(j,i,k)
282 : enddo
283 : enddo
284 : enddo
285 0 : close(222)
286 : endif !um_format
287 :
288 0 : call timestop("wann_read_umatrix")
289 :
290 0 : END SUBROUTINE wann_read_umatrix
291 : END MODULE m_wann_read_umatrix
|