Octopus
wannier.F90
Go to the documentation of this file.
1!! Copyright (C) A Buccheri, J Reimann. 2026
2!!
3!! This Source Code Form is subject to the terms of the Mozilla Public
4!! License, v. 2.0. If a copy of the MPL was not distributed with this
5!! file, You can obtain one at https://mozilla.org/MPL/2.0/.
6!!
7
8#include "global.h"
9
12program wannier
13 use, intrinsic :: iso_fortran_env
14
16 use debug_oct_m
19 use fft_oct_m
20 use global_oct_m
21 use grid_oct_m
22 use io_oct_m
24 use ions_oct_m
26 use mpi_oct_m
30 use parser_oct_m
36 use unit_oct_m
41
42 ! External Wannier 90 library
43 use w90_library
44 use w90_library_extra, only:overlaps
45
46 implicit none
47
48 type(namespace_t) :: namespace
49 type(electrons_t), pointer :: sys
50 type(wannier_opts_t) :: wannier_opts
51 integer, allocatable :: band_index(:)
52
53 type(lib_common_type) :: w90main
54
55 complex(real64), allocatable :: m_matrix(:, :, :, :)
56 complex(real64), allocatable :: u_matrix(:, :, :)
57 complex(real64), allocatable :: u_matrix_opt(:, :, :)
58 real(real64), allocatable :: eigenvals(:,:)
59
60 integer, allocatable :: nnkp(:,:), gkpb(:,:,:)
61 logical, allocatable :: exclude_list(:)
62 integer :: nn
63
64 integer :: idir, ii, itemp, ik, ist
65
66 integer :: ierr
67 real(real64) :: factor(3)
68
69 ! Initialise all Octopus singletons
70 call global_init()
71 call parser_init()
72 call messages_init()
73 call io_init()
74 namespace = namespace_t("Wannier90Lib", parent = global_namespace)
75 call profiling_init(namespace)
76 call fft_all_init(namespace)
77 call unit_system_init(namespace)
78
79 ! Initialise a system object
80 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
81 sys => electrons_t(global_namespace, mpi_world, int(option__calculationmode__dummy, int32))
82
83 if (sys%space%dim /= 3) then ! see issue #1466
84 call messages_not_implemented("oct-wannier90lib: only 3D is supported")
85 end if
86
87 ! Sanity checks
88 if (sys%kpoints%use_symmetries) then
89 message(1) = 'oct-wannier90: k-points symmetries are not allowed'
90 call messages_fatal(1)
91 end if
92 if (sys%kpoints%use_time_reversal) then
93 message(1) = 'oct-wannier90: time-reversal symmetry is not allowed'
94 call messages_fatal(1)
95 end if
96 if (sys%kpoints%reduced%nshifts > 1) then
97 message(1) = 'oct-wannier90: Wannier90 does not allow for multiple shifts of the k-point grid'
98 call messages_fatal(1)
99 end if
100 if (sys%st%d%nspin == 4) then
101 call messages_not_implemented('oct-wannier90: Spinors currently not supported', namespace)
102 end if
103
104 ! "Correct" rlattice/klattice for Wannier90 (3D periodicity assumed here).
105 ! Note: this is how the lattice is scaled in src/utils/wannier90_interface.F90
106 factor = m_one
107 do idir = sys%space%periodic_dim+1, sys%space%dim
108 factor(idir) = m_two * sys%gr%box%bounding_box_l(idir)
109 end do
110 call sys%ions%latt%scale(factor)
111
112 ! Parsing all input related to wannier90
113 call wannier_opts%parse_oct(namespace)
114 call wannier_opts%parse_win()
115
116 ! Initialize wannier90 lib object
117 call wannier90lib_init_w90main(namespace, sys%ions, sys%kpoints, sys%st, &
118 wannier_opts, w90main, stdout, stderr, ierr)
119 if (ierr /= 0) then
120 write(message(1),'(a,i0)') 'Error initializing wannier90 library: ', ierr
121 call messages_fatal(1)
122 end if
123
124 ! Check for consistency of spin polarization
125 if (sys%st%d%ispin /= unpolarized) then
126 call messages_experimental("oct-wannier90 with SpinComponents /= unpolarized")
127 elseif(sys%st%d%ispin == unpolarized .and. w90main%wvfn_read%spin_channel == 2) then
128 w90main%wvfn_read%spin_channel = 1
129 message(1) = "octopus was run with SpinComponents == unpolarized. Ignoring wannier90 spin input variable"
130 call messages_warning(1)
131 end if
132
133 ! Load octopus states
134 call wannier_calc_load_restart(global_namespace, sys%mc, sys%space, sys%st, sys%gr, sys%kpoints, wannier_opts%restart_from)
135
136 ! Convert the exclude bands to exclude list
137 safe_allocate(band_index(1:sys%st%nst))
138 safe_allocate(exclude_list(1:sys%st%nst))
139 exclude_list(1:sys%st%nst) = .false.
140
141 if (allocated(w90main%exclude_bands)) then
142 do itemp = 1, size(w90main%exclude_bands)
143 exclude_list(w90main%exclude_bands(itemp)) = .true.
144 end do
145 end if
146
147 itemp = 0
148 do ii = 1, sys%st%nst
149 if (exclude_list(ii)) cycle
150 itemp = itemp + 1
151 band_index(ii) = itemp
152 end do
153
154 ! Initial setup (previously in nnkp file)
155 call w90_get_nn(w90main, nn, stdout, stderr, ierr)
156 safe_allocate(nnkp(1:w90main%num_kpts, 1:nn))
157 safe_allocate(gkpb(1:3, 1:w90main%num_kpts, 1:nn))
158 call w90_get_nnkp(w90main, nnkp, stdout, stderr, ierr)
159 call w90_get_gkpb(w90main, gkpb, stdout, stderr, ierr)
160
161 ! Prepare computation of quantities required by Wannier90
162 safe_allocate(m_matrix(1:wannier_opts%num_bands, 1:wannier_opts%num_bands, 1:nn, 1:w90main%num_kpts))
163 safe_allocate(u_matrix_opt(1:wannier_opts%num_bands, 1:wannier_opts%num_wann, 1:w90main%num_kpts))
164 safe_allocate(eigenvals(1:wannier_opts%num_bands, 1:w90main%num_kpts))
165 safe_allocate(u_matrix(1:wannier_opts%num_wann, 1:wannier_opts%num_wann, 1:w90main%num_kpts))
166 call w90_set_m_local(w90main, m_matrix)
167 call w90_set_u_opt(w90main, u_matrix_opt)
168 call w90_set_u_matrix(w90main, u_matrix)
169
170 ! Get eigenvalues from octopus, shape is (num_bands, nk)
171 do ik = 1, w90main%num_kpts
172 do ist = 1, sys%st%nst
173 if (exclude_list(ist)) cycle
174 if (sys%st%d%ispin /= spin_polarized) then
175 eigenvals(band_index(ist), ik) = units_from_atomic(unit_ev, sys%st%eigenval(ist, ik))
176 else
177 eigenvals(band_index(ist), ik) = units_from_atomic(unit_ev, sys%st%eigenval(ist, (ik-1)*2+w90main%wvfn_read%spin_channel))
178 end if
179 end do
180 end do
181 call w90_set_eigval(w90main, eigenvals)
182
183 ! Compute or read quantities required by wannier90
184 if (wannier_opts%calc_overlaps) then
185 call wannier90lib_create_wannier90_mmn(w90main, exclude_list, band_index, sys%gr, sys%st, sys%ions, m_matrix)
186 call wannier90lib_create_wannier90_amn(w90main, wannier_opts, exclude_list, band_index, sys%space, &
187 sys%gr, sys%ions%latt, sys%st, sys%kpoints, u_matrix_opt)
188 ! TODO: Revisit spin-polarised functionality once spin-unpolarised works
189 ! call create_wannier90_spn(sys%gr, sys%st)
190 ! TODO: Add write functions here and remove dump_inputs variable in wannier90
191 else
192 call overlaps(w90main, stdout, stderr, ierr)
193 if (ierr /= 0) then
194 write(message(1), '(a,i0)') 'Error while reading overlaps with wannier90: ', ierr
195 call messages_fatal(1)
196 end if
197 end if
198
199 ! Perform wannierization
200 if (wannier_opts%wannierize) then
201 call w90_disentangle(w90main, stdout, stderr, ierr)
202 if (ierr /= 0) then
203 write(message(1), '(a,i0)') 'Error in w90_disentangle: ', ierr
204 call messages_fatal(1)
205 end if
206 call w90_project_overlap(w90main, stdout, stderr, ierr)
207
209 if (ierr /= 0) then
210 write(message(1), '(a,i0)') 'Error in w90_project_overlap: ', ierr
211 call messages_fatal(1)
212 end if
213 call w90_wannierise(w90main, stdout, stderr, ierr)
214 if (ierr /= 0) then
215 write(message(1), '(a,i0)') 'Error in w90_wannierise: ', ierr
216 call messages_fatal(1)
217 end if
218 ! Avoid plotting from Wannier90
219 w90main%w90_calculation%wannier_plot = .false.
220 ! Output U matrices, hr.dat and more output files for plotting with octopus
221 call w90_plot(w90main, stdout, stderr, ierr)
222 if (ierr /= 0) then
223 write(message(1), '(a,i0)') 'Error in w90_plot: ', ierr
224 call messages_fatal(1)
225 end if
226 end if
227 ! if (wan_opts%plot) then
228 ! Here we generate wannier functions
229 ! Currently this is only meaningful after performing the wannierization
230 ! TODO: Check functionality to read U matrix from file
231 ! generate_wannier_states(space, mesh, ions, st, kpoints)
232 ! TODO: Port function to generate wannier states from old utility
233 ! end if
234
235 ! Free memory of variables in this scope
236 safe_deallocate_a(m_matrix)
237 safe_deallocate_a(u_matrix)
238 safe_deallocate_a(u_matrix_opt)
239 safe_deallocate_a(eigenvals)
240 safe_deallocate_a(exclude_list)
241 safe_deallocate_a(band_index)
242 safe_deallocate_a(nnkp)
243 safe_deallocate_a(gkpb)
244 safe_deallocate_p(sys)
245
246 ! End all singletons
247 call fft_all_end()
248 call io_end()
249 call profiling_end(namespace)
250 call messages_end()
251 call parser_end()
252 call global_end()
253
254end program wannier
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
integer, parameter, public unpolarized
Parameters...
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
subroutine, public global_init(communicator)
Initialise Octopus.
Definition: global.F90:375
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_end()
Definition: io.F90:271
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
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 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_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
This module defines routines to write information about states.
This module handles reading and writing restart information for the states_elec_t.
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.
subroutine, public unit_system_init(namespace)
type(unit_t), public unit_ev
For output energies in eV.
Interface module to Wannier 90 library.
subroutine, public wannier90lib_create_wannier90_amn(w90main, inp_options, exclude_list, band_index, space, mesh, latt, st, kpoints, projection)
Calculate wannier90 Projection Matrix.
subroutine, public wannier90lib_create_wannier90_mmn(w90main, exclude_list, band_index, mesh, st, ions, overlap)
Kohn-Sham State Overlap Matrix.
subroutine, public wannier90lib_init_w90main(namespace, ions, kpoints, st, inp_options, w90main, stdout, stderr, ierr)
Initialize wannier90 library data.
Wannier90 related calculations.
subroutine, public wannier_calc_load_restart(namespace, mc, space, st, gr, kpoints, restart_states)
Load Octopus restart data from disk.
Wannier options module.
Class describing the electron system.
Definition: electrons.F90:223
int true(void)
program wannier
Construct maximally-localised Wannier functions through API calls to Wannier90 library.
Definition: wannier.F90:107