LCOV - code coverage report
Current view: top level - src - qs_gspace_mixing.F (source / functions) Hit Total Coverage
Test: CP2K Regtests (git:b8e0b09) Lines: 400 583 68.6 %
Date: 2024-08-31 06:31:37 Functions: 4 5 80.0 %

          Line data    Source code
       1             : !--------------------------------------------------------------------------------------------------!
       2             : !   CP2K: A general program to perform molecular dynamics simulations                              !
       3             : !   Copyright 2000-2024 CP2K developers group <https://cp2k.org>                                   !
       4             : !                                                                                                  !
       5             : !   SPDX-License-Identifier: GPL-2.0-or-later                                                      !
       6             : !--------------------------------------------------------------------------------------------------!
       7             : 
       8             : ! **************************************************************************************************
       9             : MODULE qs_gspace_mixing
      10             : 
      11             :    USE cp_control_types,                ONLY: dft_control_type
      12             :    USE kinds,                           ONLY: dp
      13             :    USE mathlib,                         ONLY: invert_matrix
      14             :    USE message_passing,                 ONLY: mp_para_env_type
      15             :    USE pw_env_types,                    ONLY: pw_env_get,&
      16             :                                               pw_env_type
      17             :    USE pw_methods,                      ONLY: pw_axpy,&
      18             :                                               pw_copy,&
      19             :                                               pw_integrate_function,&
      20             :                                               pw_scale,&
      21             :                                               pw_transfer,&
      22             :                                               pw_zero
      23             :    USE pw_pool_types,                   ONLY: pw_pool_type
      24             :    USE pw_types,                        ONLY: pw_c1d_gs_type,&
      25             :                                               pw_r3d_rs_type
      26             :    USE qs_density_mixing_types,         ONLY: broyden_mixing_nr,&
      27             :                                               gspace_mixing_nr,&
      28             :                                               mixing_storage_type,&
      29             :                                               multisecant_mixing_nr,&
      30             :                                               pulay_mixing_nr
      31             :    USE qs_environment_types,            ONLY: get_qs_env,&
      32             :                                               qs_environment_type
      33             :    USE qs_rho_atom_types,               ONLY: rho_atom_type
      34             :    USE qs_rho_types,                    ONLY: qs_rho_get,&
      35             :                                               qs_rho_type
      36             : #include "./base/base_uses.f90"
      37             : 
      38             :    IMPLICIT NONE
      39             : 
      40             :    PRIVATE
      41             : 
      42             :    CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_gspace_mixing'
      43             : 
      44             :    PUBLIC :: gspace_mixing
      45             : 
      46             : CONTAINS
      47             : 
      48             : ! **************************************************************************************************
      49             : !> \brief  Driver for the g-space mixing, calls the proper routine given the
      50             : !>requested method
      51             : !> \param qs_env ...
      52             : !> \param mixing_method ...
      53             : !> \param mixing_store ...
      54             : !> \param rho ...
      55             : !> \param para_env ...
      56             : !> \param iter_count ...
      57             : !> \par History
      58             : !>      05.2009
      59             : !>      02.2015 changed input to make g-space mixing available in linear scaling
      60             : !>              SCF [Patrick Seewald]
      61             : !> \author MI
      62             : ! **************************************************************************************************
      63        1524 :    SUBROUTINE gspace_mixing(qs_env, mixing_method, mixing_store, rho, para_env, iter_count)
      64             :       TYPE(qs_environment_type), POINTER                 :: qs_env
      65             :       INTEGER                                            :: mixing_method
      66             :       TYPE(mixing_storage_type), POINTER                 :: mixing_store
      67             :       TYPE(qs_rho_type), POINTER                         :: rho
      68             :       TYPE(mp_para_env_type), POINTER                    :: para_env
      69             :       INTEGER                                            :: iter_count
      70             : 
      71             :       CHARACTER(len=*), PARAMETER                        :: routineN = 'gspace_mixing'
      72             : 
      73             :       INTEGER                                            :: handle, iatom, ig, ispin, natom, ng, &
      74             :                                                             nimg, nspin
      75             :       LOGICAL                                            :: gapw
      76             :       REAL(dp)                                           :: alpha
      77        1524 :       REAL(KIND=dp), DIMENSION(:), POINTER               :: tot_rho_r
      78             :       TYPE(dft_control_type), POINTER                    :: dft_control
      79             :       TYPE(pw_c1d_gs_type)                               :: rho_tmp
      80        1524 :       TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER        :: rho_g
      81             :       TYPE(pw_env_type), POINTER                         :: pw_env
      82             :       TYPE(pw_pool_type), POINTER                        :: auxbas_pw_pool
      83        1524 :       TYPE(pw_r3d_rs_type), DIMENSION(:), POINTER        :: rho_r
      84        1524 :       TYPE(rho_atom_type), DIMENSION(:), POINTER         :: rho_atom
      85             : 
      86        1524 :       CALL timeset(routineN, handle)
      87             : 
      88        1524 :       CALL get_qs_env(qs_env, dft_control=dft_control, pw_env=pw_env)
      89        1524 :       IF (.NOT. (dft_control%qs_control%dftb .OR. dft_control%qs_control%xtb .OR. &
      90             :                  dft_control%qs_control%semi_empirical)) THEN
      91             : 
      92        1028 :          CPASSERT(ASSOCIATED(mixing_store))
      93        1028 :          NULLIFY (auxbas_pw_pool, rho_atom, rho_g, rho_r, tot_rho_r)
      94        1028 :          CALL qs_rho_get(rho, rho_g=rho_g, rho_r=rho_r, tot_rho_r=tot_rho_r)
      95             : 
      96        1028 :          nspin = SIZE(rho_g, 1)
      97        1028 :          nimg = dft_control%nimages
      98        1028 :          ng = SIZE(rho_g(1)%pw_grid%gsq)
      99        1028 :          CPASSERT(ng == SIZE(mixing_store%rhoin(1)%cc))
     100             : 
     101        1028 :          alpha = mixing_store%alpha
     102        1028 :          gapw = dft_control%qs_control%gapw
     103             : 
     104        1028 :          IF (nspin == 2) THEN
     105         114 :             CALL pw_env_get(pw_env=pw_env, auxbas_pw_pool=auxbas_pw_pool)
     106         114 :             CALL auxbas_pw_pool%create_pw(rho_tmp)
     107         114 :             CALL pw_zero(rho_tmp)
     108         114 :             CALL pw_copy(rho_g(1), rho_tmp)
     109         114 :             CALL pw_axpy(rho_g(2), rho_g(1), 1.0_dp)
     110         114 :             CALL pw_axpy(rho_tmp, rho_g(2), -1.0_dp)
     111         114 :             CALL pw_scale(rho_g(2), -1.0_dp)
     112             :          END IF
     113             : 
     114        1028 :          IF (iter_count + 1 <= mixing_store%nskip_mixing) THEN
     115             :             ! skip mixing
     116           8 :             DO ispin = 1, nspin
     117       27656 :                DO ig = 1, ng
     118       27652 :                   mixing_store%rhoin(ispin)%cc(ig) = rho_g(ispin)%array(ig)
     119             :                END DO
     120             :             END DO
     121           4 :             IF (mixing_store%gmix_p .AND. gapw) THEN
     122           0 :                CALL get_qs_env(qs_env=qs_env, rho_atom_set=rho_atom)
     123           0 :                natom = SIZE(rho_atom)
     124           0 :                DO ispin = 1, nspin
     125           0 :                   DO iatom = 1, natom
     126           0 :                      IF (mixing_store%paw(iatom)) THEN
     127           0 :                         mixing_store%cpc_h_in(iatom, ispin)%r_coef = rho_atom(iatom)%cpc_h(ispin)%r_coef
     128           0 :                         mixing_store%cpc_s_in(iatom, ispin)%r_coef = rho_atom(iatom)%cpc_s(ispin)%r_coef
     129             :                      END IF
     130             :                   END DO
     131             :                END DO
     132             :             END IF
     133             : 
     134           4 :             mixing_store%iter_method = "NoMix"
     135           4 :             IF (nspin == 2) THEN
     136           0 :                CALL pw_axpy(rho_g(2), rho_g(1), 1.0_dp)
     137           0 :                CALL pw_scale(rho_g(1), 0.5_dp)
     138           0 :                CALL pw_axpy(rho_g(1), rho_g(2), -1.0_dp)
     139           0 :                CALL pw_scale(rho_g(2), -1.0_dp)
     140           0 :                CALL auxbas_pw_pool%give_back_pw(rho_tmp)
     141             :             END IF
     142           4 :             CALL timestop(handle)
     143           4 :             RETURN
     144             :          END IF
     145             : 
     146        1024 :          IF ((iter_count + 1 - mixing_store%nskip_mixing) <= mixing_store%n_simple_mix) THEN
     147           6 :             CALL gmix_potential_only(qs_env, mixing_store, rho)
     148           6 :             mixing_store%iter_method = "Kerker"
     149        1018 :          ELSEIF (mixing_method == gspace_mixing_nr) THEN
     150          40 :             CALL gmix_potential_only(qs_env, mixing_store, rho)
     151          40 :             mixing_store%iter_method = "Kerker"
     152         978 :          ELSEIF (mixing_method == pulay_mixing_nr) THEN
     153         184 :             CALL pulay_mixing(qs_env, mixing_store, rho, para_env)
     154         184 :             mixing_store%iter_method = "Pulay"
     155         794 :          ELSEIF (mixing_method == broyden_mixing_nr) THEN
     156         794 :             CALL broyden_mixing(qs_env, mixing_store, rho, para_env)
     157         794 :             mixing_store%iter_method = "Broy."
     158           0 :          ELSEIF (mixing_method == multisecant_mixing_nr) THEN
     159           0 :             CPASSERT(.NOT. gapw)
     160           0 :             CALL multisecant_mixing(mixing_store, rho, para_env)
     161           0 :             mixing_store%iter_method = "MSec."
     162             :          END IF
     163             : 
     164        1024 :          IF (nspin == 2) THEN
     165         114 :             CALL pw_axpy(rho_g(2), rho_g(1), 1.0_dp)
     166         114 :             CALL pw_scale(rho_g(1), 0.5_dp)
     167         114 :             CALL pw_axpy(rho_g(1), rho_g(2), -1.0_dp)
     168         114 :             CALL pw_scale(rho_g(2), -1.0_dp)
     169         114 :             CALL auxbas_pw_pool%give_back_pw(rho_tmp)
     170             :          END IF
     171             : 
     172        2162 :          DO ispin = 1, nspin
     173        1138 :             CALL pw_transfer(rho_g(ispin), rho_r(ispin))
     174        2162 :             tot_rho_r(ispin) = pw_integrate_function(rho_r(ispin), isign=-1)
     175             :          END DO
     176             :       END IF
     177             : 
     178        1520 :       CALL timestop(handle)
     179             : 
     180        1524 :    END SUBROUTINE gspace_mixing
     181             : 
     182             : ! **************************************************************************************************
     183             : !> \brief G-space mixing performed via the Kerker damping on the density on the grid
     184             : !>        thus affecting only the caluclation of the potential, but not the denisity matrix
     185             : !> \param qs_env ...
     186             : !> \param mixing_store ...
     187             : !> \param rho ...
     188             : !> \par History
     189             : !>      02.2009
     190             : !> \author MI
     191             : ! **************************************************************************************************
     192             : 
     193          92 :    SUBROUTINE gmix_potential_only(qs_env, mixing_store, rho)
     194             : 
     195             :       TYPE(qs_environment_type), POINTER                 :: qs_env
     196             :       TYPE(mixing_storage_type), POINTER                 :: mixing_store
     197             :       TYPE(qs_rho_type), POINTER                         :: rho
     198             : 
     199             :       CHARACTER(len=*), PARAMETER :: routineN = 'gmix_potential_only'
     200             : 
     201          46 :       COMPLEX(dp), DIMENSION(:), POINTER                 :: cc_new
     202             :       INTEGER                                            :: handle, iatom, ig, ispin, natom, ng, &
     203             :                                                             nspin
     204             :       LOGICAL                                            :: gapw
     205             :       REAL(dp)                                           :: alpha, f_mix
     206             :       TYPE(dft_control_type), POINTER                    :: dft_control
     207          46 :       TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER        :: rho_g
     208          46 :       TYPE(rho_atom_type), DIMENSION(:), POINTER         :: rho_atom
     209             : 
     210           0 :       CPASSERT(ASSOCIATED(mixing_store%rhoin))
     211          46 :       CPASSERT(ASSOCIATED(mixing_store%kerker_factor))
     212             : 
     213          46 :       CALL timeset(routineN, handle)
     214          46 :       NULLIFY (cc_new, dft_control, rho_atom, rho_g)
     215             : 
     216          46 :       CALL get_qs_env(qs_env, dft_control=dft_control)
     217          46 :       CALL qs_rho_get(rho, rho_g=rho_g)
     218             : 
     219          46 :       nspin = SIZE(rho_g, 1)
     220          46 :       ng = SIZE(rho_g(1)%pw_grid%gsq)
     221             : 
     222          46 :       gapw = dft_control%qs_control%gapw
     223          46 :       alpha = mixing_store%alpha
     224             : 
     225          92 :       DO ispin = 1, nspin
     226          46 :          cc_new => rho_g(ispin)%array
     227      580654 :          DO ig = 1, mixing_store%ig_max ! ng
     228      580608 :             f_mix = mixing_store%alpha*mixing_store%kerker_factor(ig)
     229      580608 :             cc_new(ig) = (1.0_dp - f_mix)*mixing_store%rhoin(ispin)%cc(ig) + f_mix*cc_new(ig)
     230      580654 :             mixing_store%rhoin(ispin)%cc(ig) = cc_new(ig)
     231             :          END DO
     232          92 :          DO ig = mixing_store%ig_max + 1, ng
     233           0 :             f_mix = mixing_store%alpha
     234           0 :             cc_new(ig) = (1.0_dp - f_mix)*mixing_store%rhoin(ispin)%cc(ig) + f_mix*cc_new(ig)
     235          46 :             mixing_store%rhoin(ispin)%cc(ig) = cc_new(ig)
     236             :          END DO
     237             : 
     238             :       END DO
     239             : 
     240          46 :       IF (mixing_store%gmix_p .AND. gapw) THEN
     241             :          CALL get_qs_env(qs_env=qs_env, &
     242           8 :                          rho_atom_set=rho_atom)
     243           8 :          natom = SIZE(rho_atom)
     244          16 :          DO ispin = 1, nspin
     245          80 :             DO iatom = 1, natom
     246          72 :                IF (mixing_store%paw(iatom)) THEN
     247             :                   rho_atom(iatom)%cpc_h(ispin)%r_coef = alpha*rho_atom(iatom)%cpc_h(ispin)%r_coef + &
     248        7408 :                                                         mixing_store%cpc_h_in(iatom, ispin)%r_coef*(1._dp - alpha)
     249             :                   rho_atom(iatom)%cpc_s(ispin)%r_coef = alpha*rho_atom(iatom)%cpc_s(ispin)%r_coef + &
     250        7408 :                                                         mixing_store%cpc_s_in(iatom, ispin)%r_coef*(1._dp - alpha)
     251        7408 :                   mixing_store%cpc_h_in(iatom, ispin)%r_coef = rho_atom(iatom)%cpc_h(ispin)%r_coef
     252        7408 :                   mixing_store%cpc_s_in(iatom, ispin)%r_coef = rho_atom(iatom)%cpc_s(ispin)%r_coef
     253             :                END IF
     254             :             END DO
     255             :          END DO
     256             :       END IF
     257             : 
     258          46 :       CALL timestop(handle)
     259             : 
     260          46 :    END SUBROUTINE gmix_potential_only
     261             : 
     262             : ! **************************************************************************************************
     263             : !> \brief Pulay Mixing using as metrics for the residual the Kerer damping factor
     264             : !>        The mixing is applied directly on the density in reciprocal space,
     265             : !>        therefore it affects the potentials
     266             : !>        on the grid but not the density matrix
     267             : !> \param qs_env ...
     268             : !> \param mixing_store ...
     269             : !> \param rho ...
     270             : !> \param para_env ...
     271             : !> \par History
     272             : !>      03.2009
     273             : !> \author MI
     274             : ! **************************************************************************************************
     275             : 
     276         184 :    SUBROUTINE pulay_mixing(qs_env, mixing_store, rho, para_env)
     277             : 
     278             :       TYPE(qs_environment_type), POINTER                 :: qs_env
     279             :       TYPE(mixing_storage_type), POINTER                 :: mixing_store
     280             :       TYPE(qs_rho_type), POINTER                         :: rho
     281             :       TYPE(mp_para_env_type), POINTER                    :: para_env
     282             : 
     283             :       CHARACTER(len=*), PARAMETER                        :: routineN = 'pulay_mixing'
     284             : 
     285         184 :       COMPLEX(dp), ALLOCATABLE, DIMENSION(:)             :: cc_mix
     286             :       INTEGER                                            :: handle, i, iatom, ib, ibb, ig, ispin, &
     287             :                                                             jb, n1, n2, natom, nb, nb1, nbuffer, &
     288             :                                                             ng, nspin
     289             :       LOGICAL                                            :: gapw
     290             :       REAL(dp)                                           :: alpha_kerker, alpha_pulay, beta, f_mix, &
     291             :                                                             inv_err, norm, norm_c_inv, res_norm, &
     292             :                                                             vol
     293         184 :       REAL(dp), ALLOCATABLE, DIMENSION(:)                :: alpha_c, ev
     294         184 :       REAL(dp), ALLOCATABLE, DIMENSION(:, :)             :: a, b, c, c_inv, cpc_h_mix, cpc_s_mix
     295             :       TYPE(dft_control_type), POINTER                    :: dft_control
     296         184 :       TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER        :: rho_g
     297         184 :       TYPE(rho_atom_type), DIMENSION(:), POINTER         :: rho_atom
     298             : 
     299           0 :       CPASSERT(ASSOCIATED(mixing_store%res_buffer))
     300         184 :       CPASSERT(ASSOCIATED(mixing_store%rhoin_buffer))
     301         184 :       NULLIFY (dft_control, rho_atom, rho_g)
     302         184 :       CALL timeset(routineN, handle)
     303             : 
     304         184 :       CALL get_qs_env(qs_env, dft_control=dft_control)
     305         184 :       CALL qs_rho_get(rho, rho_g=rho_g)
     306         184 :       nspin = SIZE(rho_g, 1)
     307         184 :       ng = SIZE(mixing_store%res_buffer(1, 1)%cc)
     308         184 :       vol = rho_g(1)%pw_grid%vol
     309             : 
     310         184 :       alpha_kerker = mixing_store%alpha
     311         184 :       beta = mixing_store%pulay_beta
     312         184 :       alpha_pulay = mixing_store%pulay_alpha
     313         184 :       nbuffer = mixing_store%nbuffer
     314         184 :       gapw = dft_control%qs_control%gapw
     315         184 :       IF (gapw) THEN
     316          12 :          CALL get_qs_env(qs_env=qs_env, rho_atom_set=rho_atom)
     317          12 :          natom = SIZE(rho_atom)
     318             :       END IF
     319             : 
     320         552 :       ALLOCATE (cc_mix(ng))
     321             : 
     322         184 :       IF (nspin == 2) THEN
     323          14 :          CALL pw_axpy(rho_g(1), rho_g(2), 1.0_dp, -1.0_dp)
     324       96782 :          DO ig = 1, ng
     325       96782 :             mixing_store%rhoin(2)%cc(ig) = mixing_store%rhoin(1)%cc(ig) - mixing_store%rhoin(2)%cc(ig)
     326             :          END DO
     327          14 :          IF (gapw .AND. mixing_store%gmix_p) THEN
     328           0 :             DO iatom = 1, natom
     329           0 :                IF (mixing_store%paw(iatom)) THEN
     330           0 :                   rho_atom(iatom)%cpc_h(2)%r_coef = rho_atom(iatom)%cpc_h(1)%r_coef - rho_atom(iatom)%cpc_h(2)%r_coef
     331           0 :                   rho_atom(iatom)%cpc_s(2)%r_coef = rho_atom(iatom)%cpc_s(1)%r_coef - rho_atom(iatom)%cpc_s(2)%r_coef
     332             :                   mixing_store%cpc_h_in(iatom, 2)%r_coef = mixing_store%cpc_h_in(iatom, 1)%r_coef - &
     333           0 :                                                            mixing_store%cpc_h_in(iatom, 2)%r_coef
     334             :                   mixing_store%cpc_s_in(iatom, 2)%r_coef = mixing_store%cpc_s_in(iatom, 1)%r_coef - &
     335           0 :                                                            mixing_store%cpc_s_in(iatom, 2)%r_coef
     336             :                END IF
     337             :             END DO
     338             :          END IF
     339             :       END IF
     340             : 
     341         382 :       DO ispin = 1, nspin
     342         198 :          ib = MODULO(mixing_store%ncall_p(ispin), nbuffer) + 1
     343             : 
     344         198 :          mixing_store%ncall_p(ispin) = mixing_store%ncall_p(ispin) + 1
     345         198 :          nb = MIN(mixing_store%ncall_p(ispin), nbuffer)
     346         198 :          ibb = MODULO(mixing_store%ncall_p(ispin), nbuffer) + 1
     347             : 
     348     8641926 :          mixing_store%res_buffer(ib, ispin)%cc(:) = CMPLX(0._dp, 0._dp, KIND=dp)
     349         198 :          norm = 0.0_dp
     350     1306566 :          IF (nb == 1) mixing_store%rhoin_buffer(1, ispin)%cc = mixing_store%rhoin(ispin)%cc
     351         198 :          res_norm = 0.0_dp
     352     8641926 :          DO ig = 1, ng
     353     8641728 :             f_mix = mixing_store%kerker_factor(ig)
     354             :             mixing_store%res_buffer(ib, ispin)%cc(ig) = f_mix*(rho_g(ispin)%array(ig) - &
     355     8641728 :                                                                mixing_store%rhoin_buffer(ib, ispin)%cc(ig))
     356             :             res_norm = res_norm + &
     357             :                        REAL(mixing_store%res_buffer(ib, ispin)%cc(ig), dp)*REAL(mixing_store%res_buffer(ib, ispin)%cc(ig), dp) + &
     358     8641926 :                        AIMAG(mixing_store%res_buffer(ib, ispin)%cc(ig))*AIMAG(mixing_store%res_buffer(ib, ispin)%cc(ig))
     359             :          END DO
     360             : 
     361         198 :          CALL para_env%sum(res_norm)
     362         198 :          res_norm = SQRT(res_norm)
     363             : 
     364         198 :          IF (res_norm < 1.E-14_dp) THEN
     365           0 :             mixing_store%ncall_p(ispin) = mixing_store%ncall_p(ispin) - 1
     366             :          ELSE
     367         198 :             nb1 = nb + 1
     368         792 :             ALLOCATE (a(nb1, nb1))
     369        5202 :             a = 0.0_dp
     370         594 :             ALLOCATE (b(nb1, nb1))
     371        5202 :             b = 0.0_dp
     372         792 :             ALLOCATE (c(nb, nb))
     373        3518 :             c = 0.0_dp
     374         594 :             ALLOCATE (c_inv(nb, nb))
     375        3518 :             c_inv = 0.0_dp
     376         594 :             ALLOCATE (alpha_c(nb))
     377         842 :             alpha_c = 0.0_dp
     378         198 :             norm_c_inv = 0.0_dp
     379         594 :             ALLOCATE (ev(nb1))
     380        1040 :             ev = 0.0_dp
     381             : 
     382             :          END IF
     383             : 
     384         842 :          DO jb = 1, nb
     385         644 :             mixing_store%pulay_matrix(jb, ib, ispin) = 0.0_dp
     386    33413252 :             DO ig = 1, ng
     387    33412608 :                f_mix = mixing_store%special_metric(ig)
     388             :                mixing_store%pulay_matrix(jb, ib, ispin) = mixing_store%pulay_matrix(jb, ib, ispin) + &
     389             :                                                           f_mix*(REAL(mixing_store%res_buffer(jb, ispin)%cc(ig), dp) &
     390             :                                                                  *REAL(mixing_store%res_buffer(ib, ispin)%cc(ig), dp) + &
     391             :                                                                  AIMAG(mixing_store%res_buffer(jb, ispin)%cc(ig))* &
     392    33413252 :                                                                  AIMAG(mixing_store%res_buffer(ib, ispin)%cc(ig)))
     393             :             END DO
     394         644 :             CALL para_env%sum(mixing_store%pulay_matrix(jb, ib, ispin))
     395         644 :             mixing_store%pulay_matrix(jb, ib, ispin) = mixing_store%pulay_matrix(jb, ib, ispin)*vol*vol
     396         842 :             mixing_store%pulay_matrix(ib, jb, ispin) = mixing_store%pulay_matrix(jb, ib, ispin)
     397             :          END DO
     398             : 
     399         198 :          IF (nb == 1 .OR. res_norm < 1.E-14_dp) THEN
     400     1306406 :             DO ig = 1, ng
     401     1306368 :                f_mix = alpha_kerker*mixing_store%kerker_factor(ig)
     402             :                cc_mix(ig) = rho_g(ispin)%array(ig) - &
     403     1306368 :                             mixing_store%rhoin_buffer(ib, ispin)%cc(ig)
     404             :                rho_g(ispin)%array(ig) = f_mix*cc_mix(ig) + &
     405     1306368 :                                         mixing_store%rhoin_buffer(ib, ispin)%cc(ig)
     406     1306406 :                mixing_store%rhoin_buffer(ibb, ispin)%cc(ig) = rho_g(ispin)%array(ig)
     407             :             END DO
     408          38 :             IF (mixing_store%gmix_p) THEN
     409           2 :                IF (gapw) THEN
     410          18 :                   DO iatom = 1, natom
     411          18 :                      IF (mixing_store%paw(iatom)) THEN
     412             :                         mixing_store%cpc_h_res_buffer(ib, iatom, ispin)%r_coef = rho_atom(iatom)%cpc_h(ispin)%r_coef - &
     413        1850 :                                                                                  mixing_store%cpc_h_in(iatom, ispin)%r_coef
     414             :                         mixing_store%cpc_s_res_buffer(ib, iatom, ispin)%r_coef = rho_atom(iatom)%cpc_s(ispin)%r_coef - &
     415        1850 :                                                                                  mixing_store%cpc_s_in(iatom, ispin)%r_coef
     416             : 
     417             :                         rho_atom(iatom)%cpc_h(ispin)%r_coef = alpha_kerker*rho_atom(iatom)%cpc_h(ispin)%r_coef + &
     418        1850 :                                                               (1.0_dp - alpha_kerker)*mixing_store%cpc_h_in(iatom, ispin)%r_coef
     419             :                         rho_atom(iatom)%cpc_s(ispin)%r_coef = alpha_kerker*rho_atom(iatom)%cpc_s(ispin)%r_coef + &
     420        1850 :                                                               (1.0_dp - alpha_kerker)*mixing_store%cpc_s_in(iatom, ispin)%r_coef
     421             : 
     422         926 :                         mixing_store%cpc_h_in_buffer(ib, iatom, ispin)%r_coef = mixing_store%cpc_h_in(iatom, ispin)%r_coef
     423         926 :                         mixing_store%cpc_s_in_buffer(ib, iatom, ispin)%r_coef = mixing_store%cpc_s_in(iatom, ispin)%r_coef
     424             : 
     425        1850 :                         mixing_store%cpc_h_in_buffer(ibb, iatom, ispin)%r_coef = rho_atom(iatom)%cpc_h(ispin)%r_coef
     426        1850 :                         mixing_store%cpc_s_in_buffer(ibb, iatom, ispin)%r_coef = rho_atom(iatom)%cpc_s(ispin)%r_coef
     427             :                      END IF
     428             :                   END DO
     429             :                END IF
     430             :             END IF
     431             :          ELSE
     432             : 
     433        3404 :             b(1:nb, 1:nb) = mixing_store%pulay_matrix(1:nb, 1:nb, ispin)
     434        3404 :             c(1:nb, 1:nb) = b(1:nb, 1:nb)
     435         766 :             b(nb1, 1:nb) = -1.0_dp
     436         766 :             b(1:nb, nb1) = -1.0_dp
     437         160 :             b(nb1, nb1) = 0.0_dp
     438             : 
     439         160 :             CALL invert_matrix(c, c_inv, inv_err, improve=.TRUE.)
     440         766 :             alpha_c = 0.0_dp
     441         766 :             DO i = 1, nb
     442        3404 :                DO jb = 1, nb
     443        2638 :                   alpha_c(i) = alpha_c(i) + c_inv(jb, i)
     444        3244 :                   norm_c_inv = norm_c_inv + c_inv(jb, i)
     445             :                END DO
     446             :             END DO
     447         766 :             alpha_c(1:nb) = alpha_c(1:nb)/norm_c_inv
     448     7335520 :             cc_mix = CMPLX(0._dp, 0._dp, KIND=dp)
     449         766 :             DO jb = 1, nb
     450    32107006 :                DO ig = 1, ng
     451             :                   cc_mix(ig) = cc_mix(ig) + &
     452             :                                alpha_c(jb)*(mixing_store%rhoin_buffer(jb, ispin)%cc(ig) + &
     453    32106846 :                                             mixing_store%pulay_beta*mixing_store%res_buffer(jb, ispin)%cc(ig))
     454             :                END DO
     455             :             END DO
     456     7335520 :             mixing_store%rhoin_buffer(ibb, ispin)%cc = CMPLX(0._dp, 0._dp, KIND=dp)
     457         160 :             IF (alpha_pulay > 0.0_dp) THEN
     458      233290 :                DO ig = 1, ng
     459      233280 :                   f_mix = alpha_pulay*mixing_store%kerker_factor(ig)
     460             :                   rho_g(ispin)%array(ig) = f_mix*rho_g(ispin)%array(ig) + &
     461      233280 :                                            (1.0_dp - f_mix)*cc_mix(ig)
     462      233290 :                   mixing_store%rhoin_buffer(ibb, ispin)%cc(ig) = rho_g(ispin)%array(ig)
     463             :                END DO
     464             :             ELSE
     465     7102230 :                DO ig = 1, ng
     466     7102080 :                   rho_g(ispin)%array(ig) = cc_mix(ig)
     467     7102230 :                   mixing_store%rhoin_buffer(ibb, ispin)%cc(ig) = rho_g(ispin)%array(ig)
     468             :                END DO
     469             :             END IF
     470             : 
     471         160 :             IF (mixing_store%gmix_p .AND. gapw) THEN
     472          10 :                CALL get_qs_env(qs_env=qs_env, rho_atom_set=rho_atom)
     473          90 :                DO iatom = 1, natom
     474          90 :                   IF (mixing_store%paw(iatom)) THEN
     475          10 :                      n1 = SIZE(rho_atom(iatom)%cpc_h(ispin)%r_coef, 1)
     476          10 :                      n2 = SIZE(rho_atom(iatom)%cpc_h(ispin)%r_coef, 2)
     477          40 :                      ALLOCATE (cpc_h_mix(n1, n2))
     478          30 :                      ALLOCATE (cpc_s_mix(n1, n2))
     479             : 
     480             :                      mixing_store%cpc_h_res_buffer(ib, iatom, ispin)%r_coef = rho_atom(iatom)%cpc_h(ispin)%r_coef - &
     481        9250 :                                                                               mixing_store%cpc_h_in_buffer(ib, iatom, ispin)%r_coef
     482             :                      mixing_store%cpc_s_res_buffer(ib, iatom, ispin)%r_coef = rho_atom(iatom)%cpc_s(ispin)%r_coef - &
     483        9250 :                                                                               mixing_store%cpc_s_in_buffer(ib, iatom, ispin)%r_coef
     484        4630 :                      cpc_h_mix = 0.0_dp
     485        4630 :                      cpc_s_mix = 0.0_dp
     486          48 :                      DO jb = 1, nb
     487             :                         cpc_h_mix(:, :) = cpc_h_mix(:, :) + &
     488             :                                           alpha_c(jb)*mixing_store%cpc_h_in_buffer(jb, iatom, ispin)%r_coef(:, :) + &
     489       17594 :                                           alpha_c(jb)*beta*mixing_store%cpc_h_res_buffer(jb, iatom, ispin)%r_coef(:, :)
     490             :                         cpc_s_mix(:, :) = cpc_s_mix(:, :) + &
     491             :                                           alpha_c(jb)*mixing_store%cpc_s_in_buffer(jb, iatom, ispin)%r_coef(:, :) + &
     492       17604 :                                           alpha_c(jb)*beta*mixing_store%cpc_s_res_buffer(jb, iatom, ispin)%r_coef(:, :)
     493             :                      END DO
     494             :                      rho_atom(iatom)%cpc_h(ispin)%r_coef = alpha_pulay*rho_atom(iatom)%cpc_h(ispin)%r_coef + &
     495        4630 :                                                            (1.0_dp - alpha_pulay)*cpc_h_mix
     496             :                      rho_atom(iatom)%cpc_s(ispin)%r_coef = alpha_pulay*rho_atom(iatom)%cpc_s(ispin)%r_coef + &
     497        4630 :                                                            (1.0_dp - alpha_pulay)*cpc_s_mix
     498        9250 :                      mixing_store%cpc_h_in_buffer(ibb, iatom, ispin)%r_coef = rho_atom(iatom)%cpc_h(ispin)%r_coef
     499        9250 :                      mixing_store%cpc_s_in_buffer(ibb, iatom, ispin)%r_coef = rho_atom(iatom)%cpc_s(ispin)%r_coef
     500          10 :                      DEALLOCATE (cpc_h_mix)
     501          10 :                      DEALLOCATE (cpc_s_mix)
     502             :                   END IF
     503             :                END DO
     504             :             END IF
     505             :          END IF ! nb==1
     506             : 
     507         382 :          IF (res_norm > 1.E-14_dp) THEN
     508         198 :             DEALLOCATE (a)
     509         198 :             DEALLOCATE (b)
     510         198 :             DEALLOCATE (ev)
     511         198 :             DEALLOCATE (c, c_inv, alpha_c)
     512             :          END IF
     513             : 
     514             :       END DO ! ispin
     515             : 
     516         184 :       DEALLOCATE (cc_mix)
     517             : 
     518         184 :       IF (nspin == 2) THEN
     519          14 :          CALL pw_axpy(rho_g(1), rho_g(2), 1.0_dp, -1.0_dp)
     520       96782 :          DO ig = 1, ng
     521       96782 :             mixing_store%rhoin(2)%cc(ig) = mixing_store%rhoin(1)%cc(ig) - mixing_store%rhoin(2)%cc(ig)
     522             :          END DO
     523          14 :          IF (gapw .AND. mixing_store%gmix_p) THEN
     524           0 :             DO iatom = 1, natom
     525           0 :                IF (mixing_store%paw(iatom)) THEN
     526           0 :                   rho_atom(iatom)%cpc_h(2)%r_coef = rho_atom(iatom)%cpc_h(1)%r_coef - rho_atom(iatom)%cpc_h(2)%r_coef
     527           0 :                   rho_atom(iatom)%cpc_s(2)%r_coef = rho_atom(iatom)%cpc_s(1)%r_coef - rho_atom(iatom)%cpc_s(2)%r_coef
     528             :                   mixing_store%cpc_h_in(iatom, 2)%r_coef = mixing_store%cpc_h_in(iatom, 1)%r_coef - &
     529           0 :                                                            mixing_store%cpc_h_in(iatom, 2)%r_coef
     530             :                   mixing_store%cpc_s_in(iatom, 2)%r_coef = mixing_store%cpc_s_in(iatom, 1)%r_coef - &
     531           0 :                                                            mixing_store%cpc_s_in(iatom, 2)%r_coef
     532             :                END IF
     533             :             END DO
     534             :          END IF
     535             : 
     536             :       END IF
     537             : 
     538         184 :       CALL timestop(handle)
     539             : 
     540         552 :    END SUBROUTINE pulay_mixing
     541             : 
     542             : ! **************************************************************************************************
     543             : !> \brief Broyden Mixing using as metrics for the residual the Kerer damping factor
     544             : !>        The mixing is applied directly on the density expansion in reciprocal space
     545             : !> \param qs_env ...
     546             : !> \param mixing_store ...
     547             : !> \param rho ...
     548             : !> \param para_env ...
     549             : !> \par History
     550             : !>      03.2009
     551             : !> \author MI
     552             : ! **************************************************************************************************
     553             : 
     554         794 :    SUBROUTINE broyden_mixing(qs_env, mixing_store, rho, para_env)
     555             : 
     556             :       TYPE(qs_environment_type), POINTER                 :: qs_env
     557             :       TYPE(mixing_storage_type), POINTER                 :: mixing_store
     558             :       TYPE(qs_rho_type), POINTER                         :: rho
     559             :       TYPE(mp_para_env_type), POINTER                    :: para_env
     560             : 
     561             :       CHARACTER(len=*), PARAMETER                        :: routineN = 'broyden_mixing'
     562             : 
     563             :       COMPLEX(dp)                                        :: cc_mix
     564         794 :       COMPLEX(dp), ALLOCATABLE, DIMENSION(:)             :: res_rho
     565             :       INTEGER                                            :: handle, i, iatom, ib, ig, ispin, j, jb, &
     566             :                                                             kb, n1, n2, natom, nb, nbuffer, ng, &
     567             :                                                             nspin
     568             :       LOGICAL                                            :: gapw, only_kerker
     569             :       REAL(dp)                                           :: alpha, dcpc_h_res, dcpc_s_res, &
     570             :                                                             delta_norm, f_mix, inv_err, res_norm, &
     571             :                                                             rho_norm, valh, vals, w0
     572         794 :       REAL(dp), ALLOCATABLE, DIMENSION(:)                :: c, g
     573         794 :       REAL(dp), ALLOCATABLE, DIMENSION(:, :)             :: a, b
     574             :       TYPE(dft_control_type), POINTER                    :: dft_control
     575         794 :       TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER        :: rho_g
     576         794 :       TYPE(rho_atom_type), DIMENSION(:), POINTER         :: rho_atom
     577             : 
     578           0 :       CPASSERT(ASSOCIATED(mixing_store%res_buffer))
     579         794 :       CPASSERT(ASSOCIATED(mixing_store%rhoin))
     580         794 :       CPASSERT(ASSOCIATED(mixing_store%rhoin_old))
     581         794 :       CPASSERT(ASSOCIATED(mixing_store%drho_buffer))
     582         794 :       NULLIFY (dft_control, rho_atom, rho_g)
     583         794 :       CALL timeset(routineN, handle)
     584             : 
     585         794 :       CALL get_qs_env(qs_env, dft_control=dft_control)
     586         794 :       CALL qs_rho_get(rho, rho_g=rho_g)
     587             : 
     588         794 :       nspin = SIZE(rho_g, 1)
     589         794 :       ng = SIZE(mixing_store%res_buffer(1, 1)%cc)
     590             : 
     591         794 :       alpha = mixing_store%alpha
     592         794 :       w0 = mixing_store%broy_w0
     593         794 :       nbuffer = mixing_store%nbuffer
     594         794 :       gapw = dft_control%qs_control%gapw
     595             : 
     596        2382 :       ALLOCATE (res_rho(ng))
     597             : 
     598         794 :       mixing_store%ncall = mixing_store%ncall + 1
     599         794 :       IF (mixing_store%ncall == 1) THEN
     600             :          only_kerker = .TRUE.
     601             :       ELSE
     602         650 :          only_kerker = .FALSE.
     603         650 :          nb = MIN(mixing_store%ncall - 1, nbuffer)
     604         650 :          ib = MODULO(mixing_store%ncall - 2, nbuffer) + 1
     605             :       END IF
     606             : 
     607         794 :       IF (gapw) THEN
     608         102 :          CALL get_qs_env(qs_env=qs_env, rho_atom_set=rho_atom)
     609         102 :          natom = SIZE(rho_atom)
     610             :       ELSE
     611             :          natom = 0
     612             :       END IF
     613             : 
     614         794 :       IF (nspin == 2) THEN
     615         100 :          CALL pw_axpy(rho_g(1), rho_g(2), 1.0_dp, -1.0_dp)
     616     9659620 :          DO ig = 1, ng
     617     9659620 :             mixing_store%rhoin(2)%cc(ig) = mixing_store%rhoin(1)%cc(ig) - mixing_store%rhoin(2)%cc(ig)
     618             :          END DO
     619         100 :          IF (gapw .AND. mixing_store%gmix_p) THEN
     620          90 :             DO iatom = 1, natom
     621          90 :                IF (mixing_store%paw(iatom)) THEN
     622       60560 :                   rho_atom(iatom)%cpc_h(2)%r_coef = rho_atom(iatom)%cpc_h(1)%r_coef - rho_atom(iatom)%cpc_h(2)%r_coef
     623       60560 :                   rho_atom(iatom)%cpc_s(2)%r_coef = rho_atom(iatom)%cpc_s(1)%r_coef - rho_atom(iatom)%cpc_s(2)%r_coef
     624             :                   mixing_store%cpc_h_in(iatom, 2)%r_coef = mixing_store%cpc_h_in(iatom, 1)%r_coef - &
     625       60560 :                                                            mixing_store%cpc_h_in(iatom, 2)%r_coef
     626             :                   mixing_store%cpc_s_in(iatom, 2)%r_coef = mixing_store%cpc_s_in(iatom, 1)%r_coef - &
     627       60560 :                                                            mixing_store%cpc_s_in(iatom, 2)%r_coef
     628             :                END IF
     629             :             END DO
     630             :          END IF
     631             :       END IF
     632             : 
     633        1688 :       DO ispin = 1, nspin
     634    68680020 :          res_rho = CMPLX(0.0_dp, 0.0_dp, KIND=dp)
     635    68680020 :          DO ig = 1, ng
     636    68680020 :             res_rho(ig) = rho_g(ispin)%array(ig) - mixing_store%rhoin(ispin)%cc(ig)
     637             :          END DO
     638             : 
     639        1688 :          IF (only_kerker) THEN
     640     8257236 :             DO ig = 1, ng
     641     8257062 :                mixing_store%last_res(ispin)%cc(ig) = res_rho(ig)
     642     8257062 :                f_mix = alpha*mixing_store%kerker_factor(ig)
     643     8257062 :                rho_g(ispin)%array(ig) = mixing_store%rhoin(ispin)%cc(ig) + f_mix*res_rho(ig)
     644     8257062 :                mixing_store%rhoin_old(ispin)%cc(ig) = mixing_store%rhoin(ispin)%cc(ig)
     645     8257236 :                mixing_store%rhoin(ispin)%cc(ig) = rho_g(ispin)%array(ig)
     646             :             END DO
     647             : 
     648         174 :             IF (mixing_store%gmix_p) THEN
     649          24 :                IF (gapw) THEN
     650         216 :                   DO iatom = 1, natom
     651         216 :                      IF (mixing_store%paw(iatom)) THEN
     652             :                         mixing_store%cpc_h_lastres(iatom, ispin)%r_coef = rho_atom(iatom)%cpc_h(ispin)%r_coef - &
     653      245780 :                                                                           mixing_store%cpc_h_in(iatom, ispin)%r_coef
     654             :                         mixing_store%cpc_s_lastres(iatom, ispin)%r_coef = rho_atom(iatom)%cpc_s(ispin)%r_coef - &
     655      245780 :                                                                           mixing_store%cpc_s_in(iatom, ispin)%r_coef
     656             : 
     657             :                         rho_atom(iatom)%cpc_h(ispin)%r_coef = alpha*rho_atom(iatom)%cpc_h(ispin)%r_coef + &
     658      245780 :                                                               mixing_store%cpc_h_in(iatom, ispin)%r_coef*(1._dp - alpha)
     659             :                         rho_atom(iatom)%cpc_s(ispin)%r_coef = alpha*rho_atom(iatom)%cpc_s(ispin)%r_coef + &
     660      245780 :                                                               mixing_store%cpc_s_in(iatom, ispin)%r_coef*(1._dp - alpha)
     661             : 
     662      122972 :                         mixing_store%cpc_h_old(iatom, ispin)%r_coef = mixing_store%cpc_h_in(iatom, ispin)%r_coef
     663      122972 :                         mixing_store%cpc_s_old(iatom, ispin)%r_coef = mixing_store%cpc_s_in(iatom, ispin)%r_coef
     664      245780 :                         mixing_store%cpc_h_in(iatom, ispin)%r_coef = rho_atom(iatom)%cpc_h(ispin)%r_coef
     665      245780 :                         mixing_store%cpc_s_in(iatom, ispin)%r_coef = rho_atom(iatom)%cpc_s(ispin)%r_coef
     666             :                      END IF
     667             :                   END DO
     668             :                END IF
     669             :             END IF
     670             :          ELSE
     671             : 
     672        2880 :             ALLOCATE (a(nb, nb))
     673       20148 :             a = 0.0_dp
     674        2160 :             ALLOCATE (b(nb, nb))
     675       20148 :             b = 0.0_dp
     676        2160 :             ALLOCATE (c(nb))
     677        3528 :             c = 0.0_dp
     678        1440 :             ALLOCATE (g(nb))
     679        3528 :             g = 0.0_dp
     680             : 
     681         720 :             delta_norm = 0.0_dp
     682         720 :             res_norm = 0.0_dp
     683         720 :             rho_norm = 0.0_dp
     684    60422784 :             DO ig = 1, ng
     685    60422064 :                mixing_store%res_buffer(ib, ispin)%cc(ig) = res_rho(ig) - mixing_store%last_res(ispin)%cc(ig)
     686    60422064 :                mixing_store%last_res(ispin)%cc(ig) = res_rho(ig)
     687             :                res_norm = res_norm + &
     688             :                           REAL(res_rho(ig), dp)*REAL(res_rho(ig), dp) + &
     689    60422064 :                           AIMAG(res_rho(ig))*AIMAG(res_rho(ig))
     690             :                delta_norm = delta_norm + &
     691             :                             REAL(mixing_store%res_buffer(ib, ispin)%cc(ig), dp)* &
     692             :                             REAL(mixing_store%res_buffer(ib, ispin)%cc(ig), dp) + &
     693             :                             AIMAG(mixing_store%res_buffer(ib, ispin)%cc(ig))* &
     694    60422064 :                             AIMAG(mixing_store%res_buffer(ib, ispin)%cc(ig))
     695             :                rho_norm = rho_norm + &
     696             :                           REAL(rho_g(ispin)%array(ig), dp)*REAL(rho_g(ispin)%array(ig), dp) + &
     697    60422784 :                           AIMAG(rho_g(ispin)%array(ig))*AIMAG(rho_g(ispin)%array(ig))
     698             :             END DO
     699    60422784 :             DO ig = 1, ng
     700             :                mixing_store%drho_buffer(ib, ispin)%cc(ig) = &
     701             :                   mixing_store%rhoin(ispin)%cc(ig) - &
     702    60422784 :                   mixing_store%rhoin_old(ispin)%cc(ig)
     703             :             END DO
     704         720 :             CALL para_env%sum(delta_norm)
     705         720 :             delta_norm = SQRT(delta_norm)
     706         720 :             CALL para_env%sum(res_norm)
     707         720 :             res_norm = SQRT(res_norm)
     708         720 :             CALL para_env%sum(rho_norm)
     709         720 :             rho_norm = SQRT(rho_norm)
     710             : 
     711         720 :             IF (res_norm > 1.E-14_dp) THEN
     712    60422784 :                mixing_store%res_buffer(ib, ispin)%cc(:) = mixing_store%res_buffer(ib, ispin)%cc(:)/delta_norm
     713    60422784 :                mixing_store%drho_buffer(ib, ispin)%cc(:) = mixing_store%drho_buffer(ib, ispin)%cc(:)/delta_norm
     714             : 
     715         720 :                IF (mixing_store%gmix_p .AND. gapw) THEN
     716         252 :                   DO iatom = 1, natom
     717         252 :                      IF (mixing_store%paw(iatom)) THEN
     718          28 :                         n1 = SIZE(mixing_store%cpc_s_in(iatom, ispin)%r_coef, 1)
     719          28 :                         n2 = SIZE(mixing_store%cpc_s_in(iatom, ispin)%r_coef, 2)
     720         616 :                         DO i = 1, n2
     721       12964 :                            DO j = 1, n1
     722             :                               mixing_store%dcpc_h_in(ib, iatom, ispin)%r_coef(j, i) = &
     723             :                                  (mixing_store%cpc_h_in(iatom, ispin)%r_coef(j, i) - &
     724       12348 :                                   mixing_store%cpc_h_old(iatom, ispin)%r_coef(j, i))/delta_norm
     725             :                               dcpc_h_res = ((rho_atom(iatom)%cpc_h(ispin)%r_coef(j, i) - &
     726             :                                              mixing_store%cpc_h_in(iatom, ispin)%r_coef(j, i)) - &
     727       12348 :                                             mixing_store%cpc_h_lastres(iatom, ispin)%r_coef(j, i))/delta_norm
     728             :                               mixing_store%cpc_h_lastres(iatom, ispin)%r_coef(j, i) = rho_atom(iatom)%cpc_h(ispin)%r_coef(j, i) - &
     729       12348 :                                                                                     mixing_store%cpc_h_in(iatom, ispin)%r_coef(j, i)
     730             : 
     731             :                               mixing_store%dcpc_s_in(ib, iatom, ispin)%r_coef(j, i) = &
     732             :                                  (mixing_store%cpc_s_in(iatom, ispin)%r_coef(j, i) - &
     733       12348 :                                   mixing_store%cpc_s_old(iatom, ispin)%r_coef(j, i))/delta_norm
     734             :                               dcpc_s_res = ((rho_atom(iatom)%cpc_s(ispin)%r_coef(j, i) - &
     735             :                                              mixing_store%cpc_s_in(iatom, ispin)%r_coef(j, i)) - &
     736       12348 :                                             mixing_store%cpc_s_lastres(iatom, ispin)%r_coef(j, i))/delta_norm
     737             :                               mixing_store%cpc_s_lastres(iatom, ispin)%r_coef(j, i) = rho_atom(iatom)%cpc_s(ispin)%r_coef(j, i) - &
     738       12348 :                                                                                     mixing_store%cpc_s_in(iatom, ispin)%r_coef(j, i)
     739             : 
     740             :                               mixing_store%dcpc_h_in(ib, iatom, ispin)%r_coef(j, i) = &
     741             :                                  alpha*dcpc_h_res + &
     742       12348 :                                  mixing_store%dcpc_h_in(ib, iatom, ispin)%r_coef(j, i)
     743             :                               mixing_store%dcpc_s_in(ib, iatom, ispin)%r_coef(j, i) = &
     744             :                                  alpha*dcpc_s_res + &
     745       12936 :                                  mixing_store%dcpc_s_in(ib, iatom, ispin)%r_coef(j, i)
     746             :                            END DO
     747             :                         END DO
     748             :                      END IF
     749             :                   END DO
     750             :                END IF
     751             : 
     752       20148 :                a(:, :) = 0.0_dp
     753    60422784 :                DO ig = 1, ng
     754    60422064 :                   f_mix = alpha*mixing_store%kerker_factor(ig)
     755             :                   mixing_store%drho_buffer(ib, ispin)%cc(ig) = &
     756             :                      f_mix*mixing_store%res_buffer(ib, ispin)%cc(ig) + &
     757    60422784 :                      mixing_store%drho_buffer(ib, ispin)%cc(ig)
     758             :                END DO
     759        3528 :                DO jb = 1, nb
     760       13242 :                   DO kb = jb, nb
     761  1475240898 :                      DO ig = 1, mixing_store%ig_max !ng
     762             :                         a(kb, jb) = a(kb, jb) + mixing_store%p_metric(ig)*( &
     763             :                                     REAL(mixing_store%res_buffer(jb, ispin)%cc(ig), dp)* &
     764             :                                     REAL(mixing_store%res_buffer(kb, ispin)%cc(ig), dp) + &
     765             :                                     AIMAG(mixing_store%res_buffer(jb, ispin)%cc(ig))* &
     766  1475240898 :                                     AIMAG(mixing_store%res_buffer(kb, ispin)%cc(ig)))
     767             :                      END DO
     768       12522 :                      a(jb, kb) = a(kb, jb)
     769             :                   END DO
     770             :                END DO
     771         720 :                CALL para_env%sum(a)
     772             : 
     773        3528 :                C = 0.0_dp
     774        3528 :                DO jb = 1, nb
     775        2808 :                   a(jb, jb) = w0 + a(jb, jb)
     776   326542848 :                   DO ig = 1, mixing_store%ig_max !ng
     777             :                      c(jb) = c(jb) + mixing_store%p_metric(ig)*( &
     778             :                              REAL(mixing_store%res_buffer(jb, ispin)%cc(ig), dp)*REAL(res_rho(ig), dp) + &
     779   326542128 :                              AIMAG(mixing_store%res_buffer(jb, ispin)%cc(ig))*AIMAG(res_rho(ig)))
     780             :                   END DO
     781             :                END DO
     782         720 :                CALL para_env%sum(c)
     783         720 :                CALL invert_matrix(a, b, inv_err)
     784             : 
     785         720 :                CALL dgemv('T', nb, nb, 1.0_dp, B, nb, C, 1, 0.0_dp, G, 1)
     786             :             ELSE
     787           0 :                G = 0.0_dp
     788             :             END IF
     789             : 
     790    60422784 :             DO ig = 1, ng
     791    60422064 :                cc_mix = CMPLX(0.0_dp, 0.0_dp, kind=dp)
     792   386961384 :                DO jb = 1, nb
     793   386961384 :                   cc_mix = cc_mix - G(jb)*mixing_store%drho_buffer(jb, ispin)%cc(ig)
     794             :                END DO
     795    60422064 :                f_mix = alpha*mixing_store%kerker_factor(ig)
     796             : 
     797    60422064 :                IF (res_norm > 1.E-14_dp) rho_g(ispin)%array(ig) = mixing_store%rhoin(ispin)%cc(ig) + &
     798    60422064 :                                                                   f_mix*res_rho(ig) + cc_mix
     799    60422064 :                mixing_store%rhoin_old(ispin)%cc(ig) = mixing_store%rhoin(ispin)%cc(ig)
     800    60422784 :                mixing_store%rhoin(ispin)%cc(ig) = rho_g(ispin)%array(ig)
     801             :             END DO
     802             : 
     803         720 :             IF (mixing_store%gmix_p) THEN
     804          28 :                IF (gapw) THEN
     805         252 :                   DO iatom = 1, natom
     806         252 :                      IF (mixing_store%paw(iatom)) THEN
     807          28 :                         n1 = SIZE(mixing_store%cpc_s_in(iatom, ispin)%r_coef, 1)
     808          28 :                         n2 = SIZE(mixing_store%cpc_s_in(iatom, ispin)%r_coef, 2)
     809         616 :                         DO i = 1, n2
     810       12964 :                            DO j = 1, n1
     811       12348 :                               valh = 0.0_dp
     812       12348 :                               vals = 0.0_dp
     813       61740 :                               DO jb = 1, nb
     814       49392 :                                  valh = valh - G(jb)*mixing_store%dcpc_h_in(jb, iatom, ispin)%r_coef(j, i)
     815       61740 :                                  vals = vals - G(jb)*mixing_store%dcpc_s_in(jb, iatom, ispin)%r_coef(j, i)
     816             :                               END DO
     817       12936 :                               IF (res_norm > 1.E-14_dp) THEN
     818             :                                  rho_atom(iatom)%cpc_h(ispin)%r_coef(j, i) = &
     819             :                                     alpha*rho_atom(iatom)%cpc_h(ispin)%r_coef(j, i) + &
     820       12348 :                                     mixing_store%cpc_h_in(iatom, ispin)%r_coef(j, i)*(1._dp - alpha) + valh
     821             :                                  rho_atom(iatom)%cpc_s(ispin)%r_coef(j, i) = &
     822             :                                     alpha*rho_atom(iatom)%cpc_s(ispin)%r_coef(j, i) + &
     823       12348 :                                     mixing_store%cpc_s_in(iatom, ispin)%r_coef(j, i)*(1._dp - alpha) + vals
     824             :                               END IF
     825             :                            END DO
     826             :                         END DO
     827             : 
     828       12964 :                         mixing_store%cpc_h_old(iatom, ispin)%r_coef = mixing_store%cpc_h_in(iatom, ispin)%r_coef
     829       12964 :                         mixing_store%cpc_s_old(iatom, ispin)%r_coef = mixing_store%cpc_s_in(iatom, ispin)%r_coef
     830       25900 :                         mixing_store%cpc_h_in(iatom, ispin)%r_coef = rho_atom(iatom)%cpc_h(ispin)%r_coef
     831       25900 :                         mixing_store%cpc_s_in(iatom, ispin)%r_coef = rho_atom(iatom)%cpc_s(ispin)%r_coef
     832             :                      END IF
     833             :                   END DO
     834             :                END IF
     835             :             END IF
     836             : 
     837         720 :             DEALLOCATE (a, b, c, g)
     838             :          END IF
     839             : 
     840             :       END DO ! ispin
     841         794 :       IF (nspin == 2) THEN
     842         100 :          CALL pw_axpy(rho_g(1), rho_g(2), 1.0_dp, -1.0_dp)
     843     9659620 :          DO ig = 1, ng
     844     9659620 :             mixing_store%rhoin(2)%cc(ig) = mixing_store%rhoin(1)%cc(ig) - mixing_store%rhoin(2)%cc(ig)
     845             :          END DO
     846         100 :          IF (gapw .AND. mixing_store%gmix_p) THEN
     847          90 :             DO iatom = 1, natom
     848          90 :                IF (mixing_store%paw(iatom)) THEN
     849       60560 :                   rho_atom(iatom)%cpc_h(2)%r_coef = rho_atom(iatom)%cpc_h(1)%r_coef - rho_atom(iatom)%cpc_h(2)%r_coef
     850       60560 :                   rho_atom(iatom)%cpc_s(2)%r_coef = rho_atom(iatom)%cpc_s(1)%r_coef - rho_atom(iatom)%cpc_s(2)%r_coef
     851             :                   mixing_store%cpc_h_in(iatom, 2)%r_coef = mixing_store%cpc_h_in(iatom, 1)%r_coef - &
     852       60560 :                                                            mixing_store%cpc_h_in(iatom, 2)%r_coef
     853             :                   mixing_store%cpc_s_in(iatom, 2)%r_coef = mixing_store%cpc_s_in(iatom, 1)%r_coef - &
     854       60560 :                                                            mixing_store%cpc_s_in(iatom, 2)%r_coef
     855             :                END IF
     856             :             END DO
     857             :          END IF
     858             : 
     859             :       END IF
     860             : 
     861         794 :       DEALLOCATE (res_rho)
     862             : 
     863         794 :       CALL timestop(handle)
     864             : 
     865        2382 :    END SUBROUTINE broyden_mixing
     866             : 
     867             : ! **************************************************************************************************
     868             : !> \brief Multisecant scheme to perform charge density Mixing
     869             : !>        as preconditioner we use the Kerker damping factor
     870             : !>        The mixing is applied directly on the density in reciprocal space,
     871             : !>        therefore it affects the potentials
     872             : !>        on the grid but not the density matrix
     873             : !> \param mixing_store ...
     874             : !> \param rho ...
     875             : !> \param para_env ...
     876             : !> \par History
     877             : !>      05.2009
     878             : !> \author MI
     879             : ! **************************************************************************************************
     880           0 :    SUBROUTINE multisecant_mixing(mixing_store, rho, para_env)
     881             : 
     882             :       TYPE(mixing_storage_type), POINTER                 :: mixing_store
     883             :       TYPE(qs_rho_type), POINTER                         :: rho
     884             :       TYPE(mp_para_env_type), POINTER                    :: para_env
     885             : 
     886             :       CHARACTER(len=*), PARAMETER :: routineN = 'multisecant_mixing'
     887             :       COMPLEX(KIND=dp), PARAMETER                        :: cmone = (-1.0_dp, 0.0_dp), &
     888             :                                                             cone = (1.0_dp, 0.0_dp), &
     889             :                                                             czero = (0.0_dp, 0.0_dp)
     890             : 
     891             :       COMPLEX(dp)                                        :: saa, yaa
     892           0 :       COMPLEX(dp), ALLOCATABLE, DIMENSION(:)             :: gn_global, pgn, res_matrix_global, &
     893           0 :                                                             tmp_vec, ugn
     894           0 :       COMPLEX(dp), ALLOCATABLE, DIMENSION(:, :)          :: a_matrix, res_matrix, sa, step_matrix, ya
     895           0 :       COMPLEX(dp), DIMENSION(:), POINTER                 :: gn
     896             :       INTEGER                                            :: handle, ib, ib_next, ib_prev, ig, &
     897             :                                                             ig_global, iig, ispin, jb, kb, nb, &
     898             :                                                             nbuffer, ng, ng_global, nspin
     899             :       LOGICAL                                            :: use_zgemm, use_zgemm_rev
     900             :       REAL(dp) :: alpha, f_mix, gn_norm, gn_norm_old, inv_err, n_low, n_up, pgn_norm, prec, &
     901             :          r_step, reg_par, sigma_max, sigma_tilde, step_size
     902           0 :       REAL(dp), ALLOCATABLE, DIMENSION(:)                :: norm_res, norm_res_low, norm_res_up
     903           0 :       REAL(dp), ALLOCATABLE, DIMENSION(:, :)             :: b_matrix, binv_matrix
     904           0 :       REAL(dp), DIMENSION(:), POINTER                    :: g2
     905             :       REAL(dp), SAVE                                     :: sigma_old = 1.0_dp
     906           0 :       TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER        :: rho_g
     907             : 
     908           0 :       CPASSERT(ASSOCIATED(mixing_store))
     909             : 
     910           0 :       CALL timeset(routineN, handle)
     911             : 
     912           0 :       NULLIFY (gn, rho_g)
     913             : 
     914           0 :       use_zgemm = .FALSE.
     915           0 :       use_zgemm_rev = .TRUE.
     916             : 
     917             :       ! prepare the parameters
     918           0 :       CALL qs_rho_get(rho, rho_g=rho_g)
     919             : 
     920           0 :       nspin = SIZE(rho_g, 1)
     921             :       ! not implemented for large grids.
     922           0 :       CPASSERT(rho_g(1)%pw_grid%ngpts < HUGE(ng_global))
     923           0 :       ng_global = INT(rho_g(1)%pw_grid%ngpts)
     924           0 :       ng = SIZE(mixing_store%rhoin_buffer(1, 1)%cc)
     925           0 :       alpha = mixing_store%alpha
     926             : 
     927           0 :       sigma_max = mixing_store%sigma_max
     928           0 :       reg_par = mixing_store%reg_par
     929           0 :       r_step = mixing_store%r_step
     930           0 :       nbuffer = mixing_store%nbuffer
     931             : 
     932             :       ! determine the step number, and multisecant iteration
     933           0 :       nb = MIN(mixing_store%ncall, nbuffer - 1)
     934           0 :       ib = MODULO(mixing_store%ncall, nbuffer) + 1
     935           0 :       IF (mixing_store%ncall > 0) THEN
     936           0 :          ib_prev = MODULO(mixing_store%ncall - 1, nbuffer) + 1
     937             :       ELSE
     938             :          ib_prev = 0
     939             :       END IF
     940           0 :       mixing_store%ncall = mixing_store%ncall + 1
     941           0 :       ib_next = MODULO(mixing_store%ncall, nbuffer) + 1
     942             : 
     943             :       ! compute the residual gn and its norm gn_norm
     944           0 :       DO ispin = 1, nspin
     945           0 :          gn => mixing_store%res_buffer(ib, ispin)%cc
     946           0 :          gn_norm = 0.0_dp
     947           0 :          DO ig = 1, ng
     948           0 :             gn(ig) = (rho_g(ispin)%array(ig) - mixing_store%rhoin_buffer(ib, ispin)%cc(ig))
     949             :             gn_norm = gn_norm + &
     950           0 :                       REAL(gn(ig), dp)*REAL(gn(ig), dp) + AIMAG(gn(ig))*AIMAG(gn(ig))
     951             :          END DO
     952           0 :          CALL para_env%sum(gn_norm)
     953           0 :          gn_norm = SQRT(gn_norm)
     954           0 :          mixing_store%norm_res_buffer(ib, ispin) = gn_norm
     955             :       END DO
     956             : 
     957           0 :       IF (nb == 0) THEN
     958             :          !simple mixing
     959           0 :          DO ispin = 1, nspin
     960           0 :             DO ig = 1, ng
     961           0 :                f_mix = alpha*mixing_store%kerker_factor(ig)
     962             :                rho_g(ispin)%array(ig) = mixing_store%rhoin_buffer(1, ispin)%cc(ig) + &
     963           0 :                                         f_mix*mixing_store%res_buffer(1, ispin)%cc(ig)
     964           0 :                mixing_store%rhoin_buffer(ib_next, ispin)%cc(ig) = rho_g(ispin)%array(ig)
     965             :             END DO
     966             :          END DO
     967           0 :          CALL timestop(handle)
     968           0 :          RETURN
     969             :       END IF
     970             : 
     971             :       ! allocate temporary arrays
     972             :       ! step_matrix  S ngxnb
     973           0 :       ALLOCATE (step_matrix(ng, nb))
     974             :       ! res_matrix Y  ngxnb
     975           0 :       ALLOCATE (res_matrix(ng, nb))
     976             :       ! matrix A  nbxnb
     977           0 :       ALLOCATE (a_matrix(nb, ng_global))
     978             :       ! PSI nb vector of norms
     979           0 :       ALLOCATE (norm_res(nb))
     980           0 :       ALLOCATE (norm_res_low(nb))
     981           0 :       ALLOCATE (norm_res_up(nb))
     982             : 
     983             :       ! matrix B   nbxnb
     984           0 :       ALLOCATE (b_matrix(nb, nb))
     985             :       ! matrix B_inv   nbxnb
     986           0 :       ALLOCATE (binv_matrix(nb, nb))
     987             : 
     988           0 :       ALLOCATE (gn_global(ng_global))
     989           0 :       ALLOCATE (res_matrix_global(ng_global))
     990             :       IF (use_zgemm) THEN
     991             :          ALLOCATE (sa(ng, ng_global))
     992             :          ALLOCATE (ya(ng, ng_global))
     993             :       END IF
     994             :       IF (use_zgemm_rev) THEN
     995           0 :          ALLOCATE (tmp_vec(nb))
     996             :       END IF
     997           0 :       ALLOCATE (pgn(ng))
     998           0 :       ALLOCATE (ugn(ng))
     999             : 
    1000           0 :       DO ispin = 1, nspin
    1001             :          ! generate the global vector with the present residual
    1002           0 :          gn => mixing_store%res_buffer(ib, ispin)%cc
    1003           0 :          gn_global = CMPLX(0.0_dp, 0.0_dp, KIND=dp)
    1004           0 :          DO ig = 1, ng
    1005           0 :             ig_global = mixing_store%ig_global_index(ig)
    1006           0 :             gn_global(ig_global) = gn(ig)
    1007             :          END DO
    1008           0 :          CALL para_env%sum(gn_global)
    1009             : 
    1010             :          ! compute steps (matrix S) and residual differences (matrix Y)
    1011             :          ! with respect to the present
    1012             :          ! step and the present residual (use stored rho_in and res_buffer)
    1013             : 
    1014             :          ! These quantities are pre-conditioned by means of Kerker factor multipliccation
    1015           0 :          g2 => rho_g(1)%pw_grid%gsq
    1016             : 
    1017           0 :          DO jb = 1, nb + 1
    1018           0 :             IF (jb < ib) THEN
    1019             :                kb = jb
    1020           0 :             ELSEIF (jb > ib) THEN
    1021           0 :                kb = jb - 1
    1022             :             ELSE
    1023             :                CYCLE
    1024             :             END IF
    1025           0 :             norm_res(kb) = 0.0_dp
    1026           0 :             norm_res_low(kb) = 0.0_dp
    1027           0 :             norm_res_up(kb) = 0.0_dp
    1028           0 :             DO ig = 1, ng
    1029             : 
    1030           0 :                prec = 1.0_dp
    1031             : 
    1032             :                step_matrix(ig, kb) = prec*(mixing_store%rhoin_buffer(jb, ispin)%cc(ig) - &
    1033           0 :                                            mixing_store%rhoin_buffer(ib, ispin)%cc(ig))
    1034             :                res_matrix(ig, kb) = (mixing_store%res_buffer(jb, ispin)%cc(ig) - &
    1035           0 :                                      mixing_store%res_buffer(ib, ispin)%cc(ig))
    1036             :                norm_res(kb) = norm_res(kb) + REAL(res_matrix(ig, kb), dp)*REAL(res_matrix(ig, kb), dp) + &
    1037           0 :                               AIMAG(res_matrix(ig, kb))*AIMAG(res_matrix(ig, kb))
    1038           0 :                IF (g2(ig) < 4.0_dp) THEN
    1039             :                   norm_res_low(kb) = norm_res_low(kb) + &
    1040             :                                      REAL(res_matrix(ig, kb), dp)*REAL(res_matrix(ig, kb), dp) + &
    1041           0 :                                      AIMAG(res_matrix(ig, kb))*AIMAG(res_matrix(ig, kb))
    1042             :                ELSE
    1043             :                   norm_res_up(kb) = norm_res_up(kb) + &
    1044             :                                     REAL(res_matrix(ig, kb), dp)*REAL(res_matrix(ig, kb), dp) + &
    1045           0 :                                     AIMAG(res_matrix(ig, kb))*AIMAG(res_matrix(ig, kb))
    1046             :                END IF
    1047           0 :                res_matrix(ig, kb) = prec*res_matrix(ig, kb)
    1048             :             END DO
    1049             :          END DO !jb
    1050             : 
    1051             :          ! normalize each column of S and Y => Snorm Ynorm
    1052           0 :          CALL para_env%sum(norm_res)
    1053           0 :          CALL para_env%sum(norm_res_up)
    1054           0 :          CALL para_env%sum(norm_res_low)
    1055           0 :          norm_res(1:nb) = 1.0_dp/SQRT(norm_res(1:nb))
    1056             :          n_low = 0.0_dp
    1057             :          n_up = 0.0_dp
    1058           0 :          DO jb = 1, nb
    1059           0 :             n_low = n_low + norm_res_low(jb)/norm_res(jb)
    1060           0 :             n_up = n_up + norm_res_up(jb)/norm_res(jb)
    1061             :          END DO
    1062           0 :          DO ig = 1, ng
    1063           0 :             IF (g2(ig) > 4.0_dp) THEN
    1064           0 :                step_matrix(ig, 1:nb) = step_matrix(ig, 1:nb)*SQRT(n_low/n_up)
    1065           0 :                res_matrix(ig, 1:nb) = res_matrix(ig, 1:nb)*SQRT(n_low/n_up)
    1066             :             END IF
    1067             :          END DO
    1068           0 :          DO kb = 1, nb
    1069           0 :             step_matrix(1:ng, kb) = step_matrix(1:ng, kb)*norm_res(kb)
    1070           0 :             res_matrix(1:ng, kb) = res_matrix(1:ng, kb)*norm_res(kb)
    1071             :          END DO
    1072             : 
    1073             :          ! compute A as [(Ynorm^t Ynorm) + (alpha I)]^(-1) Ynorm^t
    1074             :          ! compute B
    1075           0 :          DO jb = 1, nb
    1076           0 :             DO kb = 1, nb
    1077           0 :                b_matrix(kb, jb) = 0.0_dp
    1078           0 :                DO ig = 1, ng
    1079             :                   ! it is assumed that summing over all G vector gives a real, because
    1080             :                   !  y(-G,kb) = (y(G,kb))*
    1081           0 :                   b_matrix(kb, jb) = b_matrix(kb, jb) + REAL(res_matrix(ig, kb)*res_matrix(ig, jb), dp)
    1082             :                END DO
    1083             :             END DO
    1084             :          END DO
    1085             : 
    1086           0 :          CALL para_env%sum(b_matrix)
    1087           0 :          DO jb = 1, nb
    1088           0 :             b_matrix(jb, jb) = b_matrix(jb, jb) + reg_par
    1089             :          END DO
    1090             :          ! invert B
    1091           0 :          CALL invert_matrix(b_matrix, binv_matrix, inv_err)
    1092             : 
    1093             :          ! A = Binv Ynorm^t
    1094           0 :          a_matrix = CMPLX(0.0_dp, 0.0_dp, KIND=dp)
    1095           0 :          DO kb = 1, nb
    1096           0 :             res_matrix_global = CMPLX(0.0_dp, 0.0_dp, KIND=dp)
    1097           0 :             DO ig = 1, ng
    1098           0 :                ig_global = mixing_store%ig_global_index(ig)
    1099           0 :                res_matrix_global(ig_global) = res_matrix(ig, kb)
    1100             :             END DO
    1101           0 :             CALL para_env%sum(res_matrix_global)
    1102             : 
    1103           0 :             DO jb = 1, nb
    1104           0 :                DO ig = 1, ng
    1105           0 :                   ig_global = mixing_store%ig_global_index(ig)
    1106             :                   a_matrix(jb, ig_global) = a_matrix(jb, ig_global) + &
    1107           0 :                                             binv_matrix(jb, kb)*res_matrix_global(ig_global)
    1108             :                END DO
    1109             :             END DO
    1110             :          END DO
    1111           0 :          CALL para_env%sum(a_matrix)
    1112             : 
    1113             :          ! compute the two components of gn that will be used to update rho
    1114           0 :          gn => mixing_store%res_buffer(ib, ispin)%cc
    1115           0 :          pgn_norm = 0.0_dp
    1116             : 
    1117             :          IF (use_zgemm) THEN
    1118             : 
    1119             :             CALL zgemm("N", "N", ng, ng_global, nb, cmone, step_matrix(1, 1), ng, &
    1120             :                        a_matrix(1, 1), nb, czero, sa(1, 1), ng)
    1121             :             CALL zgemm("N", "N", ng, ng_global, nb, cmone, res_matrix(1, 1), ng, &
    1122             :                        a_matrix(1, 1), nb, czero, ya(1, 1), ng)
    1123             :             DO ig = 1, ng
    1124             :                ig_global = mixing_store%ig_global_index(ig)
    1125             :                ya(ig, ig_global) = ya(ig, ig_global) + CMPLX(1.0_dp, 0.0_dp, KIND=dp)
    1126             :             END DO
    1127             : 
    1128             :             CALL zgemv("N", ng, ng_global, cone, sa(1, 1), &
    1129             :                        ng, gn_global(1), 1, czero, pgn(1), 1)
    1130             :             CALL zgemv("N", ng, ng_global, cone, ya(1, 1), &
    1131             :                        ng, gn_global(1), 1, czero, ugn(1), 1)
    1132             : 
    1133             :             DO ig = 1, ng
    1134             :                pgn_norm = pgn_norm + REAL(pgn(ig), dp)*REAL(pgn(ig), dp) + &
    1135             :                           AIMAG(pgn(ig))*AIMAG(pgn(ig))
    1136             :             END DO
    1137             :             CALL para_env%sum(pgn_norm)
    1138             :          ELSEIF (use_zgemm_rev) THEN
    1139             : 
    1140             :             CALL zgemv("N", nb, ng_global, cone, a_matrix(1, 1), &
    1141           0 :                        nb, gn_global(1), 1, czero, tmp_vec(1), 1)
    1142             : 
    1143             :             CALL zgemv("N", ng, nb, cmone, step_matrix(1, 1), ng, &
    1144           0 :                        tmp_vec(1), 1, czero, pgn(1), 1)
    1145             : 
    1146             :             CALL zgemv("N", ng, nb, cmone, res_matrix(1, 1), ng, &
    1147           0 :                        tmp_vec(1), 1, czero, ugn(1), 1)
    1148             : 
    1149           0 :             DO ig = 1, ng
    1150             :                pgn_norm = pgn_norm + REAL(pgn(ig), dp)*REAL(pgn(ig), dp) + &
    1151           0 :                           AIMAG(pgn(ig))*AIMAG(pgn(ig))
    1152           0 :                ugn(ig) = ugn(ig) + gn(ig)
    1153             :             END DO
    1154           0 :             CALL para_env%sum(pgn_norm)
    1155             : 
    1156             :          ELSE
    1157             :             DO ig = 1, ng
    1158             :                pgn(ig) = CMPLX(0.0_dp, 0.0_dp, KIND=dp)
    1159             :                ugn(ig) = CMPLX(0.0_dp, 0.0_dp, KIND=dp)
    1160             :                ig_global = mixing_store%ig_global_index(ig)
    1161             :                DO iig = 1, ng_global
    1162             :                   saa = CMPLX(0.0_dp, 0.0_dp, KIND=dp)
    1163             :                   yaa = CMPLX(0.0_dp, 0.0_dp, KIND=dp)
    1164             : 
    1165             :                   IF (ig_global == iig) yaa = CMPLX(1.0_dp, 0.0_dp, KIND=dp)
    1166             : 
    1167             :                   DO jb = 1, nb
    1168             :                      saa = saa - step_matrix(ig, jb)*a_matrix(jb, iig)
    1169             :                      yaa = yaa - res_matrix(ig, jb)*a_matrix(jb, iig)
    1170             :                   END DO
    1171             :                   pgn(ig) = pgn(ig) + saa*gn_global(iig)
    1172             :                   ugn(ig) = ugn(ig) + yaa*gn_global(iig)
    1173             :                END DO
    1174             :             END DO
    1175             :             DO ig = 1, ng
    1176             :                pgn_norm = pgn_norm + REAL(pgn(ig), dp)*REAL(pgn(ig), dp) + &
    1177             :                           AIMAG(pgn(ig))*AIMAG(pgn(ig))
    1178             :             END DO
    1179             :             CALL para_env%sum(pgn_norm)
    1180             :          END IF
    1181             : 
    1182           0 :          gn_norm = mixing_store%norm_res_buffer(ib, ispin)
    1183           0 :          gn_norm_old = mixing_store%norm_res_buffer(ib_prev, ispin)
    1184             :          IF (ib_prev /= 0) THEN
    1185           0 :             sigma_tilde = sigma_old*MAX(0.5_dp, MIN(2.0_dp, gn_norm_old/gn_norm))
    1186             :          ELSE
    1187             :             sigma_tilde = 0.5_dp
    1188             :          END IF
    1189           0 :          sigma_tilde = 0.1_dp
    1190             :          ! Step size for the unpredicted component
    1191           0 :          step_size = MIN(sigma_tilde, r_step*pgn_norm/gn_norm, sigma_max)
    1192           0 :          sigma_old = step_size
    1193             : 
    1194             :          ! update the density
    1195           0 :          DO ig = 1, ng
    1196           0 :             prec = mixing_store%kerker_factor(ig)
    1197             :             rho_g(ispin)%array(ig) = mixing_store%rhoin_buffer(ib, ispin)%cc(ig) &
    1198           0 :                                      - prec*step_size*ugn(ig) + prec*pgn(ig) ! - 0.1_dp * prec* gn(ig)
    1199           0 :             mixing_store%rhoin_buffer(ib_next, ispin)%cc(ig) = rho_g(ispin)%array(ig)
    1200             :          END DO
    1201             : 
    1202             :       END DO ! ispin
    1203             : 
    1204             :       ! Deallocate  temporary arrays
    1205           0 :       DEALLOCATE (step_matrix, res_matrix)
    1206           0 :       DEALLOCATE (norm_res)
    1207           0 :       DEALLOCATE (norm_res_up)
    1208           0 :       DEALLOCATE (norm_res_low)
    1209           0 :       DEALLOCATE (a_matrix, b_matrix, binv_matrix)
    1210           0 :       DEALLOCATE (ugn, pgn)
    1211             :       IF (use_zgemm) THEN
    1212             :          DEALLOCATE (sa, ya)
    1213             :       END IF
    1214             :       IF (use_zgemm_rev) THEN
    1215           0 :          DEALLOCATE (tmp_vec)
    1216             :       END IF
    1217           0 :       DEALLOCATE (gn_global, res_matrix_global)
    1218             : 
    1219           0 :       CALL timestop(handle)
    1220             : 
    1221           0 :    END SUBROUTINE multisecant_mixing
    1222             : 
    1223             : END MODULE qs_gspace_mixing

Generated by: LCOV version 1.15