Skip to content

Commit

Permalink
add test for write_rmn
Browse files Browse the repository at this point in the history
  • Loading branch information
Jerome Jackson authored and JeromeCCP9 committed Jan 14, 2025
1 parent bbffae2 commit b43628d
Show file tree
Hide file tree
Showing 8 changed files with 11,381 additions and 35 deletions.
76 changes: 42 additions & 34 deletions src/plot.F90
Original file line number Diff line number Diff line change
Expand Up @@ -129,7 +129,6 @@ subroutine plot_main(atom_data, band_plot, dis_manifold, fermi_energy_list, ferm

if (my_node_id == 0) on_root = .true.

! ! write extra info regarding omega_invariant
if (output_file%svd_omega) then
call plot_svd_omega_i(num_wann, num_kpts, kmesh_info, m_matrix, print_output, timer, &
dist_k, error, comm, stdout)
Expand All @@ -138,59 +137,65 @@ subroutine plot_main(atom_data, band_plot, dis_manifold, fermi_energy_list, ferm

call utility_recip_lattice_base(real_lattice, recip_lattice, volume)

if (w90_calculation%bands_plot .or. w90_calculation%fermi_surface_plot .or. &
output_file%write_hr .or. output_file%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
! setup RS cluster for calculations that need it
if (output_file%write_hr .or. &
output_file%write_r2mn .or. &
output_file%write_rmn .or. &
output_file%write_tb .or. &
w90_calculation%bands_plot .or. &
w90_calculation%fermi_surface_plot) then

call hamiltonian_setup(ham_logical, print_output, ws_region, w90_calculation, ham_k, ham_r, &
real_lattice, wannier_centres_translated, irvec, mp_grid, ndegen, &
num_kpts, num_wann, nrpts, rpt_origin, band_plot%mode, stdout, &
timer, error, transport_mode, comm)
if (allocated(error)) return
endif

! setup RS Hamilton eqn for calculations that need it
if (output_file%write_hr .or. &
output_file%write_tb .or. &
w90_calculation%bands_plot .or. &
w90_calculation%fermi_surface_plot) 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) then
write (stdout, '(1x,a)') '!!!! Kpoint grid does not include Gamma. '// &
' Interpolation may be incorrect. !!!!'
endif

! Transform Hamiltonian to WF basis
call hamiltonian_get_hr(atom_data, dis_manifold, ham_logical, real_space_ham, print_output, &
ham_k, ham_r, u_matrix, u_matrix_opt, eigval, kpt_latt, &
real_lattice, wannier_data%centres, wannier_centres_translated, &
irvec, shift_vec, nrpts, num_bands, num_kpts, num_wann, &
have_disentangled, stdout, timer, error, lsitesymmetry, comm)
if (allocated(error)) return

bands_num_spec_points = 0

if (allocated(kpoint_path%labels)) bands_num_spec_points = size(kpoint_path%labels)
endif

if (on_root) then
if (print_output%timing_level > 0) call io_stopwatch_start('plot: main', timer)

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

! if (w90_calculation%bands_plot) then
! call plot_interpolate_bands(mp_grid, real_lattice, band_plot, kpoint_path, &
! real_space_ham, ws_region, print_output, recip_lattice, &
! num_wann, wannier_data, ham_r, irvec, ndegen, nrpts, &
! wannier_centres_translated, ws_distance, &
! bands_num_spec_points, stdout, seedname, timer, error, &
! comm)
! if (allocated(error)) return
! endif

if (w90_calculation%fermi_surface_plot) then
call plot_fermi_surface(fermi_energy_list, recip_lattice, fermi_surface_plot, num_wann, &
ham_r, irvec, ndegen, nrpts, print_output%timing_level, stdout, &
Expand All @@ -211,7 +216,7 @@ subroutine plot_main(atom_data, band_plot, dis_manifold, fermi_energy_list, ferm
if (allocated(error)) return
endif

if (output_file%write_hr .or. output_file%write_rmn .or. output_file%write_tb) then
if (output_file%write_hr .or. output_file%write_tb) then
if (.not. ws_distance%done) then
call ws_translate_dist(ws_distance, ws_region, num_wann, &
wannier_data%centres, real_lattice, mp_grid, nrpts, irvec, &
Expand Down Expand Up @@ -286,6 +291,9 @@ subroutine plot_main(atom_data, band_plot, dis_manifold, fermi_energy_list, ferm
endif

if (w90_calculation%bands_plot) then
bands_num_spec_points = 0
if (allocated(kpoint_path%labels)) bands_num_spec_points = size(kpoint_path%labels)

call plot_interpolate_bands(mp_grid, real_lattice, band_plot, kpoint_path, &
real_space_ham, ws_region, print_output, recip_lattice, &
num_wann, wannier_data, ham_r, irvec, ndegen, nrpts, &
Expand Down Expand Up @@ -2401,11 +2409,11 @@ subroutine plot_write_rmn(kmesh_info, m_matrix, kpt_latt, irvec, nrpts, num_kpts
type(w90_error_type), allocatable, intent(out) :: error
type(w90_comm_type), intent(in) :: comm

integer, intent(inout) :: nrpts
integer, intent(inout) :: irvec(:, :)
integer, intent(in) :: num_wann
integer, intent(in) :: num_kpts
integer, intent(in) :: dist_k(:) ! MPI k-point distribution
integer, intent(in) :: nrpts
integer, intent(in) :: irvec(:, :)
integer, intent(in) :: num_wann
integer, intent(in) :: num_kpts
integer, intent(in) :: dist_k(:) ! MPI k-point distribution
real(kind=dp), intent(in) :: kpt_latt(:, :)
complex(kind=dp), intent(in) :: m_matrix(:, :, :, :)
character(len=50), intent(in) :: seedname
Expand Down
2 changes: 1 addition & 1 deletion src/ws_distance.F90
Original file line number Diff line number Diff line change
Expand Up @@ -80,7 +80,7 @@ subroutine ws_translate_dist(ws_distance, ws_region, num_wann, wannier_centres,
integer, intent(in) :: mp_grid(3)
integer, intent(in) :: num_wann
integer, intent(in) :: nrpts
integer, intent(in) :: irvec(3, nrpts)
integer, intent(in) :: irvec(:, :)

real(kind=dp), intent(in) :: real_lattice(3, 3)
real(kind=dp), intent(in) :: wannier_centres(:, :)
Expand Down
6 changes: 6 additions & 0 deletions test-suite/tests/jobconfig
Original file line number Diff line number Diff line change
Expand Up @@ -465,6 +465,12 @@ program = WANNIER90_DISENTANGLEMENT_SAWFS
inputs_args = ('H3S.win', '')
output = H3S.wout

# Diamond, valence states - test of the R_mn printout
[testw90_rmn]
program = WANNIER90_RMN_OK
inputs_args = ('diamond.win', '')
output = diamond_r.dat

# Test of k.p coefficients for GaAs
[testpostw90_gaas_kdotp/]
program = POSTW90_KDOTP_OK
Expand Down
Loading

0 comments on commit b43628d

Please sign in to comment.