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