Octopus
wannier90_interface.F90
Go to the documentation of this file.
1!! Copyright (C) 2017-2019 H. Huebener, N. Tancogne-Dejean
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use batch_oct_m
24 use comm_oct_m
26 use cube_oct_m
28 use debug_oct_m
30 use fft_oct_m
31 use global_oct_m
32 use grid_oct_m
33 use io_oct_m
36 use ions_oct_m
37 use, intrinsic :: iso_fortran_env
42 use loct_oct_m
44 use mesh_oct_m
47 use mpi_oct_m
51 use parser_oct_m
56 use space_oct_m
57 use string_oct_m
62 use types_oct_m
63 use unit_oct_m
65 use utils_oct_m
68
69 implicit none
70
71 integer :: w90_what, w90_mode, w90_what_default
72
73 integer :: ierr
74 integer :: dim, idim
75 integer :: ii, nik, iter, nst
76
77 type(restart_t) :: restart
78 type(electrons_t), pointer :: sys
79 logical :: w90_spinors, scdm_proj, w90_scdm
80 integer :: w90_nntot, w90_num_bands, w90_num_kpts ! w90 input parameters
81 integer, allocatable :: w90_nnk_list(:,:) !
82 character(len=80) :: w90_prefix ! w90 input file prefix
83 integer :: w90_num_wann ! input paramter
84 real(real64), allocatable :: w90_proj_centers(:,:) ! projections centers
85 integer, allocatable :: w90_proj_lmr(:,:) ! definitions for real valued Y_lm*R_r
86 integer :: w90_nproj ! number of such projections
87 integer, allocatable :: w90_spin_proj_component(:) ! up/down flag
88 real(real64), allocatable :: w90_spin_proj_axis(:,:) ! spin axis (not implemented)
89 integer :: w90_num_exclude
90 logical, allocatable :: exclude_list(:) ! list of excluded bands
91 integer, allocatable :: band_index(:) ! band index after exclusion
92 logical :: read_td_states
93 integer :: w90_spin_channel
94
95 ! scdm variables
96 integer, allocatable :: jpvt(:)
97 complex(real64), allocatable :: uk(:,:,:) ! SCDM-Wannier gauge matrices U(k)
98 complex(real64), allocatable :: psi(:,:)
99 complex(real64), allocatable :: chi(:,:), chi_diag(:,:),chi2(:,:)
100 real(real64), allocatable :: chi_eigenval(:), occ_temp(:)
101 real(real64) :: scdm_mu, scdm_sigma, smear, kvec(3), factor(3)
102 integer :: ist, jst, ik, idir
103
104 integer(int64) :: how
105
106 call global_init()
107 call parser_init()
108
109 call messages_init()
110 call io_init()
111
113
115
118
119 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
120 sys => electrons_t(global_namespace, generate_epot=.false.)
121 call sys%init_parallelization(mpi_world)
122
123 !%Variable Wannier90Prefix
124 !%Type string
125 !%Default w90
126 !%Section Utilities::oct-wannier90
127 !%Description
128 !% Prefix for wannier90 files
129 !%End
130 call parse_variable(global_namespace, 'Wannier90Prefix', 'w90', w90_prefix)
131 if (w90_prefix == 'w90') then
132 message(1) = "oct-wannier90: the prefix is set by default to w90"
133 call messages_info(1)
134 end if
135
136 !%Variable Wannier90Mode
137 !%Type integer
138 !%Default 0
139 !%Section Utilities::oct-wannier90
140 !%Description
141 !% Specifies which stage of the Wannier90 interface to use
142 !%Option none 0
143 !% Nothing is done.
144 !%Option w90_setup 1
145 !% Writes parts of the wannier90 input file <tt>w90_prefix.win</tt> corresponding to
146 !% the octopus inp file. Importantly it generates the correct form of Monkhorst-Pack mesh
147 !% written to the file w90_kpoints that has to be used in a gs calculation of Octopus by
148 !% as <tt> include w90_kpoints </tt> instead of the <tt>%KpointsGrid</tt> block.
149 !%Option w90_output 2
150 !% Generates the relevant files for a wannier90 run, specified by the variable <tt>W90_interface_files</tt>.
151 !% This needs files previously generated
152 !% by <tt>wannier90.x -pp w90 </tt>
153 !%Option w90_wannier 3
154 !% Parse the output of wannier90 to generate the Wannier states on the real-space grid.
155 !% The states will be written in the folder wannier. By default, the states are written as
156 !% binary files, similar to the Kohn-Sham states.
157 !%
158 !% Not implemented for spinor states.
159 !%End
160 call parse_variable(global_namespace, 'Wannier90Mode', 0, w90_mode)
161
162 if (w90_mode == 0) then
163 message(1) = "Wannier90Mode must be set to a value different from 0."
164 call messages_fatal(1)
165 end if
166
167 !%Variable Wannier90Files
168 !%Type flag
169 !%Default w90_mmn + w90_amn + w90_eig
170 !%Section Utilities::oct-wannier90
171 !%Description
172 !% Specifies which files to generate.
173 !% Example: <tt>w90_mmn + w90_unk</tt>
174 !%Option w90_mmn bit(1)
175 !% (see Wannier90 documentation)
176 !%Option w90_unk bit(2)
177 !% (see Wannier90 documentation)
178 !%Option w90_amn bit(3)
179 !% (see Wannier90 documentation)
180 !%Option w90_eig bit(4)
181 !% Eigenvalues. See Wannier90 documentation for more details.
182 !%End
183 w90_what_default = option__wannier90files__w90_mmn + option__wannier90files__w90_amn + option__wannier90files__w90_eig
184 call parse_variable(global_namespace, 'Wannier90Files', w90_what_default, w90_what)
185
186 !%Variable Wannier90UseTD
187 !%Type logical
188 !%Default no
189 !%Section Utilities::oct-wannier90
190 !%Description
191 !% By default oct-wannier90 uses the ground-state states to compute the necessary information.
192 !% By setting this variable to yes, oct-wannier90 will use the TD states instead.
193 !%End
194 call parse_variable(global_namespace, 'Wannier90UseTD', .false., read_td_states)
195
196 !%Variable Wannier90UseSCDM
197 !%Type logical
198 !%Default no
199 !%Section Utilities::oct-wannier90
200 !%Description
201 !% By default oct-wannier90 uses the projection method to generate the .amn file.
202 !% By setting this variable to yes, oct-wannier90 will use SCDM method instead.
203 !%End
204 call parse_variable(global_namespace, 'Wannier90UseSCDM', .false., w90_scdm)
205 if (w90_scdm) then
206 !%Variable SCDMsigma
207 !%Type float
208 !%Default 0.2
209 !%Section Utilities::oct-wannier90
210 !%Description
211 !% Broadening of SCDM smearing function.
212 !%End
213 call parse_variable(global_namespace, 'SCDMsigma', 0.2_real64, scdm_sigma)
214
215 !%Variable SCDMmu
216 !%Type float
217 !%Section Utilities::oct-wannier90
218 !%Description
219 !% Energy range up to which states are considered for SCDM.
220 !%End
221 call parse_variable(global_namespace, 'SCDMmu', m_huge, scdm_mu)
222 end if
223
224 if (sys%kpoints%use_symmetries) then
225 message(1) = 'oct-wannier90: k-points symmetries are not allowed'
226 call messages_fatal(1)
227 end if
228 if (sys%kpoints%use_time_reversal) then
229 message(1) = 'oct-wannier90: time-reversal symmetry is not allowed'
230 call messages_fatal(1)
231 end if
232 if (sys%kpoints%reduced%nshifts > 1) then
233 message(1) = 'oct-wannier90: Wannier90 does not allow for multiple shifts of the k-point grid'
234 call messages_fatal(1)
235 end if
236
237
238 if (sys%st%d%ispin /= unpolarized) then
239 call messages_experimental("oct-wannier90 with SpinComponnents /= unpolarized")
240 end if
241
242 w90_spinors = .false.
243 w90_spin_channel = 1
244
245 ! "Correct" rlattice/klattice for Wannier90 (3D periodicity assumed here).
246 factor = m_one
247 do idir = sys%space%periodic_dim+1, sys%space%dim
248 factor(idir) = m_two * sys%gr%box%bounding_box_l(idir)
249 end do
250 call sys%ions%latt%scale(factor)
251
252 ! create setup files
253 select case (w90_mode)
254 case (option__wannier90mode__w90_setup)
255 call wannier90_setup(sys%ions, sys%kpoints, sys%space)
256
257 ! load states and calculate interface files
258 case (option__wannier90mode__w90_output)
259 call wannier90_output()
260
261 case (option__wannier90mode__w90_wannier)
263
264 ! normal interface run
265 call states_elec_allocate_wfns(sys%st, sys%gr, wfs_type = type_cmplx, skip=exclude_list)
266 if (read_td_states) then
268 sys%mc, ierr, sys%gr)
269 else
271 sys%mc, ierr, sys%gr)
272 end if
273
274 if (ierr == 0) then
275 call states_elec_look(restart, nik, dim, nst, ierr)
276 if (dim == sys%st%d%dim .and. nik == sys%kpoints%reduced%npoints .and. nst == sys%st%nst) then
277 call states_elec_load(restart, global_namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
278 fixed_occ=.true., ierr=ierr, iter=iter, label = ": wannier90", skip=exclude_list)
279 else
280 write(message(1),'(a)') 'Restart structure not commensurate.'
281 call messages_fatal(1)
282 end if
283 end if
284 call restart_end(restart)
285
286 call generate_wannier_states(sys%space, sys%gr, sys%ions, sys%st, sys%kpoints)
287 case default
288 message(1) = "Wannier90Mode is set to an unsupported value."
289 call messages_fatal(1)
290 end select
291
292 safe_deallocate_a(exclude_list)
293 safe_deallocate_a(band_index)
294 safe_deallocate_a(w90_nnk_list)
295 safe_deallocate_a(w90_proj_centers)
296 safe_deallocate_a(w90_proj_lmr)
297
298 safe_deallocate_p(sys)
299 call fft_all_end()
300 call io_end()
302 call messages_end()
303 call parser_end()
304 call global_end()
305
306contains
307
308 ! --------------------------------------------------------------------------
309 subroutine wannier90_setup(ions, kpoints, space)
310 type(ions_t), intent(in) :: ions
311 type(kpoints_t), intent(in) :: kpoints
312 class(space_t), intent(in) :: space
313
314 character(len=80) :: filename
315 integer :: w90_win, ia, axis(3), npath
316
317 push_sub(wannier90_setup)
318
319 assert(space%dim == 3)
320
321 ! open win file
322 filename = trim(adjustl(w90_prefix)) //'.win'
323 w90_win = io_open(trim(filename), global_namespace, action='write')
324
325 write(w90_win,'(a)') '# this file has been created by the Octopus wannier90 utility'
326 write(w90_win,'(a)') ' '
327
328 ! write direct lattice vectors (in angstrom)
329 write(w90_win,'(a)') 'begin unit_cell_cart'
330 write(w90_win,'(a)') 'Ang'
331 do idim = 1,3
332 write(w90_win,'(f13.8,f13.8,f13.8)') units_from_atomic(unit_angstrom, ions%latt%rlattice(1:3,idim))
333 end do
334 write(w90_win,'(a)') 'end unit_cell_cart'
335 write(w90_win,'(a)') ' '
336
337 write(w90_win,'(a)') 'begin atoms_frac'
338 do ia = 1, ions%natoms
339 write(w90_win,'(a,2x,f13.8,f13.8,f13.8)') trim(ions%atom(ia)%label), ions%latt%cart_to_red(ions%pos(:, ia))
340 end do
341 write(w90_win,'(a)') 'end atoms_frac'
342 write(w90_win,'(a)') ' '
343
344 ! This is a default value. In practice, one should use projections
345 write(w90_win,'(a)') 'use_bloch_phases = .true.'
346 write(w90_win,'(a)') ' '
347
348 write(w90_win,'(a10,i4)') 'num_bands ', sys%st%nst
349 write(w90_win,'(a9,i4)') 'num_wann ', sys%st%nst
350 write(w90_win,'(a)') ' '
351
352 if (sys%st%d%ispin == spinors) then
353 write(w90_win,'(a)') 'spinors = .true.'
354 end if
355 if (sys%st%d%ispin == spin_polarized) then
356 write(w90_win, '(a)') 'spin = up'
357 end if
358
359 ! This is for convenience. This is needed for plotting the Wannier states, if requested.
360 write(w90_win,'(a)') 'write_u_matrices = .true.'
361 write(w90_win,'(a)') 'write_xyz = .true.'
362 write(w90_win,'(a)') ' '
363
364 if (kpoints%reduced%npoints == 1) then
365 write(w90_win,'(a)') 'gamma_only = .true.'
366 write(w90_win,'(a)') ' '
367 else
368 if (.not. parse_is_defined(global_namespace, 'KPointsGrid')) then
369 message(1) = 'oct-wannier90: need Monkhorst-Pack grid. Please specify %KPointsGrid'
370 call messages_fatal(1)
371 end if
372
373 ! In case the user used also a k-point path, we ignore it
374 npath = kpoints%nkpt_in_path()
375
376 axis(1:3) = kpoints%nik_axis(1:3)
377 assert(product(kpoints%nik_axis(1:3)) == kpoints%reduced%npoints - npath)
378
379 write(w90_win,'(a8,i4,i4,i4)') 'mp_grid =', axis(1:3)
380 write(w90_win,'(a)') ' '
381 write(w90_win,'(a)') 'begin kpoints '
382 ! Put a minus sign here for the wrong convention in Octopus
383
384 do ii = 1, kpoints%reduced%npoints-npath
385 write(w90_win,'(f13.8,f13.8,f13.8)') - kpoints%reduced%red_point(1:3,ii)
386 end do
387 write(w90_win,'(a)') 'end kpoints '
388 end if
389
390 call io_close(w90_win)
391
392 pop_sub(wannier90_setup)
393
394 end subroutine wannier90_setup
395
396 ! --------------------------------------------------------------------------
397 subroutine wannier90_output()
398 integer :: ik_real
399
400 push_sub(wannier90_output)
401
403
404 ! normal interface run
405 call states_elec_allocate_wfns(sys%st, sys%gr, wfs_type = type_cmplx, skip=exclude_list)
406 if (read_td_states) then
408 sys%mc, ierr, sys%gr)
409 else
411 sys%mc, ierr, sys%gr)
412 end if
413
414 if (ierr == 0) then
415 call states_elec_look(restart, nik, dim, nst, ierr)
416 if (sys%st%d%ispin == spin_polarized) then
417 nik = nik / 2
418 end if
419 if (dim == sys%st%d%dim .and. nik == sys%kpoints%reduced%npoints .and. nst == sys%st%nst) then
420 call states_elec_load(restart, global_namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
421 fixed_occ=.true., ierr=ierr, iter=iter, label = ": wannier90", skip=exclude_list)
422 else
423 write(message(1),'(a)') 'Restart structure not commensurate.'
424 call messages_fatal(1)
425 end if
426 end if
427 call restart_end(restart)
428
429 if (w90_scdm) then
430 nik = w90_num_kpts
431 safe_allocate(jpvt(1:sys%gr%np_global*sys%st%d%dim))
432 safe_allocate(psi(1:sys%gr%np, 1:sys%st%d%dim))
433 safe_allocate(occ_temp(1:w90_num_bands))
434
435 ! smear the states at gamma
436 do ist = 1, w90_num_bands
437 occ_temp(ist)= sys%st%occ(ist, 1)
438 sys%st%occ(ist, 1)=m_half*loct_erfc((sys%st%eigenval(ist, 1)-scdm_mu) / scdm_sigma)
439 end do
440
441 call zstates_elec_rrqr_decomposition(sys%st, sys%namespace, sys%gr, w90_num_bands, .true., 1, jpvt)
442
443 ! reset occupations at gamma
444 do ist = 1, w90_num_bands
445 sys%st%occ(ist, 1) = occ_temp(ist)
446 end do
447
448 safe_allocate(uk(1:w90_num_bands, 1:w90_num_bands, 1:nik))
449
450 ! auxiliary arrays for scdm procedure
451 safe_allocate(chi(1:w90_num_bands, 1:w90_num_bands))
452 safe_allocate(chi_diag(1:w90_num_bands, 1:w90_num_bands))
453 safe_allocate(chi2(1:w90_num_bands, 1:w90_num_bands))
454 safe_allocate(chi_eigenval(1:w90_num_bands))
455
456 chi(1:w90_num_bands, 1:w90_num_bands) = m_zero
457
458 do ik = 1, nik
459 kvec(:) = sys%kpoints%reduced%point(:, ik)
460
461 if (sys%st%d%ispin == spin_polarized) then
462 ik_real = (ik-1)*2 + w90_spin_channel
463 else
464 ik_real = ik
465 end if
466
467 do ist = 1, w90_num_bands
468 call states_elec_get_state(sys%st, sys%gr, ist, ik_real, psi)
469 smear=m_half * loct_erfc((sys%st%eigenval(ist, ik_real) - scdm_mu) / scdm_sigma)
470 ! NOTE: here check for domain parallelization
471 do jst = 1, w90_num_bands
472 chi(ist, jst) = smear * conjg(psi(jpvt(jst), 1)) &
473 * exp(m_zi * dot_product(sys%gr%x(jpvt(jst), 1:3), kvec(1:3)))
474 end do
475 end do
476
477 ! loewdin orhtogonalization of chi.chi
478 ! this can also be done with SVD, which might be more stable!?
479 chi_diag = matmul(conjg(transpose(chi)), chi)
480 call lalg_eigensolve(w90_num_bands, chi_diag, chi_eigenval)
481 chi2 = conjg(transpose(chi_diag))
482
483 !we need the eigenvalues to be >0
484 if (any(chi_eigenval(:) .lt. m_zero)) then
485 message(1) = 'SCDM Wannierization failed because chi matrix is'
486 message(2) = 'ill conditioned. Try increasingin scdm_sigma and/or'
487 message(3) = 'change scdm_mu.'
488 call messages_fatal(3)
489 end if
491 do ist = 1, w90_num_bands
492 chi_eigenval(ist) = m_one / sqrt(chi_eigenval(ist))
493 chi2(ist, 1:w90_num_bands) = chi_eigenval(ist) * chi2(ist, 1:w90_num_bands)
494 end do
495 ! the loewdin result would be: matmul(chi_diag,chi2)
496 ! to get the wannier gauge U(k) we multiply this with the original chi
497 uk(:,:,ik) = matmul(chi, matmul(chi_diag,chi2))
498
499 end do
500
501 safe_deallocate_a(chi)
502 safe_deallocate_a(psi)
503 safe_deallocate_a(chi_diag)
504 safe_deallocate_a(chi2)
505 safe_deallocate_a(chi_eigenval)
506 safe_deallocate_a(jpvt)
507 safe_deallocate_a(occ_temp)
508 end if
509
510 ! ---- actual interface work ----------
511 if (bitand(w90_what, option__wannier90files__w90_mmn) /= 0) then
512 call create_wannier90_mmn(sys%gr, sys%st)
513 end if
514
515 if (bitand(w90_what, option__wannier90files__w90_unk) /= 0) then
516 call write_unk(sys%space, sys%gr, sys%st)
517 end if
518
519 if (bitand(w90_what, option__wannier90files__w90_amn) /= 0) then
520 call create_wannier90_amn(sys%space, sys%gr, sys%ions%latt, sys%st, sys%kpoints)
521 end if
522
523 if (bitand(w90_what, option__wannier90files__w90_eig) /= 0) then
525 end if
526
527 safe_deallocate_a(uk)
528 safe_deallocate_a(w90_spin_proj_component)
529 safe_deallocate_a(w90_spin_proj_axis)
530
531 pop_sub(wannier90_output)
532 end subroutine wannier90_output
533
534 ! --------------------------------------------------------------------------
535 subroutine read_wannier90_files()
536 integer :: w90_nnkp, itemp, dummyint, io
537 character(len=80) :: filename, dummy, dummy1, dummy2, line
538 logical :: exist, parse_is_ok
539 real(real64) :: dummyr(7)
540
541 push_sub(read_wannier90_files)
542
543 w90_num_kpts = product(sys%kpoints%nik_axis(1:3))
544 w90_num_exclude = 0
545
546 ! open nnkp file
547 filename = trim(adjustl(w90_prefix)) //'.nnkp'
548
549 message(1) = "oct-wannier90: Parsing "//filename
550 call messages_info(1)
551
552 inquire(file=filename,exist=exist)
553 if (.not. exist) then
554 message(1) = 'oct-wannier90: Cannot find specified Wannier90 nnkp file.'
555 write(message(2),'(a)') 'Please run wannier90.x -pp '// trim(adjustl(w90_prefix)) // ' first.'
556 call messages_fatal(2)
557 end if
558
559 parse_is_ok = .false.
560
561 ! check number of k-points
562 w90_nnkp = io_open(trim(filename), global_namespace, action='read')
563 do
564 read(w90_nnkp, *, iostat=io) dummy, dummy1
565 if (io == iostat_end) exit
566
567 if (dummy == 'begin' .and. dummy1 == 'kpoints') then
568 read(w90_nnkp,*) itemp
569 if (itemp /= w90_num_kpts) then
570 message(1) = 'oct-wannier90: wannier90 setup seems to have been done with a different number of k-points.'
571 call messages_fatal(1)
572 else
573 parse_is_ok = .true.
574 exit
575 end if
576 end if
577 end do
578 call io_close(w90_nnkp)
579
580 if (.not. parse_is_ok) then
581 message(1) = 'oct-wannier90: Did not find the kpoints block in nnkp file'
582 call messages_fatal(1)
583 end if
584 parse_is_ok = .false.
585
586 ! read from nnkp file
587 ! find the nnkpts block
588 w90_nnkp = io_open(trim(filename), global_namespace, action='read', position='rewind')
589 do
590 read(w90_nnkp, *, iostat=io) dummy, dummy1
591 if (io == iostat_end) exit !End of file
592
593 if (dummy == 'begin' .and. dummy1 == 'nnkpts') then
594 read(w90_nnkp,*) w90_nntot
595 safe_allocate(w90_nnk_list(1:5, 1:w90_num_kpts * w90_nntot))
596 do ii = 1, w90_num_kpts * w90_nntot
597 read(w90_nnkp,*) w90_nnk_list(1:5, ii)
598 end do
599 !make sure we are at the end of the block
600 read(w90_nnkp,*) dummy
601 if (dummy /= 'end') then
602 message(1) = 'oct-wannier90: There dont seem to be enough k-points in nnkpts file to.'
603 call messages_fatal(1)
604 end if
605 parse_is_ok = .true.
606 exit
607 end if
608 end do
609
610 if (.not. parse_is_ok) then
611 message(1) = 'oct-wannier90: Did not find nnkpts block in nnkp file'
612 call messages_fatal(1)
613 end if
614
615 ! read from nnkp file
616 ! find the exclude block
617 safe_allocate(exclude_list(1:sys%st%nst))
618 !By default we use all the bands
619 exclude_list(1:sys%st%nst) = .false.
620 rewind(w90_nnkp)
621 do
622 read(w90_nnkp, *, iostat=io) dummy, dummy1
623 if (io == iostat_end) exit !End of file
624 if (dummy == 'begin' .and. dummy1 == 'exclude_bands') then
625 read(w90_nnkp, *) w90_num_exclude
626 do ii = 1, w90_num_exclude
627 read(w90_nnkp, *) itemp
628 if(itemp > sys%st%nst) then
629 message(1) = 'oct-wannier90: The exclude_bands list contains a state index higher than the number of states.'
630 call messages_fatal(1)
631 end if
632 exclude_list(itemp) = .true.
633 end do
634 !make sure we are at the end of the block
635 read(w90_nnkp, *) dummy
636 if (dummy /= 'end') then
637 message(1) = 'oct-wannier90: There dont seem to be enough bands in exclude_bands list.'
638 call messages_fatal(1)
639 end if
640 exit
641 end if
642 end do
643 call io_close(w90_nnkp)
644
645 !We get the number of bands
646 w90_num_bands = sys%st%nst - w90_num_exclude
647
648 safe_allocate(band_index(1:sys%st%nst))
649 itemp = 0
650 do ii = 1, sys%st%nst
651 if (exclude_list(ii)) cycle
652 itemp = itemp + 1
653 band_index(ii) = itemp
654 end do
655
656 if (bitand(w90_what, option__wannier90files__w90_amn) /= 0 &
657 .or. w90_mode == option__wannier90mode__w90_wannier ) then
658 ! parse file again for definitions of projections
659 w90_nnkp = io_open(trim(filename), global_namespace, action='read', position='rewind')
660
661 do
662 read(w90_nnkp, *, iostat=io) dummy, dummy1
663 if (io == iostat_end) then !End of file
664 message(1) = 'oct-wannier90: Did not find projections block in w90.nnkp file'
665 call messages_fatal(1)
666 end if
667
668 if (dummy == 'begin' .and. (dummy1 == 'projections' .or. dummy1 == 'spinor_projections')) then
669
670 if (dummy1 == 'spinor_projections') then
671 w90_spinors = .true.
672 if (sys%st%d%ispin /= spinors) then
673 message(1) = 'oct-wannier90: Spinor = .true. is only valid with spinors wavefunctions.'
674 call messages_fatal(1)
675 end if
676
677 message(1) = 'oct-wannier90: Spinor interface incomplete. Note there is no quantization axis implemented'
678 call messages_warning(1)
679 else
680 if (sys%st%d%ispin == spinors) then
681 message(1) = 'oct-wannier90: Octopus has spinors wavefunctions but spinor_projections is not defined.'
682 message(2) = 'oct-wannier90: Please check the input file for wannier 90.'
683 call messages_fatal(2)
684 end if
685 end if
686
687 read(w90_nnkp, *) w90_nproj
688 ! num_wann is given in w90.win, not double checked here
689 ! I assume that the wannier90.x -pp run has checked this
690 w90_num_wann = w90_nproj
691 ! In case of no projections, we use the number of bands
692 if(w90_nproj == 0) w90_num_wann = w90_num_bands
693
694 safe_allocate(w90_proj_centers(1:3, 1:w90_nproj))
695 safe_allocate(w90_proj_lmr(1:w90_nproj, 1:3))
696 if (w90_spinors) then
697 safe_allocate(w90_spin_proj_component(1:w90_nproj))
698 end if
699 if (w90_spinors) then
700 safe_allocate(w90_spin_proj_axis(1:w90_nproj, 1:3))
701 end if
702
703 do ii = 1, w90_nproj
704 read(w90_nnkp, *) w90_proj_centers(1:3, ii), w90_proj_lmr(ii, 1:3)
705 ! skip a line for now
706 read(w90_nnkp, *) dummyr(1:7)
707 if (w90_spinors) then
708 read(w90_nnkp, *) w90_spin_proj_component(ii), w90_spin_proj_axis(ii, 1:3)
709 ! use octopus spindim conventions
710 if (w90_spin_proj_component(ii) == -1) w90_spin_proj_component(ii) = 2
711 end if
712 end do
713 !make sure we are at the end of the block
714 read(w90_nnkp, *) dummy
715 if (dummy /= 'end') then
716 message(1) = 'oct-wannier90: There dont seem to be enough projections in nnkpts file to.'
717 call messages_fatal(1)
718 end if
719 exit
720 end if
721 end do
722
723 ! look for auto_projection block
724 scdm_proj = .false.
725 do
726 read(w90_nnkp, *, iostat=io) dummy, dummy1
727 if (io == iostat_end) exit !End of file
728
729 if (dummy == 'begin' .and. dummy1 == 'auto_projections') then
730 scdm_proj = .true.
731 read(w90_nnkp, *) w90_nproj
732 w90_num_wann = w90_nproj
733
734 if (.not. w90_scdm) then
735 message(1) = 'oct-wannier90: Found auto_projections block. Currently the only implemented automatic way'
736 message(2) = 'oct-wannier90: to compute projections is the SCDM method.'
737 message(3) = 'oct-wannier90: Please set Wannier90UseSCDM = yes in the inp file.'
738 call messages_fatal(3)
739 end if
740
741 if (w90_nproj /= w90_num_bands) then
742 message(1) = 'oct-wannier90: In auto_projections block first row needs to be equal to num_bands.'
743 call messages_fatal(1)
744 end if
745 read(w90_nnkp, *) dummyint
746 if (dummyint /= 0) then
747 message(1) = 'oct-wannier90: The second row in auto_projections has to be 0, per Wannier90 documentation.'
748 call messages_fatal(1)
749 end if
750 end if
751 end do
752 call io_close(w90_nnkp)
753
754 end if
755
756 message(1) = "oct-wannier90: Finished parsing "//filename
757 call messages_info(1)
758
759 ! Look extra variables variable
760 ! open win file
761 filename = trim(adjustl(w90_prefix)) //'.win'
762 message(1) = "oct-wannier90: Parsing "//filename
763 call messages_info(1)
764 w90_nnkp = io_open(trim(filename), global_namespace, action='read', position='rewind')
765 do
766 read(w90_nnkp, fmt='(a)', iostat=io) line
767 if (io == iostat_end) exit !End of file
768 if (index(line, '=') > 0) then
769 read(line, *, iostat=io) dummy, dummy2, dummy1
770 else
771 read(line, *, iostat=io) dummy, dummy1
772 end if
773
774 !Spin
775 if (dummy == 'spin') then
776 if (sys%st%d%ispin /= spin_polarized) then
777 message(1) = 'oct-wannier90: The variable spin is set for a non spin-polarized calculation.'
778 call messages_fatal(1)
779 end if
780
781 if (dummy1 == 'up') then
782 w90_spin_channel = 1
783 else if (dummy1 == 'down') then
784 w90_spin_channel = 2
785 else
786 message(1) = 'oct-wannier90: Error parsing the variable spin.'
787 call messages_fatal(1)
788 end if
789 end if
790 end do
791 call io_close(w90_nnkp)
792
793 if (sys%st%d%ispin == spin_polarized) then
794 write(message(1), '(a,i1)') 'oct-wannier90: Using spin channel ', w90_spin_channel
795 call messages_info(1)
796 end if
797
798 message(1) = "oct-wannier90: Finished parsing "//filename
799 call messages_info(1)
800
801 pop_sub(read_wannier90_files)
802
803 end subroutine read_wannier90_files
804
805 ! --------------------------------------------------------------------------
806 subroutine create_wannier90_mmn(mesh, st)
807 class(mesh_t), intent(in) :: mesh
808 type(states_elec_t), target, intent(in) :: st
809
810 integer :: ist, jst, ik, ip, w90_mmn, iknn, idim, ibind
811 real(real64) :: Gcart(3)
812 integer :: G(3)
813 character(len=80) :: filename
814 complex(real64), allocatable :: overlap(:)
815 complex(real64), allocatable :: psim(:,:), psin(:,:), phase(:)
816 type(wfs_elec_t), pointer :: batch
817 integer :: inode, node_fr, node_to
818 type(mpi_request) :: send_req
819
820 push_sub(create_wannier90_mmn)
821
822 call profiling_in("W90_MMN")
823
824 if (st%parallel_in_states) then
825 call messages_not_implemented("w90_mmn output with states parallelization")
826 end if
827
828 message(1) = "Info: Computing the overlap matrix"
829 call messages_info(1)
830
831
832 filename = './'// trim(adjustl(w90_prefix))//'.mmn'
833 w90_mmn = io_open(trim(filename), global_namespace, action='write')
834
835 ! write header
836 if (mpi_grp_is_root(mpi_world)) then
837 write(w90_mmn,*) 'Created by oct-wannier90'
838 write(w90_mmn,*) w90_num_bands, w90_num_kpts, w90_nntot
839 end if
840
841 safe_allocate(psim(1:mesh%np, 1:st%d%dim))
842 safe_allocate(psin(1:mesh%np, 1:st%d%dim))
843 safe_allocate(phase(1:mesh%np))
844 safe_allocate(overlap(1:w90_num_bands))
845
846
847 ! loop over the pairs specified in the nnkp file (read before in init)
848 do ii = 1, w90_num_kpts * w90_nntot
849 ik = w90_nnk_list(1, ii)
850 iknn = w90_nnk_list(2, ii)
851 g(1:3) = w90_nnk_list(3:5, ii)
852 if (mpi_grp_is_root(mpi_world)) write(w90_mmn, '(I10,2x,I10,2x,I3,2x,I3,2x,I3)') ik, iknn, g
853
854 ! For spin-polarized calculations, we select the right k-point
855 if (sys%st%d%ispin == spin_polarized) then
856 ik = (ik-1)*2 + w90_spin_channel
857 iknn = (iknn-1)*2 + w90_spin_channel
858 end if
859
860
861 ! Only treat local k-points
862 if(ik >= st%d%kpt%start .and. ik <= st%d%kpt%end) then
863
864 ! Wannier90 treats everything fully periodic
865 ! Conversion is done with the 3D "correct" klattice
866 call kpoints_to_absolute(sys%ions%latt, real(g, real64) , gcart)
867
868 ! Phase that gives u_{n,k+G}(r) from u_{nk}(r)
869 ! Here the minus sign of Octopus cancels with minus sign of the input G
870 ! (ik and iknn correspond in the list to minus the k-points in Octopus)
871 if (any(g /= 0)) then
872 do ip = 1, mesh%np
873 phase(ip) = exp(-m_zi*dot_product(mesh%x(ip,1:3), gcart(1:3)))
874 end do
875 end if
876
877 end if
878
879 ! Loop over distributed bands
880 do jst = 1, st%nst
881 if (exclude_list(jst)) cycle
882
883 ! Communication for the local states
884 if ( .not. st%d%kpt%parallel .and. .not. st%parallel_in_states) then
885 call states_elec_get_state(st, mesh, jst, iknn, psin)
886 else
887 node_fr = -1
888 node_to = -1
889 do inode = 0, st%d%kpt%mpi_grp%size-1
890 if(iknn >= st%st_kpt_task(inode,3) .and. iknn <= st%st_kpt_task(inode,4)) then
891 node_fr = inode
892 end if
893 if(ik >= st%st_kpt_task(inode,3) .and. ik <= st%st_kpt_task(inode,4)) then
894 node_to = inode
895 end if
896 end do
897 assert(node_fr > -1)
898 assert(node_to > -1)
900 send_req = mpi_request_null
901 ! We have locally the k-point
902 if (state_kpt_is_local(st, jst, iknn)) then
903 call states_elec_get_state(st, mesh, jst, iknn, psin)
904 ! We send it only if we don`t want to use it locally
905 if(node_to /= st%d%kpt%mpi_grp%rank) then
906 call st%d%kpt%mpi_grp%isend(psin, mesh%np*st%d%dim, mpi_double_complex, node_to, send_req)
907 end if
908 end if
909 ! We receive the desired state, only if it is not a local one
910 if(node_to == st%d%kpt%mpi_grp%rank .and. node_to /= node_fr) then
911 call st%d%kpt%mpi_grp%recv(psin, mesh%np*st%d%dim, mpi_double_complex, node_fr)
912 end if
913 if (send_req /= mpi_request_null) then
914 call st%d%kpt%mpi_grp%wait(send_req)
915 end if
916 end if
917
918 overlap = m_zero
919
920 if(ik >= st%d%kpt%start .and. ik <= st%d%kpt%end) then
921
922 ! Do not apply the phase if the phase factor is null
923 if (any(g /= 0)) then
924 ! add phase
925 do idim = 1, st%d%dim
926 do ip = 1, mesh%np
927 psin(ip, idim) = psin(ip, idim) * phase(ip)
928 end do
929 end do
930 end if
931
932
933 ! See Eq. (25) in PRB 56, 12847 (1997)
934 ! Loop over local k-points
935 do ist = 1, st%nst
936 if (exclude_list(ist)) cycle
937
938 batch => st%group%psib(st%group%iblock(ist), ik)
939
940 select case (batch%status())
941 case (batch_not_packed)
942 overlap(band_index(ist)) = m_z0
943 do idim = 1, st%d%dim
944 ibind = batch%inv_index((/ist, idim/))
945 overlap(band_index(ist)) = overlap(band_index(ist)) + &
946 zmf_dotp(mesh, batch%zff_linear(:, ibind), psin(:,idim), reduce = .false.)
947 end do
948 !Not properly done at the moment
950 call states_elec_get_state(st, mesh, ist, ik, psim)
951 overlap(band_index(ist)) = zmf_dotp(mesh, st%d%dim, psim, psin, reduce = .false.)
952 end select
953 end do !ist
954 end if
955
956 call profiling_in("W90_MMN_REDUCE")
957 call mesh%allreduce(overlap)
958 call profiling_out("W90_MMN_REDUCE")
959
960 if(st%d%kpt%parallel) then
961 call comm_allreduce(st%d%kpt%mpi_grp, overlap)
962 end if
963
964 ! write to W90 file
965 if (mpi_grp_is_root(mpi_world)) then
966 do ist = 1, st%nst
967 if (exclude_list(ist)) cycle
968 write(w90_mmn,'(e18.10,2x,e18.10)') overlap(band_index(ist))
969 end do
970 end if
971
972 end do !jst
973 end do
974
975 call io_close(w90_mmn)
976
977 safe_deallocate_a(psim)
978 safe_deallocate_a(psin)
979 safe_deallocate_a(phase)
980 safe_deallocate_a(overlap)
981
982 call profiling_out("W90_MMN")
983
984 pop_sub(create_wannier90_mmn)
985
986 end subroutine create_wannier90_mmn
987
988 ! --------------------------------------------------------------------------
989 subroutine create_wannier90_eig()
990 integer :: ist, ik, w90_eig
991 character(len=80) :: filename
992
993 push_sub(create_wannier90_eig)
994
995 if (sys%st%parallel_in_states) then
996 call messages_not_implemented("w90_eig output with states parallelization")
997 end if
998
999 if (mpi_grp_is_root(mpi_world)) then
1000 filename = './'//trim(adjustl(w90_prefix))//'.eig'
1001 w90_eig = io_open(trim(filename), global_namespace, action='write')
1002 do ik = 1, w90_num_kpts
1003 do ist = 1, sys%st%nst
1004 if (exclude_list(ist)) cycle
1005 if (sys%st%d%ispin /= spin_polarized) then
1006 write(w90_eig,'(I5,2x,I8,2x,e18.10)') band_index(ist), ik, &
1007 units_from_atomic(unit_ev, sys%st%eigenval(ist, ik))
1008 else
1009 write(w90_eig,'(I5,2x,I8,2x,e18.10)') band_index(ist), ik, &
1010 units_from_atomic(unit_ev, sys%st%eigenval(ist, (ik-1)*2+w90_spin_channel))
1011 end if
1012 end do
1013 end do
1014
1015 call io_close(w90_eig)
1016 end if
1017
1018 pop_sub(create_wannier90_eig)
1019 end subroutine create_wannier90_eig
1020
1021 ! --------------------------------------------------------------------------
1022 subroutine write_unk(space, mesh, st)
1023 class(space_t), intent(in) :: space
1024 class(mesh_t), intent(in) :: mesh
1025 type(states_elec_t), intent(in) :: st
1026
1027 integer :: ist, ik, unk_file, ispin
1028 integer :: ix, iy, iz
1029 character(len=80) :: filename
1030 complex(real64), allocatable :: psi(:)
1031 type(cube_t) :: cube
1032 type(cube_function_t) :: cf
1033
1034 push_sub(write_unk)
1035
1036 if (st%d%kpt%parallel) then
1037 call messages_not_implemented("w90_unk output with k-point parallelization")
1038 end if
1039
1040 if (sys%gr%parallel_in_domains) then
1041 call messages_not_implemented("w90_unk output with domain parallelization")
1042 end if
1043
1044 if (st%parallel_in_states) then
1045 call messages_not_implemented("w90_unk output with states parallelization")
1046 end if
1047
1048 call messages_experimental("Wannier90Files = w90_unk")
1049
1050
1051 safe_allocate(psi(1:mesh%np))
1052
1053 call cube_init(cube, mesh%idx%ll, global_namespace, space, mesh%spacing, &
1054 mesh%coord_system, need_partition=.not.mesh%parallel_in_domains)
1055 call cube_init_cube_map(cube, mesh)
1056
1057 call zcube_function_alloc_rs(cube, cf)
1058
1059 do ik = 1, w90_num_kpts
1060 do ispin = 1, st%d%dim
1061 if (mpi_grp_is_root(mpi_world)) then
1062 write(filename, '(a,i5.5,a1,i1)') './UNK', ik,'.', ispin
1063 unk_file = io_open(trim(filename), global_namespace, action='write', form='unformatted')
1064 ! write header
1065 write(unk_file) mesh%idx%ll(1:mesh%idx%dim), ik, w90_num_bands
1066 end if
1067
1068 ! states
1069 do ist = 1, st%nst
1070 if (exclude_list(ist)) cycle
1071
1072 if (sys%st%d%ispin /= spin_polarized) then
1073 call states_elec_get_state(st, mesh, ispin, ist, ik, psi)
1074 else
1075 call states_elec_get_state(st, mesh, ispin, ist, (ik-1)*2+w90_spin_channel, psi)
1076 end if
1077
1078 ! put the density in the cube
1079 ! Note: At the moment this does not work for domain parallelization
1080 if (cube%parallel_in_domains) then
1081 assert(.not. cube%parallel_in_domains)
1082 else
1083 call zmesh_to_cube(mesh, psi, cube, cf)
1084 end if
1085
1086 if (mpi_grp_is_root(mpi_world)) then
1087 write(unk_file) (((cf%zrs(ix,iy,iz), ix=1,cube%rs_n_global(1)), iy=1,cube%rs_n_global(2)), iz=1,cube%rs_n_global(3))
1088 end if
1089 end do
1090 if (mpi_grp_is_root(mpi_world)) call io_close(unk_file)
1091 end do
1092 end do
1093
1094 call dcube_function_free_rs(cube, cf)
1095 call cube_end(cube)
1096
1097 safe_deallocate_a(psi)
1098
1099 pop_sub(write_unk)
1100
1101 end subroutine write_unk
1102
1103 ! --------------------------------------------------------------------------
1104 subroutine create_wannier90_amn(space, mesh, latt, st, kpoints)
1105 class(space_t), intent(in) :: space
1106 class(mesh_t), intent(in) :: mesh
1107 type(lattice_vectors_t), intent(in) :: latt
1108 type(states_elec_t), intent(in) :: st
1109 type(kpoints_t), intent(in) :: kpoints
1110
1111 integer :: ist, ik, w90_amn, idim, iw, ip, ik_real
1112 real(real64) :: center(3), kpoint(3), threshold
1113 character(len=80) :: filename
1114 complex(real64), allocatable :: psi(:,:), phase(:), projection(:)
1115 real(real64), allocatable :: ylm(:)
1116 type(orbitalset_t), allocatable :: orbitals(:)
1117
1118 push_sub(create_wannier90_amn)
1119 call profiling_in("W90_AMN")
1120
1121 if (st%parallel_in_states) then
1122 call messages_not_implemented("w90_amn output with states parallelization")
1123 end if
1124
1125 filename = './'// trim(adjustl(w90_prefix))//'.amn'
1126 w90_amn = io_open(trim(filename), global_namespace, action='write')
1127
1128 ! write header
1129 if (mpi_grp_is_root(mpi_world)) then
1130 write(w90_amn,*) 'Created by oct-wannier90'
1131 write(w90_amn,*) w90_num_bands, w90_num_kpts, w90_num_wann
1132 end if
1133
1134 if (scdm_proj) then
1135
1136 message(1) = "Info: Writing projections obtained from SCDM."
1137 call messages_info(1)
1138
1139 if (st%d%kpt%parallel) then
1140 call messages_not_implemented("w90_amn output with k-point parallelization")
1141 end if
1142
1143
1144 do ik = 1, w90_num_kpts
1145 do ist = 1, st%nst
1146 if (exclude_list(ist)) cycle
1147 if (mpi_grp_is_root(mpi_world)) then
1148 do iw = 1, w90_nproj
1149 write (w90_amn,'(I5,2x,I5,2x,I5,2x,e18.10,2x,e18.10)') band_index(ist), iw, ik, uk(band_index(ist),iw,ik)
1150 end do
1151 end if
1152 end do !ist
1153 end do! ik
1154
1155 else
1156
1157 message(1) = "Info: Computing the projection matrix"
1158 call messages_info(1)
1159
1160 !We use the variable AOThreshold to deterine the threshold on the radii of the atomic orbitals
1161 call parse_variable(global_namespace, 'AOThreshold', 0.001_real64, threshold)
1162
1163 safe_allocate(orbitals(1:w90_nproj))
1164 ! precompute orbitals
1165 do iw=1, w90_nproj
1166 call orbitalset_init(orbitals(iw))
1167
1168 orbitals(iw)%norbs = 1
1169 orbitals(iw)%ndim = 1
1170 orbitals(iw)%radius = -log(threshold)
1171 orbitals(iw)%use_submesh = .false.
1172
1173 ! cartesian coordinate of orbital center
1174 center(1:3) = latt%red_to_cart(w90_proj_centers(1:3, iw))
1175 call submesh_init(orbitals(iw)%sphere, space, mesh, latt, center, orbitals(iw)%radius)
1176
1177 ! get dorb as submesh points
1178 safe_allocate(ylm(1:orbitals(iw)%sphere%np))
1179 ! (this is a routine from pwscf)
1180 call ylm_wannier(ylm, w90_proj_lmr(iw,1), w90_proj_lmr(iw,2), &
1181 orbitals(iw)%sphere%r, orbitals(iw)%sphere%rel_x, orbitals(iw)%sphere%np)
1182
1183 ! apply radial function
1184 if (w90_proj_lmr(iw,3) == 1) then
1185 do ip = 1,orbitals(iw)%sphere%np
1186 ylm(ip) = ylm(ip)*m_two*exp(-orbitals(iw)%sphere%r(ip))
1187 end do
1188 else
1189 call messages_not_implemented("oct-wannier90: r/=1 for the radial part")
1190 end if
1191
1192 safe_allocate(orbitals(iw)%zorb(1:orbitals(iw)%sphere%np, 1, 1))
1193 orbitals(iw)%zorb(1:orbitals(iw)%sphere%np, 1, 1) = ylm(1:orbitals(iw)%sphere%np)
1194 safe_deallocate_a(ylm)
1195
1196 safe_allocate(orbitals(iw)%phase(1:orbitals(iw)%sphere%np, st%d%kpt%start:st%d%kpt%end))
1197 orbitals(iw)%phase(:,:) = m_z0
1198 safe_allocate(orbitals(iw)%eorb_mesh(1:mesh%np, 1, 1, st%d%kpt%start:st%d%kpt%end))
1199 orbitals(iw)%eorb_mesh(:,:,:,:) = m_z0
1200
1201 call orbitalset_update_phase(orbitals(iw), space%dim, st%d%kpt, kpoints, st%d%ispin == spin_polarized, &
1202 kpt_max = w90_num_kpts)
1203
1204 end do
1205
1206 safe_allocate(psi(1:mesh%np, 1:st%d%dim))
1207 safe_allocate(phase(1:mesh%np))
1208 safe_allocate(projection(1:w90_nproj))
1209
1210 do ik = 1, w90_num_kpts
1211 kpoint(1:space%dim) = kpoints%get_point(ik)
1212 do ip = 1, mesh%np
1213 phase(ip) = exp(-m_zi* sum(mesh%x(ip, 1:space%dim) * kpoint(1:space%dim)))
1214 end do
1215
1216 !For spin-polarized calculations, we select the right k-point
1217 if (st%d%ispin == spin_polarized) then
1218 ik_real = (ik-1)*2 + w90_spin_channel
1219 else
1220 ik_real = ik
1221 end if
1222
1223
1224 do ist = 1, st%nst
1225 if (exclude_list(ist)) cycle
1226
1227 projection = m_zero
1228
1229 if(ik_real >= st%d%kpt%start .and. ik_real <= st%d%kpt%end) then
1230 call states_elec_get_state(st, mesh, ist, ik_real, psi)
1231
1232 do idim = 1, st%d%dim
1233 !The minus sign is here is for the wrong convention of Octopus
1234 do ip = 1, mesh%np
1235 psi(ip, idim) = psi(ip, idim)*phase(ip)
1236 end do
1237 end do
1238
1239 do iw = 1, w90_nproj
1240 idim = 1
1241 if (w90_spinors) idim = w90_spin_proj_component(iw)
1242
1243 !At the moment the orbitals do not depend on idim
1244 !The idim index for eorb_mesh would be for a spin-resolved orbital like j=1/2
1245 projection(iw) = zmf_dotp(mesh, psi(1:mesh%np,idim), &
1246 orbitals(iw)%eorb_mesh(1:mesh%np,1,1,ik_real), reduce = .false.)
1247 end do
1248
1249 call profiling_in("W90_AMN_REDUCE")
1250 call mesh%allreduce(projection)
1251 call profiling_out("W90_AMN_REDUCE")
1252 end if
1253
1254 if(st%d%kpt%parallel) then
1255 call comm_allreduce(st%d%kpt%mpi_grp, projection)
1256 end if
1257
1258 if (mpi_grp_is_root(mpi_world)) then
1259 do iw = 1, w90_nproj
1260 write (w90_amn,'(I5,2x,I5,2x,I5,2x,e18.10,2x,e18.10)') band_index(ist), iw, ik, projection(iw)
1261 end do
1262 end if
1263 end do !ist
1264 end do !ik
1265
1266 safe_deallocate_a(psi)
1267 safe_deallocate_a(phase)
1268 safe_deallocate_a(projection)
1269
1270 do iw = 1, w90_nproj
1271 call orbitalset_end(orbitals(iw))
1272 end do
1273 safe_deallocate_a(orbitals)
1274 end if
1275
1276 call io_close(w90_amn)
1277
1278 call profiling_out("W90_AMN")
1279
1280 pop_sub(create_wannier90_amn)
1281
1282 end subroutine create_wannier90_amn
1283
1284 ! --------------------------------------------------------------------------
1285 subroutine generate_wannier_states(space, mesh, ions, st, kpoints)
1286 class(space_t), intent(in) :: space
1287 class(mesh_t), intent(in) :: mesh
1288 type(ions_t), intent(in) :: ions
1289 type(states_elec_t), intent(in) :: st
1290 type(kpoints_t), intent(in) :: kpoints
1291
1292 integer :: w90_u_mat, w90_xyz, nwann, nik
1293 integer :: ik, iw, iw2, ip, ipmax, rankmax
1294 real(real64), allocatable :: centers(:,:)
1295 complex(real64), allocatable :: Umnk(:,:,:)
1296 complex(real64), allocatable :: zwn(:), psi(:,:), phase(:)
1297 character(len=MAX_PATH_LEN) :: fname
1298 real(real64) :: kpoint(3), wmod, wmodmax, xx(space%dim)
1299 character(len=2) :: dum
1300 logical :: exist
1301 type(unit_t) :: fn_unit
1302 complex(real64) :: scal
1303
1304 push_sub(generate_wannier_states)
1305
1306 assert(st%d%ispin /= spinors)
1307
1308 message(1) = "oct-wannier90: Constructing the Wannier states from the U matrix."
1309 call messages_info(1)
1310
1311 inquire(file=trim(trim(adjustl(w90_prefix))//'_centres.xyz'),exist=exist)
1312 if (.not. exist) then
1313 message(1) = 'oct-wannier90: Cannot find the Wannier90 file seedname_centres.xyz.'
1314 write(message(2),'(a)') 'Please run wannier90.x with "write_xyz=.true." in '// trim(adjustl(w90_prefix)) // '.'
1315 call messages_fatal(2)
1316 end if
1317
1318 w90_xyz = io_open(trim(trim(adjustl(w90_prefix))//'_centres.xyz'), global_namespace, action='read')
1319
1320 safe_allocate(centers(1:3, 1:w90_num_wann))
1321 !Skip two lines
1322 read(w90_xyz, *)
1323 read(w90_xyz, *)
1324 do iw = 1, w90_num_wann
1325 read(w90_xyz, *) dum, centers(1:3, iw)
1326 ! Wannier90 outputs the coordinates in angstrom
1327 centers(1:3, iw) = units_to_atomic(unit_angstrom, centers(1:3, iw))
1328 end do
1329 call io_close(w90_xyz)
1330
1331
1332 inquire(file=trim(trim(adjustl(w90_prefix))//'_u_dis.mat'),exist=exist)
1333 if (exist) then
1334 message(1) = 'oct-wannier90: Calculation of Wannier states with disentanglement is not yet supported.'
1335 call messages_fatal(1)
1336 end if
1337
1338 inquire(file=trim(trim(adjustl(w90_prefix))//'_u.mat'),exist=exist)
1339 if (.not. exist) then
1340 message(1) = 'oct-wannier90: Cannot find the Wannier90 seedname_u.mat file.'
1341 write(message(2),'(a)') 'Please run wannier90.x with "write_u_matrices=.true." in '// trim(adjustl(w90_prefix)) // '.'
1342 call messages_fatal(2)
1343 end if
1344 w90_u_mat = io_open(trim(trim(adjustl(w90_prefix))//'_u.mat'), global_namespace, action='read')
1345
1346
1347 !Skip one line
1348 read(w90_u_mat, *)
1349 !Read num_kpts, num_wann, num_wann for consistency check
1350 read(w90_u_mat, *) nik, nwann, nwann
1351 if (nik /= w90_num_kpts .or. nwann /= w90_num_wann) then
1352 print *, w90_num_wann, w90_num_kpts, nik, nwann
1353 message(1) = "The file contains U matrices is inconsistent with the .win file."
1354 call messages_fatal(1)
1355 end if
1356
1357 safe_allocate(umnk(1:w90_num_wann, 1:w90_num_wann, 1:w90_num_kpts))
1358
1359 do ik = 1, w90_num_kpts
1360 !Skip one line
1361 read(w90_u_mat, *)
1362 !Skip one line (k-point coordinate)
1363 read(w90_u_mat, *)
1364 read(w90_u_mat, '(f15.10,sp,f15.10)') ((umnk(iw, iw2, ik), iw=1, w90_num_wann), iw2=1, w90_num_wann)
1365 end do
1366
1367 call io_close(w90_u_mat)
1368
1369 !We read the output format for the Wannier states
1370 call parse_variable(global_namespace, 'OutputFormat', 0, how)
1371
1372 call io_mkdir('wannier', global_namespace)
1373
1374 !Computing the Wannier states in the primitive cell, from the U matrices
1375 safe_allocate(zwn(1:mesh%np))
1376 safe_allocate(psi(1:mesh%np, 1:st%d%dim))
1377 safe_allocate(phase(1:mesh%np))
1379 do iw = 1, w90_num_wann
1380
1381 zwn(:) = m_z0
1382
1383 do ik = 1, w90_num_kpts
1384 kpoint(1:space%dim) = kpoints%get_point(ik, absolute_coordinates=.true.)
1385
1386 ! We compute the Wannier orbital on a grid centered around the Wannier function
1387 do ip = 1, mesh%np
1388 xx = mesh%x(ip, 1:space%dim)-centers(1:space%dim, iw)
1389 xx = ions%latt%fold_into_cell(xx)
1390 phase(ip) = exp(-m_zi* sum( xx * kpoint(1:space%dim)))
1391 end do
1392
1393 do iw2 = 1, st%nst
1394 if (exclude_list(iw2)) cycle
1395
1396 if (st%d%ispin /= spin_polarized) then
1397 call states_elec_get_state(st, mesh, iw2, ik, psi)
1398 else
1399 call states_elec_get_state(st, mesh, iw2, (ik-1)*2+w90_spin_channel, psi)
1400 end if
1401
1402 ! The minus sign is here is for the wrong convention of Octopus
1403 do ip = 1, mesh%np
1404 zwn(ip) = zwn(ip) + umnk(band_index(iw2), iw, ik) * psi(ip, 1) * phase(ip)
1405 end do
1406 end do!ik
1407 end do!iw2
1408
1409 ! Following what Wannier90 is doing, we fix the global phase by setting the max to be real
1410 ! We also normalize to the number of k-point at this step
1411 ipmax = 0
1412 wmodmax = m_zero
1413 do ip = 1, mesh%np
1414 wmod = real(zwn(ip)*conjg(zwn(ip)), real64)
1415 if (wmod > wmodmax) then
1416 ipmax = ip
1417 wmodmax = wmod
1418 end if
1419 end do
1420 scal = sqrt(wmodmax)/zwn(ipmax)/w90_num_kpts
1421 call mesh_minmaxloc(mesh, wmodmax, rankmax, mpi_maxloc)
1422 call mesh%mpi_grp%bcast(scal, 1, mpi_double_complex, rankmax)
1423 call lalg_scal(mesh%np, scal, zwn)
1424
1425 ! Output the Wannier states
1426 fn_unit = sqrt(units_out%length**(-space%dim))
1427 write(fname, '(a,i0)') 'wannier-', iw
1428 call zio_function_output(how, "wannier", trim(fname), global_namespace, space, mesh, &
1429 zwn, fn_unit, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
1430
1431 ! Checking the ratio imag/real
1432 wmodmax = m_zero
1433 do ip = 1, mesh%np
1434 if(abs(real(zwn(ip), real64)) >= 1e-2_real64) then
1435 wmodmax = max(wmodmax, abs(aimag(zwn(ip)))/abs(real(zwn(ip), real64)))
1436 end if
1437 end do
1438 call mesh_minmaxloc(mesh, wmodmax, rankmax, mpi_maxloc)
1439
1440 write(message(1), '(a,i0,a,f11.6)') 'oct-wannier90: Wannier function ', iw, ' Max. Im/Re Ratio = ', wmodmax
1441 call messages_info(1)
1442 end do
1443
1444 safe_deallocate_a(umnk)
1445 safe_deallocate_a(zwn)
1446 safe_deallocate_a(psi)
1447 safe_deallocate_a(phase)
1448 safe_deallocate_a(centers)
1449
1451 end subroutine generate_wannier_states
1452
1453end program wannier90_interface
1454
1455!! Local Variables:
1456!! mode: f90
1457!! coding: utf-8
1458!! End:
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
This module implements batches of mesh functions.
Definition: batch.F90:133
integer, parameter, public batch_not_packed
functions are stored in CPU memory, unpacked order
Definition: batch.F90:276
integer, parameter, public batch_device_packed
functions are stored in device memory in packed order
Definition: batch.F90:276
integer, parameter, public batch_packed
functions are stored in CPU memory, in transposed (packed) order
Definition: batch.F90:276
This module handles the calculation mode.
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
integer, parameter, public p_strategy_states
parallelization in states
subroutine, public zmesh_to_cube(mesh, mf, cube, cf)
Convert a function from the mesh to the cube.
subroutine, public dcube_function_free_rs(cube, cf)
Deallocates the real space grid.
subroutine, public zcube_function_alloc_rs(cube, cf, in_device, force_alloc)
Allocates locally the real space grid, if PFFT library is not used. Otherwise, it assigns the PFFT re...
subroutine, public cube_init(cube, nn, namespace, space, spacing, coord_system, fft_type, fft_library, dont_optimize, nn_out, mpi_grp, need_partition, tp_enlarge, blocksize)
Definition: cube.F90:202
subroutine, public cube_end(cube)
Definition: cube.F90:378
subroutine, public cube_init_cube_map(cube, mesh)
Definition: cube.F90:823
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:118
subroutine, public fft_all_init(namespace)
initialize the table
Definition: fft.F90:278
subroutine, public fft_all_end()
delete all plans
Definition: fft.F90:391
real(real64), parameter, public m_two
Definition: global.F90:190
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:382
real(real64), parameter, public m_huge
Definition: global.F90:206
real(real64), parameter, public m_zero
Definition: global.F90:188
complex(real64), parameter, public m_z0
Definition: global.F90:198
complex(real64), parameter, public m_zi
Definition: global.F90:202
subroutine, public global_init(communicator)
Initialise Octopus.
Definition: global.F90:325
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
This module implements the underlying real-space grid.
Definition: grid.F90:117
Definition: io.F90:114
subroutine, public io_init(defaults)
If the argument defaults is present and set to true, then the routine will not try to read anything f...
Definition: io.F90:161
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_end()
Definition: io.F90:258
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
subroutine, public kpoints_to_absolute(latt, kin, kout)
Definition: kpoints.F90:1031
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_end()
Definition: messages.F90:278
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
subroutine, public messages_init(output_dir)
Definition: messages.F90:225
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1086
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
type(namespace_t), public global_namespace
Definition: namespace.F90:132
subroutine, public orbitalset_init(this)
Definition: orbitalset.F90:206
subroutine, public orbitalset_end(this)
Definition: orbitalset.F90:232
subroutine, public orbitalset_update_phase(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
Build the phase correction to the global phase in case the orbital crosses the border of the simulato...
Definition: orbitalset.F90:282
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
subroutine, public parser_init()
Initialise the Octopus parser.
Definition: parser.F90:450
subroutine, public parser_end()
End the Octopus parser.
Definition: parser.F90:481
subroutine, public profiling_end(namespace)
Definition: profiling.F90:413
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
Definition: profiling.F90:255
subroutine, public restart_module_init(namespace)
Definition: restart.F90:308
integer, parameter, public restart_gs
Definition: restart.F90:200
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
Definition: restart.F90:516
integer, parameter, public restart_td
Definition: restart.F90:200
integer, parameter, public restart_type_load
Definition: restart.F90:245
subroutine, public restart_end(restart)
Definition: restart.F90:722
subroutine, public zstates_elec_rrqr_decomposition(st, namespace, mesh, nst, root, ik, jpvt)
Perform RRQR on the transpose states stored in the states object and return the pivot vector.
This module defines routines to write information about states.
logical function, public state_kpt_is_local(st, ist, ik)
check whether a given state (ist, ik) is on the local node
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public states_elec_look(restart, nik, dim, nst, ierr)
Reads the 'states' file in the restart directory, and finds out the nik, dim, and nst contained in it...
This module handles reading and writing restart information for the states_elec_t.
subroutine, public states_elec_load(restart, namespace, space, st, mesh, kpoints, fixed_occ, ierr, iter, lr, lowest_missing, label, verbose, skip)
returns in ierr: <0 => Fatal error, or nothing read =0 => read all wavefunctions >0 => could only rea...
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:283
type(type_t), parameter, public type_cmplx
Definition: types.F90:134
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_t), public unit_angstrom
For XYZ files.
subroutine, public unit_system_init(namespace)
type(unit_t), public unit_ev
For output energies in eV.
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118
subroutine, public ylm_wannier(ylm, l, mr, rr, xx, nr)
Class describing the electron system.
Definition: electrons.F90:218
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Definition: wfs_elec.F90:138
int true(void)
subroutine create_wannier90_eig()
subroutine read_wannier90_files()
subroutine create_wannier90_mmn(mesh, st)
subroutine wannier90_setup(ions, kpoints, space)
subroutine write_unk(space, mesh, st)
subroutine generate_wannier_states(space, mesh, ions, st, kpoints)
program wannier90_interface
subroutine create_wannier90_amn(space, mesh, latt, st, kpoints)
subroutine wannier90_output()