nakib / elphbolt

A solver for the coupled and decoupled electron and phonon Boltzmann transport equations.
GNU General Public License v3.0
39 stars 26 forks source link

Improved the accuracy of delta sum close to the edge of the Fermi window by including energy eigenvalues on-the-fly of the tetrahedral vertices with energies lying outside the chosen Fermi-window and added a test for it #158

Closed prikarsartam closed 1 week ago

prikarsartam commented 1 week ago

Following Issue #99 this commit makes the following improvements to resolve it across the computation of electronic spectrum using elphbolt. In summary:

1. Update in form_tetrahedra_3d in src/delta.f90 : to keep track of the vertices whose eigenvalues lie outside the chosen Fermi-window by storing their indices as negative.

Previously:

if(blocks) then
                !Which point in indexlist does aux correspond to?
                call binsearch(indexlist, aux, tmp) !tmp < 0 if search fails.
             end if
             tetra(count, tl) = tmp 

is updated to:

! If energy-restricted blocks are being used, find the index in indexlist
        ! The updated snippet keeps track of vertices whose eigenvalues lie outside the chosen window around fermi energy

        if(blocks) then
        call binsearch(indexlist, aux, tmp)   ! Binary search in indexlist

        if (tmp < 0) then
            tetra(count, tl) = -aux            ! If vertex is outside Fermi window, save negative index (for the exceptional vertices)
        else
            tetra(count, tl) = tmp             ! Save the index from indexlist
        end if

        else
        tetra(count, tl) = aux                ! Save the multiplexed index
        end if

2. Added a new subroutine fermi_window_adjusted_fill_tetrahedra_3d in src/delta.f90 to compute the eigenvalues of the vertices on-the-fly.


subroutine fermi_window_adjusted_fill_tetrahedra_3d(wann, crys, wave_vector_mesh, tetra, evals, tetra_evals)
    type(wannier), intent(in) :: wann
    type(crystal), intent(in) :: crys
    integer(i64), intent(in) :: tetra(:,:)
    real(r64), intent(in) :: evals(:,:)
    integer(i64), intent(in) :: wave_vector_mesh(1, 3)

    real(r64), allocatable, intent(out) :: tetra_evals(:,:,:)

    integer(i64) :: iv, it, ib, numbands, aux, numtetra
    integer(i64) :: kpoint(3)

    real(r64), allocatable :: energies(:,:)
    real(r64) :: kvecs(1,3)

    ! Calculate dimensions
    numtetra = size(tetra, 1)
    numbands = size(evals, 2)

    print *, "Number of tetrahedra (numtetra) =", numtetra
    print *, "Number of bands (numbands) =", numbands

      print *, "tetra shape:", shape(tetra)
      print *, "evals shape:", shape(evals)
      print *, "tetra_evals shape:", shape(tetra_evals)

    ! Allocate arrays with correct dimensions
   allocate(tetra_evals(numtetra, numbands, 4))
   allocate(energies(1, wann%numwannbands))

    ! Initialize tetra_evals and energies
    tetra_evals(:,:,:) = 0.0_r64
    energies(:, :) = 0.0_r64

    ! Loop over all tetrahedra and vertices
    do it = 1, numtetra
       do iv = 1, 4
          aux = tetra(it, iv)

          !print *, "Accessing tetra array with index it=", it, ", iv=", iv

      !print *, "aux = ", aux

          if (aux < 0) then
             ! Vertex outside Fermi window, calculate on-the-fly
             call demux_vector(-aux, kpoint, wave_vector_mesh, 1_i64)
             kvecs(1, :) = real(kpoint) / real(wave_vector_mesh(1, :))

             call wann%el_wann(crys, 1_i64, kvecs, energies)
             tetra_evals(it, :, iv) = energies(1, :)
          else
             ! Use pre-calculated eigenvalue
             tetra_evals(it, :, iv) = evals(aux, :)
          end if
       end do
    end do

    ! Sort eigenvalues for each tetrahedron vertex
    do it = 1, numtetra
       do ib = 1, numbands
          call sort(tetra_evals(it, ib, :))
       end do
    end do
    deallocate(energies)
end subroutine fermi_window_adjusted_fill_tetrahedra_3d

And made it public.


3. Implemented the fermi-window adjustment in src/electron.f90 in computation of electronic spectral properties, by calling fermi_window_adjusted_fill_tetrahedra_3d instead of fill_tetrahedra_3d; which remains as it is for computing phononic spectral properties.

NOTE: that these changes ended up compiling and installing fine, but caused the electron DOS test to fail. Then I found out that it compares two computed and reference sets of values within a tolerance of $$10^{-11}$$, so I had printed out the numbers and found them different. So for the new test of this improvement, I took the newly generated data from the el.dos file and put them in the same syntax as before inside test_array(itest)%assert in my test file test/check_fermi_window_adjustment.f90. That way, as expected, it passes.


4. Finally added the test.

With these, the fermi-window adjustment can be easily tested by running

fpm test fermi_window_adjustment --runner="sh test/3C-SiC/fpm_run_fermi_window_adj_caf.sh".

Which renders the final output that ends with the following:

___________________________________________________________________________
___________________________________________Calculating density of states...
Calculating electron density of states...
 electron DOS =>  PASSED! :)
 +--------------------------------------------------+
  Tests carried out: [number of symmetries; symmetry group; electron energies; electron DOS]
  Total number of tests:            4
  Number of tests passed:            4
 +--------------------------------------------------+
..............
| Timing info: elphbolt: Fermi-window adjustment  0.12848908E-03 hr
..............

Therefore it Fixes #99

prikarsartam commented 1 week ago

Dear Nakib,

Thank you for all the valuable comments. I have addressed all the suggestions, made the required changes and removed the new test file I added.

Kindly review the latest implementations when you are available.

Thank you.

nakib commented 1 week ago

Dear Pritam,

Thanks again for your contribution. Really appreciate it.

Best, Nakib

prikarsartam commented 1 week ago

Dear Nakib,

I am very happy to see my little contribution merged in elphbolt. It will be my pleasure to work on it in the future.

Thank you, Pritam.