13 use,
intrinsic :: iso_fortran_env
44 use w90_library_extra,
only:overlaps
48 type(namespace_t) :: namespace
49 type(electrons_t),
pointer :: sys
50 type(wannier_opts_t) :: wannier_opts
51 integer,
allocatable :: band_index(:)
53 type(lib_common_type) :: w90main
55 complex(real64),
allocatable :: m_matrix(:, :, :, :)
56 complex(real64),
allocatable :: u_matrix(:, :, :)
57 complex(real64),
allocatable :: u_matrix_opt(:, :, :)
58 real(real64),
allocatable :: eigenvals(:,:)
60 integer,
allocatable :: nnkp(:,:), gkpb(:,:,:)
61 logical,
allocatable :: exclude_list(:)
64 integer :: idir, ii, itemp, ik, ist
67 real(real64) :: factor(3)
83 if (sys%space%dim /= 3)
then
88 if (sys%kpoints%use_symmetries)
then
89 message(1) =
'oct-wannier90: k-points symmetries are not allowed'
92 if (sys%kpoints%use_time_reversal)
then
93 message(1) =
'oct-wannier90: time-reversal symmetry is not allowed'
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'
100 if (sys%st%d%nspin == 4)
then
107 do idir = sys%space%periodic_dim+1, sys%space%dim
108 factor(idir) =
m_two * sys%gr%box%bounding_box_l(idir)
110 call sys%ions%latt%scale(factor)
113 call wannier_opts%parse_oct(namespace)
114 call wannier_opts%parse_win()
118 wannier_opts, w90main, stdout, stderr, ierr)
120 write(
message(1),
'(a,i0)')
'Error initializing wannier90 library: ', ierr
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"
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.
141 if (
allocated(w90main%exclude_bands))
then
142 do itemp = 1,
size(w90main%exclude_bands)
143 exclude_list(w90main%exclude_bands(itemp)) = .
true.
148 do ii = 1, sys%st%nst
149 if (exclude_list(ii)) cycle
151 band_index(ii) = itemp
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)
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)
171 do ik = 1, w90main%num_kpts
172 do ist = 1, sys%st%nst
173 if (exclude_list(ist)) cycle
177 eigenvals(band_index(ist), ik) =
units_from_atomic(
unit_ev, sys%st%eigenval(ist, (ik-1)*2+w90main%wvfn_read%spin_channel))
181 call w90_set_eigval(w90main, eigenvals)
184 if (wannier_opts%calc_overlaps)
then
187 sys%gr, sys%ions%latt, sys%st, sys%kpoints, u_matrix_opt)
192 call overlaps(w90main, stdout, stderr, ierr)
194 write(
message(1),
'(a,i0)')
'Error while reading overlaps with wannier90: ', ierr
200 if (wannier_opts%wannierize)
then
201 call w90_disentangle(w90main, stdout, stderr, ierr)
203 write(
message(1),
'(a,i0)')
'Error in w90_disentangle: ', ierr
206 call w90_project_overlap(w90main, stdout, stderr, ierr)
210 write(
message(1),
'(a,i0)')
'Error in w90_project_overlap: ', ierr
213 call w90_wannierise(w90main, stdout, stderr, ierr)
215 write(
message(1),
'(a,i0)')
'Error in w90_wannierise: ', ierr
219 w90main%w90_calculation%wannier_plot = .false.
221 call w90_plot(w90main, stdout, stderr, ierr)
223 write(
message(1),
'(a,i0)')
'Error in w90_plot: ', ierr
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)
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 ...
subroutine, public fft_all_init(namespace)
initialize the table
subroutine, public fft_all_end()
delete all plans
real(real64), parameter, public m_two
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
subroutine, public global_init(communicator)
Initialise Octopus.
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
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...
subroutine, public io_end()
subroutine, public messages_end()
subroutine, public messages_not_implemented(feature, namespace)
subroutine, public messages_init(output_dir)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_experimental(name, namespace)
type(mpi_grp_t), public mpi_world
This module handles the communicators for the various parallelization strategies.
type(namespace_t), public global_namespace
subroutine, public parser_init()
Initialise the Octopus parser.
subroutine, public parser_end()
End the Octopus parser.
subroutine, public profiling_end(namespace)
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
integer, parameter, public restart_gs
integer, parameter, public restart_td
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.
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.
Class describing the electron system.
program wannier
Construct maximally-localised Wannier functions through API calls to Wannier90 library.