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
72 integer(int64) :: w90_what_default
73
74 integer :: ierr
75 integer :: dim, idim
76 integer :: ii, nik, iter, nst
77
78 type(restart_t) :: restart
79 type(electrons_t), pointer :: sys
80 logical :: w90_spinors, scdm_proj, w90_scdm
81 integer :: w90_nntot, w90_num_bands, w90_num_kpts ! w90 input parameters
82 integer, allocatable :: w90_nnk_list(:,:) !
83 character(len=80) :: w90_prefix ! w90 input file prefix
84 integer :: w90_num_wann ! input paramter
85 real(real64), allocatable :: w90_proj_centers(:,:) ! projections centers
86 integer, allocatable :: w90_proj_lmr(:,:) ! definitions for real valued Y_lm*R_r
87 integer :: w90_nproj ! number of such projections
88 integer, allocatable :: w90_spin_proj_component(:) ! up/down flag
89 real(real64), allocatable :: w90_spin_proj_axis(:,:) ! spin axis (not implemented)
90 integer :: w90_num_exclude
91 logical, allocatable :: exclude_list(:) ! list of excluded bands
92 integer, allocatable :: band_index(:) ! band index after exclusion
93 logical :: read_td_states
94 integer :: w90_spin_channel
95 logical :: w90_bloch_sums = .false.
96
97 ! scdm variables
98 integer, allocatable :: jpvt(:)
99 complex(real64), allocatable :: uk(:,:,:) ! SCDM-Wannier gauge matrices U(k)
100 complex(real64), allocatable :: psi(:,:)
101 complex(real64), allocatable :: chi(:,:), chi_diag(:,:),chi2(:,:)
102 real(real64), allocatable :: chi_eigenval(:), occ_temp(:)
103 real(real64) :: scdm_mu, scdm_sigma, smear, kvec(3), factor(3)
104 integer :: ist, jst, ik, idir, sender
105
106 integer(int64) :: how
107
108 call global_init()
109 call parser_init()
110
111 call messages_init()
112 call io_init()
113
115
118
119 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
120 sys => electrons_t(global_namespace, int(option__calculationmode__dummy, int32))
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 !%Option w90_spn bit(5)
183 !% Spin. See Wannier90 documentation for more details.
184 !%End
185 w90_what_default = option__wannier90files__w90_mmn + option__wannier90files__w90_amn + option__wannier90files__w90_eig
186 if (sys%st%d%ispin == spinors) w90_what_default = w90_what_default + option__wannier90files__w90_spn
187 call parse_variable(global_namespace, 'Wannier90Files', w90_what_default, w90_what)
188
189 !%Variable Wannier90UseTD
190 !%Type logical
191 !%Default no
192 !%Section Utilities::oct-wannier90
193 !%Description
194 !% By default oct-wannier90 uses the ground-state states to compute the necessary information.
195 !% By setting this variable to yes, oct-wannier90 will use the TD states instead.
196 !%End
197 call parse_variable(global_namespace, 'Wannier90UseTD', .false., read_td_states)
198
199 !%Variable Wannier90UseSCDM
200 !%Type logical
201 !%Default no
202 !%Section Utilities::oct-wannier90
203 !%Description
204 !% By default oct-wannier90 uses the projection method to generate the .amn file.
205 !% By setting this variable to yes, oct-wannier90 will use SCDM method instead.
206 !%End
207 call parse_variable(global_namespace, 'Wannier90UseSCDM', .false., w90_scdm)
208 if (w90_scdm) then
209 !%Variable SCDMsigma
210 !%Type float
211 !%Default 0.2
212 !%Section Utilities::oct-wannier90
213 !%Description
214 !% Broadening of SCDM smearing function.
215 !%End
216 call parse_variable(global_namespace, 'SCDMsigma', 0.2_real64, scdm_sigma)
217
218 !%Variable SCDMmu
219 !%Type float
220 !%Section Utilities::oct-wannier90
221 !%Description
222 !% Energy range up to which states are considered for SCDM.
223 !%End
224 call parse_variable(global_namespace, 'SCDMmu', m_huge, scdm_mu)
225 end if
226
227 if (sys%kpoints%use_symmetries) then
228 message(1) = 'oct-wannier90: k-points symmetries are not allowed'
229 call messages_fatal(1)
230 end if
231 if (sys%kpoints%use_time_reversal) then
232 message(1) = 'oct-wannier90: time-reversal symmetry is not allowed'
233 call messages_fatal(1)
234 end if
235 if (sys%kpoints%reduced%nshifts > 1) then
236 message(1) = 'oct-wannier90: Wannier90 does not allow for multiple shifts of the k-point grid'
237 call messages_fatal(1)
238 end if
239
240 if (sys%st%d%ispin /= unpolarized) then
241 call messages_experimental("oct-wannier90 with SpinComponnents /= unpolarized")
242 end if
243
244 w90_spinors = .false.
245
246 !%Variable Wannier90SpinChannel
247 !%Type integer
248 !%Section Utilities::oct-wannier90
249 !%Description
250 !% Spin channel used for the Wannierization
251 !%End
252 call parse_variable(global_namespace, 'Wannier90SpinChannel', 1, w90_spin_channel)
253
254 ! "Correct" rlattice/klattice for Wannier90 (3D periodicity assumed here).
255 factor = m_one
256 do idir = sys%space%periodic_dim+1, sys%space%dim
257 factor(idir) = m_two * sys%gr%box%bounding_box_l(idir)
258 end do
259 call sys%ions%latt%scale(factor)
260
261 ! create setup files
262 select case (w90_mode)
263 case (option__wannier90mode__w90_setup)
264 call wannier90_setup(sys%ions, sys%kpoints, sys%space)
265
266 ! load states and calculate interface files
267 case (option__wannier90mode__w90_output)
268 call wannier90_output()
269
270 case (option__wannier90mode__w90_wannier)
271 !%Variable Wannier90ComputeBlochSums
272 !%Type logical
273 !%Default no
274 !%Section Utilities::oct-wannier90
275 !%Description
276 !% Only for Wannier90Mode = w90_wannier
277 !% By setting this variable to yes, oct-wannier90 will also output the Bloch sums
278 !% of the Wannier orbitals. This creates one file per k-point.
279 !%End
280 call parse_variable(global_namespace, 'Wannier90ComputeBlochSums', .false., w90_bloch_sums)
281
283
284 ! normal interface run
285 call states_elec_allocate_wfns(sys%st, sys%gr, wfs_type = type_cmplx, skip=exclude_list)
286 if (read_td_states) then
287 call restart%init(global_namespace, restart_td, restart_type_load, &
288 sys%mc, ierr, sys%gr)
289 else
290 call restart%init(global_namespace, restart_gs, restart_type_load, &
291 sys%mc, ierr, sys%gr)
292 end if
293
294 if (ierr == 0) then
295 call states_elec_look(restart, nik, dim, nst, ierr)
296 if (sys%st%d%ispin == spin_polarized) then
297 nik = nik / 2
298 end if
299 if (dim == sys%st%d%dim .and. nik == sys%kpoints%reduced%npoints .and. nst >= sys%st%nst) then
300 call states_elec_load(restart, global_namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
301 ierr, iter, label = ": wannier90", skip=exclude_list)
302 else
303 write(message(1),'(a)') 'Restart structure not commensurate.'
304 call messages_fatal(1)
305 end if
306 end if
307 call restart%end()
308
309 call generate_wannier_states(sys%space, sys%gr, sys%ions, sys%st, sys%kpoints)
310 case default
311 message(1) = "Wannier90Mode is set to an unsupported value."
312 call messages_fatal(1)
313 end select
314
315 safe_deallocate_a(exclude_list)
316 safe_deallocate_a(band_index)
317 safe_deallocate_a(w90_nnk_list)
318 safe_deallocate_a(w90_proj_centers)
319 safe_deallocate_a(w90_proj_lmr)
320
321 safe_deallocate_p(sys)
322 call fft_all_end()
323 call io_end()
325 call messages_end()
326 call parser_end()
327 call global_end()
328
329contains
330
331 ! --------------------------------------------------------------------------
332 subroutine wannier90_setup(ions, kpoints, space)
333 type(ions_t), intent(in) :: ions
334 type(kpoints_t), intent(in) :: kpoints
335 class(space_t), intent(in) :: space
336
337 character(len=80) :: filename
338 integer :: w90_win, ia, axis(3), npath
339
340 push_sub(wannier90_setup)
341
342 assert(space%dim == 3)
343
344 ! open win file
345 filename = trim(adjustl(w90_prefix)) //'.win'
346 w90_win = io_open(trim(filename), global_namespace, action='write')
347
348 write(w90_win,'(a)') '# this file has been created by the Octopus wannier90 utility'
349 write(w90_win,'(a)') ' '
350
351 ! write direct lattice vectors (in angstrom)
352 write(w90_win,'(a)') 'begin unit_cell_cart'
353 write(w90_win,'(a)') 'Ang'
354 do idim = 1,3
355 write(w90_win,'(f13.8,f13.8,f13.8)') units_from_atomic(unit_angstrom, ions%latt%rlattice(1:3,idim))
356 end do
357 write(w90_win,'(a)') 'end unit_cell_cart'
358 write(w90_win,'(a)') ' '
359
360 write(w90_win,'(a)') 'begin atoms_frac'
361 do ia = 1, ions%natoms
362 write(w90_win,'(a,2x,f13.8,f13.8,f13.8)') trim(ions%atom(ia)%label), ions%latt%cart_to_red(ions%pos(:, ia))
363 end do
364 write(w90_win,'(a)') 'end atoms_frac'
365 write(w90_win,'(a)') ' '
366
367 ! This is a default value. In practice, one should use projections
368 write(w90_win,'(a)') 'use_bloch_phases = .true.'
369 write(w90_win,'(a)') ' '
370
371 write(w90_win,'(a10,i4)') 'num_bands ', sys%st%nst
372 write(w90_win,'(a9,i4)') 'num_wann ', sys%st%nst
373 write(w90_win,'(a)') ' '
374
375 if (sys%st%d%ispin == spinors) then
376 write(w90_win,'(a)') 'spinors = .true.'
377 end if
378 if (sys%st%d%ispin == spin_polarized) then
379 if (w90_spin_channel == 1) then
380 write(w90_win, '(a)') 'spin = up'
381 else if (w90_spin_channel == 2) then
382 write(w90_win, '(a)') 'spin = down'
383 else
384 message(1) = 'Wannier90SpinChannel value is invalid.'
385 call messages_fatal(1)
386 end if
387 end if
388
389 ! This is for convenience. This is needed for plotting the Wannier states, if requested.
390 write(w90_win,'(a)') 'write_u_matrices = .true.'
391 write(w90_win,'(a)') 'write_xyz = .true.'
392 write(w90_win,'(a)') ' '
393
394 if (kpoints%reduced%npoints == 1) then
395 write(w90_win,'(a)') 'gamma_only = .true.'
396 write(w90_win,'(a)') ' '
397 else
398 if (.not. parse_is_defined(global_namespace, 'KPointsGrid')) then
399 message(1) = 'oct-wannier90: need Monkhorst-Pack grid. Please specify %KPointsGrid'
400 call messages_fatal(1)
401 end if
402
403 ! In case the user used also a k-point path, we ignore it
404 npath = kpoints%nkpt_in_path()
405
406 axis(1:3) = kpoints%nik_axis(1:3)
407 assert(product(kpoints%nik_axis(1:3)) == kpoints%reduced%npoints - npath)
408
409 write(w90_win,'(a8,i4,i4,i4)') 'mp_grid =', axis(1:3)
410 write(w90_win,'(a)') ' '
411 write(w90_win,'(a)') 'begin kpoints '
412 ! Put a minus sign here for the wrong convention in Octopus
413
414 do ii = 1, kpoints%reduced%npoints-npath
415 write(w90_win,'(f13.8,f13.8,f13.8)') - kpoints%reduced%red_point(1:3,ii)
416 end do
417 write(w90_win,'(a)') 'end kpoints '
418 end if
419
420 call io_close(w90_win)
421
422 pop_sub(wannier90_setup)
423
424 end subroutine wannier90_setup
425
426 ! --------------------------------------------------------------------------
427 subroutine wannier90_output()
428 integer :: ik_real
429
430 push_sub(wannier90_output)
431
433
434 ! normal interface run
435 call states_elec_allocate_wfns(sys%st, sys%gr, wfs_type = type_cmplx, skip=exclude_list)
436 if (read_td_states) then
437 call restart%init(global_namespace, restart_td, restart_type_load, &
438 sys%mc, ierr, sys%gr)
439 else
440 call restart%init(global_namespace, restart_gs, restart_type_load, &
441 sys%mc, ierr, sys%gr)
442 end if
443
444 if (ierr == 0) then
445 call states_elec_look(restart, nik, dim, nst, ierr)
446 if (sys%st%d%ispin == spin_polarized) then
447 nik = nik / 2
448 end if
449 if (dim == sys%st%d%dim .and. nik == sys%kpoints%reduced%npoints .and. nst >= sys%st%nst) then
450 call states_elec_load(restart, global_namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
451 ierr, iter, label = ": wannier90", skip=exclude_list)
452 else
453 write(message(1),'(a)') 'Restart structure not commensurate.'
454 call messages_fatal(1)
455 end if
456 end if
457 call restart%end()
458
459 if (w90_scdm) then
460 nik = w90_num_kpts
461 safe_allocate(jpvt(1:sys%gr%np_global*sys%st%d%dim))
462 safe_allocate(psi(1:sys%gr%np, 1:sys%st%d%dim))
463 safe_allocate(occ_temp(1:w90_num_bands))
464
465 ! smear the states at gamma
466 do ist = 1, w90_num_bands
467 occ_temp(ist)= sys%st%occ(ist, 1)
468 sys%st%occ(ist, 1)=m_half*erfc((sys%st%eigenval(ist, 1)-scdm_mu) / scdm_sigma)
469 end do
470
471 call zstates_elec_rrqr_decomposition(sys%st, sys%namespace, sys%gr, w90_num_bands, .true., 1, jpvt)
472
473 ! reset occupations at gamma
474 do ist = 1, w90_num_bands
475 sys%st%occ(ist, 1) = occ_temp(ist)
476 end do
477
478 safe_allocate(uk(1:w90_num_bands, 1:w90_num_bands, 1:nik))
479
480 ! auxiliary arrays for scdm procedure
481 safe_allocate(chi(1:w90_num_bands, 1:w90_num_bands))
482 safe_allocate(chi_diag(1:w90_num_bands, 1:w90_num_bands))
483 safe_allocate(chi2(1:w90_num_bands, 1:w90_num_bands))
484 safe_allocate(chi_eigenval(1:w90_num_bands))
485
486 chi(1:w90_num_bands, 1:w90_num_bands) = m_zero
487
488 do ik = 1, nik
489 kvec(:) = sys%kpoints%reduced%point(:, ik)
490
491 if (sys%st%d%ispin == spin_polarized) then
492 ik_real = (ik-1)*2 + w90_spin_channel
493 else
494 ik_real = ik
495 end if
496
497 do ist = 1, w90_num_bands
498 sender = -1
499 if (state_kpt_is_local(sys%st, ist, ik_real)) then
500 call states_elec_get_state(sys%st, sys%gr, ist, ik_real, psi)
501 sender = sys%st%d%kpt%mpi_grp%rank
502 endif
503 call sys%st%d%kpt%mpi_grp%allreduce_inplace(sender, 1, mpi_integer, mpi_max)
504 assert(sender >= 0)
505 call sys%st%d%kpt%mpi_grp%bcast(psi, sys%gr%np*sys%st%d%dim, mpi_double_complex, sender)
506 smear=m_half * erfc((sys%st%eigenval(ist, ik_real) - scdm_mu) / scdm_sigma)
507 ! NOTE: here check for domain parallelization
508 do jst = 1, w90_num_bands
509 chi(ist, jst) = smear * conjg(psi(jpvt(jst), 1)) &
510 * exp(m_zi * dot_product(sys%gr%x(1:3, jpvt(jst)), kvec(1:3)))
511 end do
512 end do
513
514 ! loewdin orhtogonalization of chi.chi
515 ! this can also be done with SVD, which might be more stable!?
516 chi_diag = matmul(conjg(transpose(chi)), chi)
517 call lalg_eigensolve(w90_num_bands, chi_diag, chi_eigenval)
518 chi2 = conjg(transpose(chi_diag))
519
520 !we need the eigenvalues to be >0
521 if (any(chi_eigenval(:) .lt. m_zero)) then
522 message(1) = 'SCDM Wannierization failed because chi matrix is'
523 message(2) = 'ill conditioned. Try increasingin scdm_sigma and/or'
524 message(3) = 'change scdm_mu.'
525 call messages_fatal(3)
526 end if
527
528 do ist = 1, w90_num_bands
529 chi_eigenval(ist) = m_one / sqrt(chi_eigenval(ist))
530 chi2(ist, 1:w90_num_bands) = chi_eigenval(ist) * chi2(ist, 1:w90_num_bands)
531 end do
532 ! the loewdin result would be: matmul(chi_diag,chi2)
533 ! to get the wannier gauge U(k) we multiply this with the original chi
534 uk(:,:,ik) = matmul(chi, matmul(chi_diag,chi2))
535
536 end do
537
538 safe_deallocate_a(chi)
539 safe_deallocate_a(psi)
540 safe_deallocate_a(chi_diag)
541 safe_deallocate_a(chi2)
542 safe_deallocate_a(chi_eigenval)
543 safe_deallocate_a(jpvt)
544 safe_deallocate_a(psi)
545 safe_deallocate_a(occ_temp)
546 end if
547
548 ! ---- actual interface work ----------
549 if (bitand(w90_what, option__wannier90files__w90_mmn) /= 0) then
550 call create_wannier90_mmn(sys%gr, sys%st)
551 end if
552
553 if (bitand(w90_what, option__wannier90files__w90_unk) /= 0) then
554 call write_unk(sys%space, sys%gr, sys%st, formatted=.false.)
555 end if
556
557 if (bitand(w90_what, option__wannier90files__w90_amn) /= 0) then
558 call create_wannier90_amn(sys%space, sys%gr, sys%ions%latt, sys%st, sys%kpoints)
559 end if
560
561 if (bitand(w90_what, option__wannier90files__w90_eig) /= 0) then
563 end if
564
565 if (bitand(w90_what, option__wannier90files__w90_spn) /= 0) then
566 call create_wannier90_spn(sys%gr, sys%st)
567 end if
568
569 safe_deallocate_a(uk)
570 safe_deallocate_a(w90_spin_proj_component)
571 safe_deallocate_a(w90_spin_proj_axis)
572
573 pop_sub(wannier90_output)
574 end subroutine wannier90_output
575
576 ! --------------------------------------------------------------------------
577 subroutine read_wannier90_files()
578 integer :: w90_nnkp, itemp, dummyint, io, spin_channel_win
579 character(len=80) :: filename, dummy, dummy1, dummy2, line
580 logical :: exist, parse_is_ok
581 real(real64) :: dummyr(7)
582
583 push_sub(read_wannier90_files)
584
585 w90_num_kpts = product(sys%kpoints%nik_axis(1:3))
586 assert(w90_num_kpts == sys%st%nik)
587
588
589 w90_num_exclude = 0
590
591 ! open nnkp file
592 filename = trim(adjustl(w90_prefix)) //'.nnkp'
593
594 message(1) = "oct-wannier90: Parsing "//filename
595 call messages_info(1)
596
597 inquire(file=filename,exist=exist)
598 if (.not. exist) then
599 message(1) = 'oct-wannier90: Cannot find specified Wannier90 nnkp file.'
600 write(message(2),'(a)') 'Please run wannier90.x -pp '// trim(adjustl(w90_prefix)) // ' first.'
601 call messages_fatal(2)
602 end if
603
604 parse_is_ok = .false.
605
606 ! check number of k-points
607 w90_nnkp = io_open(trim(filename), global_namespace, action='read')
608 do
609 read(w90_nnkp, *, iostat=io) dummy, dummy1
610 if (io == iostat_end) exit
611
612 if (dummy == 'begin' .and. dummy1 == 'kpoints') then
613 read(w90_nnkp,*) itemp
614 if (itemp /= w90_num_kpts) then
615 message(1) = 'oct-wannier90: wannier90 setup seems to have been done with a different number of k-points.'
616 call messages_fatal(1)
617 else
618 parse_is_ok = .true.
619 exit
620 end if
621 end if
622 end do
623 call io_close(w90_nnkp)
624
625 if (.not. parse_is_ok) then
626 message(1) = 'oct-wannier90: Did not find the kpoints block in nnkp file'
627 call messages_fatal(1)
628 end if
629 parse_is_ok = .false.
630
631 ! read from nnkp file
632 ! find the nnkpts block
633 w90_nnkp = io_open(trim(filename), global_namespace, action='read', position='rewind')
634 do
635 read(w90_nnkp, *, iostat=io) dummy, dummy1
636 if (io == iostat_end) exit !End of file
637
638 if (dummy == 'begin' .and. dummy1 == 'nnkpts') then
639 read(w90_nnkp,*) w90_nntot
640 safe_allocate(w90_nnk_list(1:5, 1:w90_num_kpts * w90_nntot))
641 do ii = 1, w90_num_kpts * w90_nntot
642 read(w90_nnkp,*) w90_nnk_list(1:5, ii)
643 end do
644 !make sure we are at the end of the block
645 read(w90_nnkp,*) dummy
646 if (dummy /= 'end') then
647 message(1) = 'oct-wannier90: There dont seem to be enough k-points in nnkpts file to.'
648 call messages_fatal(1)
649 end if
650 parse_is_ok = .true.
651 exit
652 end if
653 end do
654
655 if (.not. parse_is_ok) then
656 message(1) = 'oct-wannier90: Did not find nnkpts block in nnkp file'
657 call messages_fatal(1)
658 end if
659
660 ! read from nnkp file
661 ! find the exclude block
662 safe_allocate(exclude_list(1:sys%st%nst))
663 !By default we use all the bands
664 exclude_list(1:sys%st%nst) = .false.
665 rewind(w90_nnkp)
666 do
667 read(w90_nnkp, *, iostat=io) dummy, dummy1
668 if (io == iostat_end) exit !End of file
669 if (dummy == 'begin' .and. dummy1 == 'exclude_bands') then
670 read(w90_nnkp, *) w90_num_exclude
671 do ii = 1, w90_num_exclude
672 read(w90_nnkp, *) itemp
673 if(itemp > sys%st%nst) then
674 message(1) = 'oct-wannier90: The exclude_bands list contains a state index higher than the number of states.'
675 call messages_fatal(1)
676 end if
677 exclude_list(itemp) = .true.
678 end do
679 !make sure we are at the end of the block
680 read(w90_nnkp, *) dummy
681 if (dummy /= 'end') then
682 message(1) = 'oct-wannier90: There dont seem to be enough bands in exclude_bands list.'
683 call messages_fatal(1)
684 end if
685 exit
686 end if
687 end do
688 call io_close(w90_nnkp)
689
690 !We get the number of bands
691 w90_num_bands = sys%st%nst - w90_num_exclude
692
693 safe_allocate(band_index(1:sys%st%nst))
694 itemp = 0
695 do ii = 1, sys%st%nst
696 if (exclude_list(ii)) cycle
697 itemp = itemp + 1
698 band_index(ii) = itemp
699 end do
700
701 if (bitand(w90_what, option__wannier90files__w90_amn) /= 0 &
702 .or. w90_mode == option__wannier90mode__w90_wannier ) then
703 ! parse file again for definitions of projections
704 w90_nnkp = io_open(trim(filename), global_namespace, action='read', position='rewind')
705
706 do
707 read(w90_nnkp, *, iostat=io) dummy, dummy1
708 if (io == iostat_end) then !End of file
709 message(1) = 'oct-wannier90: Did not find projections block in w90.nnkp file'
710 call messages_fatal(1)
711 end if
712
713 if (dummy == 'begin' .and. (dummy1 == 'projections' .or. dummy1 == 'spinor_projections')) then
714
715 if (dummy1 == 'spinor_projections') then
716 w90_spinors = .true.
717 if (sys%st%d%ispin /= spinors) then
718 message(1) = 'oct-wannier90: Spinor = .true. is only valid with spinors wavefunctions.'
719 call messages_fatal(1)
720 end if
721
722 message(1) = 'oct-wannier90: Spinor interface incomplete. Note there is no quantization axis implemented'
723 call messages_warning(1)
724 else
725 if (sys%st%d%ispin == spinors) then
726 message(1) = 'oct-wannier90: Octopus has spinors wavefunctions but spinor_projections is not defined.'
727 message(2) = 'oct-wannier90: Please check the input file for wannier 90.'
728 call messages_fatal(2)
729 end if
730 end if
731
732 read(w90_nnkp, *) w90_nproj
733 ! num_wann is given in w90.win, not double checked here
734 ! I assume that the wannier90.x -pp run has checked this
735 w90_num_wann = w90_nproj
736 ! In case of no projections, we use the number of bands
737 if(w90_nproj == 0) w90_num_wann = w90_num_bands
738
739 safe_allocate(w90_proj_centers(1:3, 1:w90_nproj))
740 safe_allocate(w90_proj_lmr(1:w90_nproj, 1:3))
741 if (w90_spinors) then
742 safe_allocate(w90_spin_proj_component(1:w90_nproj))
743 end if
744 if (w90_spinors) then
745 safe_allocate(w90_spin_proj_axis(1:w90_nproj, 1:3))
746 end if
747
748 do ii = 1, w90_nproj
749 read(w90_nnkp, *) w90_proj_centers(1:3, ii), w90_proj_lmr(ii, 1:3)
750 ! skip a line for now
751 read(w90_nnkp, *) dummyr(1:7)
752 if (w90_spinors) then
753 read(w90_nnkp, *) w90_spin_proj_component(ii), w90_spin_proj_axis(ii, 1:3)
754 ! use octopus spindim conventions
755 if (w90_spin_proj_component(ii) == -1) w90_spin_proj_component(ii) = 2
756 end if
757 end do
758 !make sure we are at the end of the block
759 read(w90_nnkp, *) dummy
760 if (dummy /= 'end') then
761 message(1) = 'oct-wannier90: There dont seem to be enough projections in nnkpts file to.'
762 call messages_fatal(1)
763 end if
764 exit
765 end if
766 end do
767
768 ! look for auto_projection block
769 scdm_proj = .false.
770 do
771 read(w90_nnkp, *, iostat=io) dummy, dummy1
772 if (io == iostat_end) exit !End of file
773
774 if (dummy == 'begin' .and. dummy1 == 'auto_projections') then
775 scdm_proj = .true.
776 read(w90_nnkp, *) w90_nproj
777 w90_num_wann = w90_nproj
778
779 if (.not. w90_scdm) then
780 message(1) = 'oct-wannier90: Found auto_projections block. Currently the only implemented automatic way'
781 message(2) = 'oct-wannier90: to compute projections is the SCDM method.'
782 message(3) = 'oct-wannier90: Please set Wannier90UseSCDM = yes in the inp file.'
783 call messages_fatal(3)
784 end if
785
786 if (w90_nproj /= w90_num_bands) then
787 message(1) = 'oct-wannier90: In auto_projections block first row needs to be equal to num_bands.'
788 call messages_fatal(1)
789 end if
790 read(w90_nnkp, *) dummyint
791 if (dummyint /= 0) then
792 message(1) = 'oct-wannier90: The second row in auto_projections has to be 0, per Wannier90 documentation.'
793 call messages_fatal(1)
794 end if
795 end if
796 end do
797 call io_close(w90_nnkp)
798
799 end if
800
801 message(1) = "oct-wannier90: Finished parsing "//filename
802 call messages_info(1)
803
804 ! Look extra variables variable
805 ! open win file
806 filename = trim(adjustl(w90_prefix)) //'.win'
807 message(1) = "oct-wannier90: Parsing "//filename
808 call messages_info(1)
809 w90_nnkp = io_open(trim(filename), global_namespace, action='read', position='rewind')
810 do
811 read(w90_nnkp, fmt='(a)', iostat=io) line
812 if (io == iostat_end) exit !End of file
813 if (index(line, '=') > 0) then
814 read(line, *, iostat=io) dummy, dummy2, dummy1
815 else
816 read(line, *, iostat=io) dummy, dummy1
817 end if
818
819 !Spin
820 if (dummy == 'spin') then
821 if (sys%st%d%ispin /= spin_polarized) then
822 message(1) = 'oct-wannier90: The variable spin is set for a non spin-polarized calculation.'
823 call messages_fatal(1)
824 end if
825
826 if (dummy1 == 'up') then
827 spin_channel_win = 1
828 else if (dummy1 == 'down') then
829 spin_channel_win = 2
830 else
831 message(1) = 'oct-wannier90: Error parsing the variable spin.'
832 call messages_fatal(1)
833 end if
834 if (spin_channel_win /= w90_spin_channel) then
835 message(1) = 'spin polarization input from .win and Wannier90SpinChannel do not agree.'
836 call messages_fatal(1)
837 end if
838 end if
839 end do
840 call io_close(w90_nnkp)
841
842 if (sys%st%d%ispin == spin_polarized) then
843 write(message(1), '(a,i1)') 'oct-wannier90: Using spin channel ', w90_spin_channel
844 call messages_info(1)
845 end if
846
847 message(1) = "oct-wannier90: Finished parsing "//filename
848 call messages_info(1)
849
850 pop_sub(read_wannier90_files)
851
852 end subroutine read_wannier90_files
853
854 ! --------------------------------------------------------------------------
855 subroutine create_wannier90_mmn(mesh, st)
856 class(mesh_t), intent(in) :: mesh
857 type(states_elec_t), target, intent(in) :: st
858
859 integer :: ist, jst, ik, ip, w90_mmn, iknn, idim, ibind
860 real(real64) :: Gcart(3)
861 integer :: G(3)
862 character(len=80) :: filename
863 complex(real64), allocatable :: overlap(:)
864 complex(real64), allocatable :: psim(:,:), psin(:,:), phase(:)
865 type(wfs_elec_t), pointer :: batch
866 integer :: inode, node_fr, node_to
867 type(mpi_request) :: send_req
868
869 push_sub(create_wannier90_mmn)
870
871 call profiling_in("W90_MMN")
872
873 if (st%parallel_in_states) then
874 call messages_not_implemented("w90_mmn output with states parallelization")
875 end if
876
877 message(1) = "Info: Computing the overlap matrix"
878 call messages_info(1)
879
880
881 filename = './'// trim(adjustl(w90_prefix))//'.mmn'
882 w90_mmn = io_open(trim(filename), global_namespace, action='write')
883
884 ! write header
885 if (st%system_grp%is_root()) then
886 write(w90_mmn,*) 'Created by oct-wannier90'
887 write(w90_mmn,*) w90_num_bands, w90_num_kpts, w90_nntot
888 end if
889
890 safe_allocate(psim(1:mesh%np, 1:st%d%dim))
891 safe_allocate(psin(1:mesh%np, 1:st%d%dim))
892 safe_allocate(phase(1:mesh%np))
893 safe_allocate(overlap(1:w90_num_bands))
894
895
896 ! loop over the pairs specified in the nnkp file (read before in init)
897 do ii = 1, w90_num_kpts * w90_nntot
898 ik = w90_nnk_list(1, ii)
899 iknn = w90_nnk_list(2, ii)
900 g(1:3) = w90_nnk_list(3:5, ii)
901 if (st%system_grp%is_root()) write(w90_mmn, '(I10,2x,I10,2x,I3,2x,I3,2x,I3)') ik, iknn, g
902 ! For mixed periodicity we remove the G vectors along the non-periodic directions
903 ! as this correspond to the infinite vacuum limit
904 g(sys%space%periodic_dim+1:sys%space%dim) = 0
905
906 ! For spin-polarized calculations, we select the right k-point
907 if (sys%st%d%ispin == spin_polarized) then
908 ik = (ik-1)*2 + w90_spin_channel
909 iknn = (iknn-1)*2 + w90_spin_channel
910 end if
911
912
913 ! Only treat local k-points
914 if(ik >= st%d%kpt%start .and. ik <= st%d%kpt%end) then
915
916 ! Wannier90 treats everything fully periodic
917 ! Conversion is done with the 3D "correct" klattice
918 call kpoints_to_absolute(sys%ions%latt, real(g, real64) , gcart)
919
920 ! Phase that gives u_{n,k+G}(r) from u_{nk}(r)
921 ! Here the minus sign of Octopus cancels with minus sign of the input G
922 ! (ik and iknn correspond in the list to minus the k-points in Octopus)
923 if (any(g /= 0)) then
924 !$omp parallel do
925 do ip = 1, mesh%np
926 phase(ip) = exp(-m_zi*dot_product(mesh%x(1:3, ip), gcart(1:3)))
927 end do
928 end if
929
930 end if
931
932 ! Loop over distributed bands
933 ! Slow index in M_{mn}^{k,b}, so \psi_n
934 do jst = 1, st%nst
935 if (exclude_list(jst)) cycle
936
937 ! Communication for the local states
938 if ( .not. st%d%kpt%parallel .and. .not. st%parallel_in_states) then
939 call states_elec_get_state(st, mesh, jst, iknn, psin)
940 else
941 node_fr = -1
942 node_to = -1
943 do inode = 0, st%d%kpt%mpi_grp%size-1
944 if(iknn >= st%st_kpt_task(inode,3) .and. iknn <= st%st_kpt_task(inode,4)) then
945 node_fr = inode
946 end if
947 if(ik >= st%st_kpt_task(inode,3) .and. ik <= st%st_kpt_task(inode,4)) then
948 node_to = inode
949 end if
950 end do
951 assert(node_fr > -1)
952 assert(node_to > -1)
953
954 send_req = mpi_request_null
955 ! We have locally the k-point
956 if (state_kpt_is_local(st, jst, iknn)) then
957 call states_elec_get_state(st, mesh, jst, iknn, psin)
958 ! We send it only if we don`t want to use it locally
959 if(node_to /= st%d%kpt%mpi_grp%rank) then
960 call st%d%kpt%mpi_grp%isend(psin, mesh%np*st%d%dim, mpi_double_complex, node_to, send_req)
961 end if
962 end if
963 ! We receive the desired state, only if it is not a local one
964 if(node_to == st%d%kpt%mpi_grp%rank .and. node_to /= node_fr) then
965 call st%d%kpt%mpi_grp%recv(psin, mesh%np*st%d%dim, mpi_double_complex, node_fr)
966 end if
967 if (send_req /= mpi_request_null) then
968 call st%d%kpt%mpi_grp%wait(send_req)
969 end if
970 end if
971
972 overlap = m_zero
973
974 if(ik >= st%d%kpt%start .and. ik <= st%d%kpt%end) then
975
976 ! Do not apply the phase if the phase factor is null
977 if (any(g /= 0)) then
978 ! add phase
979 !$omp parallel
980 do idim = 1, st%d%dim
981 !$omp do
982 do ip = 1, mesh%np
983 psin(ip, idim) = psin(ip, idim) * phase(ip)
984 end do
985 end do
986 !$omp end parallel
987 end if
988
989
990 ! See Eq. (25) in PRB 56, 12847 (1997)
991 ! Loop over local k-points
992 ! Fast index in M_{mn}^{k,b}, so m
993 do ist = 1, st%nst
994 if (exclude_list(ist)) cycle
995
996 batch => st%group%psib(st%group%iblock(ist), ik)
997
998 select case (batch%status())
999 case (batch_not_packed)
1000 overlap(band_index(ist)) = m_z0
1001 do idim = 1, st%d%dim
1002 ibind = batch%inv_index((/ist, idim/))
1003 overlap(band_index(ist)) = overlap(band_index(ist)) + &
1004 zmf_dotp(mesh, batch%zff_linear(:, ibind), psin(:,idim), reduce = .false.)
1005 end do
1006 !Not properly done at the moment
1008 call states_elec_get_state(st, mesh, ist, ik, psim)
1009 overlap(band_index(ist)) = zmf_dotp(mesh, st%d%dim, psim, psin, reduce = .false.)
1010 end select
1011 end do !ist
1012 end if
1013
1014 call profiling_in("W90_MMN_REDUCE")
1015 call mesh%allreduce(overlap)
1016 call profiling_out("W90_MMN_REDUCE")
1017
1018 if(st%d%kpt%parallel) then
1019 call comm_allreduce(st%d%kpt%mpi_grp, overlap)
1020 end if
1021
1022 ! write to W90 file
1023 if (st%system_grp%is_root()) then
1024 do ist = 1, st%nst
1025 if (exclude_list(ist)) cycle
1026 write(w90_mmn,'(e18.10,2x,e18.10)') overlap(band_index(ist))
1027 end do
1028 end if
1029
1030 end do !jst
1031 end do
1032
1033 call io_close(w90_mmn)
1034
1035 safe_deallocate_a(psim)
1036 safe_deallocate_a(psin)
1037 safe_deallocate_a(phase)
1038 safe_deallocate_a(overlap)
1039
1040 call profiling_out("W90_MMN")
1041
1042 pop_sub(create_wannier90_mmn)
1043
1044 end subroutine create_wannier90_mmn
1045
1046 ! --------------------------------------------------------------------------
1047 subroutine create_wannier90_eig()
1048 integer :: ist, ik, w90_eig
1049 character(len=80) :: filename
1050
1051 push_sub(create_wannier90_eig)
1052
1053 if (sys%st%parallel_in_states) then
1054 call messages_not_implemented("w90_eig output with states parallelization")
1055 end if
1056
1057 if (sys%st%system_grp%is_root()) then
1058 filename = './'//trim(adjustl(w90_prefix))//'.eig'
1059 w90_eig = io_open(trim(filename), global_namespace, action='write')
1060 do ik = 1, w90_num_kpts
1061 do ist = 1, sys%st%nst
1062 if (exclude_list(ist)) cycle
1063 if (sys%st%d%ispin /= spin_polarized) then
1064 write(w90_eig,'(I5,2x,I8,2x,e18.10)') band_index(ist), ik, &
1065 units_from_atomic(unit_ev, sys%st%eigenval(ist, ik))
1066 else
1067 write(w90_eig,'(I5,2x,I8,2x,e18.10)') band_index(ist), ik, &
1068 units_from_atomic(unit_ev, sys%st%eigenval(ist, (ik-1)*2+w90_spin_channel))
1069 end if
1070 end do
1071 end do
1072
1073 call io_close(w90_eig)
1074 end if
1075
1076 pop_sub(create_wannier90_eig)
1077 end subroutine create_wannier90_eig
1078
1079 ! --------------------------------------------------------------------------
1080 subroutine write_unk(space, mesh, st, formatted)
1081 class(space_t), intent(in) :: space
1082 class(mesh_t), intent(in) :: mesh
1083 type(states_elec_t), intent(in) :: st
1084 logical, intent(in) :: formatted
1085
1086 integer :: ist, ik, unk_file, ispin
1087 integer :: ix, iy, iz
1088 real(real64) :: w_real, w_imag
1089 character(len=80) :: filename
1090 complex(real64), allocatable :: psi(:)
1091 type(cube_t) :: cube
1092 type(cube_function_t) :: cf
1093
1094 push_sub(write_unk)
1095
1096 if (st%d%kpt%parallel) then
1097 call messages_not_implemented("w90_unk output with k-point parallelization")
1098 end if
1099
1100 if (sys%gr%parallel_in_domains) then
1101 call messages_not_implemented("w90_unk output with domain parallelization")
1102 end if
1103
1104 if (st%parallel_in_states) then
1105 call messages_not_implemented("w90_unk output with states parallelization")
1106 end if
1107
1108 call messages_experimental("Wannier90Files = w90_unk")
1109
1110
1111 safe_allocate(psi(1:mesh%np))
1112
1113 call cube_init(cube, mesh%idx%ll, global_namespace, space, mesh%spacing, &
1114 mesh%coord_system, need_partition=.not.mesh%parallel_in_domains)
1115 call cube_init_cube_map(cube, mesh)
1116
1117 call zcube_function_alloc_rs(cube, cf)
1118
1119 do ik = 1, w90_num_kpts
1120 do ispin = 1, st%d%dim
1121 if (st%system_grp%is_root()) then
1122 write(filename, '(a,i5.5,a1,i1)') './UNK', ik,'.', ispin
1123 ! write header
1124 if (formatted) then
1125 unk_file = io_open(trim(filename), global_namespace, action='write', form='formatted')
1126 write(unk_file, *) mesh%idx%ll(1:mesh%idx%dim), ik, w90_num_bands
1127 else
1128 unk_file = io_open(trim(filename), global_namespace, action='write', form='unformatted')
1129 write(unk_file) mesh%idx%ll(1:mesh%idx%dim), ik, w90_num_bands
1130 end if
1131 end if
1132
1133 ! states
1134 do ist = 1, st%nst
1135 if (exclude_list(ist)) cycle
1136
1137 if (sys%st%d%ispin /= spin_polarized) then
1138 call states_elec_get_state(st, mesh, ispin, ist, ik, psi)
1139 else
1140 call states_elec_get_state(st, mesh, ispin, ist, (ik-1)*2+w90_spin_channel, psi)
1141 end if
1143 ! put the density in the cube
1144 ! Note: At the moment this does not work for domain parallelization
1145 assert(.not. cube%parallel_in_domains)
1146 call zmesh_to_cube(mesh, psi, cube, cf)
1147
1148 if (st%system_grp%is_root()) then
1149 if (formatted) then
1150 do iz=1,cube%rs_n_global(3)
1151 do iy = 1,cube%rs_n_global(2)
1152 do ix = 1,cube%rs_n_global(1)
1153 w_real = real(cf%zrs(ix,iy,iz), real64)
1154 w_imag = aimag(cf%zrs(ix,iy,iz))
1155 write (unk_file, *) w_real, w_imag
1156 end do
1157 end do
1158 end do
1159 else
1160 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))
1161 end if
1162 end if
1163 end do
1164 if (st%system_grp%is_root()) call io_close(unk_file)
1165 end do
1166 end do
1167
1168 call zcube_function_free_rs(cube, cf)
1169 call cube_end(cube)
1170
1171 safe_deallocate_a(psi)
1172
1173 pop_sub(write_unk)
1174
1175 end subroutine write_unk
1176
1177 ! --------------------------------------------------------------------------
1178 subroutine create_wannier90_amn(space, mesh, latt, st, kpoints)
1179 class(space_t), intent(in) :: space
1180 class(mesh_t), intent(in) :: mesh
1181 type(lattice_vectors_t), intent(in) :: latt
1182 type(states_elec_t), intent(in) :: st
1183 type(kpoints_t), intent(in) :: kpoints
1184
1185 integer :: ist, ik, w90_amn, idim, iw, ip, ik_real
1186 real(real64) :: center(3), kpoint(3), threshold
1187 character(len=80) :: filename
1188 complex(real64), allocatable :: psi(:,:), phase(:), projection(:)
1189 real(real64), allocatable :: ylm(:)
1190 type(orbitalset_t), allocatable :: orbitals(:)
1191
1192 push_sub(create_wannier90_amn)
1193 call profiling_in("W90_AMN")
1194
1195 if (st%parallel_in_states) then
1196 call messages_not_implemented("w90_amn output with states parallelization")
1197 end if
1198
1199 filename = './'// trim(adjustl(w90_prefix))//'.amn'
1200 w90_amn = io_open(trim(filename), global_namespace, action='write')
1201
1202 ! write header
1203 if (st%system_grp%is_root()) then
1204 write(w90_amn,*) 'Created by oct-wannier90'
1205 write(w90_amn,*) w90_num_bands, w90_num_kpts, w90_num_wann
1206 end if
1207
1208 if (scdm_proj) then
1209
1210 message(1) = "Info: Writing projections obtained from SCDM."
1211 call messages_info(1)
1212
1213
1214 do ik = 1, w90_num_kpts
1215 do ist = 1, st%nst
1216 if (exclude_list(ist)) cycle
1217 if (st%system_grp%is_root()) then
1218 do iw = 1, w90_nproj
1219 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)
1220 end do
1221 end if
1222 end do !ist
1223 end do! ik
1224
1225 else
1226
1227 message(1) = "Info: Computing the projection matrix"
1228 call messages_info(1)
1229
1230 !We use the variable AOThreshold to deterine the threshold on the radii of the atomic orbitals
1231 call parse_variable(global_namespace, 'AOThreshold', 0.001_real64, threshold)
1232
1233 safe_allocate(orbitals(1:w90_nproj))
1234 ! precompute orbitals
1235 do iw=1, w90_nproj
1236 call orbitalset_init(orbitals(iw))
1237
1238 orbitals(iw)%norbs = 1
1239 orbitals(iw)%ndim = 1
1240 orbitals(iw)%radius = -log(threshold)
1241 orbitals(iw)%use_submesh = .false.
1242
1243 ! cartesian coordinate of orbital center
1244 center(1:3) = latt%red_to_cart(w90_proj_centers(1:3, iw))
1245 call submesh_init(orbitals(iw)%sphere, space, mesh, latt, center, orbitals(iw)%radius)
1246
1247 ! get dorb as submesh points
1248 safe_allocate(ylm(1:orbitals(iw)%sphere%np))
1249 ! (this is a routine from pwscf)
1250 call ylm_wannier(ylm, w90_proj_lmr(iw,1), w90_proj_lmr(iw,2), &
1251 orbitals(iw)%sphere%r, orbitals(iw)%sphere%rel_x, orbitals(iw)%sphere%np)
1252
1253 ! apply radial function
1254 if (w90_proj_lmr(iw,3) == 1) then
1255 do ip = 1,orbitals(iw)%sphere%np
1256 ylm(ip) = ylm(ip)*m_two*exp(-orbitals(iw)%sphere%r(ip))
1257 end do
1258 else
1259 call messages_not_implemented("oct-wannier90: r/=1 for the radial part")
1260 end if
1261
1262 safe_allocate(orbitals(iw)%zorb(1:orbitals(iw)%sphere%np, 1, 1))
1263 orbitals(iw)%zorb(1:orbitals(iw)%sphere%np, 1, 1) = ylm(1:orbitals(iw)%sphere%np)
1264 safe_deallocate_a(ylm)
1265
1266 safe_allocate(orbitals(iw)%phase(1:orbitals(iw)%sphere%np, st%d%kpt%start:st%d%kpt%end))
1267 orbitals(iw)%phase(:,:) = m_z0
1268 safe_allocate(orbitals(iw)%eorb_mesh(1:mesh%np, 1, 1, st%d%kpt%start:st%d%kpt%end))
1269 orbitals(iw)%eorb_mesh(:,:,:,:) = m_z0
1270
1271 call orbitalset_update_phase(orbitals(iw), space%dim, st%d%kpt, kpoints, st%d%ispin == spin_polarized, &
1272 kpt_max = w90_num_kpts)
1274 end do
1275
1276 safe_allocate(psi(1:mesh%np, 1:st%d%dim))
1277 safe_allocate(phase(1:mesh%np))
1278 safe_allocate(projection(1:w90_nproj))
1279
1280 do ik = 1, w90_num_kpts
1281 kpoint(1:space%dim) = kpoints%get_point(ik)
1282 !$omp parallel do
1283 do ip = 1, mesh%np
1284 phase(ip) = exp(-m_zi* sum(mesh%x(1:space%dim, ip) * kpoint(1:space%dim)))
1285 end do
1286
1287 !For spin-polarized calculations, we select the right k-point
1288 if (st%d%ispin == spin_polarized) then
1289 ik_real = (ik-1)*2 + w90_spin_channel
1290 else
1291 ik_real = ik
1292 end if
1293
1294
1295 do ist = 1, st%nst
1296 if (exclude_list(ist)) cycle
1297
1298 projection = m_zero
1299
1300 if(ik_real >= st%d%kpt%start .and. ik_real <= st%d%kpt%end) then
1301 call states_elec_get_state(st, mesh, ist, ik_real, psi)
1302
1303 !$omp parallel
1304 do idim = 1, st%d%dim
1305 !The minus sign is here is for the wrong convention of Octopus
1306 !$omp do
1307 do ip = 1, mesh%np
1308 psi(ip, idim) = psi(ip, idim)*phase(ip)
1309 end do
1310 end do
1311 !$omp end parallel
1312
1313 do iw = 1, w90_nproj
1314 idim = 1
1315 if (w90_spinors) idim = w90_spin_proj_component(iw)
1316
1317 !At the moment the orbitals do not depend on idim
1318 !The idim index for eorb_mesh would be for a spin-resolved orbital like j=1/2
1319 projection(iw) = zmf_dotp(mesh, psi(1:mesh%np,idim), &
1320 orbitals(iw)%eorb_mesh(1:mesh%np,1,1,ik_real), reduce = .false.)
1321 end do
1322
1323 call profiling_in("W90_AMN_REDUCE")
1324 call mesh%allreduce(projection)
1325 call profiling_out("W90_AMN_REDUCE")
1326 end if
1327
1328 if(st%d%kpt%parallel) then
1329 call comm_allreduce(st%d%kpt%mpi_grp, projection)
1330 end if
1331
1332 if (st%system_grp%is_root()) then
1333 do iw = 1, w90_nproj
1334 write (w90_amn,'(I5,2x,I5,2x,I5,2x,e18.10,2x,e18.10)') band_index(ist), iw, ik, projection(iw)
1335 end do
1336 end if
1337 end do !ist
1338 end do !ik
1339
1340 safe_deallocate_a(psi)
1341 safe_deallocate_a(phase)
1342 safe_deallocate_a(projection)
1343
1344 do iw = 1, w90_nproj
1345 call orbitalset_end(orbitals(iw))
1346 end do
1347 safe_deallocate_a(orbitals)
1348 end if
1349
1350 call io_close(w90_amn)
1351
1352 call profiling_out("W90_AMN")
1353
1354 pop_sub(create_wannier90_amn)
1355
1356 end subroutine create_wannier90_amn
1357
1358 ! --------------------------------------------------------------------------
1362 subroutine create_wannier90_spn(mesh, st)
1363 class(mesh_t), intent(in) :: mesh
1364 type(states_elec_t), target, intent(in) :: st
1365
1366 integer :: ist, jst, ik, w90_spn, counter
1367 character(len=80) :: filename
1368 complex(real64), allocatable :: spin(:,:,:)
1369 complex(real64), allocatable :: psim(:,:), psin(:,:)
1370 complex(real64) :: dot_upup, dot_updown, dot_downup, dot_downdown
1371
1372 push_sub(create_wannier90_spn)
1373 call profiling_in("W90_SPN")
1374
1375 assert(st%d%ispin == spinors)
1376
1377 if (st%parallel_in_states) then
1378 call messages_not_implemented("w90_spn output with states parallelization")
1379 end if
1380
1381 message(1) = "Info: Computing the spin file"
1382 call messages_info(1)
1383
1384 filename = './'// trim(adjustl(w90_prefix))//'.spn'
1385 w90_spn = io_open(trim(filename), global_namespace, action='write')
1386
1387 ! write header
1388 if (st%system_grp%is_root()) then
1389 write(w90_spn,*) 'Created by oct-wannier90'
1390 write(w90_spn,*) w90_num_bands, w90_num_kpts
1391 end if
1392
1393 safe_allocate(psim(1:mesh%np, 1:st%d%dim))
1394 safe_allocate(psin(1:mesh%np, 1:st%d%dim))
1395 safe_allocate(spin(1:3, 1:(w90_num_bands*(w90_num_bands+1))/2, 1:w90_num_kpts))
1396 spin = m_zero
1397
1398 ! loop over the pairs specified in the nnkp file (read before in init)
1399 do ik = st%d%kpt%start, st%d%kpt%end
1400 counter = 0
1401 do jst = 1, st%nst
1402 if (exclude_list(jst)) cycle
1403
1404 call states_elec_get_state(st, mesh, jst, ik, psim)
1405 do ist = 1, jst
1406 if (exclude_list(ist)) cycle
1407
1408 counter = counter + 1
1409
1410 call states_elec_get_state(st, mesh, ist, ik, psin)
1411
1412 dot_upup = zmf_dotp(mesh, psin(:, 1), psim(:, 1), reduce = .false.)
1413 dot_downdown = zmf_dotp(mesh, psin(:, 2), psim(:, 2), reduce = .false.)
1414 dot_updown = zmf_dotp(mesh, psin(:, 1), psim(:, 2), reduce = .false.)
1415 dot_downup = zmf_dotp(mesh, psin(:, 2), psim(:, 1), reduce = .false.)
1416
1417 spin(1, counter, ik) = dot_updown + dot_downup
1418 spin(2, counter, ik) = -m_zi * dot_updown + m_zi * dot_downup
1419 spin(3, counter, ik) = dot_upup - dot_downdown
1420 end do !ist
1421 end do
1422 end do
1423
1424 call profiling_in("W90_SPN_REDUCE")
1425 call mesh%allreduce(spin)
1426
1427 if(st%d%kpt%parallel) then
1428 call comm_allreduce(st%d%kpt%mpi_grp, spin)
1429 end if
1430 call profiling_out("W90_SPN_REDUCE")
1431
1432 ! write to W90 file
1433 if (st%system_grp%is_root()) then
1434 do ik = 1, w90_num_kpts
1435 counter = 0
1436 do jst = 1, st%nst
1437 if (exclude_list(jst)) cycle
1438
1439 do ist = 1, jst
1440 if (exclude_list(ist)) cycle
1441
1442 counter = counter + 1
1443 write(w90_spn, '(e18.10,2x,e18.10)') spin(1, counter, ik)
1444 write(w90_spn, '(e18.10,2x,e18.10)') spin(2, counter, ik)
1445 write(w90_spn, '(e18.10,2x,e18.10)') spin(3, counter, ik)
1446 end do
1447 end do
1448 end do !ik
1449 end if
1450
1451 call io_close(w90_spn)
1452
1453 safe_deallocate_a(psim)
1454 safe_deallocate_a(psin)
1455 safe_deallocate_a(spin)
1456
1457 call profiling_out("W90_SPN")
1458
1459 pop_sub(create_wannier90_spn)
1460 end subroutine create_wannier90_spn
1461
1462
1463 ! --------------------------------------------------------------------------
1464 subroutine generate_wannier_states(space, mesh, ions, st, kpoints)
1465 class(space_t), intent(in) :: space
1466 class(mesh_t), intent(in) :: mesh
1467 type(ions_t), intent(in) :: ions
1468 type(states_elec_t), intent(in) :: st
1469 type(kpoints_t), intent(in) :: kpoints
1470
1471 integer :: w90_u_mat, w90_xyz, nwann, nik
1472 integer :: ik, iw, iw2, ip, ipmax, rankmax, idmmax
1473 real(real64), allocatable :: centers(:,:), supercell_centers(:,:), new_centers(:,:)
1474 complex(real64), allocatable :: Umnk(:,:,:)
1475 complex(real64), allocatable :: zwn(:,:,:), psi(:,:), phase(:,:), phase_bloch(:), zwn_bloch(:,:,:)
1476 character(len=MAX_PATH_LEN) :: fname
1477 real(real64) :: kpoint(3), wmod, wmodmax, xx(space%dim)
1478 character(len=2) :: dum
1479 logical :: exist
1480 type(unit_t) :: fn_unit
1481 complex(real64) :: scal
1482 type(block_t) :: blk
1483 integer :: supercell(space%dim), ii, jj, kk, ncols, Nreplica, irep, irepmax
1484 real(real64) :: min_image_displ(3), offset(space%dim)
1485 integer :: jw
1486 integer, allocatable :: parent(:)
1487 real(real64), parameter :: tol_cluster = 0.25_real64 ! In Bohr
1488
1489
1490 push_sub(generate_wannier_states)
1491
1492 message(1) = "oct-wannier90: Constructing the Wannier states from the U matrix."
1493 call messages_info(1)
1494
1495 inquire(file=trim(trim(adjustl(w90_prefix))//'_centres.xyz'),exist=exist)
1496 if (.not. exist) then
1497 message(1) = 'oct-wannier90: Cannot find the Wannier90 file seedname_centres.xyz.'
1498 write(message(2),'(a)') 'Please run wannier90.x with "write_xyz=.true." in '// trim(adjustl(w90_prefix)) // '.'
1499 call messages_fatal(2)
1500 end if
1501
1502 w90_xyz = io_open(trim(trim(adjustl(w90_prefix))//'_centres.xyz'), global_namespace, action='read')
1503
1504 safe_allocate(centers(1:3, 1:w90_num_wann))
1505 !Skip two lines
1506 read(w90_xyz, *)
1507 read(w90_xyz, *)
1508 do iw = 1, w90_num_wann
1509 read(w90_xyz, *) dum, centers(1:3, iw)
1510 ! Wannier90 outputs the coordinates in angstrom
1511 centers(1:3, iw) = units_to_atomic(unit_angstrom, centers(1:3, iw))
1512 end do
1513 call io_close(w90_xyz)
1514
1515 ! In order to find "more consistent" Wannier centers, the centers are clustered in groups,
1516 ! Using periodic boundary condition to find the centers that are the closest
1517 safe_allocate(parent(1:w90_num_wann))
1518 do iw = 1, w90_num_wann
1519 parent(iw) = iw
1520 do jw = 1, iw-1
1521 min_image_displ(:) = ions%latt%cart_to_red(centers(:, iw) - centers(:, jw))
1522 min_image_displ(1:space%periodic_dim) = min_image_displ(1:space%periodic_dim) - nint(min_image_displ(1:space%periodic_dim))
1523 min_image_displ = ions%latt%red_to_cart(min_image_displ)
1524 if (norm2(min_image_displ) < tol_cluster) then
1525 parent(iw) = parent(jw)
1526 exit
1527 end if
1528 end do
1529 end do
1530
1531 message(1) = "Info : Clustering of the Wannier centers"
1532 call messages_info(1)
1533
1534 safe_allocate(new_centers(1:3, 1:w90_num_wann))
1535 do iw = 1, w90_num_wann
1536 write(message(1), '(a,i0,a,3(f7.3,a))') 'Info : Original Wannier center ', &
1537 iw, ' (', centers(1, iw), ', ', centers(2, iw), ', ', centers(3, iw), ')'
1538
1539 if (parent(iw) == iw) then
1540 ! We are computing |w_{n,0}>, so we need to take only a shift to get to the central cell
1541 new_centers(1:3, iw) = ions%latt%fold_into_cell(centers(1:3, iw))
1542 else
1543 ! This center has a closeby periodic parent. We then move it to be as close as possible to its "parent"
1544 ! Thanks to the loop above, this should always be inside of very close to a parent inside the unit cell
1545 min_image_displ(:) = ions%latt%cart_to_red(centers(:, iw) - new_centers(:, parent(iw)))
1546 min_image_displ(1:space%periodic_dim) = real(nint(min_image_displ(1:space%periodic_dim)), real64)
1547 min_image_displ = ions%latt%red_to_cart(min_image_displ)
1548 new_centers(1:3, iw) = centers(1:3, iw) - min_image_displ
1549 end if
1550
1551 write(message(2), '(a,i0,a,3(f7.3,a))') 'Info : New Wannier center ', &
1552 iw, ' (', new_centers(1, iw), ', ', new_centers(2, iw), ', ', new_centers(3, iw), ')'
1553 write(message(3), '(a)') ''
1554 call messages_info(3)
1555 end do
1556
1557 ! Getting the offset, for the correct phase factor
1558 centers = new_centers - centers
1559 safe_deallocate_a(parent)
1560 safe_deallocate_a(new_centers)
1561
1562 inquire(file=trim(trim(adjustl(w90_prefix))//'_u_dis.mat'),exist=exist)
1563 if (exist) then
1564 message(1) = 'oct-wannier90: Calculation of Wannier states with disentanglement is not yet supported.'
1565 call messages_fatal(1)
1566 end if
1567
1568 inquire(file=trim(trim(adjustl(w90_prefix))//'_u.mat'),exist=exist)
1569 if (.not. exist) then
1570 message(1) = 'oct-wannier90: Cannot find the Wannier90 seedname_u.mat file.'
1571 write(message(2),'(a)') 'Please run wannier90.x with "write_u_matrices=.true." in '// trim(adjustl(w90_prefix)) // '.'
1572 call messages_fatal(2)
1573 end if
1574 w90_u_mat = io_open(trim(trim(adjustl(w90_prefix))//'_u.mat'), global_namespace, action='read')
1575
1576
1577 !Skip one line
1578 read(w90_u_mat, *)
1579 !Read num_kpts, num_wann, num_wann for consistency check
1580 read(w90_u_mat, *) nik, nwann, nwann
1581 if (nik /= w90_num_kpts .or. nwann /= w90_num_wann) then
1582 message(1) = "The file contains U matrices is inconsistent with the .win file."
1583 call messages_fatal(1)
1584 end if
1585
1586 safe_allocate(umnk(1:w90_num_wann, 1:w90_num_wann, 1:w90_num_kpts))
1587
1588 do ik = 1, w90_num_kpts
1589 !Skip one line
1590 read(w90_u_mat, *)
1591 !Skip one line (k-point coordinate)
1592 read(w90_u_mat, *)
1593 read(w90_u_mat, '(f15.10,sp,f15.10)') ((umnk(iw, iw2, ik), iw=1, w90_num_wann), iw2=1, w90_num_wann)
1594 end do
1595
1596 call io_close(w90_u_mat)
1597
1598 !We read the output format for the Wannier states
1599 call parse_variable(global_namespace, 'OutputFormat', 0, how)
1600 if (how == 0) then
1601 message(1) = "OutputFormat must be specified for outputing Wannier functions."
1602 call messages_fatal(1)
1603 end if
1604
1605 !%Variable Wannier90Supercell
1606 !%Type block
1607 !%Section Utilities::oct-wannier90
1608 !%Description
1609 !% This block allows to specify the size of the supercell used to compute the Wannier functions
1610 !% If not specified, the code uses a default 1x1x1 cell, i.e., it plots the Wannier90 in the primitive cell.
1611 !%End
1612 if (parse_is_defined(sys%namespace, 'Wannier90Supercell')) then
1613 if (parse_block(sys%namespace, 'Wannier90Supercell', blk) == 0) then
1614 ncols = parse_block_cols(blk, 0)
1615 if (ncols /= space%dim) then
1616 write(message(1),'(a,i3,a,i3)') 'Wannier90Supercell has ', ncols, ' columns but must have ', sys%space%dim
1617 call messages_fatal(1, namespace=sys%namespace)
1618 end if
1619 do ii = 1, space%dim
1620 call parse_block_integer(blk, 0, ii - 1, supercell(ii))
1621 end do
1622
1623 call parse_block_end(blk)
1624 end if
1625 else
1626 supercell = 1
1627 end if
1628
1629 nreplica = product(supercell)
1630
1631 ! The center of each replica of the unit cell
1632 safe_allocate(supercell_centers(1:space%dim, 1:nreplica))
1633 offset(1:space%dim) = -floor((real(supercell(1:space%dim), real64) - m_one) / m_two)
1634 irep = 1
1635 do ii = 0, supercell(1)-1
1636 do jj = 0, supercell(2)-1
1637 do kk = 0, supercell(3)-1
1638 supercell_centers(:, irep) = ions%latt%red_to_cart(offset &
1639 + [real(ii, real64), real(jj, real64), real(kk, real64)])
1640 irep = irep + 1
1641 end do
1642 end do
1643 end do
1644
1645
1646
1647 call io_mkdir('wannier', global_namespace)
1648
1649 !Computing the Wannier states in the primitive cell, from the U matrices
1650 safe_allocate(zwn(1:mesh%np, 1:nreplica, 1:st%d%dim))
1651 safe_allocate(psi(1:mesh%np, 1:st%d%dim))
1652 safe_allocate(phase(1:mesh%np, 1:nreplica))
1653 safe_allocate(phase_bloch(1:mesh%np))
1654 if (w90_bloch_sums) then
1655 safe_allocate(zwn_bloch(1:mesh%np, 1:st%d%dim, st%d%kpt%start:st%d%kpt%end))
1656 end if
1657
1658 do iw = 1, w90_num_wann
1659
1660 zwn(:, :, :) = m_z0
1661
1662 do ik = 1, w90_num_kpts
1663
1664 if (.not. (ik >= st%d%kpt%start .and. ik <= st%d%kpt%end)) cycle
1665
1666 if (w90_bloch_sums) zwn_bloch(:,:,ik) = m_z0
1667
1668 kpoint(1:space%dim) = kpoints%get_point(ik, absolute_coordinates=.true.)
1669
1670 ! We compute the Wannier orbital on a grid centered around the Wannier function
1671 ! The minus sign is here is for the wrong convention of Octopus
1672 do irep = 1, nreplica
1673 !$omp parallel do
1674 do ip = 1, mesh%np
1675 xx = mesh%x(1:space%dim, ip)-centers(1:space%dim, iw) + supercell_centers(1:space%dim, irep)
1676 phase(ip, irep) = exp(-m_zi* sum( xx * kpoint(1:space%dim)))
1677 end do
1678 end do
1679 if (w90_bloch_sums) then
1680 do ip = 1, mesh%np
1681 xx = mesh%x(1:space%dim, ip)
1682 phase_bloch(ip) = exp(-m_zi* sum( xx * kpoint(1:space%dim)))
1683 end do
1684 end if
1685
1686 do iw2 = 1, st%nst
1687 if (exclude_list(iw2)) cycle
1688
1689 if (st%d%ispin /= spin_polarized) then
1690 call states_elec_get_state(st, mesh, iw2, ik, psi)
1691 else
1692 call states_elec_get_state(st, mesh, iw2, (ik-1)*2+w90_spin_channel, psi)
1693 end if
1694
1695 !$omp parallel
1696 do idim = 1, st%d%dim
1697 do irep = 1, nreplica
1698 !$omp do
1699 do ip = 1, mesh%np
1700 zwn(ip, irep, idim) = zwn(ip, irep, idim) + umnk(band_index(iw2), iw, ik) * psi(ip, idim) * phase(ip, irep)
1701 end do
1702 end do
1703 end do
1704 !$omp end parallel
1705
1706 ! See Eq. 19 in arXiv:0605539
1707 if (w90_bloch_sums) then
1708 !$omp parallel
1709 do idim = 1, st%d%dim
1710 !$omp do
1711 do ip = 1, mesh%np
1712 zwn_bloch(ip, idim, ik) = zwn_bloch(ip, idim, ik) &
1713 + umnk(band_index(iw2), iw, ik) * psi(ip, idim) * phase_bloch(ip)
1714 end do
1715 end do
1716 !$omp end parallel
1717 end if
1718 end do!iw2
1719 end do!ik
1720
1721 if(st%d%kpt%parallel) then
1722 call comm_allreduce(st%d%kpt%mpi_grp, zwn)
1723 end if
1724
1725 ! Following what Wannier90 is doing, we fix the global phase by setting the max to be real
1726 ! We also normalize to the number of k-point at this step
1727 if (sys%st%d%ispin /= spinors) then
1728 ipmax = 0
1729 irepmax = 0
1730 wmodmax = m_zero
1731 do irep = 1, nreplica
1732 do ip = 1, mesh%np
1733 wmod = real(zwn(ip, irep, 1)*conjg(zwn(ip, irep, 1)), real64)
1734 if (wmod > wmodmax) then
1735 ipmax = ip
1736 irepmax = irep
1737 wmodmax = wmod
1738 end if
1739 end do
1740 end do
1741 scal = sqrt(wmodmax)/zwn(ipmax, irepmax, 1)/w90_num_kpts
1742 call mesh_minmaxloc(mesh, wmodmax, rankmax, mpi_maxloc)
1743 call mesh%mpi_grp%bcast(scal, 1, mpi_double_complex, rankmax)
1744 call lalg_scal(mesh%np, nreplica, scal, zwn(:,:,1))
1745
1746 if (w90_bloch_sums) then
1747 do ik = st%d%kpt%start, st%d%kpt%end
1748 ipmax = 0
1749 idmmax= 0
1750 wmodmax = m_zero
1751 do idim = 1, st%d%dim
1752 do ip = 1, mesh%np
1753 wmod = real(zwn_bloch(ip, idim, ik)*conjg(zwn_bloch(ip, idim, ik)), real64)
1754 if (wmod > wmodmax) then
1755 ipmax = ip
1756 wmodmax = wmod
1757 idmmax = idim
1758 end if
1759 end do
1760 end do
1761 scal = sqrt(wmodmax)/zwn_bloch(ipmax, idmmax, ik)/w90_num_kpts
1762 call mesh_minmaxloc(mesh, wmodmax, rankmax, mpi_maxloc)
1763 call mesh%mpi_grp%bcast(scal, 1, mpi_double_complex, rankmax)
1764 call lalg_scal(mesh%np, st%d%dim, scal, zwn_bloch(:,:,ik))
1765 end do
1766 end if
1767 end if
1768
1769 ! Output the Wannier states
1770 fn_unit = sqrt(units_out%length**(-space%dim))
1771 do idim = 1, st%d%dim
1772 if (st%d%ispin == spinors) then
1773 write(fname, '(a,i3.3,a4,i1)') 'wannier-', iw, '-isp', idim
1774 else
1775 write(fname, '(a,i3.3,a4,i1)') 'wannier-', iw
1776 end if
1777 if (any(supercell_centers > 1)) then
1778 call io_function_output_supercell(how, "wannier", trim(fname), mesh, &
1779 space, ions%latt, zwn(:, :, idim), supercell_centers, supercell, fn_unit, ierr, global_namespace, &
1780 pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
1781
1782 else
1783 call zio_function_output(how, "wannier", trim(fname), global_namespace, space, mesh, &
1784 zwn(:, 1, idim), fn_unit, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
1785 end if
1786 end do
1787
1788 ! Output Bloch sums of Wannier states
1789 if (w90_bloch_sums) then
1790 do ik = st%d%kpt%start, st%d%kpt%end
1791 do idim = 1, st%d%dim
1792 if (st%d%ispin == spinors) then
1793 write(fname, '(a,i3.3,a4,i1,a3,i5.5)') 'wannier_bloch-', iw, '-isp', idim, '-ik', ik
1794 else
1795 write(fname, '(a,i3.3,a4,i1,a3,i5.5)') 'wannier_bloch-', iw, '-ik', ik
1796 end if
1797 call zio_function_output(how, "wannier", trim(fname), global_namespace, space, mesh, &
1798 zwn_bloch(:, idim, ik), fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
1799 end do
1800 end do
1801 end if
1802
1803 ! Checking the ratio imag/real
1804 if (sys%st%d%ispin /= spinors) then
1805 wmodmax = m_zero
1806 do irep = 1, nreplica
1807 do ip = 1, mesh%np
1808 if(abs(real(zwn(ip, irep, 1), real64)) >= 1e-2_real64) then
1809 wmodmax = max(wmodmax, abs(aimag(zwn(ip, irep, 1)))/abs(real(zwn(ip, irep, 1), real64)))
1810 end if
1811 end do
1812 end do
1813 call mesh_minmaxloc(mesh, wmodmax, rankmax, mpi_maxloc)
1814
1815 write(message(1), '(a,i4,a,f11.6)') 'oct-wannier90: Wannier function ', iw, ' Max. Im/Re Ratio = ', wmodmax
1816 call messages_info(1)
1817 else
1818 write(message(1), '(a,i4)') 'oct-wannier90: Wannier function done ', iw
1819 call messages_info(1)
1820 end if
1821 end do
1822
1823 safe_deallocate_a(umnk)
1824 safe_deallocate_a(zwn)
1825 safe_deallocate_a(psi)
1826 safe_deallocate_a(phase)
1827 safe_deallocate_a(phase_bloch)
1828 safe_deallocate_a(centers)
1829 safe_deallocate_a(supercell_centers)
1830
1832 end subroutine generate_wannier_states
1833
1834end program wannier90_interface
1835
1836!! Local Variables:
1837!! mode: f90
1838!! coding: utf-8
1839!! End:
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double floor(double __x) __attribute__((__nothrow__
This module implements batches of mesh functions.
Definition: batch.F90:135
integer, parameter, public batch_not_packed
functions are stored in CPU memory, unpacked order
Definition: batch.F90:286
integer, parameter, public batch_device_packed
functions are stored in device memory in packed order
Definition: batch.F90:286
integer, parameter, public batch_packed
functions are stored in CPU memory, in transposed (packed) order
Definition: batch.F90:286
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 zcube_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:204
subroutine, public cube_end(cube)
Definition: cube.F90:387
subroutine, public cube_init_cube_map(cube, mesh)
Definition: cube.F90:824
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:120
subroutine, public fft_all_init(namespace)
initialize the table
Definition: fft.F90:269
subroutine, public fft_all_end()
delete all plans
Definition: fft.F90:380
real(real64), parameter, public m_two
Definition: global.F90:193
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:458
real(real64), parameter, public m_huge
Definition: global.F90:209
real(real64), parameter, public m_zero
Definition: global.F90:191
complex(real64), parameter, public m_z0
Definition: global.F90:201
complex(real64), parameter, public m_zi
Definition: global.F90:205
subroutine, public global_init(communicator)
Initialise Octopus.
Definition: global.F90:365
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:119
Definition: io.F90:116
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:165
subroutine, public io_close(iunit, grp)
Definition: io.F90:467
subroutine, public io_end()
Definition: io.F90:271
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:402
subroutine, public kpoints_to_absolute(latt, kin, kout)
Definition: kpoints.F90:1137
System information (time, memory, sysname)
Definition: loct.F90:117
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_end()
Definition: messages.F90:272
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
subroutine, public messages_init(output_dir)
Definition: messages.F90:219
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:524
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:409
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1039
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:593
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
type(namespace_t), public global_namespace
Definition: namespace.F90:135
subroutine, public orbitalset_init(this)
Definition: orbitalset.F90:210
subroutine, public orbitalset_end(this)
Definition: orbitalset.F90:236
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:286
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:455
subroutine, public parser_init()
Initialise the Octopus parser.
Definition: parser.F90:402
subroutine, public parser_end()
End the Octopus parser.
Definition: parser.F90:434
subroutine, public profiling_end(namespace)
Definition: profiling.F90:415
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
Definition: profiling.F90:257
integer, parameter, public restart_gs
Definition: restart.F90:156
integer, parameter, public restart_td
Definition: restart.F90:156
integer, parameter, public restart_type_load
Definition: restart.F90:183
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, 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:227
type(type_t), public type_cmplx
Definition: types.F90:136
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
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:120
subroutine, public ylm_wannier(ylm, l, mr, rr, xx, nr)
Class describing the electron system.
Definition: electrons.F90:220
Describes mesh distribution to nodes.
Definition: mesh.F90:187
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Definition: wfs_elec.F90:141
int true(void)
subroutine create_wannier90_spn(mesh, st)
Write the spn file containing .
subroutine create_wannier90_eig()
subroutine read_wannier90_files()
subroutine write_unk(space, mesh, st, formatted)
subroutine create_wannier90_mmn(mesh, st)
subroutine wannier90_setup(ions, kpoints, space)
subroutine generate_wannier_states(space, mesh, ions, st, kpoints)
program wannier90_interface
subroutine create_wannier90_amn(space, mesh, latt, st, kpoints)
subroutine wannier90_output()