Octopus
wannier90lib_interface.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#include "global.h"
8
13
15 use debug_oct_m
16 use global_oct_m
17 use ions_oct_m, only: ions_t
18 use io_oct_m, only: io_open, io_close
19 use kpoints_oct_m, only: kpoints_t
21 use mesh_oct_m, only: mesh_t
23 use mpi_oct_m, only: mpi_request
27 use space_oct_m, only: space_t
32 only: &
37 use wfs_elec_oct_m, only: wfs_elec_t
38
39 ! External library
40 ! Wannier90 lib type and routines are not dressed with the preprocessor guard as they
41 ! are limited to this scope, and this module should only get included in compilation
42 ! if Wannier90 lib is linked to
43 use w90_library, only: lib_common_type, &
44 w90_set_option, &
45 w90_input_reader, &
46 w90_input_setopt, &
47 w90_set_comm
48
49 implicit none
50 private
51 public :: wannier90lib_init_w90main, &
54
55contains
56
62 subroutine wannier90lib_init_w90main(namespace, ions, kpoints, st, &
63 inp_options, w90main, stdout, stderr, ierr)
64 type(namespace_t), intent(in) :: namespace
65 type(ions_t), intent(in) :: ions
66 type(kpoints_t), intent(in) :: kpoints
67 type(states_elec_t), intent(in) :: st
68 type(wannier_opts_t), intent(in) :: inp_options
69 type(lib_common_type), intent(inout) :: w90main
70 integer, intent(inout) :: stdout, stderr, ierr
71
72#ifdef HAVE_MPI
73 integer, allocatable :: distk(:)
74 integer :: ik, iq, nk
75#endif
76
78
79 ! Set required options
80 ! Note: num_bands does not include excluded bands!
81 ! The order is nst > num_bands > num_wann
82 ! TODO: Correctly handle this depending on inputs
83 ! TODO: Distinguish everywhere between nst and num_bands
84 call w90_set_option(w90main, 'kpoints', -kpoints%reduced%red_point(1:3,:)) ! minus sign for wrong octopus convention
85 call w90_set_option(w90main, 'mp_grid', kpoints%nik_axis(1:3))
86 call w90_set_option(w90main, 'num_bands', inp_options%num_bands)
87 call w90_set_option(w90main, 'num_wann', inp_options%num_wann)
88 call w90_set_option(w90main, 'unit_cell_cart', units_from_atomic(unit_angstrom, ions%latt%rlattice(1:3,:)))
89
90#ifdef HAVE_MPI
91 nk = product(kpoints%nik_axis(1:3))
92 safe_allocate(distk(1:nk))
93 if (st%d%kpt%parallel) then
94 ! pass octopus communicator to wannier90 library
95 call w90_set_comm(w90main, st%d%kpt%mpi_grp%comm%MPI_VAL)
96 ! Feed kpoint distribution to wannier90 library
97 distk = -1
98 ! (k, spin)
99 do iq = 1, st%nik
100 ! k-only
101 ik = st%d%get_kpoint_index(iq)
102 if (distk(ik) == -1) then
103 distk(ik) = st%d%kpt%node(iq)
104 else if (distk(ik) /= st%d%kpt%node(iq)) then
105 message(1) = 'Error: Wannier90 k-point parallelization: same k-point' // &
106 'assigned to different nodes.'
107 call messages_fatal(1, namespace=namespace)
108 end if
109 end do
110 else
111 ! octopus default communicator is not acceptable in that case
112 call w90_set_comm(w90main, st%system_grp%comm%MPI_VAL)
113 ! No parallelization, assign all k-points to node 0
114 distk = 0
115 endif
116 call w90_set_option(w90main, 'distk', distk)
117 safe_deallocate_a(distk)
118#endif
119
120 ! Set atom data for projection parsing
121 call wannier90lib_set_atom_data(w90main, ions)
122
123 ! Parse projections from input file
124 call wannier90lib_set_projections(w90main, inp_options)
125
126 ! Option to dump files
127 if (inp_options%dump_inputs) then
128 call w90_set_option(w90main, "dump_inputs", .true.)
129 end if
130
131 ! Initialize library type with seedname
132 call w90_input_setopt(w90main, inp_options%prefix, stdout, stderr, ierr)
133
134 ! Read optional options from file
135 ! Note: Mandatory options are only read from Octopus's input.
136 ! If any mandatory options are set in `.win`, they will be ignored triggering a warning.
137 call wannier90lib_warn_inputs(inp_options)
138 call w90_input_reader(w90main, stdout, stderr, ierr)
139 if (ierr /= 0) then
140 write(message(1), '(a,i0)') 'Error in w90_input_reader: ', ierr
141 call messages_fatal(1)
142 end if
143
145 end subroutine wannier90lib_init_w90main
146
147
153 subroutine wannier90lib_set_atom_data(w90main, ions)
154 type(lib_common_type), intent(inout) :: w90main
155 type(ions_t), intent(in) :: ions
156
157 integer :: ia, n_atoms
158 real(real64), allocatable :: atom_pos_frac(:,:)
159 character(len=2), allocatable :: atom_labels(:)
160
162
163 n_atoms = ions%natoms
164
165 ! Allocate arrays for atom data
166 safe_allocate(atom_pos_frac(1:3, 1:n_atoms))
167 safe_allocate(atom_labels(1:n_atoms))
168
169 ! Extract atomic positions and labels from ions
170 do ia = 1, n_atoms
171 ! Get Cartesian positions from octopus
172 ! Convert to fractional (reduced) coordinates for wannier90
173 atom_pos_frac(1:3, ia) = ions%latt%cart_to_red(ions%pos(1:3, ia))
174 ! Get atomic labels
175 atom_labels(ia) = trim(ions%atom(ia)%label)
176 end do
177
178 ! Set atom positions in fractional coordinates
179 call w90_set_option(w90main, 'atoms_frac', atom_pos_frac)
180 ! Set symbols
181 call w90_set_option(w90main, 'symbols', atom_labels)
182
183 safe_deallocate_a(atom_pos_frac)
184 safe_deallocate_a(atom_labels)
185
187 end subroutine wannier90lib_set_atom_data
188
189
194 subroutine wannier90lib_set_projections(w90main, inp_options)
195 type(lib_common_type), intent(inout) :: w90main
196 type(wannier_opts_t), intent(in) :: inp_options
197
198 character(len=80) :: filename, line
199 integer :: w90_win, io
200 logical :: exist
201
203
204 ! Correctly handle projections by parsing them from .win file
205 filename = trim(adjustl(inp_options%prefix)) //'.win'
206
207 message(1) = "oct-wannier90: Parsing projections from "//filename
208 call messages_info(1)
209
210 inquire(file=filename,exist=exist)
211 if (.not. exist) then
212 message(1) = 'oct-wannier90: Cannot find specified Wannier90 win file.'
213 call messages_fatal(1)
214 end if
215
216 ! read projections and set them in w90main
217 w90_win = io_open(trim(filename), action='read')
218 do
219 read(w90_win, '(A)', iostat=io) line
220 if (is_iostat_end(io)) exit
221
222 if (trim(line) == 'begin projections') then
223 do
224 read(w90_win, '(A)') line
225 if (trim(line) == 'end projections') exit
226 call w90_set_option(w90main, 'projections', trim(line))
227 end do
228 end if
229 end do
230 call io_close(w90_win)
231
233 end subroutine wannier90lib_set_projections
234
240 subroutine wannier90lib_warn_inputs(inp_options)
241 type(wannier_opts_t), intent(in) :: inp_options
242
243 character(len=80) :: filename, line, dummy, dummy1
244 integer :: w90_win, io
245 logical :: exist
246
248
249 filename = trim(adjustl(inp_options%prefix)) //'.win'
250
251 message(1) = "oct-wannier90: Parsing "//filename
252 call messages_info(1)
253
254 inquire(file=filename,exist=exist)
255 if (.not. exist) then
256 message(1) = 'oct-wannier90: Cannot find specified Wannier90 win file.'
257 call messages_fatal(1)
258 end if
260 ! warn which options are not considered from .win file
261 w90_win = io_open(trim(filename), action='read')
262 do
263 read(w90_win, '(A)', iostat=io) line
264 if (is_iostat_end(io)) exit
265
266 if (trim(line) == 'begin kpoints') then
267 message(1) = 'oct-wannier90: Warning: kpoints specified in .win file' // &
268 ' are ignored. Using kpoints from octopus input instead.'
269 call messages_warning(1)
270 else if (trim(line) == 'begin atoms_frac') then
271 message(1) = 'oct-wannier90: Warning: atomic positions specified in .win file' // &
272 ' are ignored. Using positions from octopus input instead.'
273 call messages_warning(1)
274 else if (trim(line) == 'begin unit_cell_cart') then
275 message(1) = 'oct-wannier90: Warning: unit cell vectors specified in .win file' // &
276 ' are ignored. Using unit cell vectors from octopus input instead.'
277 call messages_warning(1)
278 end if
279 if (index(line, '=') <= 0) then
280 read(line, *, iostat=io) dummy, dummy1
281 if(dummy == 'mp_grid') then
282 message(1) = 'oct-wannier90: Warning: mp_grid specified in .win file' // &
283 ' is ignored. Using k-point mesh from octopus input instead.'
284 call messages_warning(1)
285 end if
286 end if
287 end do
288 call io_close(w90_win)
289
291 end subroutine wannier90lib_warn_inputs
292
296 subroutine wannier90lib_create_wannier90_amn(w90main, inp_options, exclude_list, &
297 band_index, space, mesh, latt, st, kpoints, projection)
298 type(lib_common_type), intent(in) :: w90main
299 type(wannier_opts_t), intent(in) :: inp_options
300 class(space_t), intent(in) :: space
301 class(mesh_t), intent(in) :: mesh
302 type(lattice_vectors_t), intent(in) :: latt
303 type(states_elec_t), intent(in) :: st
304 type(kpoints_t), intent(in) :: kpoints
305 logical, intent(in) :: exclude_list(:)
306 integer, intent(in) :: band_index(:)
307
308 complex(real64), contiguous, intent(out) :: projection(:, :, :)
309 type(wannier_calc_proj_t), allocatable :: proj_input_oct(:)
310
311 integer :: iw
312
314 call profiling_in("W90_AMN")
315
316 ! Note: in w90main, proj_input contains the projections from input file,
317 ! proj contains the projections that are actually used in wannier90
318 ! We use proj_input since that is also what is written to file in -pp mode
319
320
321 safe_allocate(proj_input_oct(1:w90main%num_proj))
322 do iw = 1, w90main%num_proj
323 proj_input_oct(iw)%l = w90main%proj_input(iw)%l
324 proj_input_oct(iw)%m = w90main%proj_input(iw)%m
325 proj_input_oct(iw)%s = w90main%proj_input(iw)%s
326 proj_input_oct(iw)%radial = w90main%proj_input(iw)%radial
327 proj_input_oct(iw)%site = w90main%proj_input(iw)%site
328 proj_input_oct(iw)%s_qaxis = w90main%proj_input(iw)%s_qaxis
329 proj_input_oct(iw)%z = w90main%proj_input(iw)%z
330 proj_input_oct(iw)%x = w90main%proj_input(iw)%x
331 proj_input_oct(iw)%zona = w90main%proj_input(iw)%zona
332 end do
333
334 call wannier_calc_create_amn(inp_options, exclude_list, &
335 band_index, space, mesh, latt, st, kpoints, &
336 w90main%wvfn_read%spin_channel, w90main%num_proj, &
337 proj_input_oct, projection)
338
339 safe_deallocate_a(proj_input_oct)
340
341 call profiling_out("W90_AMN")
342
344
346
350 subroutine wannier90lib_create_wannier90_mmn(w90main, exclude_list, band_index, mesh, st, ions, overlap)
351 type(lib_common_type), intent(in) :: w90main
352 logical, intent(in) :: exclude_list(:)
353 integer, intent(in) :: band_index(:)
354 class(mesh_t), intent(in) :: mesh
355 type(states_elec_t), target, intent(in) :: st
356 type(ions_t), intent(in) :: ions
357
358 complex(real64), contiguous, intent(out) :: overlap(:, :, :, :)
359
362 call wannier_calc_create_mmn(exclude_list, band_index, &
363 mesh, st, ions, w90main%wvfn_read%spin_channel, &
364 w90main%kmesh_info%nnlist, w90main%kmesh_info%nncell, &
365 overlap)
366
368
370
371
Definition: io.F90:116
subroutine, public io_close(iunit, grp)
Definition: io.F90:467
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:402
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
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_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
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
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.
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 wannier90lib_set_atom_data(w90main, ions)
Set atom data in wannier90 library.
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.
subroutine wannier90lib_set_projections(w90main, inp_options)
Parse projections.
subroutine wannier90lib_warn_inputs(inp_options)
Check for options in .win file.
Wannier90 related calculations.
subroutine, public wannier_calc_create_amn(wan_opts, exclude_list, band_index, space, mesh, latt, st, kpoints, spin_channel, num_proj, proj_input, projection)
Calculation of Wannier90 Projection Matrix.
subroutine, public wannier_calc_create_mmn(exclude_list, band_index, mesh, st, ions, spin_channel, nnlist, nncell, overlap)
Kohn-Sham State Overlap Matrix.
Wannier options module.
Describes mesh distribution to nodes.
Definition: mesh.F90:187
The states_elec_t class contains all electronic wave functions.
Mocks the projection type from wannier90.
Wannier related options.
batches of electronic states
Definition: wfs_elec.F90:141
int true(void)