This routine should be called first from a code calling the library mode to setup all the variables. NOTE! The library mode currently works ONLY in serial (when called from a parallel code, make sure to run it only on 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_tot | |||
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 | |||
logical, | intent(in) | :: | spinors_loc | |||
integer, | intent(out) | :: | nntot_loc | |||
integer, | intent(out), | dimension(num_kpts_loc, num_nnmax) | :: | nnlist_loc | ||
integer, | intent(out), | dimension(3, num_kpts_loc, num_nnmax) | :: | nncell_loc | ||
integer, | intent(out) | :: | num_bands_loc | |||
integer, | intent(out) | :: | num_wann_loc | |||
real(kind=dp), | intent(out), | dimension(3, num_bands_tot) | :: | proj_site_loc | ||
integer, | intent(out), | dimension(num_bands_tot) | :: | proj_l_loc | ||
integer, | intent(out), | dimension(num_bands_tot) | :: | proj_m_loc | ||
integer, | intent(out), | dimension(num_bands_tot) | :: | proj_radial_loc | ||
real(kind=dp), | intent(out), | dimension(3, num_bands_tot) | :: | proj_z_loc | ||
real(kind=dp), | intent(out), | dimension(3, num_bands_tot) | :: | proj_x_loc | ||
real(kind=dp), | intent(out), | dimension(num_bands_tot) | :: | proj_zona_loc | ||
integer, | intent(out), | dimension(num_bands_tot) | :: | exclude_bands_loc | ||
integer, | intent(out), | optional | dimension(num_bands_tot) | :: | proj_s_loc | |
real(kind=dp), | intent(out), | optional | dimension(3, num_bands_tot) | :: | proj_s_qaxis_loc |
subroutine wannier_setup(seed__name, mp_grid_loc, num_kpts_loc, &
real_lattice_loc, recip_lattice_loc, kpt_latt_loc, num_bands_tot, &
num_atoms_loc, atom_symbols_loc, atoms_cart_loc, gamma_only_loc, spinors_loc, &
nntot_loc, nnlist_loc, nncell_loc, num_bands_loc, num_wann_loc, &
proj_site_loc, proj_l_loc, proj_m_loc, proj_radial_loc, proj_z_loc, &
proj_x_loc, proj_zona_loc, exclude_bands_loc, proj_s_loc, proj_s_qaxis_loc)
!! This routine should be called first from a code calling the library
!! mode to setup all the variables.
!! NOTE! The library mode currently works ONLY in serial (when called from
!! a parallel code, make sure to run it only on 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_kmesh
use w90_comms, only: comms_setup_vars
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_tot
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
logical, intent(in) :: spinors_loc
integer, intent(out) :: nntot_loc
integer, dimension(num_kpts_loc, num_nnmax), intent(out) :: nnlist_loc
integer, dimension(3, num_kpts_loc, num_nnmax), intent(out) :: nncell_loc
integer, intent(out) :: num_bands_loc
integer, intent(out) :: num_wann_loc
real(kind=dp), dimension(3, num_bands_tot), intent(out) :: proj_site_loc
integer, dimension(num_bands_tot), intent(out) :: proj_l_loc
integer, dimension(num_bands_tot), intent(out) :: proj_m_loc
integer, dimension(num_bands_tot), intent(out) :: proj_radial_loc
real(kind=dp), dimension(3, num_bands_tot), intent(out) :: proj_z_loc
real(kind=dp), dimension(3, num_bands_tot), intent(out) :: proj_x_loc
real(kind=dp), dimension(num_bands_tot), intent(out) :: proj_zona_loc
integer, dimension(num_bands_tot), intent(out) :: exclude_bands_loc
integer, dimension(num_bands_tot), optional, intent(out) :: proj_s_loc
real(kind=dp), dimension(3, num_bands_tot), optional, intent(out) :: proj_s_qaxis_loc
real(kind=dp) time0, time1, time2
character(len=9) :: stat, pos, cdate, ctime
integer :: ierr
logical :: wout_found
time0 = io_time()
call comms_setup_vars
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 param_write_header()
write (stdout, '(/a/)') ' Wannier90 is running in LIBRARY MODE'
write (stdout, '(a/)') ' Setting up k-point neighbours...'
! copy local data into module variables
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
num_atoms = num_atoms_loc
call param_lib_set_atoms(atom_symbols_loc, atoms_cart_loc)
gamma_only = gamma_only_loc
spinors = spinors_loc
! GP: at this point we don't know yet the number of excluded bands...
num_bands = num_bands_tot
library_param_read_first_pass = .true.
call param_read()
! Following calls will all NOT be first_pass, and I need to pass
! directly num_bands, that is already set internally now to num_bands = num_bands_tot - num_exclude_bands
library_param_read_first_pass = .false.
! set cell_volume as it is written to output in param_write
cell_volume = real_lattice(1, 1)*(real_lattice(2, 2)*real_lattice(3, 3) - real_lattice(3, 2)*real_lattice(2, 3)) + &
real_lattice(1, 2)*(real_lattice(2, 3)*real_lattice(3, 1) - real_lattice(3, 3)*real_lattice(2, 1)) + &
real_lattice(1, 3)*(real_lattice(2, 1)*real_lattice(3, 2) - real_lattice(3, 1)*real_lattice(2, 2))
call param_write()
time1 = io_time()
write (stdout, '(1x,a25,f11.3,a)') 'Time to read parameters ', time1 - time0, ' (sec)'
if (.not. explicit_nnkpts) call kmesh_get()
! Now we zero all of the local output data, then copy in the data
! from the parameters module
nntot_loc = 0
nnlist_loc = 0
nncell_loc = 0
proj_site_loc = 0.0_dp
proj_l_loc = 0
proj_m_loc = 0
proj_z_loc = 0.0_dp
proj_x_loc = 0.0_dp
proj_radial_loc = 0
proj_zona_loc = 0.0_dp
exclude_bands_loc = 0
nntot_loc = nntot
nnlist_loc(:, 1:nntot) = nnlist(:, 1:nntot)
nncell_loc(:, :, 1:nntot) = nncell(:, :, 1:nntot)
num_bands_loc = num_bands
num_wann_loc = num_wann
if (allocated(proj_site)) then
proj_site_loc(:, 1:num_proj) = proj_site(:, 1:num_proj)
proj_l_loc(1:num_proj) = proj_l(1:num_proj)
proj_m_loc(1:num_proj) = proj_m(1:num_proj)
proj_z_loc(:, 1:num_proj) = proj_z(:, 1:num_proj)
proj_x_loc(:, 1:num_proj) = proj_x(:, 1:num_proj)
proj_radial_loc(1:num_proj) = proj_radial(1:num_proj)
proj_zona_loc(1:num_proj) = proj_zona(1:num_proj)
if (allocated(proj_s) .and. present(proj_s_loc) .and. present(proj_s_qaxis_loc)) then
proj_s_loc(1:num_proj) = proj_s(1:num_proj)
proj_s_qaxis_loc(:, 1:num_proj) = proj_s_qaxis(:, 1:num_proj)
end if
endif
if (allocated(exclude_bands)) then
exclude_bands_loc(1:num_exclude_bands) = exclude_bands(1:num_exclude_bands)
end if
if (postproc_setup) then
call kmesh_write()
write (stdout, '(1x,a25,f11.3,a)') 'Time to write kmesh ', io_time(), ' (sec)'
write (stdout, '(/a)') ' '//trim(seedname)//'.nnkp written.'
endif
call kmesh_dealloc()
call param_dealloc()
write (stdout, '(1x,a25,f11.3,a)') 'Time to write kmesh ', io_time(), ' (sec)'
write (stdout, '(/a/)') ' Finished setting up k-point neighbours.'
call io_date(cdate, ctime)
write (stdout, '(2a)') ' Exiting wannier_setup in wannier90 ', ctime
close (stdout)
end subroutine wannier_setup