This routine should be called after wannier_setup from a code calling the library mode to actually run the Wannier code.
NOTE! The library mode currently works ONLY in serial. When called from an external code, wannier90 needs to be compiled in sequential and wannier_run called with 1 MPI process.
For more information, check a (minimal) example of how it can be used in the folder test-suite/library-mode-test/test_library.F90
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
character(len=*), | intent(in) | :: | seed__name | |||
integer, | intent(in), | dimension(3) | :: | mp_grid_loc | ||
integer, | intent(in) | :: | num_kpts_loc | |||
real(kind=dp), | intent(in), | dimension(3, 3) | :: | real_lattice_loc | ||
real(kind=dp), | intent(in), | dimension(3, 3) | :: | recip_lattice_loc | ||
real(kind=dp), | intent(in), | dimension(3, num_kpts_loc) | :: | kpt_latt_loc | ||
integer, | intent(in) | :: | num_bands_loc | |||
integer, | intent(in) | :: | num_wann_loc | |||
integer, | intent(in) | :: | nntot_loc | |||
integer, | intent(in) | :: | num_atoms_loc | |||
character(len=*), | intent(in), | dimension(num_atoms_loc) | :: | atom_symbols_loc | ||
real(kind=dp), | intent(in), | dimension(3, num_atoms_loc) | :: | atoms_cart_loc | ||
logical, | intent(in) | :: | gamma_only_loc | |||
complex(kind=dp), | intent(in), | dimension(num_bands_loc, num_bands_loc, nntot_loc, num_kpts_loc) | :: | M_matrix_loc | ||
complex(kind=dp), | intent(in), | dimension(num_bands_loc, num_wann_loc, num_kpts_loc) | :: | A_matrix_loc | ||
real(kind=dp), | intent(in), | dimension(num_bands_loc, num_kpts_loc) | :: | eigenvalues_loc | ||
complex(kind=dp), | intent(out), | dimension(num_wann_loc, num_wann_loc, num_kpts_loc) | :: | U_matrix_loc | ||
complex(kind=dp), | intent(out), | optional | dimension(num_bands_loc, num_wann_loc, num_kpts_loc) | :: | U_matrix_opt_loc | |
logical, | intent(out), | optional | dimension(num_bands_loc, num_kpts_loc) | :: | lwindow_loc | |
real(kind=dp), | intent(out), | optional | dimension(3, num_wann_loc) | :: | wann_centres_loc | |
real(kind=dp), | intent(out), | optional | dimension(num_wann_loc) | :: | wann_spreads_loc | |
real(kind=dp), | intent(out), | optional | dimension(3) | :: | spread_loc |
subroutine wannier_run(seed__name, mp_grid_loc, num_kpts_loc, &
real_lattice_loc, recip_lattice_loc, kpt_latt_loc, num_bands_loc, &
num_wann_loc, nntot_loc, num_atoms_loc, atom_symbols_loc, &
atoms_cart_loc, gamma_only_loc, M_matrix_loc, A_matrix_loc, eigenvalues_loc, &
U_matrix_loc, U_matrix_opt_loc, lwindow_loc, wann_centres_loc, &
wann_spreads_loc, spread_loc)
!! This routine should be called after wannier_setup from a code calling
!! the library mode to actually run the Wannier code.
!!
!! NOTE! The library mode currently works ONLY in serial.
!! When called from an external code, wannier90 needs to be compiled
!! in sequential and wannier_run called with 1 MPI process.
!!
!! For more information, check a (minimal) example of how it can be used
!! in the folder test-suite/library-mode-test/test_library.F90
use w90_constants
use w90_parameters
use w90_io
use w90_hamiltonian
use w90_kmesh
use w90_disentangle
use w90_overlap
use w90_wannierise
use w90_plot
use w90_transport
use w90_comms, only: my_node_id, num_nodes, &
comms_array_split, comms_scatterv
implicit none
character(len=*), intent(in) :: seed__name
integer, dimension(3), intent(in) :: mp_grid_loc
integer, intent(in) :: num_kpts_loc
real(kind=dp), dimension(3, 3), intent(in) :: real_lattice_loc
real(kind=dp), dimension(3, 3), intent(in) :: recip_lattice_loc
real(kind=dp), dimension(3, num_kpts_loc), intent(in) :: kpt_latt_loc
integer, intent(in) :: num_bands_loc
integer, intent(in) :: num_wann_loc
integer, intent(in) :: nntot_loc
integer, intent(in) :: num_atoms_loc
character(len=*), dimension(num_atoms_loc), intent(in) :: atom_symbols_loc
real(kind=dp), dimension(3, num_atoms_loc), intent(in) :: atoms_cart_loc
logical, intent(in) :: gamma_only_loc
complex(kind=dp), dimension(num_bands_loc, num_bands_loc, nntot_loc, num_kpts_loc), intent(in) :: M_matrix_loc
complex(kind=dp), dimension(num_bands_loc, num_wann_loc, num_kpts_loc), intent(in) :: A_matrix_loc
real(kind=dp), dimension(num_bands_loc, num_kpts_loc), intent(in) :: eigenvalues_loc
complex(kind=dp), dimension(num_wann_loc, num_wann_loc, num_kpts_loc), intent(out) :: U_matrix_loc
complex(kind=dp), dimension(num_bands_loc, num_wann_loc, num_kpts_loc), optional, intent(out) :: U_matrix_opt_loc
logical, dimension(num_bands_loc, num_kpts_loc), optional, intent(out) :: lwindow_loc
real(kind=dp), dimension(3, num_wann_loc), optional, intent(out) :: wann_centres_loc
real(kind=dp), dimension(num_wann_loc), optional, intent(out) :: wann_spreads_loc
real(kind=dp), dimension(3), optional, intent(out) :: spread_loc
real(kind=dp) time0, time1, time2
character(len=9) :: stat, pos, cdate, ctime
integer :: ierr, loop_k, loop_w
logical :: wout_found
integer :: nkp, nn, n, m
! Needed to split an array on different nodes
integer, dimension(0:num_nodes - 1) :: counts
integer, dimension(0:num_nodes - 1) :: displs
time0 = io_time()
library = .true.
! seedname="wannier"
seedname = trim(adjustl(seed__name))
inquire (file=trim(seedname)//'.wout', exist=wout_found)
if (wout_found) then
stat = 'old'
else
stat = 'replace'
endif
pos = 'append'
stdout = io_file_unit()
open (unit=stdout, file=trim(seedname)//'.wout', status=trim(stat), position=trim(pos))
call io_date(cdate, ctime)
write (stdout, '(/,2a,/)') ' Resuming Wannier90 at ', ctime
! call param_write_header
! copy local data into module variables
num_bands = num_bands_loc
mp_grid = mp_grid_loc
num_kpts = num_kpts_loc
real_lattice = real_lattice_loc
recip_lattice = recip_lattice_loc
allocate (kpt_latt(3, num_kpts), stat=ierr)
if (ierr /= 0) call io_error('Error allocating kpt_latt in wannier_setup')
kpt_latt = kpt_latt_loc
allocate (eigval(num_bands, num_kpts), stat=ierr)
if (ierr /= 0) call io_error('Error allocating eigval in wannier_setup')
eigval = eigenvalues_loc
num_atoms = num_atoms_loc
gamma_only = gamma_only_loc
call param_lib_set_atoms(atom_symbols_loc, atoms_cart_loc)
call param_read()
call param_write()
time1 = io_time()
write (stdout, '(1x,a25,f11.3,a)') 'Time to read parameters ', time1 - time0, ' (sec)'
call kmesh_get()
time2 = io_time()
write (stdout, '(1x,a25,f11.3,a)') 'Time to get kmesh ', time2 - time1, ' (sec)'
call comms_array_split(num_kpts, counts, displs)
call overlap_allocate()
if (disentanglement) then
m_matrix_orig = m_matrix_loc
a_matrix = a_matrix_loc
u_matrix_opt = cmplx_0
u_matrix = cmplx_0
else
m_matrix = m_matrix_loc
u_matrix = a_matrix_loc
endif
! IMPORTANT NOTE: _loc are variables local to this function, passed in as variables
! Instead, _local are variables local to the MPI process.
if (disentanglement) then
call comms_scatterv(m_matrix_orig_local, num_bands*num_bands*nntot*counts(my_node_id), &
m_matrix_orig, num_bands*num_bands*nntot*counts, num_bands*num_bands*nntot*displs)
else
call comms_scatterv(m_matrix_local, num_wann*num_wann*nntot*counts(my_node_id), &
m_matrix, num_wann*num_wann*nntot*counts, num_wann*num_wann*nntot*displs)
endif
!~ ! Check Mmn(k,b) is symmetric in m and n for gamma_only case
!~ if (gamma_only) call overlap_check_m_symmetry()
if (disentanglement) then
have_disentangled = .false.
call dis_main()
have_disentangled = .true.
call param_write_chkpt('postdis')
time1 = io_time()
write (stdout, '(1x,a25,f11.3,a)') 'Time to disentangle ', time1 - time2, ' (sec)'
else
if (gamma_only) then
call overlap_project_gamma()
else
call overlap_project()
endif
time1 = io_time()
write (stdout, '(1x,a25,f11.3,a)') 'Time to project overlaps ', time1 - time2, ' (sec)'
end if
if (gamma_only) then
call wann_main_gamma()
else
call wann_main()
endif
call param_write_chkpt('postwann')
time2 = io_time()
write (stdout, '(1x,a25,f11.3,a)') 'Time for wannierise ', time2 - time1, ' (sec)'
if (wannier_plot .or. bands_plot .or. fermi_surface_plot .or. write_hr) then
call plot_main()
time1 = io_time()
write (stdout, '(1x,a25,f11.3,a)') 'Time for plotting ', time1 - time2, ' (sec)'
end if
time2 = io_time()
if (transport) then
call tran_main()
time1 = io_time()
write (stdout, '(1x,a25,f11.3,a)') 'Time for transport ', time1 - time2, ' (sec)'
end if
! Now we zero all of the local output data, then copy in the data
! from the parameters module
u_matrix_loc = u_matrix
if (present(u_matrix_opt_loc) .and. present(lwindow_loc)) then
if (disentanglement) then
u_matrix_opt_loc = u_matrix_opt
lwindow_loc = lwindow
else
u_matrix_opt_loc = cmplx_0
do loop_k = 1, num_kpts
do loop_w = 1, num_wann
u_matrix_opt_loc(loop_w, loop_w, loop_k) = cmplx_1
end do
end do
lwindow_loc = .true.
end if
end if
if (present(wann_centres_loc)) wann_centres_loc = wannier_centres
if (present(wann_spreads_loc)) wann_spreads_loc = wannier_spreads
if (present(spread_loc)) then
spread_loc(1) = omega_total
spread_loc(2) = omega_invariant
spread_loc(3) = omega_tilde
endif
call hamiltonian_dealloc()
call overlap_dealloc()
call kmesh_dealloc()
call param_dealloc()
write (stdout, '(1x,a25,f11.3,a)') 'Total Execution Time ', io_time() - time0, ' (sec)'
if (timing_level > 0) call io_print_timings()
write (stdout, *)
write (stdout, '(1x,a)') 'All done: wannier90 exiting'
close (stdout)
end subroutine wannier_run