Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

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

Merged
merged 5 commits into from
Nov 1, 2024

Conversation

prikarsartam
Copy link
Contributor

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

  • Added a way to keep track of the exceptional vertices during analytical tetrahedron method.
  • Added a new subroutine to compute the eigenvalues of the tetrahedral vertices on-the-fly.
  • Updated Fermi-window adjustment across the computation of electronic spectral properties.
  • Added a test to compute the electronic density of states with Fermi-window adjusted analytical tetrahedron method.
    (fpm test fermi_window_adjustment --runner="sh test/3C-SiC/fpm_run_fermi_window_adj_caf.sh")

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.

  • Added test/check_fermi_window_adjustment.f90 similar to bte_regression.f90, except with minimum dependencies to test the computation of electronic density of states.
  • Added test/3C-SiC/fpm_run_fermi_window_adj_caf.sh similar to the runner of bte_regression test which executes cafrun -n 2 ../build/caf_*/test/fermi_window_adjustment | tee 3C_SiC_test.output.
  • Updated fpm.toml with an addition of
[[test]]
name = "fermi_window_adjustment"
source-dir="test"
main = "check_fermi_window_adjustment.f90"

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

Copy link
Owner

@nakib nakib left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Dear Pritam,

Thanks a lot for your contribution. Really appreciate it!

I left my review below. Please take a look when you have time.

src/delta.f90 Outdated Show resolved Hide resolved
src/delta.f90 Show resolved Hide resolved
src/delta.f90 Outdated Show resolved Hide resolved
src/delta.f90 Outdated Show resolved Hide resolved
src/delta.f90 Outdated Show resolved Hide resolved
src/electron.f90 Show resolved Hide resolved
test/check_fermi_window_adjustment.f90 Outdated Show resolved Hide resolved
test/check_fermi_window_adjustment.f90 Outdated Show resolved Hide resolved
test/check_fermi_window_adjustment.f90 Outdated Show resolved Hide resolved
test/check_fermi_window_adjustment.f90 Outdated Show resolved Hide resolved
@prikarsartam prikarsartam requested a review from nakib November 1, 2024 16:18
@prikarsartam
Copy link
Contributor Author

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 nakib merged commit 1e72b6f into nakib:prot Nov 1, 2024
1 check passed
@nakib
Copy link
Owner

nakib commented Nov 1, 2024

Dear Pritam,

Thanks again for your contribution. Really appreciate it.

Best,
Nakib

@prikarsartam
Copy link
Contributor Author

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.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
None yet
Projects
None yet
Development

Successfully merging this pull request may close these issues.

2 participants