LCOV - code coverage report
Current view: top level - src - dm_ls_chebyshev.F (source / functions) Hit Total Coverage
Test: CP2K Regtests (git:15a58fb) Lines: 173 174 99.4 %
Date: 2025-02-18 08:24:35 Functions: 3 3 100.0 %

          Line data    Source code
       1             : !--------------------------------------------------------------------------------------------------!
       2             : !   CP2K: A general program to perform molecular dynamics simulations                              !
       3             : !   Copyright 2000-2025 CP2K developers group <https://cp2k.org>                                   !
       4             : !                                                                                                  !
       5             : !   SPDX-License-Identifier: GPL-2.0-or-later                                                      !
       6             : !--------------------------------------------------------------------------------------------------!
       7             : 
       8             : ! **************************************************************************************************
       9             : !> \brief Routines using linear scaling chebyshev methods
      10             : !> \par History
      11             : !>       2012.10 created [Jinwoong Cha]
      12             : !> \author Jinwoong Cha
      13             : ! **************************************************************************************************
      14             : MODULE dm_ls_chebyshev
      15             :    USE arnoldi_api,                     ONLY: arnoldi_extremal
      16             :    USE cp_dbcsr_api,                    ONLY: &
      17             :         dbcsr_add, dbcsr_add_on_diag, dbcsr_copy, dbcsr_create, dbcsr_get_info, &
      18             :         dbcsr_get_occupation, dbcsr_multiply, dbcsr_release, dbcsr_scale, dbcsr_set, dbcsr_trace, &
      19             :         dbcsr_type, dbcsr_type_no_symmetry
      20             :    USE cp_dbcsr_contrib,                ONLY: dbcsr_frobenius_norm
      21             :    USE cp_log_handling,                 ONLY: cp_get_default_logger,&
      22             :                                               cp_logger_get_default_unit_nr,&
      23             :                                               cp_logger_type
      24             :    USE cp_output_handling,              ONLY: cp_p_file,&
      25             :                                               cp_print_key_finished_output,&
      26             :                                               cp_print_key_should_output,&
      27             :                                               cp_print_key_unit_nr
      28             :    USE dm_ls_scf_qs,                    ONLY: write_matrix_to_cube
      29             :    USE dm_ls_scf_types,                 ONLY: ls_scf_env_type
      30             :    USE input_section_types,             ONLY: section_get_ivals,&
      31             :                                               section_vals_val_get
      32             :    USE kinds,                           ONLY: default_string_length,&
      33             :                                               dp
      34             :    USE machine,                         ONLY: m_flush,&
      35             :                                               m_walltime
      36             :    USE mathconstants,                   ONLY: pi
      37             :    USE qs_environment_types,            ONLY: qs_environment_type
      38             : #include "./base/base_uses.f90"
      39             : 
      40             :    IMPLICIT NONE
      41             : 
      42             :    PRIVATE
      43             : 
      44             :    CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'dm_ls_chebyshev'
      45             : 
      46             :    PUBLIC :: compute_chebyshev
      47             : 
      48             : CONTAINS
      49             : 
      50             : ! **************************************************************************************************
      51             : !> \brief compute chebyshev polynomials up to order n for a given value of x
      52             : !> \param value ...
      53             : !> \param x ...
      54             : !> \param n ...
      55             : !> \par History
      56             : !>       2012.11 created [Jinwoong Cha]
      57             : !> \author Jinwoong Cha
      58             : ! **************************************************************************************************
      59     3729300 :    SUBROUTINE chebyshev_poly(value, x, n)
      60             :       REAL(KIND=dp), INTENT(OUT)                         :: value
      61             :       REAL(KIND=dp), INTENT(IN)                          :: x
      62             :       INTEGER, INTENT(IN)                                :: n
      63             : 
      64             : !polynomial values
      65             : !number of chev polynomials
      66             : 
      67     3729300 :       value = COS((n - 1)*ACOS(x))
      68             : 
      69     3729300 :    END SUBROUTINE chebyshev_poly
      70             : 
      71             : ! **************************************************************************************************
      72             : !> \brief kernel for chebyshev polynomials expansion (Jackson kernel)
      73             : !> \param value ...
      74             : !> \param n ...
      75             : !> \param nc ...
      76             : !> \par History
      77             : !>       2012.11 created [Jinwoong Cha]
      78             : !> \author Jinwoong Cha
      79             : ! **************************************************************************************************
      80        3000 :    SUBROUTINE kernel(value, n, nc)
      81             :       REAL(KIND=dp), INTENT(OUT)                         :: value
      82             :       INTEGER, INTENT(IN)                                :: n, nc
      83             : 
      84             : !kernel at n
      85             : !n-1 order of chebyshev polynomials
      86             : !number of total chebyshev polynomials
      87             : !Kernel define
      88             : 
      89             :       value = 1.0_dp/(nc + 1.0_dp)*((nc - (n - 1) + 1.0_dp)* &
      90        3000 :                                     COS(pi*(n - 1)/(nc + 1.0_dp)) + SIN(pi*(n - 1)/(nc + 1.0_dp))*1.0_dp/TAN(pi/(nc + 1.0_dp)))
      91             : 
      92        3000 :    END SUBROUTINE kernel
      93             : 
      94             : ! **************************************************************************************************
      95             : !> \brief compute properties based on chebyshev expansion
      96             : !> \param qs_env ...
      97             : !> \param ls_scf_env ...
      98             : !> \par History
      99             : !>       2012.10 created [Jinwoong Cha]
     100             : !> \author Jinwoong Cha
     101             : ! **************************************************************************************************
     102           6 :    SUBROUTINE compute_chebyshev(qs_env, ls_scf_env)
     103             :       TYPE(qs_environment_type), POINTER                 :: qs_env
     104             :       TYPE(ls_scf_env_type)                              :: ls_scf_env
     105             : 
     106             :       CHARACTER(len=*), PARAMETER                        :: routineN = 'compute_chebyshev'
     107             :       REAL(KIND=dp), PARAMETER                           :: scale_evals = 1.01_dp
     108             : 
     109             :       CHARACTER(LEN=30)                                  :: middle_name
     110             :       CHARACTER(LEN=default_string_length)               :: title
     111             :       INTEGER                                            :: handle, icheb, igrid, iinte, ispin, &
     112             :                                                             iwindow, n_gridpoint_dos, ncheb, &
     113             :                                                             ninte, Nrows, nwindow, unit_cube, &
     114             :                                                             unit_dos, unit_nr
     115             :       LOGICAL                                            :: converged, write_cubes
     116             :       REAL(KIND=dp) :: chev_T, chev_T_dos, dummy1, final, frob_matrix, initial, interval_a, &
     117             :          interval_b, max_ev, min_ev, occ, orbital_occ, summa, t1, t2
     118           6 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:)           :: chev_E, chev_Es_dos, dos, dummy2, ev1, &
     119           6 :                                                             ev2, kernel_g, mu, sev1, sev2, trace_dm
     120           6 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :)        :: aitchev_T, E_inte, gdensity, sqrt_vec
     121           6 :       REAL(KIND=dp), DIMENSION(:), POINTER               :: tmp_r
     122             :       TYPE(cp_logger_type), POINTER                      :: logger
     123             :       TYPE(dbcsr_type)                                   :: matrix_dummy1, matrix_F, matrix_tmp1, &
     124             :                                                             matrix_tmp2, matrix_tmp3
     125           6 :       TYPE(dbcsr_type), DIMENSION(:), POINTER            :: matrix_dummy2
     126             : 
     127           6 :       IF (.NOT. ls_scf_env%chebyshev%compute_chebyshev) RETURN
     128             : 
     129           6 :       CALL timeset(routineN, handle)
     130             : 
     131             :       ! get a useful output_unit
     132           6 :       logger => cp_get_default_logger()
     133           6 :       IF (logger%para_env%is_source()) THEN
     134           3 :          unit_nr = cp_logger_get_default_unit_nr(logger, local=.TRUE.)
     135             :       ELSE
     136           3 :          unit_nr = -1
     137             :       END IF
     138             : 
     139           6 :       ncheb = ls_scf_env%chebyshev%n_chebyshev
     140           6 :       ninte = 2*ncheb
     141           6 :       n_gridpoint_dos = ls_scf_env%chebyshev%n_gridpoint_dos
     142             : 
     143           6 :       write_cubes = BTEST(cp_print_key_should_output(logger%iter_info, ls_scf_env%chebyshev%print_key_cube), cp_p_file)
     144           6 :       IF (write_cubes) THEN
     145           2 :          IF (ASSOCIATED(ls_scf_env%chebyshev%min_energy)) DEALLOCATE (ls_scf_env%chebyshev%min_energy)
     146           2 :          CALL section_vals_val_get(ls_scf_env%chebyshev%print_key_cube, "MIN_ENERGY", r_vals=tmp_r)
     147           6 :          ALLOCATE (ls_scf_env%chebyshev%min_energy(SIZE(tmp_r)))
     148           8 :          ls_scf_env%chebyshev%min_energy = tmp_r
     149             : 
     150           2 :          IF (ASSOCIATED(ls_scf_env%chebyshev%max_energy)) DEALLOCATE (ls_scf_env%chebyshev%max_energy)
     151           2 :          CALL section_vals_val_get(ls_scf_env%chebyshev%print_key_cube, "MAX_ENERGY", r_vals=tmp_r)
     152           6 :          ALLOCATE (ls_scf_env%chebyshev%max_energy(SIZE(tmp_r)))
     153           8 :          ls_scf_env%chebyshev%max_energy = tmp_r
     154             : 
     155           2 :          nwindow = SIZE(ls_scf_env%chebyshev%min_energy)
     156             :       ELSE
     157             :          nwindow = 0
     158             :       END IF
     159             : 
     160          14 :       ALLOCATE (ev1(1:nwindow))
     161           8 :       ALLOCATE (ev2(1:nwindow))
     162           8 :       ALLOCATE (sev1(1:nwindow))
     163           8 :       ALLOCATE (sev2(1:nwindow))
     164           8 :       ALLOCATE (trace_dm(1:nwindow))
     165          20 :       ALLOCATE (matrix_dummy2(1:nwindow))
     166             : 
     167          12 :       DO iwindow = 1, nwindow
     168           6 :          ev1(iwindow) = ls_scf_env%chebyshev%min_energy(iwindow)
     169          12 :          ev2(iwindow) = ls_scf_env%chebyshev%max_energy(iwindow)
     170             :       END DO
     171             : 
     172           6 :       IF (unit_nr > 0) THEN
     173           3 :          WRITE (unit_nr, '()')
     174           3 :          WRITE (unit_nr, '(T2,A)') "STARTING CHEBYSHEV CALCULATION"
     175             :       END IF
     176             : 
     177             :       ! create 3 temporary matrices
     178           6 :       CALL dbcsr_create(matrix_tmp1, template=ls_scf_env%matrix_s, matrix_type=dbcsr_type_no_symmetry)
     179           6 :       CALL dbcsr_create(matrix_tmp2, template=ls_scf_env%matrix_s, matrix_type=dbcsr_type_no_symmetry)
     180           6 :       CALL dbcsr_create(matrix_tmp3, template=ls_scf_env%matrix_s, matrix_type=dbcsr_type_no_symmetry)
     181           6 :       CALL dbcsr_create(matrix_F, template=ls_scf_env%matrix_s, matrix_type=dbcsr_type_no_symmetry)
     182           6 :       CALL dbcsr_create(matrix_dummy1, template=ls_scf_env%matrix_s, matrix_type=dbcsr_type_no_symmetry)
     183             : 
     184          12 :       DO iwindow = 1, nwindow
     185             :          CALL dbcsr_create(matrix_dummy2(iwindow), template=ls_scf_env%matrix_s, &
     186          12 :                            matrix_type=dbcsr_type_no_symmetry)
     187             :       END DO
     188             : 
     189          12 :       DO ispin = 1, SIZE(ls_scf_env%matrix_ks)
     190             :          ! create matrix_F=inv(sqrt(S))*H*inv(sqrt(S))
     191             :          CALL dbcsr_multiply("N", "N", 1.0_dp, ls_scf_env%matrix_s_sqrt_inv, ls_scf_env%matrix_ks(ispin), &
     192           6 :                              0.0_dp, matrix_tmp1, filter_eps=ls_scf_env%eps_filter)
     193             :          CALL dbcsr_multiply("N", "N", 1.0_dp, matrix_tmp1, ls_scf_env%matrix_s_sqrt_inv, &
     194           6 :                              0.0_dp, matrix_F, filter_eps=ls_scf_env%eps_filter)
     195             : 
     196             :          ! find largest and smallest eigenvalues
     197             :          CALL arnoldi_extremal(matrix_F, max_ev, min_ev, converged=converged, max_iter=ls_scf_env%max_iter_lanczos, &
     198           6 :                                threshold=ls_scf_env%eps_lanczos) !Lanczos algorithm to calculate eigenvalue
     199           6 :          IF (unit_nr > 0) WRITE (unit_nr, '(T2,A,2F16.8,A,L2)') &
     200           3 :             "smallest largest eigenvalue", min_ev, max_ev, " converged ", converged
     201           6 :          IF (nwindow > 0) THEN
     202           2 :             IF (unit_nr > 0) WRITE (unit_nr, '(T2,A,1000F16.8)') "requested interval-min_energy", ev1(:)
     203           2 :             IF (unit_nr > 0) WRITE (unit_nr, '(T2,A,1000F16.8)') "requested interval-max_energy", ev2(:)
     204             :          END IF
     205           6 :          interval_a = (max_ev - min_ev)*scale_evals/2
     206           6 :          interval_b = (max_ev + min_ev)/2
     207             : 
     208          12 :          sev1(:) = (ev1(:) - interval_b)/interval_a !scaled ev1 vector
     209          12 :          sev2(:) = (ev2(:) - interval_b)/interval_a !scaled ev2 vector
     210             : 
     211             :          !chebyshev domain,pi*sqrt(1-x^2) vector construction and chebyshev polynomials for integration (for g(E))
     212          20 :          ALLOCATE (E_inte(1:ninte + 1, 1:nwindow))
     213          14 :          ALLOCATE (sqrt_vec(1:ninte + 1, 1:nwindow))
     214             : 
     215          12 :          DO iwindow = 1, nwindow
     216        4818 :             DO iinte = 1, ninte + 1
     217        4806 :                E_inte(iinte, iwindow) = sev1(iwindow) + ((sev2(iwindow) - sev1(iwindow))/ninte)*(iinte - 1)
     218        4812 :                sqrt_vec(iinte, iwindow) = pi*SQRT(1.0_dp - E_inte(iinte, iwindow)*E_inte(iinte, iwindow))
     219             :             END DO
     220             :          END DO
     221             : 
     222             :          !integral.. (identical to the coefficient for g(E))
     223             : 
     224          20 :          ALLOCATE (aitchev_T(1:ncheb, 1:nwindow)) !after intergral. =>ainte
     225             : 
     226          12 :          DO iwindow = 1, nwindow
     227        2412 :             DO icheb = 1, ncheb
     228        2400 :                CALL chebyshev_poly(initial, E_inte(1, iwindow), icheb)
     229        2400 :                CALL chebyshev_poly(final, E_inte(1, iwindow), icheb)
     230        2400 :           summa = (sev2(iwindow) - sev1(iwindow))/(2.0_dp*ninte)*(initial/sqrt_vec(1, iwindow) + final/sqrt_vec(ninte + 1, iwindow))
     231     1920000 :                DO iinte = 2, ninte
     232     1917600 :                   CALL chebyshev_poly(chev_T, E_inte(iinte, iwindow), icheb)
     233     1920000 :                   summa = summa + ((sev2(iwindow) - sev1(iwindow))/ninte)*(chev_T/sqrt_vec(iinte, iwindow))
     234             :                END DO
     235        2400 :                aitchev_T(icheb, iwindow) = summa
     236        4806 :                summa = 0
     237             :             END DO
     238             :          END DO
     239             : 
     240             :          ! scale the matrix to get evals in the interval -1,1
     241           6 :          CALL dbcsr_add_on_diag(matrix_F, -interval_b)
     242           6 :          CALL dbcsr_scale(matrix_F, 1/interval_a)
     243             : 
     244             :          ! compute chebyshev matrix recursion
     245           6 :          CALL dbcsr_get_info(matrix=matrix_F, nfullrows_total=Nrows) !get information about a matrix
     246           6 :          CALL dbcsr_set(matrix_dummy1, 0.0_dp) !empty matrix creation(for density matrix)
     247             : 
     248          12 :          DO iwindow = 1, nwindow
     249          12 :             CALL dbcsr_set(matrix_dummy2(iwindow), 0.0_dp) !empty matrix creation(for density matrix)
     250             :          END DO
     251             : 
     252          18 :          ALLOCATE (mu(1:ncheb))
     253          12 :          ALLOCATE (kernel_g(1:ncheb))
     254           6 :          CALL kernel(kernel_g(1), 1, ncheb)
     255           6 :          CALL kernel(kernel_g(2), 2, ncheb)
     256             : 
     257           6 :          CALL dbcsr_set(matrix_tmp1, 0.0_dp) !matrix creation
     258           6 :          CALL dbcsr_add_on_diag(matrix_tmp1, 1.0_dp) !add a only number to diagonal elements
     259           6 :          CALL dbcsr_trace(matrix_tmp1, trace=mu(1))
     260           6 :          CALL dbcsr_copy(matrix_tmp2, matrix_F) !make matrix_tmp2 = matrix_F
     261           6 :          CALL dbcsr_trace(matrix_tmp2, trace=mu(2))
     262             : 
     263          12 :          DO iwindow = 1, nwindow
     264           6 :             CALL dbcsr_copy(matrix_dummy1, matrix_tmp1)
     265           6 :             CALL dbcsr_copy(matrix_dummy2(iwindow), matrix_tmp2) !matrix_dummy2=
     266           6 :             CALL dbcsr_scale(matrix_dummy1, kernel_g(1)*aitchev_T(1, iwindow)) !first term of chebyshev poly(matrix)
     267           6 :             CALL dbcsr_scale(matrix_dummy2(iwindow), 2.0_dp*kernel_g(2)*aitchev_T(2, iwindow)) !second term of chebyshev poly(matrix)
     268             : 
     269          12 :             CALL dbcsr_add(matrix_dummy2(iwindow), matrix_dummy1, 1.0_dp, 1.0_dp)
     270             :          END DO
     271             : 
     272        2994 :          DO icheb = 2, ncheb - 1
     273        2988 :             t1 = m_walltime()
     274             :             CALL dbcsr_multiply("N", "N", 2.0_dp, matrix_F, matrix_tmp2, &
     275        2988 :                                 -1.0_dp, matrix_tmp1, filter_eps=ls_scf_env%eps_filter) !matrix multiplication(Recursion)
     276        2988 :             CALL dbcsr_copy(matrix_tmp3, matrix_tmp1)
     277        2988 :             CALL dbcsr_copy(matrix_tmp1, matrix_tmp2)
     278        2988 :             CALL dbcsr_copy(matrix_tmp2, matrix_tmp3)
     279        2988 :             CALL dbcsr_trace(matrix_tmp2, trace=mu(icheb + 1)) !icheb+1 th coefficient
     280        2988 :             CALL kernel(kernel_g(icheb + 1), icheb + 1, ncheb)
     281             : 
     282        5376 :             DO iwindow = 1, nwindow
     283             : 
     284        2388 :                CALL dbcsr_copy(matrix_dummy1, matrix_tmp2)
     285        2388 :                CALL dbcsr_scale(matrix_dummy1, 2.0_dp*kernel_g(icheb + 1)*aitchev_T(icheb + 1, iwindow)) !second term of chebyshev poly(matrix)
     286        2388 :                CALL dbcsr_add(matrix_dummy2(iwindow), matrix_dummy1, 1.0_dp, 1.0_dp)
     287        5376 :                CALL dbcsr_trace(matrix_dummy2(iwindow), trace=trace_dm(iwindow)) !icheb+1 th coefficient
     288             : 
     289             :             END DO
     290             : 
     291        2988 :             occ = dbcsr_get_occupation(matrix_tmp1)
     292        2988 :             t2 = m_walltime()
     293        2994 :             IF (unit_nr > 0 .AND. MOD(icheb, 20) == 0) THEN
     294          72 :                CALL m_flush(unit_nr)
     295          72 :                IF (nwindow > 0) THEN
     296             :                   WRITE (unit_nr, '(T2,A,I5,1X,A,1X,F8.3,1X,A,1X,F8.6,1X,A,1X,1000F16.8)') &
     297          19 :                      "Iter.", icheb, "time=", t2 - t1, "occ=", occ, "traces=", trace_dm(:)
     298             :                ELSE
     299             :                   WRITE (unit_nr, '(T2,A,I5,1X,A,1X,F8.3,1X,A,1X,F8.6)') &
     300          53 :                      "Iter.", icheb, "time=", t2 - t1, "occ=", occ
     301             :                END IF
     302             :             END IF
     303             :          END DO
     304             : 
     305          12 :          DO iwindow = 1, nwindow
     306           6 :             IF (SIZE(ls_scf_env%matrix_ks) == 1) THEN
     307           6 :                orbital_occ = 2.0_dp
     308             :             ELSE
     309           0 :                orbital_occ = 1.0_dp
     310             :             END IF
     311             :             CALL dbcsr_multiply("N", "N", 1.0_dp, ls_scf_env%matrix_s_sqrt_inv, matrix_dummy2(iwindow), &
     312           6 :                                 0.0_dp, matrix_tmp1, filter_eps=ls_scf_env%eps_filter)
     313             :             CALL dbcsr_multiply("N", "N", orbital_occ, matrix_tmp1, ls_scf_env%matrix_s_sqrt_inv, &
     314           6 :                                 0.0_dp, matrix_tmp2, filter_eps=ls_scf_env%eps_filter)
     315           6 :             CALL dbcsr_copy(matrix_dummy2(iwindow), matrix_tmp2)
     316             : 
     317             :             ! look at the difference with the density matrix from the ls routines
     318           6 :             IF (.FALSE.) THEN
     319             :                CALL dbcsr_copy(matrix_tmp1, matrix_tmp2)
     320             :                CALL dbcsr_add(matrix_tmp1, ls_scf_env%matrix_p(ispin), 1.0_dp, -1.0_dp) !comparison
     321             :                frob_matrix = dbcsr_frobenius_norm(matrix_tmp1)
     322             :                IF (unit_nr > 0) WRITE (unit_nr, *) "Difference between Chebyshev DM and LS DM", frob_matrix
     323             :             END IF
     324             :          END DO
     325             : 
     326             :          write_cubes = BTEST(cp_print_key_should_output(logger%iter_info, &
     327           6 :                                                         ls_scf_env%chebyshev%print_key_cube), cp_p_file)
     328          24 :          IF (write_cubes) THEN
     329           8 :             DO iwindow = 1, nwindow
     330           6 :                WRITE (middle_name, "(A,I0)") "E_DENSITY_WINDOW_", iwindow
     331           6 :                WRITE (title, "(A,1X,F16.8,1X,A,1X,F16.8)") "Energy range : ", ev1(iwindow), "to", ev2(iwindow)
     332             :                unit_cube = cp_print_key_unit_nr(logger, ls_scf_env%chebyshev%print_key_cube, &
     333             :                                                 "", extension=".cube", & !added 01/22/2012
     334           6 :                                                 middle_name=TRIM(middle_name), log_filename=.FALSE.)
     335             :                CALL write_matrix_to_cube(qs_env, ls_scf_env, matrix_dummy2(iwindow), unit_cube, title, &
     336           6 :                                          section_get_ivals(ls_scf_env%chebyshev%print_key_cube, "STRIDE"))
     337           8 :                CALL cp_print_key_finished_output(unit_cube, logger, ls_scf_env%chebyshev%print_key_cube, "")
     338             :             END DO
     339             :          END IF
     340             : 
     341             :       END DO
     342             : 
     343             :       ! Chebyshev expansion with calculated coefficient
     344             :       ! grid construction and rescaling (by J)
     345             :       unit_dos = cp_print_key_unit_nr(logger, ls_scf_env%chebyshev%print_key_dos, "", extension=".xy", &
     346           6 :                                       middle_name="DOS", log_filename=.FALSE.)
     347             : 
     348           6 :       IF (unit_dos > 0) THEN
     349           9 :          ALLOCATE (dos(1:n_gridpoint_dos))
     350          10 :          ALLOCATE (gdensity(1:n_gridpoint_dos, 1:nwindow))
     351           6 :          ALLOCATE (chev_E(1:n_gridpoint_dos))
     352           6 :          ALLOCATE (chev_Es_dos(1:n_gridpoint_dos))
     353           4 :          ALLOCATE (dummy2(1:nwindow))
     354        3103 :          DO igrid = 1, n_gridpoint_dos
     355        3100 :             chev_E(igrid) = min_ev + (igrid - 1)*(max_ev - min_ev)/(n_gridpoint_dos - 1)
     356        3103 :             chev_Es_dos(igrid) = (chev_E(igrid) - interval_b)/interval_a
     357             :          END DO
     358        3103 :          DO igrid = 1, n_gridpoint_dos
     359        9100 :             dummy1 = 0.0_dp !summation of polynomials
     360        9100 :             dummy2(:) = 0.0_dp !summation of polynomials
     361     1810000 :             DO icheb = 2, ncheb
     362     1806900 :                CALL chebyshev_poly(chev_T_dos, chev_Es_dos(igrid), icheb)
     363     1806900 :                dummy1 = dummy1 + kernel_g(icheb)*mu(icheb)*chev_T_dos
     364     4204000 :                DO iwindow = 1, nwindow
     365     4200900 :                   dummy2(iwindow) = dummy2(iwindow) + kernel_g(icheb)*aitchev_T(icheb, iwindow)*chev_T_dos
     366             :                END DO
     367             :             END DO
     368             :             dos(igrid) = 1.0_dp/(interval_a*Nrows* &
     369        3100 :                                  (pi*SQRT(1.0_dp - chev_Es_dos(igrid)*chev_Es_dos(igrid))))*(kernel_g(1)*mu(1) + 2.0_dp*dummy1)
     370        9100 :             DO iwindow = 1, nwindow
     371        9100 :                gdensity(igrid, iwindow) = kernel_g(1)*aitchev_T(1, iwindow) + 2.0_dp*dummy2(iwindow)
     372             :             END DO
     373        3103 :             WRITE (unit_dos, '(1000F16.8)') chev_E(igrid), dos(igrid), gdensity(igrid, :)
     374             :          END DO
     375           3 :          DEALLOCATE (chev_Es_dos, chev_E, dos, gdensity)
     376             :       END IF
     377           6 :       CALL cp_print_key_finished_output(unit_dos, logger, ls_scf_env%chebyshev%print_key_dos, "")
     378             : 
     379             :       ! free the matrices
     380           6 :       CALL dbcsr_release(matrix_tmp1)
     381           6 :       CALL dbcsr_release(matrix_tmp2)
     382           6 :       CALL dbcsr_release(matrix_tmp3)
     383           6 :       CALL dbcsr_release(matrix_F)
     384           6 :       CALL dbcsr_release(matrix_dummy1)
     385             : 
     386          12 :       DO iwindow = 1, nwindow
     387          12 :          CALL dbcsr_release(matrix_dummy2(iwindow))
     388             :       END DO
     389             : 
     390           6 :       DEALLOCATE (ev1, ev2, sev1, sev2, matrix_dummy2)
     391             : 
     392             :       !Need deallocation
     393           6 :       DEALLOCATE (mu, kernel_g, aitchev_T, E_inte, sqrt_vec)
     394             : 
     395           6 :       IF (unit_nr > 0) WRITE (unit_nr, '(T2,A)') "ENDING CHEBYSHEV CALCULATION"
     396             : 
     397           6 :       CALL timestop(handle)
     398             : 
     399          12 :    END SUBROUTINE compute_chebyshev
     400             : 
     401             : END MODULE dm_ls_chebyshev

Generated by: LCOV version 1.15