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
45 use mesh_oct_m
48 use mpi_oct_m
52 use parser_oct_m
57 use space_oct_m
58 use string_oct_m
63 use types_oct_m
64 use unit_oct_m
66 use utils_oct_m
69
70 implicit none
71
72 integer :: w90_what, w90_mode
73 integer(int64) :: w90_what_default
74
75 integer :: ierr
76 integer :: dim, idim
77 integer :: ii, nik, iter, nst
78
79 type(restart_t) :: restart
80 type(electrons_t), pointer :: sys
81 logical :: w90_spinors, scdm_proj, w90_scdm
82 integer :: w90_nntot, w90_num_bands, w90_num_kpts ! w90 input parameters
83 integer, allocatable :: w90_nnk_list(:,:) !
84 character(len=80) :: w90_prefix ! w90 input file prefix
85 integer :: w90_num_wann ! input paramter
86 real(real64), allocatable :: w90_proj_centers(:,:) ! projections centers
87 integer, allocatable :: w90_proj_lmr(:,:) ! definitions for real valued Y_lm*R_r
88 integer :: w90_nproj ! number of such projections
89 integer, allocatable :: w90_spin_proj_component(:) ! up/down flag
90 real(real64), allocatable :: w90_spin_proj_axis(:,:) ! spin axis (not implemented)
91 integer :: w90_num_exclude
92 logical, allocatable :: exclude_list(:) ! list of excluded bands
93 integer, allocatable :: band_index(:) ! band index after exclusion
94 logical :: read_td_states
95 integer :: w90_spin_channel
96 logical :: w90_bloch_sums = .false.
97
98 ! scdm variables
99 integer, allocatable :: jpvt(:)
100 complex(real64), allocatable :: uk(:,:,:) ! SCDM-Wannier gauge matrices U(k)
101 complex(real64), allocatable :: psi(:,:)
102 complex(real64), allocatable :: chi(:,:), chi_diag(:,:),chi2(:,:)
103 real(real64), allocatable :: chi_eigenval(:), occ_temp(:)
104 real(real64) :: scdm_mu, scdm_sigma, smear, kvec(3), factor(3)
105 integer :: ist, jst, ik, idir, sender
106
107 integer(int64) :: how
108
109 call global_init()
110 call parser_init()
111
112 call messages_init()
113 call io_init()
114
119
120 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
121 sys => electrons_t(global_namespace, mpi_world, int(option__calculationmode__dummy, int32))
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(occ_temp)
545 end if
546
547 ! ---- actual interface work ----------
548 if (bitand(w90_what, option__wannier90files__w90_mmn) /= 0) then
549 call create_wannier90_mmn(sys%gr, sys%st)
550 end if
551
552 if (bitand(w90_what, option__wannier90files__w90_unk) /= 0) then
553 call write_unk(sys%space, sys%gr, sys%st, formatted=.false.)
554 end if
555
556 if (bitand(w90_what, option__wannier90files__w90_amn) /= 0) then
557 call create_wannier90_amn(sys%space, sys%gr, sys%ions%latt, sys%st, sys%kpoints)
558 end if
559
560 if (bitand(w90_what, option__wannier90files__w90_eig) /= 0) then
562 end if
563
564 if (bitand(w90_what, option__wannier90files__w90_spn) /= 0) then
565 call create_wannier90_spn(sys%gr, sys%st)
566 end if
567
568 safe_deallocate_a(uk)
569 safe_deallocate_a(w90_spin_proj_component)
570 safe_deallocate_a(w90_spin_proj_axis)
571
572 pop_sub(wannier90_output)
573 end subroutine wannier90_output
574
575 ! --------------------------------------------------------------------------
576 subroutine read_wannier90_files()
577 integer :: w90_nnkp, itemp, dummyint, io, spin_channel_win
578 character(len=80) :: filename, dummy, dummy1, dummy2, line
579 logical :: exist, parse_is_ok
580 real(real64) :: dummyr(7)
581
582 push_sub(read_wannier90_files)
583
584 w90_num_kpts = product(sys%kpoints%nik_axis(1:3))
585 assert(w90_num_kpts == sys%st%nik)
586
587
588 w90_num_exclude = 0
589
590 ! open nnkp file
591 filename = trim(adjustl(w90_prefix)) //'.nnkp'
592
593 message(1) = "oct-wannier90: Parsing "//filename
594 call messages_info(1)
595
596 inquire(file=filename,exist=exist)
597 if (.not. exist) then
598 message(1) = 'oct-wannier90: Cannot find specified Wannier90 nnkp file.'
599 write(message(2),'(a)') 'Please run wannier90.x -pp '// trim(adjustl(w90_prefix)) // ' first.'
600 call messages_fatal(2)
601 end if
602
603 parse_is_ok = .false.
604
605 ! check number of k-points
606 w90_nnkp = io_open(trim(filename), global_namespace, action='read')
607 do
608 read(w90_nnkp, *, iostat=io) dummy, dummy1
609 if (io == iostat_end) exit
610
611 if (dummy == 'begin' .and. dummy1 == 'kpoints') then
612 read(w90_nnkp,*) itemp
613 if (itemp /= w90_num_kpts) then
614 message(1) = 'oct-wannier90: wannier90 setup seems to have been done with a different number of k-points.'
615 call messages_fatal(1)
616 else
617 parse_is_ok = .true.
618 exit
619 end if
620 end if
621 end do
622 call io_close(w90_nnkp)
623
624 if (.not. parse_is_ok) then
625 message(1) = 'oct-wannier90: Did not find the kpoints block in nnkp file'
626 call messages_fatal(1)
627 end if
628 parse_is_ok = .false.
629
630 ! read from nnkp file
631 ! find the nnkpts block
632 w90_nnkp = io_open(trim(filename), global_namespace, action='read', position='rewind')
633 do
634 read(w90_nnkp, *, iostat=io) dummy, dummy1
635 if (io == iostat_end) exit !End of file
636
637 if (dummy == 'begin' .and. dummy1 == 'nnkpts') then
638 read(w90_nnkp,*) w90_nntot
639 safe_allocate(w90_nnk_list(1:5, 1:w90_num_kpts * w90_nntot))
640 do ii = 1, w90_num_kpts * w90_nntot
641 read(w90_nnkp,*) w90_nnk_list(1:5, ii)
642 end do
643 !make sure we are at the end of the block
644 read(w90_nnkp,*) dummy
645 if (dummy /= 'end') then
646 message(1) = 'oct-wannier90: There dont seem to be enough k-points in nnkpts file to.'
647 call messages_fatal(1)
648 end if
649 parse_is_ok = .true.
650 exit
651 end if
652 end do
653
654 if (.not. parse_is_ok) then
655 message(1) = 'oct-wannier90: Did not find nnkpts block in nnkp file'
656 call messages_fatal(1)
657 end if
658
659 ! read from nnkp file
660 ! find the exclude block
661 safe_allocate(exclude_list(1:sys%st%nst))
662 !By default we use all the bands
663 exclude_list(1:sys%st%nst) = .false.
664 rewind(w90_nnkp)
665 do
666 read(w90_nnkp, *, iostat=io) dummy, dummy1
667 if (io == iostat_end) exit !End of file
668 if (dummy == 'begin' .and. dummy1 == 'exclude_bands') then
669 read(w90_nnkp, *) w90_num_exclude
670 do ii = 1, w90_num_exclude
671 read(w90_nnkp, *) itemp
672 if(itemp > sys%st%nst) then
673 message(1) = 'oct-wannier90: The exclude_bands list contains a state index higher than the number of states.'
674 call messages_fatal(1)
675 end if
676 exclude_list(itemp) = .true.
677 end do
678 !make sure we are at the end of the block
679 read(w90_nnkp, *) dummy
680 if (dummy /= 'end') then
681 message(1) = 'oct-wannier90: There dont seem to be enough bands in exclude_bands list.'
682 call messages_fatal(1)
683 end if
684 exit
685 end if
686 end do
687 call io_close(w90_nnkp)
688
689 !We get the number of bands
690 w90_num_bands = sys%st%nst - w90_num_exclude
691
692 safe_allocate(band_index(1:sys%st%nst))
693 itemp = 0
694 do ii = 1, sys%st%nst
695 if (exclude_list(ii)) cycle
696 itemp = itemp + 1
697 band_index(ii) = itemp
698 end do
699
700 if (bitand(w90_what, option__wannier90files__w90_amn) /= 0 &
701 .or. w90_mode == option__wannier90mode__w90_wannier ) then
702 ! parse file again for definitions of projections
703 w90_nnkp = io_open(trim(filename), global_namespace, action='read', position='rewind')
704
705 do
706 read(w90_nnkp, *, iostat=io) dummy, dummy1
707 if (io == iostat_end) then !End of file
708 message(1) = 'oct-wannier90: Did not find projections block in w90.nnkp file'
709 call messages_fatal(1)
710 end if
711
712 if (dummy == 'begin' .and. (dummy1 == 'projections' .or. dummy1 == 'spinor_projections')) then
713
714 if (dummy1 == 'spinor_projections') then
715 w90_spinors = .true.
716 if (sys%st%d%ispin /= spinors) then
717 message(1) = 'oct-wannier90: Spinor = .true. is only valid with spinors wavefunctions.'
718 call messages_fatal(1)
719 end if
720
721 message(1) = 'oct-wannier90: Spinor interface incomplete. Note there is no quantization axis implemented'
722 call messages_warning(1)
723 else
724 if (sys%st%d%ispin == spinors) then
725 message(1) = 'oct-wannier90: Octopus has spinors wavefunctions but spinor_projections is not defined.'
726 message(2) = 'oct-wannier90: Please check the input file for wannier 90.'
727 call messages_fatal(2)
728 end if
729 end if
730
731 read(w90_nnkp, *) w90_nproj
732 ! num_wann is given in w90.win, not double checked here
733 ! I assume that the wannier90.x -pp run has checked this
734 w90_num_wann = w90_nproj
735 ! In case of no projections, we use the number of bands
736 if(w90_nproj == 0) w90_num_wann = w90_num_bands
737
738 safe_allocate(w90_proj_centers(1:3, 1:w90_nproj))
739 safe_allocate(w90_proj_lmr(1:w90_nproj, 1:3))
740 if (w90_spinors) then
741 safe_allocate(w90_spin_proj_component(1:w90_nproj))
742 end if
743 if (w90_spinors) then
744 safe_allocate(w90_spin_proj_axis(1:w90_nproj, 1:3))
745 end if
746
747 do ii = 1, w90_nproj
748 read(w90_nnkp, *) w90_proj_centers(1:3, ii), w90_proj_lmr(ii, 1:3)
749 ! skip a line for now
750 read(w90_nnkp, *) dummyr(1:7)
751 if (w90_spinors) then
752 read(w90_nnkp, *) w90_spin_proj_component(ii), w90_spin_proj_axis(ii, 1:3)
753 ! use octopus spindim conventions
754 if (w90_spin_proj_component(ii) == -1) w90_spin_proj_component(ii) = 2
755 end if
756 end do
757 !make sure we are at the end of the block
758 read(w90_nnkp, *) dummy
759 if (dummy /= 'end') then
760 message(1) = 'oct-wannier90: There dont seem to be enough projections in nnkpts file to.'
761 call messages_fatal(1)
762 end if
763 exit
764 end if
765 end do
766
767 ! look for auto_projection block
768 scdm_proj = .false.
769 do
770 read(w90_nnkp, *, iostat=io) dummy, dummy1
771 if (io == iostat_end) exit !End of file
772
773 if (dummy == 'begin' .and. dummy1 == 'auto_projections') then
774 scdm_proj = .true.
775 read(w90_nnkp, *) w90_nproj
776 w90_num_wann = w90_nproj
777
778 if (.not. w90_scdm) then
779 message(1) = 'oct-wannier90: Found auto_projections block. Currently the only implemented automatic way'
780 message(2) = 'oct-wannier90: to compute projections is the SCDM method.'
781 message(3) = 'oct-wannier90: Please set Wannier90UseSCDM = yes in the inp file.'
782 call messages_fatal(3)
783 end if
784
785 if (w90_nproj /= w90_num_bands) then
786 message(1) = 'oct-wannier90: In auto_projections block first row needs to be equal to num_bands.'
787 call messages_fatal(1)
788 end if
789 read(w90_nnkp, *) dummyint
790 if (dummyint /= 0) then
791 message(1) = 'oct-wannier90: The second row in auto_projections has to be 0, per Wannier90 documentation.'
792 call messages_fatal(1)
793 end if
794 end if
795 end do
796 call io_close(w90_nnkp)
797
798 end if
799
800 message(1) = "oct-wannier90: Finished parsing "//filename
801 call messages_info(1)
802
803 ! Look extra variables variable
804 ! open win file
805 filename = trim(adjustl(w90_prefix)) //'.win'
806 message(1) = "oct-wannier90: Parsing "//filename
807 call messages_info(1)
808 w90_nnkp = io_open(trim(filename), global_namespace, action='read', position='rewind')
809 do
810 read(w90_nnkp, fmt='(a)', iostat=io) line
811 if (io == iostat_end) exit !End of file
812 if (index(line, '=') > 0) then
813 read(line, *, iostat=io) dummy, dummy2, dummy1
814 else
815 read(line, *, iostat=io) dummy, dummy1
816 end if
817
818 !Spin
819 if (dummy == 'spin') then
820 if (sys%st%d%ispin /= spin_polarized) then
821 message(1) = 'oct-wannier90: The variable spin is set for a non spin-polarized calculation.'
822 call messages_fatal(1)
823 end if
824
825 if (dummy1 == 'up') then
826 spin_channel_win = 1
827 else if (dummy1 == 'down') then
828 spin_channel_win = 2
829 else
830 message(1) = 'oct-wannier90: Error parsing the variable spin.'
831 call messages_fatal(1)
832 end if
833 if (spin_channel_win /= w90_spin_channel) then
834 message(1) = 'spin polarization input from .win and Wannier90SpinChannel do not agree.'
835 call messages_fatal(1)
836 end if
837 end if
838 end do
839 call io_close(w90_nnkp)
840
841 if (sys%st%d%ispin == spin_polarized) then
842 write(message(1), '(a,i1)') 'oct-wannier90: Using spin channel ', w90_spin_channel
843 call messages_info(1)
844 end if
845
846 message(1) = "oct-wannier90: Finished parsing "//filename
847 call messages_info(1)
848
849 pop_sub(read_wannier90_files)
850
851 end subroutine read_wannier90_files
852
853 ! --------------------------------------------------------------------------
854 subroutine create_wannier90_mmn(mesh, st)
855 class(mesh_t), intent(in) :: mesh
856 type(states_elec_t), target, intent(in) :: st
857
858 integer :: ist, jst, ik, ip, w90_mmn, iknn, ib
859 real(real64) :: Gcart(3)
860 integer :: G(3)
861 character(len=80) :: filename
862 complex(real64), allocatable :: overlap(:), ss(:)
863 complex(real64), allocatable :: psin(:,:), phase(:)
864 type(wfs_elec_t), pointer :: batch
865 integer :: inode, node_fr, node_to
866 type(mpi_request) :: send_req
867
868 push_sub(create_wannier90_mmn)
869
870 call profiling_in("W90_MMN")
871
872 if (st%parallel_in_states) then
873 call messages_not_implemented("w90_mmn output with states parallelization")
874 end if
875
876 message(1) = "Info: Computing the overlap matrix"
877 call messages_info(1)
878
879
880 filename = './'// trim(adjustl(w90_prefix))//'.mmn'
881 w90_mmn = io_open(trim(filename), global_namespace, action='write')
882
883 ! write header
884 if (st%system_grp%is_root()) then
885 write(w90_mmn,*) 'Created by oct-wannier90'
886 write(w90_mmn,*) w90_num_bands, w90_num_kpts, w90_nntot
887 end if
888
889 safe_allocate(psin(1:mesh%np, 1:st%d%dim))
890 safe_allocate(phase(1:mesh%np))
891 safe_allocate(overlap(1:w90_num_bands))
892 safe_allocate(ss(1:st%nst))
893
894
895 ! loop over the pairs specified in the nnkp file (read before in init)
896 do ii = 1, w90_num_kpts * w90_nntot
897 ik = w90_nnk_list(1, ii)
898 iknn = w90_nnk_list(2, ii)
899 g(1:3) = w90_nnk_list(3:5, ii)
900 if (st%system_grp%is_root()) write(w90_mmn, '(I10,2x,I10,2x,I3,2x,I3,2x,I3)') ik, iknn, g
901 ! For mixed periodicity we remove the G vectors along the non-periodic directions
902 ! as this correspond to the infinite vacuum limit
903 g(sys%space%periodic_dim+1:sys%space%dim) = 0
904
905 ! For spin-polarized calculations, we select the right k-point
906 if (sys%st%d%ispin == spin_polarized) then
907 ik = (ik-1)*2 + w90_spin_channel
908 iknn = (iknn-1)*2 + w90_spin_channel
909 end if
910
911
912 ! Only treat local k-points
913 if(ik >= st%d%kpt%start .and. ik <= st%d%kpt%end) then
914
915 ! Wannier90 treats everything fully periodic
916 ! Conversion is done with the 3D "correct" klattice
917 call kpoints_to_absolute(sys%ions%latt, real(g, real64) , gcart)
918
919 ! Phase that gives u_{n,k+G}(r) from u_{nk}(r)
920 ! Here the minus sign of Octopus cancels with minus sign of the input G
921 ! (ik and iknn correspond in the list to minus the k-points in Octopus)
922 if (any(g /= 0)) then
923 !$omp parallel do
924 do ip = 1, mesh%np
925 phase(ip) = exp(-m_zi*dot_product(mesh%x(1:3, ip), gcart(1:3)))
926 end do
927 end if
928
929 end if
930
931 ! Loop over distributed bands
932 ! Slow index in M_{mn}^{k,b}, so \psi_n
933 do jst = 1, st%nst
934 if (exclude_list(jst)) cycle
935
936 ! Communication for the local states
937 if ( .not. st%d%kpt%parallel .and. .not. st%parallel_in_states) then
938 call states_elec_get_state(st, mesh, jst, iknn, psin)
939 else
940 node_fr = -1
941 node_to = -1
942 do inode = 0, st%d%kpt%mpi_grp%size-1
943 if(iknn >= st%st_kpt_task(inode,3) .and. iknn <= st%st_kpt_task(inode,4)) then
944 node_fr = inode
945 end if
946 if(ik >= st%st_kpt_task(inode,3) .and. ik <= st%st_kpt_task(inode,4)) then
947 node_to = inode
948 end if
949 end do
950 assert(node_fr > -1)
951 assert(node_to > -1)
952
953 send_req = mpi_request_null
954 ! We have locally the k-point
955 if (state_kpt_is_local(st, jst, iknn)) then
956 call states_elec_get_state(st, mesh, jst, iknn, psin)
957 ! We send it only if we don`t want to use it locally
958 if(node_to /= st%d%kpt%mpi_grp%rank) then
959 call st%d%kpt%mpi_grp%isend(psin, mesh%np*st%d%dim, mpi_double_complex, node_to, send_req)
960 end if
961 end if
962 ! We receive the desired state, only if it is not a local one
963 if(node_to == st%d%kpt%mpi_grp%rank .and. node_to /= node_fr) then
964 call st%d%kpt%mpi_grp%recv(psin, mesh%np*st%d%dim, mpi_double_complex, node_fr)
965 end if
966 if (send_req /= mpi_request_null) then
967 call st%d%kpt%mpi_grp%wait(send_req)
968 end if
969 end if
970
971 overlap = m_zero
972
973 if(ik >= st%d%kpt%start .and. ik <= st%d%kpt%end) then
974
975 ! Do not apply the phase if the phase factor is null
976 if (any(g /= 0)) then
977 ! add phase
978 !$omp parallel
979 do idim = 1, st%d%dim
980 !$omp do
981 do ip = 1, mesh%np
982 psin(ip, idim) = psin(ip, idim) * phase(ip)
983 end do
984 end do
985 !$omp end parallel
986 end if
987
988
989 ! See Eq. (25) in PRB 56, 12847 (1997)
990 ! Fast index in M_{mn}^{k,b}, so m.
991 do ib = st%group%block_start, st%group%block_end
992 batch => st%group%psib(ib, ik)
993 if (all(exclude_list(batch%ist(1:batch%nst)))) cycle
994 call zmesh_batch_mf_dotp(mesh, batch, psin, ss(1:batch%nst), reduce = .false.)
995 do ist = 1, batch%nst
996 if (exclude_list(batch%ist(ist))) cycle
997 overlap(band_index(batch%ist(ist))) = ss(ist)
998 end do
999 end do
1000 end if
1001
1002 call profiling_in("W90_MMN_REDUCE")
1003 call mesh%allreduce(overlap)
1004 call profiling_out("W90_MMN_REDUCE")
1005
1006 if(st%d%kpt%parallel) then
1007 call comm_allreduce(st%d%kpt%mpi_grp, overlap)
1008 end if
1009
1010 ! write to W90 file
1011 if (st%system_grp%is_root()) then
1012 do ist = 1, st%nst
1013 if (exclude_list(ist)) cycle
1014 write(w90_mmn,'(e18.10,2x,e18.10)') overlap(band_index(ist))
1015 end do
1016 end if
1017
1018 end do !jst
1019 end do
1020
1021 call io_close(w90_mmn)
1022
1023 safe_deallocate_a(psin)
1024 safe_deallocate_a(phase)
1025 safe_deallocate_a(overlap)
1026 safe_deallocate_a(ss)
1027
1028 call profiling_out("W90_MMN")
1029
1030 pop_sub(create_wannier90_mmn)
1031
1032 end subroutine create_wannier90_mmn
1033
1034 ! --------------------------------------------------------------------------
1035 subroutine create_wannier90_eig()
1036 integer :: ist, ik, w90_eig
1037 character(len=80) :: filename
1038
1039 push_sub(create_wannier90_eig)
1040
1041 if (sys%st%parallel_in_states) then
1042 call messages_not_implemented("w90_eig output with states parallelization")
1043 end if
1044
1045 if (sys%st%system_grp%is_root()) then
1046 filename = './'//trim(adjustl(w90_prefix))//'.eig'
1047 w90_eig = io_open(trim(filename), global_namespace, action='write')
1048 do ik = 1, w90_num_kpts
1049 do ist = 1, sys%st%nst
1050 if (exclude_list(ist)) cycle
1051 if (sys%st%d%ispin /= spin_polarized) then
1052 write(w90_eig,'(I5,2x,I8,2x,e18.10)') band_index(ist), ik, &
1053 units_from_atomic(unit_ev, sys%st%eigenval(ist, ik))
1054 else
1055 write(w90_eig,'(I5,2x,I8,2x,e18.10)') band_index(ist), ik, &
1056 units_from_atomic(unit_ev, sys%st%eigenval(ist, (ik-1)*2+w90_spin_channel))
1057 end if
1058 end do
1059 end do
1060
1061 call io_close(w90_eig)
1062 end if
1063
1064 pop_sub(create_wannier90_eig)
1065 end subroutine create_wannier90_eig
1066
1067 ! --------------------------------------------------------------------------
1068 subroutine write_unk(space, mesh, st, formatted)
1069 class(space_t), intent(in) :: space
1070 class(mesh_t), intent(in) :: mesh
1071 type(states_elec_t), intent(in) :: st
1072 logical, intent(in) :: formatted
1073
1074 integer :: ist, ik, unk_file, ispin
1075 integer :: ix, iy, iz
1076 real(real64) :: w_real, w_imag
1077 character(len=80) :: filename
1078 complex(real64), allocatable :: psi(:)
1079 type(cube_t) :: cube
1080 type(cube_function_t) :: cf
1081
1082 push_sub(write_unk)
1083
1084 if (st%d%kpt%parallel) then
1085 call messages_not_implemented("w90_unk output with k-point parallelization")
1086 end if
1087
1088 if (sys%gr%parallel_in_domains) then
1089 call messages_not_implemented("w90_unk output with domain parallelization")
1090 end if
1091
1092 if (st%parallel_in_states) then
1093 call messages_not_implemented("w90_unk output with states parallelization")
1094 end if
1095
1096 call messages_experimental("Wannier90Files = w90_unk")
1097
1098
1099 safe_allocate(psi(1:mesh%np))
1100
1101 call cube_init(cube, mesh%idx%ll, global_namespace, space, mesh%spacing, &
1102 mesh%coord_system, need_partition=.not.mesh%parallel_in_domains)
1103 call cube_init_cube_map(cube, mesh)
1104
1105 call zcube_function_alloc_rs(cube, cf)
1106
1107 do ik = 1, w90_num_kpts
1108 do ispin = 1, st%d%dim
1109 if (st%system_grp%is_root()) then
1110 write(filename, '(a,i5.5,a1,i1)') './UNK', ik,'.', ispin
1111 ! write header
1112 if (formatted) then
1113 unk_file = io_open(trim(filename), global_namespace, action='write', form='formatted')
1114 write(unk_file, *) mesh%idx%ll(1:mesh%idx%dim), ik, w90_num_bands
1115 else
1116 unk_file = io_open(trim(filename), global_namespace, action='write', form='unformatted')
1117 write(unk_file) mesh%idx%ll(1:mesh%idx%dim), ik, w90_num_bands
1118 end if
1119 end if
1120
1121 ! states
1122 do ist = 1, st%nst
1123 if (exclude_list(ist)) cycle
1124
1125 if (sys%st%d%ispin /= spin_polarized) then
1126 call states_elec_get_state(st, mesh, ispin, ist, ik, psi)
1127 else
1128 call states_elec_get_state(st, mesh, ispin, ist, (ik-1)*2+w90_spin_channel, psi)
1129 end if
1131 ! put the density in the cube
1132 ! Note: At the moment this does not work for domain parallelization
1133 assert(.not. cube%parallel_in_domains)
1134 call zmesh_to_cube(mesh, psi, cube, cf)
1135
1136 if (st%system_grp%is_root()) then
1137 if (formatted) then
1138 do iz=1,cube%rs_n_global(3)
1139 do iy = 1,cube%rs_n_global(2)
1140 do ix = 1,cube%rs_n_global(1)
1141 w_real = real(cf%zrs(ix,iy,iz), real64)
1142 w_imag = aimag(cf%zrs(ix,iy,iz))
1143 write (unk_file, *) w_real, w_imag
1144 end do
1145 end do
1146 end do
1147 else
1148 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))
1149 end if
1150 end if
1151 end do
1152 if (st%system_grp%is_root()) call io_close(unk_file)
1153 end do
1154 end do
1155
1156 call zcube_function_free_rs(cube, cf)
1157 call cube_end(cube)
1158
1159 safe_deallocate_a(psi)
1160
1161 pop_sub(write_unk)
1162
1163 end subroutine write_unk
1164
1165 ! --------------------------------------------------------------------------
1166 subroutine create_wannier90_amn(space, mesh, latt, st, kpoints)
1167 class(space_t), intent(in) :: space
1168 class(mesh_t), intent(in) :: mesh
1169 type(lattice_vectors_t), intent(in) :: latt
1170 type(states_elec_t), intent(in) :: st
1171 type(kpoints_t), intent(in) :: kpoints
1172
1173 integer :: ist, ik, w90_amn, idim, iw, ip, ik_real
1174 real(real64) :: center(3), kpoint(3), threshold
1175 character(len=80) :: filename
1176 complex(real64), allocatable :: psi(:,:), phase(:), projection(:)
1177 real(real64), allocatable :: ylm(:)
1178 type(orbitalset_t), allocatable :: orbitals(:)
1179
1180 push_sub(create_wannier90_amn)
1181 call profiling_in("W90_AMN")
1182
1183 if (st%parallel_in_states) then
1184 call messages_not_implemented("w90_amn output with states parallelization")
1185 end if
1186
1187 filename = './'// trim(adjustl(w90_prefix))//'.amn'
1188 w90_amn = io_open(trim(filename), global_namespace, action='write')
1189
1190 ! write header
1191 if (st%system_grp%is_root()) then
1192 write(w90_amn,*) 'Created by oct-wannier90'
1193 write(w90_amn,*) w90_num_bands, w90_num_kpts, w90_num_wann
1194 end if
1195
1196 if (scdm_proj) then
1197
1198 message(1) = "Info: Writing projections obtained from SCDM."
1199 call messages_info(1)
1200
1201
1202 do ik = 1, w90_num_kpts
1203 do ist = 1, st%nst
1204 if (exclude_list(ist)) cycle
1205 if (st%system_grp%is_root()) then
1206 do iw = 1, w90_nproj
1207 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)
1208 end do
1209 end if
1210 end do !ist
1211 end do! ik
1212
1213 else
1214
1215 message(1) = "Info: Computing the projection matrix"
1216 call messages_info(1)
1217
1218 !We use the variable AOThreshold to deterine the threshold on the radii of the atomic orbitals
1219 call parse_variable(global_namespace, 'AOThreshold', 0.001_real64, threshold)
1220
1221 safe_allocate(orbitals(1:w90_nproj))
1222 ! precompute orbitals
1223 do iw=1, w90_nproj
1224 call orbitalset_init(orbitals(iw))
1225
1226 orbitals(iw)%norbs = 1
1227 orbitals(iw)%ndim = 1
1228 orbitals(iw)%radius = -log(threshold)
1229 orbitals(iw)%use_submesh = .false.
1230
1231 ! cartesian coordinate of orbital center
1232 center(1:3) = latt%red_to_cart(w90_proj_centers(1:3, iw))
1233 call submesh_init(orbitals(iw)%sphere, space, mesh, latt, center, orbitals(iw)%radius)
1234
1235 ! get dorb as submesh points
1236 safe_allocate(ylm(1:orbitals(iw)%sphere%np))
1237 ! (this is a routine from pwscf)
1238 call ylm_wannier(ylm, w90_proj_lmr(iw,1), w90_proj_lmr(iw,2), &
1239 orbitals(iw)%sphere%r, orbitals(iw)%sphere%rel_x, orbitals(iw)%sphere%np)
1240
1241 ! apply radial function
1242 if (w90_proj_lmr(iw,3) == 1) then
1243 do ip = 1,orbitals(iw)%sphere%np
1244 ylm(ip) = ylm(ip)*m_two*exp(-orbitals(iw)%sphere%r(ip))
1245 end do
1246 else
1247 call messages_not_implemented("oct-wannier90: r/=1 for the radial part")
1248 end if
1249
1250 safe_allocate(orbitals(iw)%zorb(1:orbitals(iw)%sphere%np, 1, 1))
1251 orbitals(iw)%zorb(1:orbitals(iw)%sphere%np, 1, 1) = ylm(1:orbitals(iw)%sphere%np)
1252 safe_deallocate_a(ylm)
1253
1254 safe_allocate(orbitals(iw)%phase(1:orbitals(iw)%sphere%np, st%d%kpt%start:st%d%kpt%end))
1255 orbitals(iw)%phase(:,:) = m_z0
1256 safe_allocate(orbitals(iw)%eorb_mesh(1:mesh%np, 1, 1, st%d%kpt%start:st%d%kpt%end))
1257 orbitals(iw)%eorb_mesh(:,:,:,:) = m_z0
1258
1259 call orbitalset_update_phase(orbitals(iw), space%dim, st%d%kpt, kpoints, st%d%ispin == spin_polarized, &
1260 kpt_max = w90_num_kpts)
1262 end do
1263
1264 safe_allocate(psi(1:mesh%np, 1:st%d%dim))
1265 safe_allocate(phase(1:mesh%np))
1266 safe_allocate(projection(1:w90_nproj))
1267
1268 do ik = 1, w90_num_kpts
1269 kpoint(1:space%dim) = kpoints%get_point(ik)
1270 !$omp parallel do
1271 do ip = 1, mesh%np
1272 phase(ip) = exp(-m_zi* sum(mesh%x(1:space%dim, ip) * kpoint(1:space%dim)))
1273 end do
1274
1275 !For spin-polarized calculations, we select the right k-point
1276 if (st%d%ispin == spin_polarized) then
1277 ik_real = (ik-1)*2 + w90_spin_channel
1278 else
1279 ik_real = ik
1280 end if
1281
1282
1283 do ist = 1, st%nst
1284 if (exclude_list(ist)) cycle
1285
1286 projection = m_zero
1287
1288 if(ik_real >= st%d%kpt%start .and. ik_real <= st%d%kpt%end) then
1289 call states_elec_get_state(st, mesh, ist, ik_real, psi)
1290
1291 !$omp parallel
1292 do idim = 1, st%d%dim
1293 !The minus sign is here is for the wrong convention of Octopus
1294 !$omp do
1295 do ip = 1, mesh%np
1296 psi(ip, idim) = psi(ip, idim)*phase(ip)
1297 end do
1298 end do
1299 !$omp end parallel
1300
1301 do iw = 1, w90_nproj
1302 idim = 1
1303 if (w90_spinors) idim = w90_spin_proj_component(iw)
1304
1305 !At the moment the orbitals do not depend on idim
1306 !The idim index for eorb_mesh would be for a spin-resolved orbital like j=1/2
1307 projection(iw) = zmf_dotp(mesh, psi(1:mesh%np,idim), &
1308 orbitals(iw)%eorb_mesh(1:mesh%np,1,1,ik_real), reduce = .false.)
1309 end do
1310
1311 call profiling_in("W90_AMN_REDUCE")
1312 call mesh%allreduce(projection)
1313 call profiling_out("W90_AMN_REDUCE")
1314 end if
1315
1316 if(st%d%kpt%parallel) then
1317 call comm_allreduce(st%d%kpt%mpi_grp, projection)
1318 end if
1319
1320 if (st%system_grp%is_root()) then
1321 do iw = 1, w90_nproj
1322 write (w90_amn,'(I5,2x,I5,2x,I5,2x,e18.10,2x,e18.10)') band_index(ist), iw, ik, projection(iw)
1323 end do
1324 end if
1325 end do !ist
1326 end do !ik
1327
1328 safe_deallocate_a(psi)
1329 safe_deallocate_a(phase)
1330 safe_deallocate_a(projection)
1331
1332 do iw = 1, w90_nproj
1333 call orbitalset_end(orbitals(iw))
1334 end do
1335 safe_deallocate_a(orbitals)
1336 end if
1337
1338 call io_close(w90_amn)
1339
1340 call profiling_out("W90_AMN")
1341
1342 pop_sub(create_wannier90_amn)
1343
1344 end subroutine create_wannier90_amn
1345
1346 ! --------------------------------------------------------------------------
1350 subroutine create_wannier90_spn(mesh, st)
1351 class(mesh_t), intent(in) :: mesh
1352 type(states_elec_t), target, intent(in) :: st
1353
1354 integer :: ist, jst, ik, w90_spn, counter
1355 character(len=80) :: filename
1356 complex(real64), allocatable :: spin(:,:,:)
1357 complex(real64), allocatable :: psim(:,:), psin(:,:)
1358 complex(real64) :: dot_upup, dot_updown, dot_downup, dot_downdown
1359
1360 push_sub(create_wannier90_spn)
1361 call profiling_in("W90_SPN")
1362
1363 assert(st%d%ispin == spinors)
1364
1365 if (st%parallel_in_states) then
1366 call messages_not_implemented("w90_spn output with states parallelization")
1367 end if
1368
1369 message(1) = "Info: Computing the spin file"
1370 call messages_info(1)
1371
1372 filename = './'// trim(adjustl(w90_prefix))//'.spn'
1373 w90_spn = io_open(trim(filename), global_namespace, action='write')
1374
1375 ! write header
1376 if (st%system_grp%is_root()) then
1377 write(w90_spn,*) 'Created by oct-wannier90'
1378 write(w90_spn,*) w90_num_bands, w90_num_kpts
1379 end if
1380
1381 safe_allocate(psim(1:mesh%np, 1:st%d%dim))
1382 safe_allocate(psin(1:mesh%np, 1:st%d%dim))
1383 safe_allocate(spin(1:3, 1:(w90_num_bands*(w90_num_bands+1))/2, 1:w90_num_kpts))
1384 spin = m_zero
1385
1386 ! loop over the pairs specified in the nnkp file (read before in init)
1387 do ik = st%d%kpt%start, st%d%kpt%end
1388 counter = 0
1389 do jst = 1, st%nst
1390 if (exclude_list(jst)) cycle
1391
1392 call states_elec_get_state(st, mesh, jst, ik, psim)
1393 do ist = 1, jst
1394 if (exclude_list(ist)) cycle
1395
1396 counter = counter + 1
1397
1398 call states_elec_get_state(st, mesh, ist, ik, psin)
1399
1400 dot_upup = zmf_dotp(mesh, psin(:, 1), psim(:, 1), reduce = .false.)
1401 dot_downdown = zmf_dotp(mesh, psin(:, 2), psim(:, 2), reduce = .false.)
1402 dot_updown = zmf_dotp(mesh, psin(:, 1), psim(:, 2), reduce = .false.)
1403 dot_downup = zmf_dotp(mesh, psin(:, 2), psim(:, 1), reduce = .false.)
1404
1405 spin(1, counter, ik) = dot_updown + dot_downup
1406 spin(2, counter, ik) = -m_zi * dot_updown + m_zi * dot_downup
1407 spin(3, counter, ik) = dot_upup - dot_downdown
1408 end do !ist
1409 end do
1410 end do
1411
1412 call profiling_in("W90_SPN_REDUCE")
1413 call mesh%allreduce(spin)
1414
1415 if(st%d%kpt%parallel) then
1416 call comm_allreduce(st%d%kpt%mpi_grp, spin)
1417 end if
1418 call profiling_out("W90_SPN_REDUCE")
1419
1420 ! write to W90 file
1421 if (st%system_grp%is_root()) then
1422 do ik = 1, w90_num_kpts
1423 counter = 0
1424 do jst = 1, st%nst
1425 if (exclude_list(jst)) cycle
1426
1427 do ist = 1, jst
1428 if (exclude_list(ist)) cycle
1429
1430 counter = counter + 1
1431 write(w90_spn, '(e18.10,2x,e18.10)') spin(1, counter, ik)
1432 write(w90_spn, '(e18.10,2x,e18.10)') spin(2, counter, ik)
1433 write(w90_spn, '(e18.10,2x,e18.10)') spin(3, counter, ik)
1434 end do
1435 end do
1436 end do !ik
1437 end if
1438
1439 call io_close(w90_spn)
1440
1441 safe_deallocate_a(psim)
1442 safe_deallocate_a(psin)
1443 safe_deallocate_a(spin)
1444
1445 call profiling_out("W90_SPN")
1446
1447 pop_sub(create_wannier90_spn)
1448 end subroutine create_wannier90_spn
1449
1450
1451 ! --------------------------------------------------------------------------
1452 subroutine generate_wannier_states(space, mesh, ions, st, kpoints)
1453 class(space_t), intent(in) :: space
1454 class(mesh_t), intent(in) :: mesh
1455 type(ions_t), intent(in) :: ions
1456 type(states_elec_t), intent(in) :: st
1457 type(kpoints_t), intent(in) :: kpoints
1458
1459 integer :: w90_u_mat, w90_xyz, nwann, nik
1460 integer :: ik, iw, iw2, ip, ipmax, rankmax, idmmax
1461 real(real64), allocatable :: centers(:,:), supercell_centers(:,:), new_centers(:,:)
1462 complex(real64), allocatable :: Umnk(:,:,:)
1463 complex(real64), allocatable :: zwn(:,:,:), psi(:,:), phase(:,:), phase_bloch(:), zwn_bloch(:,:,:)
1464 character(len=MAX_PATH_LEN) :: fname
1465 real(real64) :: kpoint(3), wmod, wmodmax, xx(space%dim)
1466 character(len=2) :: dum
1467 logical :: exist
1468 type(unit_t) :: fn_unit
1469 complex(real64) :: scal
1470 type(block_t) :: blk
1471 integer :: supercell(space%dim), ii, jj, kk, ncols, Nreplica, irep, irepmax
1472 real(real64) :: min_image_displ(3), offset(space%dim)
1473 integer :: jw
1474 integer, allocatable :: parent(:)
1475 real(real64), parameter :: tol_cluster = 0.25_real64 ! In Bohr
1476
1477
1478 push_sub(generate_wannier_states)
1479
1480 message(1) = "oct-wannier90: Constructing the Wannier states from the U matrix."
1481 call messages_info(1)
1482
1483 inquire(file=trim(trim(adjustl(w90_prefix))//'_centres.xyz'),exist=exist)
1484 if (.not. exist) then
1485 message(1) = 'oct-wannier90: Cannot find the Wannier90 file seedname_centres.xyz.'
1486 write(message(2),'(a)') 'Please run wannier90.x with "write_xyz=.true." in '// trim(adjustl(w90_prefix)) // '.'
1487 call messages_fatal(2)
1488 end if
1489
1490 w90_xyz = io_open(trim(trim(adjustl(w90_prefix))//'_centres.xyz'), global_namespace, action='read')
1491
1492 safe_allocate(centers(1:3, 1:w90_num_wann))
1493 !Skip two lines
1494 read(w90_xyz, *)
1495 read(w90_xyz, *)
1496 do iw = 1, w90_num_wann
1497 read(w90_xyz, *) dum, centers(1:3, iw)
1498 ! Wannier90 outputs the coordinates in angstrom
1499 centers(1:3, iw) = units_to_atomic(unit_angstrom, centers(1:3, iw))
1500 end do
1501 call io_close(w90_xyz)
1502
1503 ! In order to find "more consistent" Wannier centers, the centers are clustered in groups,
1504 ! Using periodic boundary condition to find the centers that are the closest
1505 safe_allocate(parent(1:w90_num_wann))
1506 do iw = 1, w90_num_wann
1507 parent(iw) = iw
1508 do jw = 1, iw-1
1509 min_image_displ(:) = ions%latt%cart_to_red(centers(:, iw) - centers(:, jw))
1510 min_image_displ(1:space%periodic_dim) = min_image_displ(1:space%periodic_dim) - nint(min_image_displ(1:space%periodic_dim))
1511 min_image_displ = ions%latt%red_to_cart(min_image_displ)
1512 if (norm2(min_image_displ) < tol_cluster) then
1513 parent(iw) = parent(jw)
1514 exit
1515 end if
1516 end do
1517 end do
1518
1519 message(1) = "Info : Clustering of the Wannier centers"
1520 call messages_info(1)
1521
1522 safe_allocate(new_centers(1:3, 1:w90_num_wann))
1523 do iw = 1, w90_num_wann
1524 write(message(1), '(a,i0,a,3(f7.3,a))') 'Info : Original Wannier center ', &
1525 iw, ' (', centers(1, iw), ', ', centers(2, iw), ', ', centers(3, iw), ')'
1526
1527 if (parent(iw) == iw) then
1528 ! We are computing |w_{n,0}>, so we need to take only a shift to get to the central cell
1529 new_centers(1:3, iw) = ions%latt%fold_into_cell(centers(1:3, iw))
1530 else
1531 ! This center has a closeby periodic parent. We then move it to be as close as possible to its "parent"
1532 ! Thanks to the loop above, this should always be inside of very close to a parent inside the unit cell
1533 min_image_displ(:) = ions%latt%cart_to_red(centers(:, iw) - new_centers(:, parent(iw)))
1534 min_image_displ(1:space%periodic_dim) = real(nint(min_image_displ(1:space%periodic_dim)), real64)
1535 min_image_displ = ions%latt%red_to_cart(min_image_displ)
1536 new_centers(1:3, iw) = centers(1:3, iw) - min_image_displ
1537 end if
1538
1539 write(message(2), '(a,i0,a,3(f7.3,a))') 'Info : New Wannier center ', &
1540 iw, ' (', new_centers(1, iw), ', ', new_centers(2, iw), ', ', new_centers(3, iw), ')'
1541 write(message(3), '(a)') ''
1542 call messages_info(3)
1543 end do
1544
1545 ! Getting the offset, for the correct phase factor
1546 centers = new_centers - centers
1547 safe_deallocate_a(parent)
1548 safe_deallocate_a(new_centers)
1549
1550 inquire(file=trim(trim(adjustl(w90_prefix))//'_u_dis.mat'),exist=exist)
1551 if (exist) then
1552 message(1) = 'oct-wannier90: Calculation of Wannier states with disentanglement is not yet supported.'
1553 call messages_fatal(1)
1554 end if
1555
1556 inquire(file=trim(trim(adjustl(w90_prefix))//'_u.mat'),exist=exist)
1557 if (.not. exist) then
1558 message(1) = 'oct-wannier90: Cannot find the Wannier90 seedname_u.mat file.'
1559 write(message(2),'(a)') 'Please run wannier90.x with "write_u_matrices=.true." in '// trim(adjustl(w90_prefix)) // '.'
1560 call messages_fatal(2)
1561 end if
1562 w90_u_mat = io_open(trim(trim(adjustl(w90_prefix))//'_u.mat'), global_namespace, action='read')
1563
1564
1565 !Skip one line
1566 read(w90_u_mat, *)
1567 !Read num_kpts, num_wann, num_wann for consistency check
1568 read(w90_u_mat, *) nik, nwann, nwann
1569 if (nik /= w90_num_kpts .or. nwann /= w90_num_wann) then
1570 message(1) = "The file contains U matrices is inconsistent with the .win file."
1571 call messages_fatal(1)
1572 end if
1573
1574 safe_allocate(umnk(1:w90_num_wann, 1:w90_num_wann, 1:w90_num_kpts))
1575
1576 do ik = 1, w90_num_kpts
1577 !Skip one line
1578 read(w90_u_mat, *)
1579 !Skip one line (k-point coordinate)
1580 read(w90_u_mat, *)
1581 read(w90_u_mat, '(f15.10,sp,f15.10)') ((umnk(iw, iw2, ik), iw=1, w90_num_wann), iw2=1, w90_num_wann)
1582 end do
1583
1584 call io_close(w90_u_mat)
1585
1586 !We read the output format for the Wannier states
1587 call parse_variable(global_namespace, 'OutputFormat', 0, how)
1588 if (how == 0) then
1589 message(1) = "OutputFormat must be specified for outputing Wannier functions."
1590 call messages_fatal(1)
1591 end if
1592
1593 !%Variable Wannier90Supercell
1594 !%Type block
1595 !%Section Utilities::oct-wannier90
1596 !%Description
1597 !% This block allows to specify the size of the supercell used to compute the Wannier functions
1598 !% If not specified, the code uses a default 1x1x1 cell, i.e., it plots the Wannier90 in the primitive cell.
1599 !%End
1600 if (parse_is_defined(sys%namespace, 'Wannier90Supercell')) then
1601 if (parse_block(sys%namespace, 'Wannier90Supercell', blk) == 0) then
1602 ncols = parse_block_cols(blk, 0)
1603 if (ncols /= space%dim) then
1604 write(message(1),'(a,i3,a,i3)') 'Wannier90Supercell has ', ncols, ' columns but must have ', sys%space%dim
1605 call messages_fatal(1, namespace=sys%namespace)
1606 end if
1607 do ii = 1, space%dim
1608 call parse_block_integer(blk, 0, ii - 1, supercell(ii))
1609 end do
1610
1611 call parse_block_end(blk)
1612 end if
1613 else
1614 supercell = 1
1615 end if
1616
1617 nreplica = product(supercell)
1618
1619 ! The center of each replica of the unit cell
1620 safe_allocate(supercell_centers(1:space%dim, 1:nreplica))
1621 offset(1:space%dim) = -floor((real(supercell(1:space%dim), real64) - m_one) / m_two)
1622 irep = 1
1623 do ii = 0, supercell(1)-1
1624 do jj = 0, supercell(2)-1
1625 do kk = 0, supercell(3)-1
1626 supercell_centers(:, irep) = ions%latt%red_to_cart(offset &
1627 + [real(ii, real64), real(jj, real64), real(kk, real64)])
1628 irep = irep + 1
1629 end do
1630 end do
1631 end do
1632
1633
1634
1635 call io_mkdir('wannier', global_namespace)
1636
1637 !Computing the Wannier states in the primitive cell, from the U matrices
1638 safe_allocate(zwn(1:mesh%np, 1:nreplica, 1:st%d%dim))
1639 safe_allocate(psi(1:mesh%np, 1:st%d%dim))
1640 safe_allocate(phase(1:mesh%np, 1:nreplica))
1641 safe_allocate(phase_bloch(1:mesh%np))
1642 if (w90_bloch_sums) then
1643 safe_allocate(zwn_bloch(1:mesh%np, 1:st%d%dim, st%d%kpt%start:st%d%kpt%end))
1644 end if
1645
1646 do iw = 1, w90_num_wann
1647
1648 zwn(:, :, :) = m_z0
1649
1650 do ik = 1, w90_num_kpts
1651
1652 if (.not. (ik >= st%d%kpt%start .and. ik <= st%d%kpt%end)) cycle
1653
1654 if (w90_bloch_sums) zwn_bloch(:,:,ik) = m_z0
1655
1656 kpoint(1:space%dim) = kpoints%get_point(ik, absolute_coordinates=.true.)
1657
1658 ! We compute the Wannier orbital on a grid centered around the Wannier function
1659 ! The minus sign is here is for the wrong convention of Octopus
1660 do irep = 1, nreplica
1661 !$omp parallel do
1662 do ip = 1, mesh%np
1663 xx = mesh%x(1:space%dim, ip)-centers(1:space%dim, iw) + supercell_centers(1:space%dim, irep)
1664 phase(ip, irep) = exp(-m_zi* sum( xx * kpoint(1:space%dim)))
1665 end do
1666 end do
1667 if (w90_bloch_sums) then
1668 do ip = 1, mesh%np
1669 xx = mesh%x(1:space%dim, ip)
1670 phase_bloch(ip) = exp(-m_zi* sum( xx * kpoint(1:space%dim)))
1671 end do
1672 end if
1673
1674 do iw2 = 1, st%nst
1675 if (exclude_list(iw2)) cycle
1676
1677 if (st%d%ispin /= spin_polarized) then
1678 call states_elec_get_state(st, mesh, iw2, ik, psi)
1679 else
1680 call states_elec_get_state(st, mesh, iw2, (ik-1)*2+w90_spin_channel, psi)
1681 end if
1682
1683 !$omp parallel
1684 do idim = 1, st%d%dim
1685 do irep = 1, nreplica
1686 !$omp do
1687 do ip = 1, mesh%np
1688 zwn(ip, irep, idim) = zwn(ip, irep, idim) + umnk(band_index(iw2), iw, ik) * psi(ip, idim) * phase(ip, irep)
1689 end do
1690 end do
1691 end do
1692 !$omp end parallel
1693
1694 ! See Eq. 19 in arXiv:0605539
1695 if (w90_bloch_sums) then
1696 !$omp parallel
1697 do idim = 1, st%d%dim
1698 !$omp do
1699 do ip = 1, mesh%np
1700 zwn_bloch(ip, idim, ik) = zwn_bloch(ip, idim, ik) &
1701 + umnk(band_index(iw2), iw, ik) * psi(ip, idim) * phase_bloch(ip)
1702 end do
1703 end do
1704 !$omp end parallel
1705 end if
1706 end do!iw2
1707 end do!ik
1708
1709 if(st%d%kpt%parallel) then
1710 call comm_allreduce(st%d%kpt%mpi_grp, zwn)
1711 end if
1712
1713 ! Following what Wannier90 is doing, we fix the global phase by setting the max to be real
1714 ! We also normalize to the number of k-point at this step
1715 if (sys%st%d%ispin /= spinors) then
1716 ipmax = 0
1717 irepmax = 0
1718 wmodmax = m_zero
1719 do irep = 1, nreplica
1720 do ip = 1, mesh%np
1721 wmod = real(zwn(ip, irep, 1)*conjg(zwn(ip, irep, 1)), real64)
1722 if (wmod > wmodmax) then
1723 ipmax = ip
1724 irepmax = irep
1725 wmodmax = wmod
1726 end if
1727 end do
1728 end do
1729 scal = sqrt(wmodmax)/zwn(ipmax, irepmax, 1)/w90_num_kpts
1730 call mesh_minmaxloc(mesh, wmodmax, rankmax, mpi_maxloc)
1731 call mesh%mpi_grp%bcast(scal, 1, mpi_double_complex, rankmax)
1732 call lalg_scal(mesh%np, nreplica, scal, zwn(:,:,1))
1733
1734 if (w90_bloch_sums) then
1735 do ik = st%d%kpt%start, st%d%kpt%end
1736 ipmax = 0
1737 idmmax= 0
1738 wmodmax = m_zero
1739 do idim = 1, st%d%dim
1740 do ip = 1, mesh%np
1741 wmod = real(zwn_bloch(ip, idim, ik)*conjg(zwn_bloch(ip, idim, ik)), real64)
1742 if (wmod > wmodmax) then
1743 ipmax = ip
1744 wmodmax = wmod
1745 idmmax = idim
1746 end if
1747 end do
1748 end do
1749 scal = sqrt(wmodmax)/zwn_bloch(ipmax, idmmax, ik)/w90_num_kpts
1750 call mesh_minmaxloc(mesh, wmodmax, rankmax, mpi_maxloc)
1751 call mesh%mpi_grp%bcast(scal, 1, mpi_double_complex, rankmax)
1752 call lalg_scal(mesh%np, st%d%dim, scal, zwn_bloch(:,:,ik))
1753 end do
1754 end if
1755 end if
1756
1757 ! Output the Wannier states
1758 fn_unit = sqrt(units_out%length**(-space%dim))
1759 do idim = 1, st%d%dim
1760 if (st%d%ispin == spinors) then
1761 write(fname, '(a,i3.3,a4,i1)') 'wannier-', iw, '-isp', idim
1762 else
1763 write(fname, '(a,i3.3,a4,i1)') 'wannier-', iw
1764 end if
1765 if (any(supercell_centers > 1)) then
1766 call io_function_output_supercell(how, "wannier", trim(fname), mesh, &
1767 space, ions%latt, zwn(:, :, idim), supercell_centers, supercell, fn_unit, ierr, global_namespace, &
1768 pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
1769
1770 else
1771 call zio_function_output(how, "wannier", trim(fname), global_namespace, space, mesh, &
1772 zwn(:, 1, idim), fn_unit, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
1773 end if
1774 end do
1775
1776 ! Output Bloch sums of Wannier states
1777 if (w90_bloch_sums) then
1778 do ik = st%d%kpt%start, st%d%kpt%end
1779 do idim = 1, st%d%dim
1780 if (st%d%ispin == spinors) then
1781 write(fname, '(a,i3.3,a4,i1,a3,i5.5)') 'wannier_bloch-', iw, '-isp', idim, '-ik', ik
1782 else
1783 write(fname, '(a,i3.3,a4,i1,a3,i5.5)') 'wannier_bloch-', iw, '-ik', ik
1784 end if
1785 call zio_function_output(how, "wannier", trim(fname), global_namespace, space, mesh, &
1786 zwn_bloch(:, idim, ik), fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
1787 end do
1788 end do
1789 end if
1790
1791 ! Checking the ratio imag/real
1792 if (sys%st%d%ispin /= spinors) then
1793 wmodmax = m_zero
1794 do irep = 1, nreplica
1795 do ip = 1, mesh%np
1796 if(abs(real(zwn(ip, irep, 1), real64)) >= 1e-2_real64) then
1797 wmodmax = max(wmodmax, abs(aimag(zwn(ip, irep, 1)))/abs(real(zwn(ip, irep, 1), real64)))
1798 end if
1799 end do
1800 end do
1801 call mesh_minmaxloc(mesh, wmodmax, rankmax, mpi_maxloc)
1802
1803 write(message(1), '(a,i4,a,f11.6)') 'oct-wannier90: Wannier function ', iw, ' Max. Im/Re Ratio = ', wmodmax
1804 call messages_info(1)
1805 else
1806 write(message(1), '(a,i4)') 'oct-wannier90: Wannier function done ', iw
1807 call messages_info(1)
1808 end if
1809 end do
1810
1811 safe_deallocate_a(umnk)
1812 safe_deallocate_a(zwn)
1813 safe_deallocate_a(psi)
1814 safe_deallocate_a(phase)
1815 safe_deallocate_a(phase_bloch)
1816 safe_deallocate_a(centers)
1817 safe_deallocate_a(supercell_centers)
1818
1820 end subroutine generate_wannier_states
1821
1822end program wannier90_interface
1823
1824!! Local Variables:
1825!! mode: f90
1826!! coding: utf-8
1827!! 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
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:202
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:494
real(real64), parameter, public m_huge
Definition: global.F90:218
real(real64), parameter, public m_zero
Definition: global.F90:200
complex(real64), parameter, public m_z0
Definition: global.F90:210
complex(real64), parameter, public m_zi
Definition: global.F90:214
subroutine, public global_init(communicator)
Initialise Octopus.
Definition: global.F90:375
real(real64), parameter, public m_half
Definition: global.F90:206
real(real64), parameter, public m_one
Definition: global.F90:201
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 functions over batches of mesh functions.
Definition: mesh_batch.F90:118
subroutine, public zmesh_batch_mf_dotp(mesh, aa, psi, dot, reduce, nst)
calculate the dot products between a batch and a vector of mesh functions
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:273
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1068
subroutine, public messages_init(output_dir)
Definition: messages.F90:220
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1040
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
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:463
subroutine, public parser_init()
Initialise the Octopus parser.
Definition: parser.F90:410
subroutine, public parser_end()
End the Octopus parser.
Definition: parser.F90:442
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:184
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:223
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()