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