wannier_run Subroutine

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)

Uses

  • proc~~wannier_run~~UsesGraph proc~wannier_run wannier_run module~w90_transport w90_transport proc~wannier_run->module~w90_transport module~w90_disentangle w90_disentangle proc~wannier_run->module~w90_disentangle module~w90_constants w90_constants proc~wannier_run->module~w90_constants module~w90_overlap w90_overlap proc~wannier_run->module~w90_overlap module~w90_hamiltonian w90_hamiltonian proc~wannier_run->module~w90_hamiltonian module~w90_kmesh w90_kmesh proc~wannier_run->module~w90_kmesh module~w90_io w90_io proc~wannier_run->module~w90_io module~w90_wannierise w90_wannierise proc~wannier_run->module~w90_wannierise module~w90_plot w90_plot proc~wannier_run->module~w90_plot module~w90_comms w90_comms proc~wannier_run->module~w90_comms module~w90_parameters w90_parameters proc~wannier_run->module~w90_parameters module~w90_transport->module~w90_constants module~w90_disentangle->module~w90_constants module~w90_disentangle->module~w90_io module~w90_disentangle->module~w90_comms module~w90_disentangle->module~w90_parameters module~w90_sitesym w90_sitesym module~w90_disentangle->module~w90_sitesym module~w90_overlap->module~w90_constants module~w90_overlap->module~w90_io module~w90_overlap->module~w90_comms module~w90_overlap->module~w90_parameters module~w90_hamiltonian->module~w90_constants module~w90_hamiltonian->module~w90_comms module~w90_kmesh->module~w90_constants module~w90_kmesh->module~w90_comms module~w90_kmesh->module~w90_parameters module~w90_io->module~w90_constants module~w90_wannierise->module~w90_constants module~w90_wannierise->module~w90_comms module~w90_comms->module~w90_constants module~w90_comms->module~w90_io module~w90_parameters->module~w90_constants module~w90_parameters->module~w90_io module~w90_sitesym->module~w90_constants module~w90_sitesym->module~w90_io

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

Arguments

Type IntentOptional AttributesName
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

Calls

proc~~wannier_run~~CallsGraph proc~wannier_run wannier_run proc~hamiltonian_dealloc hamiltonian_dealloc proc~wannier_run->proc~hamiltonian_dealloc proc~overlap_project_gamma overlap_project_gamma proc~wannier_run->proc~overlap_project_gamma proc~io_time io_time proc~wannier_run->proc~io_time proc~overlap_dealloc overlap_dealloc proc~wannier_run->proc~overlap_dealloc proc~overlap_allocate overlap_allocate proc~wannier_run->proc~overlap_allocate proc~param_write param_write proc~wannier_run->proc~param_write proc~io_date io_date proc~wannier_run->proc~io_date proc~param_lib_set_atoms param_lib_set_atoms proc~wannier_run->proc~param_lib_set_atoms proc~dis_main dis_main proc~wannier_run->proc~dis_main proc~wann_main_gamma wann_main_gamma proc~wannier_run->proc~wann_main_gamma proc~wann_main wann_main proc~wannier_run->proc~wann_main proc~plot_main plot_main proc~wannier_run->proc~plot_main proc~param_write_chkpt param_write_chkpt proc~wannier_run->proc~param_write_chkpt proc~param_read param_read proc~wannier_run->proc~param_read proc~kmesh_get kmesh_get proc~wannier_run->proc~kmesh_get proc~param_dealloc param_dealloc proc~wannier_run->proc~param_dealloc proc~tran_main tran_main proc~wannier_run->proc~tran_main proc~comms_array_split comms_array_split proc~wannier_run->proc~comms_array_split proc~overlap_project overlap_project proc~wannier_run->proc~overlap_project proc~kmesh_dealloc kmesh_dealloc proc~wannier_run->proc~kmesh_dealloc proc~io_file_unit io_file_unit proc~wannier_run->proc~io_file_unit proc~io_error io_error proc~overlap_project_gamma->proc~io_error dgemm dgemm proc~overlap_project_gamma->dgemm dgesvd dgesvd proc~overlap_project_gamma->dgesvd proc~overlap_allocate->proc~comms_array_split proc~utility_lowercase utility_lowercase proc~param_lib_set_atoms->proc~utility_lowercase proc~dis_main->proc~comms_array_split proc~dis_main->proc~io_file_unit proc~dis_windows dis_windows proc~dis_main->proc~dis_windows proc~dis_extract_gamma dis_extract_gamma proc~dis_main->proc~dis_extract_gamma proc~dis_extract dis_extract proc~dis_main->proc~dis_extract proc~dis_proj_froz dis_proj_froz proc~dis_main->proc~dis_proj_froz proc~dis_project dis_project proc~dis_main->proc~dis_project proc~wann_main_gamma->proc~io_time proc~wann_main_gamma->proc~param_write_chkpt proc~wann_phases wann_phases proc~wann_main_gamma->proc~wann_phases proc~wann_spread_copy wann_spread_copy proc~wann_main_gamma->proc~wann_spread_copy proc~utility_zgemm utility_zgemm proc~wann_main_gamma->proc~utility_zgemm proc~wann_check_unitarity wann_check_unitarity proc~wann_main_gamma->proc~wann_check_unitarity proc~wann_omega_gamma wann_omega_gamma proc~wann_main_gamma->proc~wann_omega_gamma proc~wann_main->proc~comms_array_split proc~wann_main->proc~io_file_unit proc~hamiltonian_get_hr hamiltonian_get_hr proc~wann_main->proc~hamiltonian_get_hr proc~wann_main->proc~io_error proc~hamiltonian_setup hamiltonian_setup proc~wann_main->proc~hamiltonian_setup proc~io_wallclocktime io_wallclocktime proc~wann_main->proc~io_wallclocktime proc~wann_domega wann_domega proc~wann_main->proc~wann_domega proc~wann_main->proc~wann_phases proc~wann_main->proc~wann_spread_copy proc~wann_omega wann_omega proc~wann_main->proc~wann_omega proc~wann_main->proc~wann_check_unitarity proc~plot_main->proc~hamiltonian_get_hr proc~plot_main->proc~hamiltonian_setup proc~ws_write_vec ws_write_vec proc~plot_main->proc~ws_write_vec proc~param_write_chkpt->proc~io_date proc~param_write_chkpt->proc~io_file_unit proc~param_read->proc~io_file_unit proc~param_in_file param_in_file proc~param_read->proc~param_in_file proc~param_get_range_vector param_get_range_vector proc~param_read->proc~param_get_range_vector proc~param_get_keyword_block param_get_keyword_block proc~param_read->proc~param_get_keyword_block proc~param_get_centre_constraints param_get_centre_constraints proc~param_read->proc~param_get_centre_constraints proc~param_read->proc~io_error proc~get_smearing_index get_smearing_index proc~param_read->proc~get_smearing_index proc~param_get_projections param_get_projections proc~param_read->proc~param_get_projections proc~get_module_kmesh get_module_kmesh proc~param_read->proc~get_module_kmesh proc~utility_metric utility_metric proc~param_read->proc~utility_metric proc~param_get_block_length param_get_block_length proc~param_read->proc~param_get_block_length proc~param_get_keyword_vector param_get_keyword_vector proc~param_read->proc~param_get_keyword_vector proc~param_get_keyword param_get_keyword proc~param_read->proc~param_get_keyword proc~param_get_keyword_kpath param_get_keyword_kpath proc~param_read->proc~param_get_keyword_kpath proc~param_uppercase param_uppercase proc~param_read->proc~param_uppercase proc~internal_set_kmesh internal_set_kmesh proc~param_read->proc~internal_set_kmesh proc~param_get_atoms param_get_atoms proc~param_read->proc~param_get_atoms proc~param_get_vector_length param_get_vector_length proc~param_read->proc~param_get_vector_length proc~kmesh_get->proc~io_error proc~kmesh_supercell_sort kmesh_supercell_sort proc~kmesh_get->proc~kmesh_supercell_sort proc~kmesh_shell_fixed kmesh_shell_fixed proc~kmesh_get->proc~kmesh_shell_fixed proc~kmesh_shell_automatic kmesh_shell_automatic proc~kmesh_get->proc~kmesh_shell_automatic proc~kmesh_shell_from_file kmesh_shell_from_file proc~kmesh_get->proc~kmesh_shell_from_file proc~tran_cut_hr_one_dim tran_cut_hr_one_dim proc~tran_main->proc~tran_cut_hr_one_dim proc~tran_lcr_2c2_sort tran_lcr_2c2_sort proc~tran_main->proc~tran_lcr_2c2_sort proc~tran_main->proc~hamiltonian_get_hr proc~tran_bulk tran_bulk proc~tran_main->proc~tran_bulk proc~tran_main->proc~hamiltonian_setup proc~tran_lcr_2c2_build_ham tran_lcr_2c2_build_ham proc~tran_main->proc~tran_lcr_2c2_build_ham proc~tran_lcr tran_lcr proc~tran_main->proc~tran_lcr proc~tran_find_integral_signatures tran_find_integral_signatures proc~tran_main->proc~tran_find_integral_signatures proc~tran_parity_enforce tran_parity_enforce proc~tran_main->proc~tran_parity_enforce proc~tran_reduce_hr tran_reduce_hr proc~tran_main->proc~tran_reduce_hr proc~tran_get_ht tran_get_ht proc~tran_main->proc~tran_get_ht proc~overlap_project->proc~comms_array_split proc~overlap_project->proc~io_error proc~tran_lcr_2c2_sort->proc~tran_cut_hr_one_dim proc~tran_lcr_2c2_sort->proc~io_error proc~tran_lcr_2c2_sort->proc~tran_reduce_hr proc~check_and_sort_similar_centres check_and_sort_similar_centres proc~tran_lcr_2c2_sort->proc~check_and_sort_similar_centres proc~group group proc~tran_lcr_2c2_sort->proc~group proc~sort sort proc~tran_lcr_2c2_sort->proc~sort proc~master_sort_and_group master_sort_and_group proc~tran_lcr_2c2_sort->proc~master_sort_and_group proc~param_in_file->proc~io_file_unit proc~param_in_file->proc~utility_lowercase proc~dis_windows->proc~io_error 10 10 proc~dis_windows->10 proc~dis_extract_gamma->proc~io_time proc~dis_extract_gamma->proc~io_error dspevx dspevx proc~dis_extract_gamma->dspevx proc~param_get_keyword_block->proc~io_error proc~param_get_centre_constraint_from_column param_get_centre_constraint_from_column proc~param_get_centre_constraints->proc~param_get_centre_constraint_from_column proc~tran_bulk->proc~io_date proc~tran_bulk->proc~io_file_unit proc~tran_read_htx tran_read_htX proc~tran_bulk->proc~tran_read_htx proc~tran_green tran_green proc~tran_bulk->proc~tran_green zgemm zgemm proc~tran_bulk->zgemm proc~tran_transfer tran_transfer proc~tran_bulk->proc~tran_transfer proc~hamiltonian_wigner_seitz hamiltonian_wigner_seitz proc~hamiltonian_setup->proc~hamiltonian_wigner_seitz proc~sitesym_symmetrize_gradient sitesym_symmetrize_gradient proc~wann_domega->proc~sitesym_symmetrize_gradient proc~utility_inv3 utility_inv3 proc~wann_phases->proc~utility_inv3 proc~internal_maxloc internal_maxloc proc~kmesh_supercell_sort->proc~internal_maxloc proc~param_get_projections->proc~io_error proc~utility_strip utility_strip proc~param_get_projections->proc~utility_strip proc~utility_string_to_coord utility_string_to_coord proc~param_get_projections->proc~utility_string_to_coord proc~tran_lcr_2c2_build_ham->proc~io_date proc~tran_lcr_2c2_build_ham->proc~io_file_unit proc~tran_lcr_2c2_build_ham->proc~tran_cut_hr_one_dim proc~tran_lcr_2c2_build_ham->proc~io_error proc~tran_lcr_2c2_build_ham->proc~tran_reduce_hr proc~ws_write_vec->proc~io_date proc~ws_write_vec->proc~io_file_unit proc~tran_lcr->proc~io_date proc~tran_lcr->proc~io_file_unit proc~tran_lcr->proc~io_error proc~tran_read_htxy tran_read_htXY proc~tran_lcr->proc~tran_read_htxy proc~tran_lcr->proc~tran_read_htx zgbsv zgbsv proc~tran_lcr->zgbsv proc~tran_lcr->proc~tran_green proc~tran_lcr->zgemm proc~tran_lcr->proc~tran_transfer proc~tran_read_htc tran_read_htC proc~tran_lcr->proc~tran_read_htc proc~get_module_kmesh->proc~internal_set_kmesh proc~kmesh_shell_fixed->proc~io_error proc~kmesh_shell_fixed->dgesvd proc~utility_zgemm->zgemm proc~kmesh_shell_automatic->proc~io_error proc~kmesh_shell_automatic->dgesvd proc~tran_find_integral_signatures->proc~io_file_unit proc~tran_find_integral_signatures->proc~io_error proc~param_get_keyword_vector->proc~io_error proc~dis_extract->proc~comms_array_split proc~dis_extract->proc~io_error proc~dis_extract->proc~io_wallclocktime interface~comms_allreduce comms_allreduce proc~dis_extract->interface~comms_allreduce zhpevx zhpevx proc~dis_extract->zhpevx proc~sitesym_symmetrize_zmatrix sitesym_symmetrize_zmatrix proc~dis_extract->proc~sitesym_symmetrize_zmatrix proc~dis_proj_froz->proc~io_error proc~tran_reduce_hr->proc~io_error proc~wann_omega->interface~comms_allreduce proc~param_get_atoms->proc~io_error proc~param_get_atoms->proc~param_get_block_length proc~dis_project->proc~io_error proc~kmesh_shell_from_file->proc~io_file_unit proc~kmesh_shell_from_file->proc~io_error proc~kmesh_shell_from_file->dgesvd proc~tran_get_ht->proc~io_date proc~tran_get_ht->proc~io_file_unit proc~wann_check_unitarity->proc~io_error proc~tran_read_htxy->proc~io_file_unit proc~tran_read_htx->proc~io_file_unit proc~check_and_sort_similar_centres->proc~io_error proc~comms_allreduce_real comms_allreduce_real interface~comms_allreduce->proc~comms_allreduce_real proc~comms_allreduce_cmplx comms_allreduce_cmplx interface~comms_allreduce->proc~comms_allreduce_cmplx proc~hamiltonian_wigner_seitz->proc~io_error proc~tran_green->proc~io_error proc~tran_green->zgemm zgesv zgesv proc~tran_green->zgesv proc~tran_transfer->proc~io_error zcopy zcopy proc~tran_transfer->zcopy zaxpy zaxpy proc~tran_transfer->zaxpy proc~tran_transfer->zgesv proc~master_sort_and_group->proc~group proc~master_sort_and_group->proc~sort proc~tran_read_htc->proc~io_file_unit

Contents

Source Code


Source Code

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