plot_main Subroutine

public subroutine plot_main()

Uses

  • proc~~plot_main~~UsesGraph proc~plot_main plot_main module~w90_ws_distance w90_ws_distance proc~plot_main->module~w90_ws_distance module~w90_constants w90_constants proc~plot_main->module~w90_constants module~w90_io w90_io proc~plot_main->module~w90_io module~w90_parameters w90_parameters proc~plot_main->module~w90_parameters module~w90_hamiltonian w90_hamiltonian proc~plot_main->module~w90_hamiltonian module~w90_ws_distance->module~w90_constants module~w90_ws_distance->module~w90_parameters module~w90_io->module~w90_constants module~w90_parameters->module~w90_constants module~w90_parameters->module~w90_io module~w90_hamiltonian->module~w90_constants module~w90_comms w90_comms module~w90_hamiltonian->module~w90_comms module~w90_comms->module~w90_constants module~w90_comms->module~w90_io

Main plotting routine

Arguments

None

Calls

proc~~plot_main~~CallsGraph proc~plot_main plot_main proc~hamiltonian_setup hamiltonian_setup proc~plot_main->proc~hamiltonian_setup proc~ws_write_vec ws_write_vec proc~plot_main->proc~ws_write_vec proc~hamiltonian_get_hr hamiltonian_get_hr proc~plot_main->proc~hamiltonian_get_hr proc~hamiltonian_wigner_seitz hamiltonian_wigner_seitz proc~hamiltonian_setup->proc~hamiltonian_wigner_seitz proc~io_date io_date proc~ws_write_vec->proc~io_date proc~io_file_unit io_file_unit proc~ws_write_vec->proc~io_file_unit proc~io_error io_error proc~hamiltonian_wigner_seitz->proc~io_error

Called by

proc~~plot_main~~CalledByGraph proc~plot_main plot_main program~wannier wannier program~wannier->proc~plot_main proc~wannier_run wannier_run proc~wannier_run->proc~plot_main

Contents

Source Code


Source Code

  subroutine plot_main()
    !! Main plotting routine
    !============================================!

    use w90_constants, only: eps6
    use w90_io, only: stdout, io_stopwatch
    use w90_parameters, only: num_kpts, bands_plot, dos_plot, &
      kpt_latt, fermi_surface_plot, &
      wannier_plot, timing_level, write_bvec, &
      write_hr, write_rmn, write_tb, write_u_matrices
    use w90_hamiltonian, only: hamiltonian_get_hr, hamiltonian_write_hr, &
      hamiltonian_setup, hamiltonian_write_rmn, &
      hamiltonian_write_tb, nrpts, irvec
    use w90_ws_distance, only: done_ws_distance, ws_translate_dist, &
      ws_write_vec

    implicit none

    integer :: nkp
    logical :: have_gamma

    if (timing_level > 0) call io_stopwatch('plot: main', 1)

    ! Print the header only if there is something to plot
    if (bands_plot .or. dos_plot .or. fermi_surface_plot .or. write_hr .or. &
        wannier_plot .or. write_u_matrices .or. write_tb) then
      write (stdout, '(1x,a)') '*---------------------------------------------------------------------------*'
      write (stdout, '(1x,a)') '|                               PLOTTING                                    |'
      write (stdout, '(1x,a)') '*---------------------------------------------------------------------------*'
      write (stdout, *)
    end if

    if (bands_plot .or. dos_plot .or. fermi_surface_plot .or. write_hr .or. &
        write_tb) then
      ! Check if the kmesh includes the gamma point
      have_gamma = .false.
      do nkp = 1, num_kpts
        if (all(abs(kpt_latt(:, nkp)) < eps6)) have_gamma = .true.
      end do
      if (.not. have_gamma) &
           write (stdout, '(1x,a)') '!!!! Kpoint grid does not include Gamma. '// &
           & ' Interpolation may be incorrect. !!!!'
      ! Transform Hamiltonian to WF basis
      !
      call hamiltonian_setup()
      !
      call hamiltonian_get_hr()
      !
      if (bands_plot) call plot_interpolate_bands
      !
      if (fermi_surface_plot) call plot_fermi_surface
      !
      if (write_hr) call hamiltonian_write_hr()
      !
      if (write_rmn) call hamiltonian_write_rmn()
      !
      if (write_tb) call hamiltonian_write_tb()
      !
      if (write_hr .or. write_rmn .or. write_tb) then
        if (.not. done_ws_distance) call ws_translate_dist(nrpts, irvec)
        call ws_write_vec(nrpts, irvec)
      end if
    end if

    if (wannier_plot) call plot_wannier

    if (write_bvec) call plot_bvec

    if (write_u_matrices) call plot_u_matrices

    if (timing_level > 0) call io_stopwatch('plot: main', 2)

  end subroutine plot_main