LCOV - code coverage report
Current view: top level - src - hartree_local_methods.F (source / functions) Coverage Total Hit
Test: CP2K Regtests (git:cccd2f3) Lines: 99.2 % 358 355
Test Date: 2026-05-06 07:07:47 Functions: 100.0 % 6 6

            Line data    Source code
       1              : !--------------------------------------------------------------------------------------------------!
       2              : !   CP2K: A general program to perform molecular dynamics simulations                              !
       3              : !   Copyright 2000-2026 CP2K developers group <https://cp2k.org>                                   !
       4              : !                                                                                                  !
       5              : !   SPDX-License-Identifier: GPL-2.0-or-later                                                      !
       6              : !--------------------------------------------------------------------------------------------------!
       7              : ! **************************************************************************************************
       8              : MODULE hartree_local_methods
       9              :    USE atomic_kind_types,               ONLY: atomic_kind_type,&
      10              :                                               get_atomic_kind
      11              :    USE basis_set_types,                 ONLY: get_gto_basis_set,&
      12              :                                               gto_basis_set_type
      13              :    USE cell_types,                      ONLY: cell_type
      14              :    USE cp_control_types,                ONLY: dft_control_type
      15              :    USE hartree_local_types,             ONLY: allocate_ecoul_1center,&
      16              :                                               ecoul_1center_type,&
      17              :                                               hartree_local_type,&
      18              :                                               set_ecoul_1c
      19              :    USE kinds,                           ONLY: dp
      20              :    USE mathconstants,                   ONLY: fourpi,&
      21              :                                               pi
      22              :    USE message_passing,                 ONLY: mp_para_env_type
      23              :    USE orbital_pointers,                ONLY: indso,&
      24              :                                               nsoset
      25              :    USE pw_env_types,                    ONLY: pw_env_get,&
      26              :                                               pw_env_type
      27              :    USE pw_poisson_types,                ONLY: pw_poisson_periodic,&
      28              :                                               pw_poisson_type
      29              :    USE qs_charges_types,                ONLY: qs_charges_type
      30              :    USE qs_cneo_methods,                 ONLY: Vh_1c_nuc_integrals,&
      31              :                                               calculate_rhoz_cneo
      32              :    USE qs_cneo_types,                   ONLY: cneo_potential_type,&
      33              :                                               rhoz_cneo_type
      34              :    USE qs_environment_types,            ONLY: get_qs_env,&
      35              :                                               qs_environment_type
      36              :    USE qs_grid_atom,                    ONLY: grid_atom_type
      37              :    USE qs_harmonics_atom,               ONLY: get_none0_cg_list,&
      38              :                                               harmonics_atom_type
      39              :    USE qs_kind_types,                   ONLY: get_qs_kind,&
      40              :                                               get_qs_kind_set,&
      41              :                                               qs_kind_type
      42              :    USE qs_local_rho_types,              ONLY: get_local_rho,&
      43              :                                               local_rho_type,&
      44              :                                               rhoz_type
      45              :    USE qs_rho0_types,                   ONLY: get_rho0_mpole,&
      46              :                                               rho0_atom_type,&
      47              :                                               rho0_mpole_type
      48              :    USE qs_rho_atom_types,               ONLY: get_rho_atom,&
      49              :                                               rho_atom_coeff,&
      50              :                                               rho_atom_type
      51              :    USE util,                            ONLY: get_limit
      52              : #include "./base/base_uses.f90"
      53              : 
      54              :    IMPLICIT NONE
      55              : 
      56              :    PRIVATE
      57              : 
      58              :    CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'hartree_local_methods'
      59              : 
      60              :    ! Public Subroutine
      61              : 
      62              :    PUBLIC :: init_coulomb_local, calculate_Vh_1center, Vh_1c_gg_integrals
      63              : 
      64              : CONTAINS
      65              : 
      66              : ! **************************************************************************************************
      67              : !> \brief ...
      68              : !> \param hartree_local ...
      69              : !> \param natom ...
      70              : ! **************************************************************************************************
      71         3026 :    SUBROUTINE init_coulomb_local(hartree_local, natom)
      72              : 
      73              :       TYPE(hartree_local_type), POINTER                  :: hartree_local
      74              :       INTEGER, INTENT(IN)                                :: natom
      75              : 
      76              :       CHARACTER(len=*), PARAMETER :: routineN = 'init_coulomb_local'
      77              : 
      78              :       INTEGER                                            :: handle
      79         3026 :       TYPE(ecoul_1center_type), DIMENSION(:), POINTER    :: ecoul_1c
      80              : 
      81         3026 :       CALL timeset(routineN, handle)
      82              : 
      83         3026 :       NULLIFY (ecoul_1c)
      84              :       !   Allocate and Initialize 1-center Potentials and Integrals
      85         3026 :       CALL allocate_ecoul_1center(ecoul_1c, natom)
      86         3026 :       hartree_local%ecoul_1c => ecoul_1c
      87              : 
      88         3026 :       CALL timestop(handle)
      89              : 
      90         3026 :    END SUBROUTINE init_coulomb_local
      91              : 
      92              : ! **************************************************************************************************
      93              : !> \brief Calculates Hartree potential for hard and soft densities (including
      94              : !>        nuclear charge and compensation charges) using numerical integration
      95              : !> \param vrad_h ...
      96              : !> \param vrad_s ...
      97              : !> \param rrad_h ...
      98              : !> \param rrad_s ...
      99              : !> \param rrad_0 ...
     100              : !> \param rrad_z ...
     101              : !> \param grid_atom ...
     102              : !> \par History
     103              : !>      05.2012 JGH refactoring
     104              : !> \author ??
     105              : ! **************************************************************************************************
     106           15 :    SUBROUTINE calculate_Vh_1center(vrad_h, vrad_s, rrad_h, rrad_s, rrad_0, rrad_z, grid_atom)
     107              : 
     108              :       REAL(dp), DIMENSION(:, :), INTENT(INOUT)           :: vrad_h, vrad_s
     109              :       TYPE(rho_atom_coeff), DIMENSION(:), INTENT(IN)     :: rrad_h, rrad_s
     110              :       REAL(dp), DIMENSION(:, :), INTENT(IN)              :: rrad_0
     111              :       REAL(dp), DIMENSION(:), INTENT(IN)                 :: rrad_z
     112              :       TYPE(grid_atom_type), POINTER                      :: grid_atom
     113              : 
     114              :       CHARACTER(len=*), PARAMETER :: routineN = 'calculate_Vh_1center'
     115              : 
     116              :       INTEGER                                            :: handle, ir, iso, ispin, l_ang, &
     117              :                                                             max_s_harm, nchannels, nr, nspins
     118              :       REAL(dp)                                           :: I1_down, I1_up, I2_down, I2_up, prefactor
     119           15 :       REAL(dp), ALLOCATABLE, DIMENSION(:)                :: rho_1, rho_2
     120           15 :       REAL(dp), DIMENSION(:), POINTER                    :: wr
     121           15 :       REAL(dp), DIMENSION(:, :), POINTER                 :: oor2l, r2l
     122              : 
     123           15 :       CALL timeset(routineN, handle)
     124              : 
     125           15 :       nr = grid_atom%nr
     126           15 :       max_s_harm = SIZE(vrad_h, 2)
     127           15 :       nspins = SIZE(rrad_h, 1)
     128           15 :       nchannels = SIZE(rrad_0, 2)
     129              : 
     130           15 :       r2l => grid_atom%rad2l
     131           15 :       oor2l => grid_atom%oorad2l
     132           15 :       wr => grid_atom%wr
     133              : 
     134           60 :       ALLOCATE (rho_1(nr), rho_2(nr))
     135              : 
     136          198 :       DO iso = 1, max_s_harm
     137         9333 :          rho_1(:) = 0.0_dp
     138         9333 :          rho_2(:) = 0.0_dp
     139          948 :          IF (iso == 1) rho_1(:) = rrad_z(:)
     140         6933 :          IF (iso <= nchannels) rho_2(:) = rrad_0(:, iso)
     141          549 :          DO ispin = 1, nspins
     142        18666 :             rho_1(:) = rho_1(:) + rrad_h(ispin)%r_coef(:, iso)
     143        18849 :             rho_2(:) = rho_2(:) + rrad_s(ispin)%r_coef(:, iso)
     144              :          END DO
     145              : 
     146          183 :          l_ang = indso(1, iso)
     147          183 :          prefactor = fourpi/(2._dp*l_ang + 1._dp)
     148              : 
     149         9333 :          rho_1(:) = rho_1(:)*wr(:)
     150         9333 :          rho_2(:) = rho_2(:)*wr(:)
     151              : 
     152          183 :          I1_up = 0.0_dp
     153          183 :          I1_down = 0.0_dp
     154          183 :          I2_up = 0.0_dp
     155          183 :          I2_down = 0.0_dp
     156              : 
     157          183 :          I1_up = r2l(nr, l_ang)*rho_1(nr)
     158          183 :          I2_up = r2l(nr, l_ang)*rho_2(nr)
     159              : 
     160         9150 :          DO ir = nr - 1, 1, -1
     161         8967 :             I1_down = I1_down + oor2l(ir, l_ang + 1)*rho_1(ir)
     162         9150 :             I2_down = I2_down + oor2l(ir, l_ang + 1)*rho_2(ir)
     163              :          END DO
     164              : 
     165              :          vrad_h(nr, iso) = vrad_h(nr, iso) + prefactor* &
     166          183 :                            (oor2l(nr, l_ang + 1)*I1_up + r2l(nr, l_ang)*I1_down)
     167              :          vrad_s(nr, iso) = vrad_s(nr, iso) + prefactor* &
     168          183 :                            (oor2l(nr, l_ang + 1)*I2_up + r2l(nr, l_ang)*I2_down)
     169              : 
     170         9165 :          DO ir = nr - 1, 1, -1
     171         8967 :             I1_up = I1_up + r2l(ir, l_ang)*rho_1(ir)
     172         8967 :             I1_down = I1_down - oor2l(ir, l_ang + 1)*rho_1(ir)
     173         8967 :             I2_up = I2_up + r2l(ir, l_ang)*rho_2(ir)
     174         8967 :             I2_down = I2_down - oor2l(ir, l_ang + 1)*rho_2(ir)
     175              : 
     176              :             vrad_h(ir, iso) = vrad_h(ir, iso) + prefactor* &
     177         8967 :                               (oor2l(ir, l_ang + 1)*I1_up + r2l(ir, l_ang)*I1_down)
     178              :             vrad_s(ir, iso) = vrad_s(ir, iso) + prefactor* &
     179         9150 :                               (oor2l(ir, l_ang + 1)*I2_up + r2l(ir, l_ang)*I2_down)
     180              : 
     181              :          END DO
     182              : 
     183              :       END DO
     184              : 
     185           15 :       DEALLOCATE (rho_1, rho_2)
     186              : 
     187           15 :       CALL timestop(handle)
     188              : 
     189           15 :    END SUBROUTINE calculate_Vh_1center
     190              : 
     191              : ! **************************************************************************************************
     192              : !> \brief Calculates one center GAPW Hartree energies and matrix elements
     193              : !>        Hartree potentials are input
     194              : !>        Takes possible background charge into account
     195              : !>        Special case for densities without core charge
     196              : !> \param qs_env ...
     197              : !> \param energy_hartree_1c ...
     198              : !> \param ecoul_1c ...
     199              : !> \param local_rho_set ...
     200              : !> \param para_env ...
     201              : !> \param tddft ...
     202              : !> \param local_rho_set_2nd ...
     203              : !> \param core_2nd ...
     204              : !> \par History
     205              : !>      05.2012 JGH refactoring
     206              : !> \author ??
     207              : ! **************************************************************************************************
     208        24550 :    SUBROUTINE Vh_1c_gg_integrals(qs_env, energy_hartree_1c, ecoul_1c, local_rho_set, para_env, tddft, local_rho_set_2nd, &
     209              :                                  core_2nd)
     210              : 
     211              :       TYPE(qs_environment_type), POINTER                 :: qs_env
     212              :       REAL(kind=dp), INTENT(out)                         :: energy_hartree_1c
     213              :       TYPE(ecoul_1center_type), DIMENSION(:), POINTER    :: ecoul_1c
     214              :       TYPE(local_rho_type), POINTER                      :: local_rho_set
     215              :       TYPE(mp_para_env_type), POINTER                    :: para_env
     216              :       LOGICAL, INTENT(IN)                                :: tddft
     217              :       TYPE(local_rho_type), OPTIONAL, POINTER            :: local_rho_set_2nd
     218              :       LOGICAL, INTENT(IN), OPTIONAL                      :: core_2nd
     219              : 
     220              :       CHARACTER(LEN=*), PARAMETER :: routineN = 'Vh_1c_gg_integrals'
     221              : 
     222              :       INTEGER :: bo(2), handle, iat, iatom, ikind, ipgf1, is1, iset1, iso, l_ang, llmax, &
     223              :          llmax_nuc, lmax0, lmax0_2nd, lmax_0, m1, max_iso, max_iso_not0, max_iso_not0_nuc, &
     224              :          max_s_harm, max_s_harm_nuc, maxl, maxl_nuc, maxso, maxso_nuc, mepos, n1, nat, nchan_0, &
     225              :          nkind, nr, nset, nset_nuc, nsotot, nsotot_nuc, nspins, num_pe
     226        24550 :       INTEGER, ALLOCATABLE, DIMENSION(:)                 :: cg_n_list
     227        24550 :       INTEGER, ALLOCATABLE, DIMENSION(:, :, :)           :: cg_list
     228        24550 :       INTEGER, DIMENSION(:), POINTER                     :: atom_list, lmax, lmax_nuc, lmin, &
     229        24550 :                                                             lmin_nuc, npgf, npgf_nuc
     230              :       LOGICAL                                            :: cneo, core_charge, l_2nd_local_rho, &
     231              :                                                             my_core_2nd, my_periodic, paw_atom
     232              :       REAL(dp)                                           :: back_ch, ecoul_1_z_cneo, factor
     233        24550 :       REAL(dp), ALLOCATABLE, DIMENSION(:)                :: gexp, sqrtwr
     234        24550 :       REAL(dp), ALLOCATABLE, DIMENSION(:, :)             :: aVh1b_00, aVh1b_00_nuc, aVh1b_hh, &
     235        24550 :                                                             aVh1b_hh_nuc, aVh1b_ss, aVh1b_ss_nuc, &
     236        24550 :                                                             g0_h_w
     237        24550 :       REAL(dp), DIMENSION(:), POINTER                    :: rrad_z, vrrad_z
     238        24550 :       REAL(dp), DIMENSION(:, :), POINTER                 :: g0_h, g0_h_2nd, gsph, gsph_nuc, rrad_0, &
     239        24550 :                                                             Vh1_h, Vh1_s, vrrad_0, zet, zet_nuc
     240        24550 :       REAL(dp), DIMENSION(:, :, :), POINTER              :: my_CG, Qlm_gg, Qlm_gg_2nd
     241        24550 :       TYPE(atomic_kind_type), DIMENSION(:), POINTER      :: atomic_kind_set
     242              :       TYPE(cell_type), POINTER                           :: cell
     243              :       TYPE(cneo_potential_type), POINTER                 :: cneo_potential
     244              :       TYPE(dft_control_type), POINTER                    :: dft_control
     245              :       TYPE(grid_atom_type), POINTER                      :: grid_atom
     246              :       TYPE(gto_basis_set_type), POINTER                  :: basis_1c, nuc_basis
     247              :       TYPE(harmonics_atom_type), POINTER                 :: harmonics
     248              :       TYPE(pw_env_type), POINTER                         :: pw_env
     249              :       TYPE(pw_poisson_type), POINTER                     :: poisson_env
     250              :       TYPE(qs_charges_type), POINTER                     :: qs_charges
     251        24550 :       TYPE(qs_kind_type), DIMENSION(:), POINTER          :: qs_kind_set
     252        24550 :       TYPE(rho0_atom_type), DIMENSION(:), POINTER        :: rho0_atom_set, rho0_atom_set_2nd
     253              :       TYPE(rho0_mpole_type), POINTER                     :: rho0_mpole, rho0_mpole_2nd
     254        24550 :       TYPE(rho_atom_type), DIMENSION(:), POINTER         :: rho_atom_set, rho_atom_set_2nd
     255              :       TYPE(rho_atom_type), POINTER                       :: rho_atom
     256        24550 :       TYPE(rhoz_cneo_type), DIMENSION(:), POINTER        :: rhoz_cneo_set
     257              :       TYPE(rhoz_cneo_type), POINTER                      :: rhoz_cneo
     258        24550 :       TYPE(rhoz_type), DIMENSION(:), POINTER             :: rhoz_set, rhoz_set_2nd
     259              : 
     260        24550 :       CALL timeset(routineN, handle)
     261              : 
     262        24550 :       NULLIFY (cell, dft_control, poisson_env, pw_env, qs_charges)
     263        24550 :       NULLIFY (atomic_kind_set, qs_kind_set, rho_atom_set, rho0_atom_set)
     264        24550 :       NULLIFY (rho0_mpole, rhoz_set)
     265        24550 :       NULLIFY (atom_list, grid_atom, harmonics)
     266        24550 :       NULLIFY (basis_1c, lmin, lmax, npgf, zet)
     267        24550 :       NULLIFY (gsph)
     268        24550 :       NULLIFY (rhoz_cneo, rhoz_cneo_set)
     269              : 
     270              :       CALL get_qs_env(qs_env=qs_env, &
     271              :                       cell=cell, dft_control=dft_control, &
     272              :                       atomic_kind_set=atomic_kind_set, &
     273              :                       qs_kind_set=qs_kind_set, &
     274        24550 :                       pw_env=pw_env, qs_charges=qs_charges)
     275              : 
     276        24550 :       CALL pw_env_get(pw_env, poisson_env=poisson_env)
     277        24550 :       my_periodic = (poisson_env%method == pw_poisson_periodic)
     278              : 
     279        24550 :       back_ch = qs_charges%background*cell%deth
     280              : 
     281              :       ! rhoz_set is not accessed in TDDFT
     282              :       CALL get_local_rho(local_rho_set, rho_atom_set, rho0_atom_set, rho0_mpole, rhoz_set, &
     283        24550 :                          rhoz_cneo_set) ! for integral space
     284              : 
     285              :       ! for forces we need a second local_rho_set
     286        24550 :       l_2nd_local_rho = .FALSE.
     287        24550 :       IF (PRESENT(local_rho_set_2nd)) THEN ! for potential
     288          352 :          l_2nd_local_rho = .TRUE.
     289          352 :          NULLIFY (rho_atom_set_2nd, rho0_atom_set_2nd, rhoz_set_2nd) ! for potential
     290          352 :          CALL get_local_rho(local_rho_set_2nd, rho_atom_set_2nd, rho0_atom_set_2nd, rho0_mpole_2nd, rhoz_set=rhoz_set_2nd)
     291              :       END IF
     292              : 
     293        24550 :       nkind = SIZE(atomic_kind_set, 1)
     294        24550 :       nspins = dft_control%nspins
     295              : 
     296        24550 :       core_charge = .NOT. tddft  ! for forces mixed version
     297        24550 :       my_core_2nd = .TRUE.
     298        24550 :       IF (PRESENT(core_2nd)) my_core_2nd = .NOT. core_2nd   ! if my_core_2nd true, include core charge
     299              : 
     300              :       ! The aim of the following code was to return immediately if the subroutine
     301              :       ! was called for triplet excited states in spin-restricted case. This check
     302              :       ! is also performed before invocation of this subroutine. It should be save
     303              :       ! to remove the optional argument 'do_triplet' from the subroutine interface.
     304              :       !IF (tddft) THEN
     305              :       !   CPASSERT(PRESENT(do_triplet))
     306              :       !   IF (nspins == 1 .AND. do_triplet) RETURN
     307              :       !END IF
     308              : 
     309        24550 :       CALL get_qs_kind_set(qs_kind_set, maxg_iso_not0=max_iso)
     310        24550 :       CALL get_rho0_mpole(rho0_mpole=rho0_mpole, lmax_0=lmax_0)
     311              : 
     312              :       !   Put to 0 the local hartree energy contribution from 1 center integrals
     313        24550 :       energy_hartree_1c = 0.0_dp
     314              :       !   Restore total quantum nuclear density to zero
     315        24550 :       qs_charges%total_rho1_hard_nuc = 0.0_dp
     316        24550 :       rho0_mpole%tot_rhoz_cneo_s = 0.0_dp
     317              : 
     318              :       !   Here starts the loop over all the atoms
     319        72964 :       DO ikind = 1, nkind
     320              : 
     321        48414 :          CALL get_atomic_kind(atomic_kind_set(ikind), atom_list=atom_list, natom=nat)
     322        48414 :          NULLIFY (cneo_potential)
     323              :          CALL get_qs_kind(qs_kind_set(ikind), &
     324              :                           grid_atom=grid_atom, &
     325              :                           harmonics=harmonics, ngrid_rad=nr, &
     326              :                           max_iso_not0=max_iso_not0, paw_atom=paw_atom, &
     327        48414 :                           cneo_potential=cneo_potential)
     328              :          CALL get_qs_kind(qs_kind_set(ikind), &
     329        48414 :                           basis_set=basis_1c, basis_type="GAPW_1C")
     330              : 
     331        48414 :          cneo = ASSOCIATED(cneo_potential)
     332        48414 :          IF (cneo .AND. tddft) &
     333            0 :             CPABORT("Electronic TDDFT with CNEO quantum nuclei is not implemented.")
     334              : 
     335        48414 :          NULLIFY (nuc_basis)
     336        48414 :          max_iso_not0_nuc = 0
     337        48414 :          IF (cneo) THEN
     338           48 :             CPASSERT(paw_atom)
     339              :             CALL get_qs_kind(qs_kind_set(ikind), &
     340           48 :                              basis_set=nuc_basis, basis_type="NUC")
     341           48 :             max_iso_not0_nuc = cneo_potential%harmonics%max_iso_not0
     342              :          END IF
     343              : 
     344       121378 :          IF (paw_atom) THEN
     345              : !===========    PAW   ===============
     346              :             CALL get_gto_basis_set(gto_basis_set=basis_1c, lmax=lmax, lmin=lmin, &
     347              :                                    maxso=maxso, npgf=npgf, maxl=maxl, &
     348        46570 :                                    nset=nset, zet=zet)
     349              : 
     350        46570 :             max_s_harm = harmonics%max_s_harm
     351        46570 :             llmax = harmonics%llmax
     352              : 
     353        46570 :             nsotot = maxso*nset
     354       186280 :             ALLOCATE (gsph(nr, nsotot))
     355       139710 :             ALLOCATE (gexp(nr))
     356       232850 :             ALLOCATE (sqrtwr(nr), g0_h_w(nr, 0:lmax_0))
     357              : 
     358       186280 :             ALLOCATE (aVh1b_hh(nsotot, nsotot))
     359       139710 :             ALLOCATE (aVh1b_ss(nsotot, nsotot))
     360       139710 :             ALLOCATE (aVh1b_00(nsotot, nsotot))
     361              : 
     362        46570 :             NULLIFY (Qlm_gg, g0_h)
     363              :             CALL get_rho0_mpole(rho0_mpole=rho0_mpole, ikind=ikind, &
     364              :                                 l0_ikind=lmax0, &
     365        46570 :                                 Qlm_gg=Qlm_gg, g0_h=g0_h) ! Qlm_gg of density
     366              : 
     367        46570 :             IF (PRESENT(local_rho_set_2nd)) THEN ! for potential
     368          716 :                NULLIFY (Qlm_gg_2nd, g0_h_2nd)
     369              :                CALL get_rho0_mpole(rho0_mpole=rho0_mpole_2nd, ikind=ikind, &
     370              :                                    l0_ikind=lmax0_2nd, &
     371          716 :                                    Qlm_gg=Qlm_gg_2nd, g0_h=g0_h_2nd) ! Qlm_gg of density
     372              :             END IF
     373        46570 :             nchan_0 = nsoset(lmax0)
     374              : 
     375        46570 :             IF (nchan_0 > MAX(max_iso_not0, max_iso_not0_nuc)) &
     376            0 :                CPABORT("channels for rho0 > # max of spherical harmonics")
     377              : 
     378        46570 :             NULLIFY (Vh1_h, Vh1_s)
     379       186280 :             ALLOCATE (Vh1_h(nr, max_iso_not0))
     380       186280 :             ALLOCATE (Vh1_s(nr, MAX(max_iso_not0, max_iso_not0_nuc, nchan_0)))
     381              : 
     382        46570 :             NULLIFY (lmax_nuc, lmin_nuc, npgf_nuc, zet_nuc, gsph_nuc)
     383        46570 :             maxso_nuc = 0
     384        46570 :             maxl_nuc = -1
     385        46570 :             nset_nuc = 0
     386        46570 :             max_s_harm_nuc = 0
     387        46570 :             llmax_nuc = -1
     388        46570 :             nsotot_nuc = 0
     389        46570 :             IF (cneo) THEN
     390              :                CALL get_gto_basis_set(gto_basis_set=nuc_basis, lmax=lmax_nuc, &
     391              :                                       lmin=lmin_nuc, maxso=maxso_nuc, npgf=npgf_nuc, &
     392           48 :                                       maxl=maxl_nuc, nset=nset_nuc, zet=zet_nuc)
     393              : 
     394           48 :                max_s_harm_nuc = cneo_potential%harmonics%max_s_harm
     395           48 :                llmax_nuc = cneo_potential%harmonics%llmax
     396           48 :                nsotot_nuc = maxso_nuc*nset_nuc
     397          192 :                ALLOCATE (gsph_nuc(nr, nsotot_nuc))
     398       198336 :                gsph_nuc = 0.0_dp
     399          192 :                ALLOCATE (aVh1b_hh_nuc(nsotot_nuc, nsotot_nuc))
     400          144 :                ALLOCATE (aVh1b_ss_nuc(nsotot_nuc, nsotot_nuc))
     401          192 :                ALLOCATE (aVh1b_00_nuc(nsotot_nuc, nsotot_nuc))
     402              :             END IF
     403              : 
     404            0 :             ALLOCATE (cg_list(2, nsoset(MAX(maxl, maxl_nuc))**2, &
     405              :                               MAX(max_s_harm, max_s_harm_nuc)), &
     406       279420 :                       cg_n_list(MAX(max_s_harm, max_s_harm_nuc)))
     407              : 
     408        46570 :             NULLIFY (rrad_z, my_CG)
     409        46570 :             my_CG => harmonics%my_CG
     410              : 
     411              :             !     set to zero temporary arrays
     412      2602390 :             sqrtwr = 0.0_dp
     413      7934024 :             g0_h_w = 0.0_dp
     414      2602390 :             gexp = 0.0_dp
     415     98927168 :             gsph = 0.0_dp
     416              : 
     417      2602390 :             sqrtwr(1:nr) = SQRT(grid_atom%wr(1:nr))
     418       176584 :             DO l_ang = 0, lmax0
     419      7213484 :                g0_h_w(1:nr, l_ang) = g0_h(1:nr, l_ang)*grid_atom%wr(1:nr)
     420              :             END DO
     421              : 
     422              :             m1 = 0
     423       168680 :             DO iset1 = 1, nset
     424       122110 :                n1 = nsoset(lmax(iset1))
     425       449884 :                DO ipgf1 = 1, npgf(iset1)
     426     18238174 :                   gexp(1:nr) = EXP(-zet(ipgf1, iset1)*grid_atom%rad2(1:nr))*sqrtwr(1:nr)
     427      1369688 :                   DO is1 = nsoset(lmin(iset1) - 1) + 1, nsoset(lmax(iset1))
     428       919804 :                      iso = is1 + (ipgf1 - 1)*n1 + m1
     429       919804 :                      l_ang = indso(1, is1)
     430    100638658 :                      gsph(1:nr, iso) = grid_atom%rad2l(1:nr, l_ang)*gexp(1:nr)
     431              :                   END DO ! is1
     432              :                END DO ! ipgf1
     433       168680 :                m1 = m1 + maxso
     434              :             END DO ! iset1
     435              : 
     436        46570 :             IF (cneo) THEN
     437              :                ! initialize nuclear pmat, cpc and e_core to zero
     438          140 :                DO iat = 1, nat
     439           92 :                   iatom = atom_list(iat)
     440        50876 :                   rhoz_cneo_set(iatom)%pmat = 0.0_dp
     441        50876 :                   rhoz_cneo_set(iatom)%cpc_h = 0.0_dp
     442        50876 :                   rhoz_cneo_set(iatom)%cpc_s = 0.0_dp
     443           92 :                   rhoz_cneo_set(iatom)%e_core = 0.0_dp
     444          140 :                   rhoz_cneo_set(iatom)%ready = .FALSE.
     445              :                END DO
     446              : 
     447              :                ! calculate nuclear gsph
     448           48 :                m1 = 0
     449          480 :                DO iset1 = 1, nset_nuc
     450          432 :                   n1 = nsoset(lmax_nuc(iset1))
     451          864 :                   DO ipgf1 = 1, npgf_nuc(iset1)
     452        22032 :                      gexp(1:nr) = EXP(-zet_nuc(ipgf1, iset1)*grid_atom%rad2(1:nr))*sqrtwr(1:nr)
     453         1968 :                      DO is1 = nsoset(lmin_nuc(iset1) - 1) + 1, nsoset(lmax_nuc(iset1))
     454         1104 :                         iso = is1 + (ipgf1 - 1)*n1 + m1
     455         1104 :                         l_ang = indso(1, is1)
     456       111936 :                         gsph_nuc(1:nr, iso) = cneo_potential%rad2l(1:nr, l_ang)*gexp(1:nr)
     457              :                      END DO ! is1
     458              :                   END DO ! ipgf1
     459          480 :                   m1 = m1 + maxso_nuc
     460              :                END DO ! iset1
     461              :             END IF
     462              : 
     463              :             !     Distribute the atoms of this kind
     464        46570 :             num_pe = para_env%num_pe
     465        46570 :             mepos = para_env%mepos
     466        46570 :             bo = get_limit(nat, num_pe, mepos)
     467              : 
     468        83478 :             DO iat = bo(1), bo(2) !1,nat
     469        36908 :                iatom = atom_list(iat)
     470        36908 :                rho_atom => rho_atom_set(iatom)
     471              : 
     472        36908 :                NULLIFY (rrad_z, vrrad_z, rrad_0, vrrad_0)
     473        36908 :                IF (core_charge .AND. .NOT. cneo) THEN
     474        30630 :                   rrad_z => rhoz_set(ikind)%r_coef ! for density
     475              :                END IF
     476        36908 :                IF (my_core_2nd .AND. .NOT. cneo) THEN
     477        30838 :                   IF (l_2nd_local_rho) THEN
     478          263 :                      vrrad_z => rhoz_set_2nd(ikind)%vr_coef ! for potential
     479              :                   ELSE
     480        30575 :                      vrrad_z => rhoz_set(ikind)%vr_coef ! for potential
     481              :                   END IF
     482              :                END IF
     483        36908 :                rrad_0 => rho0_atom_set(iatom)%rho0_rad_h%r_coef ! for density
     484        36908 :                vrrad_0 => rho0_atom_set(iatom)%vrho0_rad_h%r_coef
     485        36908 :                IF (l_2nd_local_rho) THEN
     486          526 :                   rho_atom => rho_atom_set_2nd(iatom)
     487          526 :                   vrrad_0 => rho0_atom_set_2nd(iatom)%vrho0_rad_h%r_coef ! for potential
     488              :                END IF
     489        36908 :                IF (my_periodic .AND. back_ch > 1.E-3_dp) THEN
     490          954 :                   factor = -2.0_dp*pi/3.0_dp*SQRT(fourpi)*qs_charges%background
     491              :                ELSE
     492        35954 :                   factor = 0._dp
     493              :                END IF
     494              : 
     495              :                CALL Vh_1c_atom_potential(rho_atom, vrrad_0, &
     496              :                                          grid_atom, my_core_2nd .AND. .NOT. cneo, & ! core charge for potential (2nd)
     497              :                                          vrrad_z, Vh1_h, Vh1_s, &
     498        36908 :                                          nchan_0, nspins, max_iso_not0, factor)
     499              : 
     500        36908 :                IF (l_2nd_local_rho) rho_atom => rho_atom_set(iatom) ! rho_atom for density
     501              : 
     502        36908 :                ecoul_1_z_cneo = 0.0_dp
     503        36908 :                IF (cneo) THEN
     504           46 :                   rhoz_cneo => rhoz_cneo_set(iatom)
     505              :                   ! Add the soft tail of nuclear Hartree potential to total Vh1_s first.
     506              :                   ! vrho already contains the -zeff factor.
     507         1196 :                   DO iso = 1, max_iso_not0_nuc
     508       116196 :                      Vh1_s(:, iso) = Vh1_s(:, iso) + rhoz_cneo%vrho_rad_s(:, iso)
     509              :                   END DO
     510              : 
     511              :                   ! Build nuclear 1c integrals according to Vh1_h from electronic density only
     512              :                   ! and Vh1_s from electron density and last step (or initial guess) nuclear density
     513              :                   ! Vh1_h_nuc = -Z*Vh1_h, Vh1_s_nuc = -Z*Vh1_s
     514              :                   CALL Vh_1c_nuc_integrals(rhoz_cneo, cneo_potential%zeff, &
     515              :                                            aVh1b_hh_nuc, aVh1b_ss_nuc, aVh1b_00_nuc, Vh1_h, Vh1_s, &
     516              :                                            max_iso_not0, max_iso_not0_nuc, &
     517              :                                            max_s_harm_nuc, llmax_nuc, cg_list, cg_n_list, &
     518              :                                            nset_nuc, npgf_nuc, lmin_nuc, lmax_nuc, nsotot_nuc, maxso_nuc, &
     519              :                                            nchan_0, gsph_nuc, g0_h_w, cneo_potential%harmonics%my_CG, &
     520           46 :                                            cneo_potential%Qlm_gg)
     521              : 
     522              :                   ! Solve the nuclear 1c problem
     523              :                   CALL calculate_rhoz_cneo(rhoz_cneo, cneo_potential, cg_list, cg_n_list, nset_nuc, &
     524           46 :                                            npgf_nuc, lmin_nuc, lmax_nuc, maxl_nuc, maxso_nuc)
     525              :                   ! slm_int(iso=1) = sqrt(4*Pi), slm_int(iso>1) = 0, without using Lebedev grid
     526              :                   ! when printing, nuclear density is positive, thus needing the minus sign
     527              :                   qs_charges%total_rho1_hard_nuc = qs_charges%total_rho1_hard_nuc - SQRT(fourpi) &
     528         2346 :                                                    *SUM(rhoz_cneo%rho_rad_h(:, 1)*grid_atom%wr(:))
     529              :                   rho0_mpole%tot_rhoz_cneo_s = rho0_mpole%tot_rhoz_cneo_s - SQRT(fourpi) &
     530         2346 :                                                *SUM(rhoz_cneo%rho_rad_s(:, 1)*grid_atom%wr(:))
     531              : 
     532              :                   ! Calculate the contributions to Ecoul coming from Vh1_h*rhoz_cneo_h and Vh1_s*rhoz_cneo_s.
     533              :                   ! Self-interaction of the quantum nucleus is already removed in ecoul_1_z_cneo.
     534              :                   ! rho already contains the -zeff factor.
     535          460 :                   DO iso = 1, MIN(max_iso_not0, max_iso_not0_nuc)
     536              :                      ecoul_1_z_cneo = ecoul_1_z_cneo + 0.5_dp* &
     537              :                                       (SUM(Vh1_h(:, iso)*rhoz_cneo%rho_rad_h(:, iso)*grid_atom%wr(:)) &
     538        41860 :                                        - SUM(Vh1_s(:, iso)*rhoz_cneo%rho_rad_s(:, iso)*grid_atom%wr(:)))
     539              :                   END DO
     540          782 :                   DO iso = max_iso_not0 + 1, max_iso_not0_nuc
     541              :                      ecoul_1_z_cneo = ecoul_1_z_cneo - 0.5_dp* &
     542        37582 :                                       SUM(Vh1_s(:, iso)*rhoz_cneo%rho_rad_s(:, iso)*grid_atom%wr(:))
     543              :                   END DO
     544              : 
     545              :                   ! Add nuclear Hartree potential to total Vh1_h after solving the nuclear 1c problem
     546              :                   ! to avoid nuclear self-interaction.
     547              :                   ! vrho already contains the -zeff factor.
     548              :                   ! Here the min of two max_iso_not0's is chosen, because even when the nuclear one
     549              :                   ! is larger, it is meaningless to let Vh1_h have higher angular momentum components
     550              :                   ! as Vh1_h now is only used for the electronic part.
     551          460 :                   DO iso = 1, MIN(max_iso_not0, max_iso_not0_nuc)
     552        41860 :                      Vh1_h(:, iso) = Vh1_h(:, iso) + rhoz_cneo%vrho_rad_h(:, iso)
     553              :                   END DO
     554              :                END IF
     555              : 
     556              :                CALL Vh_1c_atom_energy(energy_hartree_1c, ecoul_1c, rho_atom, rrad_0, &
     557              :                                       grid_atom, iatom, core_charge .AND. .NOT. cneo, & ! core charge for density
     558        36908 :                                       rrad_z, Vh1_h, Vh1_s, nchan_0, nspins, max_iso_not0)
     559              : 
     560        36908 :                IF (l_2nd_local_rho) rho_atom => rho_atom_set_2nd(iatom) ! rho_atom for potential (2nd)
     561              : 
     562        36908 :                IF (cneo) THEN
     563           46 :                   CALL set_ecoul_1c(ecoul_1c, iatom, ecoul_1_z=ecoul_1c(iatom)%ecoul_1_z + ecoul_1_z_cneo)
     564           46 :                   energy_hartree_1c = energy_hartree_1c + ecoul_1_z_cneo
     565              :                END IF
     566              : 
     567              :                CALL Vh_1c_atom_integrals(rho_atom, &  ! results (int_local_h and int_local_s) written on rho_atom_2nd
     568              :                                          ! int_local_h and int_local_s are used in update_ks_atom
     569              :                                          ! on int_local_h mixed core / non-core
     570              :                                          aVh1b_hh, aVh1b_ss, aVh1b_00, Vh1_h, Vh1_s, max_iso_not0, &
     571              :                                          max_s_harm, llmax, cg_list, cg_n_list, &
     572              :                                          nset, npgf, lmin, lmax, nsotot, maxso, nspins, nchan_0, gsph, &
     573        83478 :                                          g0_h_w, my_CG, Qlm_gg) ! Qlm_gg for density from local_rho_set
     574              : 
     575              :             END DO ! iat
     576              : 
     577        46570 :             DEALLOCATE (aVh1b_hh)
     578        46570 :             DEALLOCATE (aVh1b_ss)
     579        46570 :             DEALLOCATE (aVh1b_00)
     580        46570 :             IF (ALLOCATED(aVh1b_hh_nuc)) DEALLOCATE (aVh1b_hh_nuc)
     581        46570 :             IF (ALLOCATED(aVh1b_ss_nuc)) DEALLOCATE (aVh1b_ss_nuc)
     582        46570 :             IF (ALLOCATED(aVh1b_00_nuc)) DEALLOCATE (aVh1b_00_nuc)
     583        46570 :             DEALLOCATE (Vh1_h, Vh1_s)
     584        46570 :             DEALLOCATE (cg_list, cg_n_list)
     585        46570 :             DEALLOCATE (gsph)
     586        46570 :             IF (ASSOCIATED(gsph_nuc)) DEALLOCATE (gsph_nuc)
     587        46570 :             DEALLOCATE (gexp)
     588        46570 :             DEALLOCATE (sqrtwr, g0_h_w)
     589              : 
     590       139710 :             IF (cneo) THEN
     591              :                ! broadcast nuclear pmat, cpc and e_core
     592          140 :                DO iat = 1, nat
     593           92 :                   iatom = atom_list(iat)
     594       101660 :                   CALL para_env%sum(rhoz_cneo_set(iatom)%pmat)
     595       101660 :                   CALL para_env%sum(rhoz_cneo_set(iatom)%cpc_h)
     596       101660 :                   CALL para_env%sum(rhoz_cneo_set(iatom)%cpc_s)
     597           92 :                   CALL para_env%sum(rhoz_cneo_set(iatom)%e_core)
     598          140 :                   rhoz_cneo_set(iatom)%ready = .TRUE.
     599              :                END DO
     600              :             END IF
     601              :          ELSE
     602              : !===========   NO  PAW   ===============
     603              : !  This term is taken care of using the core density as in GPW
     604              :             CYCLE
     605              :          END IF ! paw
     606              :       END DO ! ikind
     607              : 
     608        24550 :       CALL para_env%sum(energy_hartree_1c)
     609        24550 :       CALL para_env%sum(qs_charges%total_rho1_hard_nuc)
     610        24550 :       CALL para_env%sum(rho0_mpole%tot_rhoz_cneo_s)
     611              : 
     612        24550 :       CALL timestop(handle)
     613              : 
     614        49100 :    END SUBROUTINE Vh_1c_gg_integrals
     615              : 
     616              : ! **************************************************************************************************
     617              : 
     618              : ! **************************************************************************************************
     619              : !> \brief ...
     620              : !> \param rho_atom ...
     621              : !> \param vrrad_0 ...
     622              : !> \param grid_atom ...
     623              : !> \param core_charge ...
     624              : !> \param vrrad_z ...
     625              : !> \param Vh1_h ...
     626              : !> \param Vh1_s ...
     627              : !> \param nchan_0 ...
     628              : !> \param nspins ...
     629              : !> \param max_iso_not0 ...
     630              : !> \param bfactor ...
     631              : ! **************************************************************************************************
     632        36908 :    SUBROUTINE Vh_1c_atom_potential(rho_atom, vrrad_0, &
     633              :                                    grid_atom, core_charge, vrrad_z, Vh1_h, Vh1_s, &
     634              :                                    nchan_0, nspins, max_iso_not0, bfactor)
     635              : 
     636              :       TYPE(rho_atom_type), POINTER                       :: rho_atom
     637              :       REAL(dp), DIMENSION(:, :), POINTER                 :: vrrad_0
     638              :       TYPE(grid_atom_type), POINTER                      :: grid_atom
     639              :       LOGICAL, INTENT(IN)                                :: core_charge
     640              :       REAL(dp), DIMENSION(:), POINTER                    :: vrrad_z
     641              :       REAL(dp), DIMENSION(:, :), POINTER                 :: Vh1_h, Vh1_s
     642              :       INTEGER, INTENT(IN)                                :: nchan_0, nspins, max_iso_not0
     643              :       REAL(dp), INTENT(IN)                               :: bfactor
     644              : 
     645              :       INTEGER                                            :: ir, iso, ispin, nr
     646        36908 :       TYPE(rho_atom_coeff), DIMENSION(:), POINTER        :: vr_h, vr_s
     647              : 
     648        36908 :       nr = grid_atom%nr
     649              : 
     650        36908 :       NULLIFY (vr_h, vr_s)
     651        36908 :       CALL get_rho_atom(rho_atom=rho_atom, vrho_rad_h=vr_h, vrho_rad_s=vr_s)
     652              : 
     653     27134972 :       Vh1_h = 0.0_dp
     654     27172508 :       Vh1_s = 0.0_dp
     655              : 
     656      3469826 :       IF (core_charge) Vh1_h(:, 1) = vrrad_z(:)
     657              : 
     658       327984 :       DO iso = 1, nchan_0
     659     31344780 :          Vh1_s(:, iso) = vrrad_0(:, iso)
     660              :       END DO
     661              : 
     662        78835 :       DO ispin = 1, nspins
     663       657122 :          DO iso = 1, max_iso_not0
     664     32014257 :             Vh1_h(:, iso) = Vh1_h(:, iso) + vr_h(ispin)%r_coef(:, iso)
     665     32056184 :             Vh1_s(:, iso) = Vh1_s(:, iso) + vr_s(ispin)%r_coef(:, iso)
     666              :          END DO
     667              :       END DO
     668              : 
     669        36908 :       IF (bfactor /= 0._dp) THEN
     670        56854 :          DO ir = 1, nr
     671        55900 :             Vh1_h(ir, 1) = Vh1_h(ir, 1) + bfactor*grid_atom%rad2(ir)*grid_atom%wr(ir)
     672        56854 :             Vh1_s(ir, 1) = Vh1_s(ir, 1) + bfactor*grid_atom%rad2(ir)*grid_atom%wr(ir)
     673              :          END DO
     674              :       END IF
     675              : 
     676        36908 :    END SUBROUTINE Vh_1c_atom_potential
     677              : 
     678              : ! **************************************************************************************************
     679              : 
     680              : ! **************************************************************************************************
     681              : !> \brief ...
     682              : !> \param energy_hartree_1c ...
     683              : !> \param ecoul_1c ...
     684              : !> \param rho_atom ...
     685              : !> \param rrad_0 ...
     686              : !> \param grid_atom ...
     687              : !> \param iatom ...
     688              : !> \param core_charge ...
     689              : !> \param rrad_z ...
     690              : !> \param Vh1_h ...
     691              : !> \param Vh1_s ...
     692              : !> \param nchan_0 ...
     693              : !> \param nspins ...
     694              : !> \param max_iso_not0 ...
     695              : ! **************************************************************************************************
     696        36908 :    SUBROUTINE Vh_1c_atom_energy(energy_hartree_1c, ecoul_1c, rho_atom, rrad_0, &
     697              :                                 grid_atom, iatom, core_charge, rrad_z, Vh1_h, Vh1_s, &
     698              :                                 nchan_0, nspins, max_iso_not0)
     699              : 
     700              :       REAL(dp), INTENT(INOUT)                            :: energy_hartree_1c
     701              :       TYPE(ecoul_1center_type), DIMENSION(:), POINTER    :: ecoul_1c
     702              :       TYPE(rho_atom_type), POINTER                       :: rho_atom
     703              :       REAL(dp), DIMENSION(:, :), POINTER                 :: rrad_0
     704              :       TYPE(grid_atom_type), POINTER                      :: grid_atom
     705              :       INTEGER, INTENT(IN)                                :: iatom
     706              :       LOGICAL, INTENT(IN)                                :: core_charge
     707              :       REAL(dp), DIMENSION(:), POINTER                    :: rrad_z
     708              :       REAL(dp), DIMENSION(:, :), POINTER                 :: Vh1_h, Vh1_s
     709              :       INTEGER, INTENT(IN)                                :: nchan_0, nspins, max_iso_not0
     710              : 
     711              :       INTEGER                                            :: iso, ispin, nr
     712              :       REAL(dp)                                           :: ecoul_1_0, ecoul_1_h, ecoul_1_s, &
     713              :                                                             ecoul_1_z
     714        36908 :       TYPE(rho_atom_coeff), DIMENSION(:), POINTER        :: r_h, r_s
     715              : 
     716        36908 :       nr = grid_atom%nr
     717              : 
     718        36908 :       NULLIFY (r_h, r_s)
     719        36908 :       CALL get_rho_atom(rho_atom=rho_atom, rho_rad_h=r_h, rho_rad_s=r_s)
     720              : 
     721              :       !       Calculate the contributions to Ecoul coming from Vh1_h*rhoz
     722        36908 :       ecoul_1_z = 0.0_dp
     723        36908 :       IF (core_charge) THEN
     724      1721270 :          ecoul_1_z = 0.5_dp*SUM(Vh1_h(:, 1)*rrad_z(:)*grid_atom%wr(:))
     725              :       END IF
     726              : 
     727              :       !       Calculate the contributions to Ecoul coming from  Vh1_s*rho0
     728        36908 :       ecoul_1_0 = 0.0_dp
     729       327984 :       DO iso = 1, nchan_0
     730     15690844 :          ecoul_1_0 = ecoul_1_0 + 0.5_dp*SUM(Vh1_s(:, iso)*rrad_0(:, iso)*grid_atom%wr(:))
     731              :       END DO
     732              : 
     733              :       !       Calculate the contributions to Ecoul coming from Vh1_h*rho1_h and Vh1_s*rho1_s
     734        36908 :       ecoul_1_s = 0.0_dp
     735        36908 :       ecoul_1_h = 0.0_dp
     736        78835 :       DO ispin = 1, nspins
     737       657122 :          DO iso = 1, max_iso_not0
     738     32014257 :             ecoul_1_s = ecoul_1_s + 0.5_dp*SUM(Vh1_s(:, iso)*r_s(ispin)%r_coef(:, iso)*grid_atom%wr(:))
     739     32056184 :             ecoul_1_h = ecoul_1_h + 0.5_dp*SUM(Vh1_h(:, iso)*r_h(ispin)%r_coef(:, iso)*grid_atom%wr(:))
     740              :          END DO
     741              :       END DO
     742              : 
     743        36908 :       CALL set_ecoul_1c(ecoul_1c, iatom, ecoul_1_z=ecoul_1_z, ecoul_1_0=ecoul_1_0)
     744        36908 :       CALL set_ecoul_1c(ecoul_1c=ecoul_1c, iatom=iatom, ecoul_1_h=ecoul_1_h, ecoul_1_s=ecoul_1_s)
     745              : 
     746        36908 :       energy_hartree_1c = energy_hartree_1c + ecoul_1_z - ecoul_1_0
     747        36908 :       energy_hartree_1c = energy_hartree_1c + ecoul_1_h - ecoul_1_s
     748              : 
     749        36908 :    END SUBROUTINE Vh_1c_atom_energy
     750              : 
     751              : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     752              : 
     753              : ! **************************************************************************************************
     754              : !> \brief ...
     755              : !> \param rho_atom ...
     756              : !> \param aVh1b_hh ...
     757              : !> \param aVh1b_ss ...
     758              : !> \param aVh1b_00 ...
     759              : !> \param Vh1_h ...
     760              : !> \param Vh1_s ...
     761              : !> \param max_iso_not0 ...
     762              : !> \param max_s_harm ...
     763              : !> \param llmax ...
     764              : !> \param cg_list ...
     765              : !> \param cg_n_list ...
     766              : !> \param nset ...
     767              : !> \param npgf ...
     768              : !> \param lmin ...
     769              : !> \param lmax ...
     770              : !> \param nsotot ...
     771              : !> \param maxso ...
     772              : !> \param nspins ...
     773              : !> \param nchan_0 ...
     774              : !> \param gsph ...
     775              : !> \param g0_h_w ...
     776              : !> \param my_CG ...
     777              : !> \param Qlm_gg ...
     778              : ! **************************************************************************************************
     779        36908 :    SUBROUTINE Vh_1c_atom_integrals(rho_atom, &
     780        36908 :                                    aVh1b_hh, aVh1b_ss, aVh1b_00, Vh1_h, Vh1_s, max_iso_not0, &
     781        36908 :                                    max_s_harm, llmax, cg_list, cg_n_list, &
     782              :                                    nset, npgf, lmin, lmax, nsotot, maxso, nspins, nchan_0, gsph, &
     783        36908 :                                    g0_h_w, my_CG, Qlm_gg)
     784              : 
     785              :       TYPE(rho_atom_type), POINTER                       :: rho_atom
     786              :       REAL(dp), DIMENSION(:, :)                          :: aVh1b_hh, aVh1b_ss, aVh1b_00
     787              :       REAL(dp), DIMENSION(:, :), POINTER                 :: Vh1_h, Vh1_s
     788              :       INTEGER, INTENT(IN)                                :: max_iso_not0, max_s_harm, llmax
     789              :       INTEGER, DIMENSION(:, :, :)                        :: cg_list
     790              :       INTEGER, DIMENSION(:)                              :: cg_n_list
     791              :       INTEGER, INTENT(IN)                                :: nset
     792              :       INTEGER, DIMENSION(:), POINTER                     :: npgf, lmin, lmax
     793              :       INTEGER, INTENT(IN)                                :: nsotot, maxso, nspins, nchan_0
     794              :       REAL(dp), DIMENSION(:, :), POINTER                 :: gsph
     795              :       REAL(dp), DIMENSION(:, 0:)                         :: g0_h_w
     796              :       REAL(dp), DIMENSION(:, :, :), POINTER              :: my_CG, Qlm_gg
     797              : 
     798              :       INTEGER                                            :: icg, ipgf1, ipgf2, ir, is1, is2, iset1, &
     799              :                                                             iset2, iso, iso1, iso2, ispin, l_ang, &
     800              :                                                             m1, m2, max_iso_not0_local, n1, n2, nr
     801              :       REAL(dp)                                           :: gVg_0, gVg_h, gVg_s
     802        36908 :       TYPE(rho_atom_coeff), DIMENSION(:), POINTER        :: int_local_h, int_local_s
     803              : 
     804        36908 :       NULLIFY (int_local_h, int_local_s)
     805              :       CALL get_rho_atom(rho_atom=rho_atom, &
     806              :                         ga_Vlocal_gb_h=int_local_h, &
     807        36908 :                         ga_Vlocal_gb_s=int_local_s)
     808              : 
     809              :       !       Calculate the integrals of the potential with 2 primitives
     810    132814350 :       aVh1b_hh = 0.0_dp
     811    132814350 :       aVh1b_ss = 0.0_dp
     812    132814350 :       aVh1b_00 = 0.0_dp
     813              : 
     814        36908 :       nr = SIZE(gsph, 1)
     815              : 
     816        36908 :       m1 = 0
     817       130525 :       DO iset1 = 1, nset
     818              :          m2 = 0
     819       404680 :          DO iset2 = 1, nset
     820              :             CALL get_none0_cg_list(my_CG, lmin(iset1), lmax(iset1), lmin(iset2), lmax(iset2), &
     821       311063 :                                    max_s_harm, llmax, cg_list, cg_n_list, max_iso_not0_local)
     822              : 
     823       311063 :             n1 = nsoset(lmax(iset1))
     824      1084025 :             DO ipgf1 = 1, npgf(iset1)
     825       772962 :                n2 = nsoset(lmax(iset2))
     826      3281891 :                DO ipgf2 = 1, npgf(iset2)
     827              :                   !               with contributions to  V1_s*rho0
     828     22714116 :                   DO iso = 1, MIN(nchan_0, max_iso_not0)
     829     20516250 :                      l_ang = indso(1, iso)
     830   1102922360 :                      gVg_0 = SUM(Vh1_s(:, iso)*g0_h_w(:, l_ang))
     831     44210251 :                      DO icg = 1, cg_n_list(iso)
     832     21496135 :                         is1 = cg_list(1, icg, iso)
     833     21496135 :                         is2 = cg_list(2, icg, iso)
     834              : 
     835     21496135 :                         iso1 = is1 + n1*(ipgf1 - 1) + m1
     836     21496135 :                         iso2 = is2 + n2*(ipgf2 - 1) + m2
     837     21496135 :                         gVg_h = 0.0_dp
     838     21496135 :                         gVg_s = 0.0_dp
     839              : 
     840   1146215845 :                         DO ir = 1, nr
     841   1124719710 :                            gVg_h = gVg_h + gsph(ir, iso1)*gsph(ir, iso2)*Vh1_h(ir, iso)
     842   1146215845 :                            gVg_s = gVg_s + gsph(ir, iso1)*gsph(ir, iso2)*Vh1_s(ir, iso)
     843              :                         END DO ! ir
     844              : 
     845     21496135 :                         aVh1b_hh(iso1, iso2) = aVh1b_hh(iso1, iso2) + gVg_h*my_CG(is1, is2, iso)
     846     21496135 :                         aVh1b_ss(iso1, iso2) = aVh1b_ss(iso1, iso2) + gVg_s*my_CG(is1, is2, iso)
     847     42012385 :                         aVh1b_00(iso1, iso2) = aVh1b_00(iso1, iso2) + gVg_0*Qlm_gg(iso1, iso2, iso)
     848              : 
     849              :                      END DO !icg
     850              :                   END DO ! iso
     851              :                   !               without contributions to  V1_s*rho0
     852     28536756 :                   DO iso = nchan_0 + 1, max_iso_not0
     853     37778359 :                      DO icg = 1, cg_n_list(iso)
     854     10014565 :                         is1 = cg_list(1, icg, iso)
     855     10014565 :                         is2 = cg_list(2, icg, iso)
     856              : 
     857     10014565 :                         iso1 = is1 + n1*(ipgf1 - 1) + m1
     858     10014565 :                         iso2 = is2 + n2*(ipgf2 - 1) + m2
     859     10014565 :                         gVg_h = 0.0_dp
     860     10014565 :                         gVg_s = 0.0_dp
     861              : 
     862    519652045 :                         DO ir = 1, nr
     863    509637480 :                            gVg_h = gVg_h + gsph(ir, iso1)*gsph(ir, iso2)*Vh1_h(ir, iso)
     864    519652045 :                            gVg_s = gVg_s + gsph(ir, iso1)*gsph(ir, iso2)*Vh1_s(ir, iso)
     865              :                         END DO ! ir
     866              : 
     867     10014565 :                         aVh1b_hh(iso1, iso2) = aVh1b_hh(iso1, iso2) + gVg_h*my_CG(is1, is2, iso)
     868     35580493 :                         aVh1b_ss(iso1, iso2) = aVh1b_ss(iso1, iso2) + gVg_s*my_CG(is1, is2, iso)
     869              : 
     870              :                      END DO !icg
     871              :                   END DO ! iso
     872              :                END DO ! ipgf2
     873              :             END DO ! ipgf1
     874       404680 :             m2 = m2 + maxso
     875              :          END DO ! iset2
     876       130525 :          m1 = m1 + maxso
     877              :       END DO !iset1
     878        78835 :       DO ispin = 1, nspins
     879        41927 :          CALL daxpy(nsotot*nsotot, 1.0_dp, aVh1b_hh, 1, int_local_h(ispin)%r_coef, 1)
     880        41927 :          CALL daxpy(nsotot*nsotot, 1.0_dp, aVh1b_ss, 1, int_local_s(ispin)%r_coef, 1)
     881        41927 :          CALL daxpy(nsotot*nsotot, -1.0_dp, aVh1b_00, 1, int_local_h(ispin)%r_coef, 1)
     882        78835 :          CALL daxpy(nsotot*nsotot, -1.0_dp, aVh1b_00, 1, int_local_s(ispin)%r_coef, 1)
     883              :       END DO ! ispin
     884              : 
     885        36908 :    END SUBROUTINE Vh_1c_atom_integrals
     886              : 
     887              : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     888              : 
     889              : END MODULE hartree_local_methods
        

Generated by: LCOV version 2.0-1