LCOV - code coverage report
Current view: top level - src - core_ppnl.F (source / functions) Coverage Total Hit
Test: CP2K Regtests (git:936074a) Lines: 97.5 % 122 119
Test Date: 2025-12-04 06:27:48 Functions: 100.0 % 1 1

            Line data    Source code
       1              : !--------------------------------------------------------------------------------------------------!
       2              : !   CP2K: A general program to perform molecular dynamics simulations                              !
       3              : !   Copyright 2000-2025 CP2K developers group <https://cp2k.org>                                   !
       4              : !                                                                                                  !
       5              : !   SPDX-License-Identifier: GPL-2.0-or-later                                                      !
       6              : !--------------------------------------------------------------------------------------------------!
       7              : ! **************************************************************************************************
       8              : !> \brief Calculation of the non-local pseudopotential contribution to the core Hamiltonian
       9              : !>         <a|V(non-local)|b> = <a|p(l,i)>*h(i,j)*<p(l,j)|b>
      10              : !> \par History
      11              : !>      - refactered from qs_core_hamiltian [Joost VandeVondele, 2008-11-01]
      12              : !>      - full rewrite [jhu, 2009-01-23]
      13              : !>      - Extended by the derivatives for DFPT [Sandra Luber, Edward Ditler, 2021]
      14              : ! **************************************************************************************************
      15              : MODULE core_ppnl
      16              :    USE ai_angmom,                       ONLY: angmom
      17              :    USE ai_overlap,                      ONLY: overlap
      18              :    USE atomic_kind_types,               ONLY: atomic_kind_type,&
      19              :                                               get_atomic_kind_set
      20              :    USE basis_set_types,                 ONLY: gto_basis_set_p_type,&
      21              :                                               gto_basis_set_type
      22              :    USE cp_dbcsr_api,                    ONLY: dbcsr_add,&
      23              :                                               dbcsr_get_block_p,&
      24              :                                               dbcsr_p_type
      25              :    USE external_potential_types,        ONLY: gth_potential_p_type,&
      26              :                                               gth_potential_type,&
      27              :                                               sgp_potential_p_type,&
      28              :                                               sgp_potential_type
      29              :    USE kinds,                           ONLY: dp,&
      30              :                                               int_8
      31              :    USE orbital_pointers,                ONLY: init_orbital_pointers,&
      32              :                                               nco,&
      33              :                                               ncoset
      34              :    USE particle_types,                  ONLY: particle_type
      35              :    USE qs_force_types,                  ONLY: qs_force_type
      36              :    USE qs_kind_types,                   ONLY: get_qs_kind,&
      37              :                                               get_qs_kind_set,&
      38              :                                               qs_kind_type
      39              :    USE qs_neighbor_list_types,          ONLY: neighbor_list_set_p_type
      40              :    USE sap_kind_types,                  ONLY: alist_type,&
      41              :                                               clist_type,&
      42              :                                               get_alist,&
      43              :                                               release_sap_int,&
      44              :                                               sap_int_type,&
      45              :                                               sap_sort
      46              :    USE virial_methods,                  ONLY: virial_pair_force
      47              :    USE virial_types,                    ONLY: virial_type
      48              : 
      49              : !$ USE OMP_LIB, ONLY: omp_lock_kind, &
      50              : !$                    omp_init_lock, omp_set_lock, &
      51              : !$                    omp_unset_lock, omp_destroy_lock
      52              : 
      53              : #include "./base/base_uses.f90"
      54              : 
      55              :    IMPLICIT NONE
      56              : 
      57              :    PRIVATE
      58              : 
      59              :    CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'core_ppnl'
      60              : 
      61              :    PUBLIC :: build_core_ppnl
      62              : 
      63              : CONTAINS
      64              : 
      65              : ! **************************************************************************************************
      66              : !> \brief ...
      67              : !> \param matrix_h ...
      68              : !> \param matrix_p ...
      69              : !> \param force ...
      70              : !> \param virial ...
      71              : !> \param calculate_forces ...
      72              : !> \param use_virial ...
      73              : !> \param nder ...
      74              : !> \param qs_kind_set ...
      75              : !> \param atomic_kind_set ...
      76              : !> \param particle_set ...
      77              : !> \param sab_orb ...
      78              : !> \param sap_ppnl ...
      79              : !> \param eps_ppnl ...
      80              : !> \param nimages ...
      81              : !> \param cell_to_index ...
      82              : !> \param basis_type ...
      83              : !> \param deltaR Weighting factors of the derivatives wrt. nuclear positions
      84              : !> \param matrix_l ...
      85              : !> \param atcore ...
      86              : ! **************************************************************************************************
      87        15018 :    SUBROUTINE build_core_ppnl(matrix_h, matrix_p, force, virial, calculate_forces, use_virial, nder, &
      88              :                               qs_kind_set, atomic_kind_set, particle_set, sab_orb, sap_ppnl, eps_ppnl, &
      89        15018 :                               nimages, cell_to_index, basis_type, deltaR, matrix_l, atcore)
      90              : 
      91              :       TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER       :: matrix_h, matrix_p
      92              :       TYPE(qs_force_type), DIMENSION(:), POINTER         :: force
      93              :       TYPE(virial_type), POINTER                         :: virial
      94              :       LOGICAL, INTENT(IN)                                :: calculate_forces
      95              :       LOGICAL                                            :: use_virial
      96              :       INTEGER                                            :: nder
      97              :       TYPE(qs_kind_type), DIMENSION(:), POINTER          :: qs_kind_set
      98              :       TYPE(atomic_kind_type), DIMENSION(:), POINTER      :: atomic_kind_set
      99              :       TYPE(particle_type), DIMENSION(:), POINTER         :: particle_set
     100              :       TYPE(neighbor_list_set_p_type), DIMENSION(:), &
     101              :          POINTER                                         :: sab_orb, sap_ppnl
     102              :       REAL(KIND=dp), INTENT(IN)                          :: eps_ppnl
     103              :       INTEGER, INTENT(IN)                                :: nimages
     104              :       INTEGER, DIMENSION(:, :, :), OPTIONAL, POINTER     :: cell_to_index
     105              :       CHARACTER(LEN=*), INTENT(IN)                       :: basis_type
     106              :       REAL(KIND=dp), DIMENSION(:, :), INTENT(IN), &
     107              :          OPTIONAL                                        :: deltaR
     108              :       TYPE(dbcsr_p_type), DIMENSION(:, :), OPTIONAL, &
     109              :          POINTER                                         :: matrix_l
     110              :       REAL(KIND=dp), DIMENSION(:), INTENT(INOUT), &
     111              :          OPTIONAL                                        :: atcore
     112              : 
     113              :       CHARACTER(LEN=*), PARAMETER                        :: routineN = 'build_core_ppnl'
     114              : 
     115              :       INTEGER :: atom_a, first_col, handle, i, i_dim, iab, iac, iatom, ib, ibc, icol, ikind, &
     116              :          ilist, img, irow, iset, j, jatom, jb, jkind, jneighbor, kac, katom, kbc, kkind, l, &
     117              :          lc_max, lc_min, ldai, ldsab, lppnl, maxco, maxder, maxl, maxlgto, maxlppnl, maxppnl, &
     118              :          maxsgf, na, natom, nb, ncoa, ncoc, nkind, nlist, nneighbor, nnl, np, nppnl, nprjc, nseta, &
     119              :          nsgfa, prjc, sgfa, slot
     120        15018 :       INTEGER, ALLOCATABLE, DIMENSION(:)                 :: atom_of_kind, kind_of
     121              :       INTEGER, DIMENSION(3)                              :: cell_b, cell_c
     122        15018 :       INTEGER, DIMENSION(:), POINTER                     :: la_max, la_min, npgfa, nprj_ppnl, &
     123        15018 :                                                             nsgf_seta
     124        15018 :       INTEGER, DIMENSION(:, :), POINTER                  :: first_sgfa
     125              :       LOGICAL                                            :: do_dR, do_gth, do_kp, do_soc, doat, &
     126              :                                                             found, ppnl_present
     127              :       REAL(KIND=dp)                                      :: atk, dac, f0, ppnl_radius
     128        15018 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:)           :: radp
     129        15018 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :)        :: sab, work
     130        15018 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :)     :: ai_work, lab, work_l
     131              :       REAL(KIND=dp), DIMENSION(1)                        :: rprjc, zetc
     132              :       REAL(KIND=dp), DIMENSION(3)                        :: fa, fb, rab, rac, rbc
     133              :       REAL(KIND=dp), DIMENSION(3, 3)                     :: pv_thread
     134              :       TYPE(gto_basis_set_type), POINTER                  :: orb_basis_set
     135        15018 :       TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER  :: basis_set
     136              :       TYPE(gth_potential_type), POINTER                  :: gth_potential
     137        15018 :       TYPE(gth_potential_p_type), DIMENSION(:), POINTER  :: gpotential
     138              :       TYPE(clist_type), POINTER                          :: clist
     139              :       TYPE(alist_type), POINTER                          :: alist_ac, alist_bc
     140        30036 :       REAL(KIND=dp), DIMENSION(SIZE(particle_set))       :: at_thread
     141        15018 :       REAL(KIND=dp), DIMENSION(:, :, :), POINTER         :: achint, acint, alkint, bchint, bcint, &
     142        15018 :                                                             blkint
     143        15018 :       REAL(KIND=dp), DIMENSION(:, :), POINTER            :: cprj, h_block, l_block_x, l_block_y, &
     144        15018 :                                                             l_block_z, p_block, r_2block, &
     145        15018 :                                                             r_3block, rpgfa, sphi_a, vprj_ppnl, &
     146        15018 :                                                             wprj_ppnl, zeta
     147        15018 :       REAL(KIND=dp), DIMENSION(:), POINTER               :: a_nl, alpha_ppnl, hprj, set_radius_a
     148        30036 :       REAL(KIND=dp), DIMENSION(3, SIZE(particle_set))    :: force_thread
     149        15018 :       TYPE(sap_int_type), DIMENSION(:), POINTER          :: sap_int
     150        15018 :       TYPE(sgp_potential_p_type), DIMENSION(:), POINTER  :: spotential
     151              :       TYPE(sgp_potential_type), POINTER                  :: sgp_potential
     152              : 
     153              : !$    INTEGER(kind=omp_lock_kind), &
     154        15018 : !$       ALLOCATABLE, DIMENSION(:) :: locks
     155              : !$    INTEGER(KIND=int_8)                                :: iatom8
     156              : !$    INTEGER                                            :: lock_num, hash
     157              : !$    INTEGER, PARAMETER                                 :: nlock = 501
     158              : 
     159              :       MARK_USED(int_8)
     160              : 
     161        15018 :       do_dR = .FALSE.
     162           72 :       IF (PRESENT(deltaR)) do_dR = .TRUE.
     163        15018 :       doat = .FALSE.
     164        15018 :       IF (PRESENT(atcore)) doat = .TRUE.
     165        15018 :       IF ((calculate_forces .OR. doat) .AND. do_dR) THEN
     166            0 :          CPABORT("core_ppl: incompatible options")
     167              :       END IF
     168              : 
     169        15018 :       IF (calculate_forces) THEN
     170         6431 :          CALL timeset(routineN//"_forces", handle)
     171              :       ELSE
     172         8587 :          CALL timeset(routineN, handle)
     173              :       END IF
     174              : 
     175        15018 :       do_soc = PRESENT(matrix_l)
     176              : 
     177        15018 :       ppnl_present = ASSOCIATED(sap_ppnl)
     178              : 
     179        15018 :       IF (ppnl_present) THEN
     180              : 
     181        15018 :          nkind = SIZE(atomic_kind_set)
     182        15018 :          natom = SIZE(particle_set)
     183              : 
     184        15018 :          do_kp = (nimages > 1)
     185              : 
     186        15018 :          IF (do_kp) THEN
     187          228 :             CPASSERT(PRESENT(cell_to_index) .AND. ASSOCIATED(cell_to_index))
     188              :          END IF
     189              : 
     190        15018 :          IF (calculate_forces .OR. doat) THEN
     191         6493 :             IF (SIZE(matrix_p, 1) == 2) THEN
     192         1982 :                DO img = 1, nimages
     193              :                   CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
     194         1292 :                                  alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
     195              :                   CALL dbcsr_add(matrix_p(2, img)%matrix, matrix_p(1, img)%matrix, &
     196         1982 :                                  alpha_scalar=-2.0_dp, beta_scalar=1.0_dp)
     197              :                END DO
     198              :             END IF
     199              :          END IF
     200              : 
     201        15018 :          maxder = ncoset(nder)
     202              : 
     203              :          CALL get_qs_kind_set(qs_kind_set, &
     204              :                               maxco=maxco, &
     205              :                               maxlgto=maxlgto, &
     206              :                               maxsgf=maxsgf, &
     207              :                               maxlppnl=maxlppnl, &
     208              :                               maxppnl=maxppnl, &
     209        15018 :                               basis_type=basis_type)
     210              : 
     211        15018 :          maxl = MAX(maxlgto, maxlppnl)
     212        15018 :          CALL init_orbital_pointers(maxl + nder + 1)
     213              : 
     214        15018 :          ldsab = MAX(maxco, ncoset(maxlppnl), maxsgf, maxppnl)
     215        15018 :          ldai = ncoset(maxl + nder + 1)
     216              : 
     217              :          ! sap_int needs to be shared as multiple threads need to access this
     218       103142 :          ALLOCATE (sap_int(nkind*nkind))
     219        73106 :          DO i = 1, nkind*nkind
     220        58088 :             NULLIFY (sap_int(i)%alist, sap_int(i)%asort, sap_int(i)%aindex)
     221        73106 :             sap_int(i)%nalist = 0
     222              :          END DO
     223              : 
     224              :          ! Set up direct access to basis and potential
     225       160026 :          ALLOCATE (basis_set(nkind), gpotential(nkind), spotential(nkind))
     226        43330 :          DO ikind = 1, nkind
     227        28312 :             CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set, basis_type=basis_type)
     228        28312 :             IF (ASSOCIATED(orb_basis_set)) THEN
     229        28312 :                basis_set(ikind)%gto_basis_set => orb_basis_set
     230              :             ELSE
     231            0 :                NULLIFY (basis_set(ikind)%gto_basis_set)
     232              :             END IF
     233        28312 :             CALL get_qs_kind(qs_kind_set(ikind), gth_potential=gth_potential, sgp_potential=sgp_potential)
     234        28312 :             NULLIFY (gpotential(ikind)%gth_potential)
     235        28312 :             NULLIFY (spotential(ikind)%sgp_potential)
     236        43330 :             IF (ASSOCIATED(gth_potential)) THEN
     237        28064 :                gpotential(ikind)%gth_potential => gth_potential
     238        28064 :                IF (do_soc .AND. (.NOT. gth_potential%soc)) THEN
     239            0 :                   CPABORT("Spin-orbit coupling selected, but GTH potential without SOC parameters provided")
     240              :                END IF
     241          248 :             ELSE IF (ASSOCIATED(sgp_potential)) THEN
     242           10 :                spotential(ikind)%sgp_potential => sgp_potential
     243              :             END IF
     244              :          END DO
     245              : 
     246              :          ! Allocate sap int
     247       788882 :          DO slot = 1, sap_ppnl(1)%nl_size
     248              : 
     249       773864 :             ikind = sap_ppnl(1)%nlist_task(slot)%ikind
     250       773864 :             kkind = sap_ppnl(1)%nlist_task(slot)%jkind
     251       773864 :             iatom = sap_ppnl(1)%nlist_task(slot)%iatom
     252       773864 :             katom = sap_ppnl(1)%nlist_task(slot)%jatom
     253       773864 :             nlist = sap_ppnl(1)%nlist_task(slot)%nlist
     254       773864 :             ilist = sap_ppnl(1)%nlist_task(slot)%ilist
     255       773864 :             nneighbor = sap_ppnl(1)%nlist_task(slot)%nnode
     256              : 
     257       773864 :             iac = ikind + nkind*(kkind - 1)
     258       773864 :             IF (.NOT. ASSOCIATED(basis_set(ikind)%gto_basis_set)) CYCLE
     259       773864 :             IF (.NOT. ASSOCIATED(gpotential(kkind)%gth_potential) .AND. &
     260              :                 .NOT. ASSOCIATED(spotential(kkind)%sgp_potential)) CYCLE
     261       773864 :             IF (.NOT. ASSOCIATED(sap_int(iac)%alist)) THEN
     262        31766 :                sap_int(iac)%a_kind = ikind
     263        31766 :                sap_int(iac)%p_kind = kkind
     264        31766 :                sap_int(iac)%nalist = nlist
     265       158650 :                ALLOCATE (sap_int(iac)%alist(nlist))
     266        95118 :                DO i = 1, nlist
     267        63352 :                   NULLIFY (sap_int(iac)%alist(i)%clist)
     268        63352 :                   sap_int(iac)%alist(i)%aatom = 0
     269        95118 :                   sap_int(iac)%alist(i)%nclist = 0
     270              :                END DO
     271              :             END IF
     272       788882 :             IF (.NOT. ASSOCIATED(sap_int(iac)%alist(ilist)%clist)) THEN
     273        63294 :                sap_int(iac)%alist(ilist)%aatom = iatom
     274        63294 :                sap_int(iac)%alist(ilist)%nclist = nneighbor
     275      1343510 :                ALLOCATE (sap_int(iac)%alist(ilist)%clist(nneighbor))
     276       837158 :                DO i = 1, nneighbor
     277       837158 :                   sap_int(iac)%alist(ilist)%clist(i)%catom = 0
     278              :                END DO
     279              :             END IF
     280              :          END DO
     281              : 
     282              :          ! Calculate the overlap integrals <a|p>
     283              : !$OMP PARALLEL &
     284              : !$OMP DEFAULT (NONE) &
     285              : !$OMP SHARED  (basis_set, gpotential, spotential, maxder, ncoset, &
     286              : !$OMP          sap_ppnl, sap_int, nkind, ldsab, ldai, nder, nco, do_soc ) &
     287              : !$OMP PRIVATE (ikind, kkind, iatom, katom, nlist, ilist, nneighbor, jneighbor, &
     288              : !$OMP          cell_c, rac, iac, first_sgfa, la_max, la_min, npgfa, nseta, nsgfa, nsgf_seta, &
     289              : !$OMP          slot, sphi_a, zeta, cprj, hprj, lppnl, nppnl, nprj_ppnl, &
     290              : !$OMP          clist, iset, ncoa, sgfa, prjc, work, work_l, sab, lab, ai_work, nprjc, &
     291              : !$OMP          ppnl_radius, ncoc, rpgfa, first_col, vprj_ppnl, wprj_ppnl, i, j, l, do_gth, &
     292              : !$OMP          set_radius_a, rprjc, dac, lc_max, lc_min, zetc, alpha_ppnl, &
     293        15018 : !$OMP          na, nb, np, nnl, a_nl, radp, i_dim, ib, jb)
     294              : 
     295              :          ALLOCATE (sab(ldsab, ldsab*maxder), work(ldsab, ldsab*maxder))
     296              :          sab = 0.0_dp
     297              :          ALLOCATE (ai_work(ldai, ldai, ncoset(nder + 1)))
     298              :          ai_work = 0.0_dp
     299              :          IF (do_soc) THEN
     300              :             ALLOCATE (lab(ldsab, ldsab, 3), work_l(ldsab, ldsab, 3))
     301              :             lab = 0.0_dp
     302              :          END IF
     303              : 
     304              : !$OMP DO SCHEDULE(GUIDED)
     305              :          DO slot = 1, sap_ppnl(1)%nl_size
     306              : 
     307              :             ikind = sap_ppnl(1)%nlist_task(slot)%ikind
     308              :             kkind = sap_ppnl(1)%nlist_task(slot)%jkind
     309              :             iatom = sap_ppnl(1)%nlist_task(slot)%iatom
     310              :             katom = sap_ppnl(1)%nlist_task(slot)%jatom
     311              :             nlist = sap_ppnl(1)%nlist_task(slot)%nlist
     312              :             ilist = sap_ppnl(1)%nlist_task(slot)%ilist
     313              :             nneighbor = sap_ppnl(1)%nlist_task(slot)%nnode
     314              :             jneighbor = sap_ppnl(1)%nlist_task(slot)%inode
     315              :             cell_c(:) = sap_ppnl(1)%nlist_task(slot)%cell(:)
     316              :             rac(1:3) = sap_ppnl(1)%nlist_task(slot)%r(1:3)
     317              : 
     318              :             iac = ikind + nkind*(kkind - 1)
     319              :             IF (.NOT. ASSOCIATED(basis_set(ikind)%gto_basis_set)) CYCLE
     320              :             ! Get definition of basis set
     321              :             first_sgfa => basis_set(ikind)%gto_basis_set%first_sgf
     322              :             la_max => basis_set(ikind)%gto_basis_set%lmax
     323              :             la_min => basis_set(ikind)%gto_basis_set%lmin
     324              :             npgfa => basis_set(ikind)%gto_basis_set%npgf
     325              :             nseta = basis_set(ikind)%gto_basis_set%nset
     326              :             nsgfa = basis_set(ikind)%gto_basis_set%nsgf
     327              :             nsgf_seta => basis_set(ikind)%gto_basis_set%nsgf_set
     328              :             rpgfa => basis_set(ikind)%gto_basis_set%pgf_radius
     329              :             set_radius_a => basis_set(ikind)%gto_basis_set%set_radius
     330              :             sphi_a => basis_set(ikind)%gto_basis_set%sphi
     331              :             zeta => basis_set(ikind)%gto_basis_set%zet
     332              :             ! Get definition of PP projectors
     333              :             IF (ASSOCIATED(gpotential(kkind)%gth_potential)) THEN
     334              :                ! GTH potential
     335              :                do_gth = .TRUE.
     336              :                alpha_ppnl => gpotential(kkind)%gth_potential%alpha_ppnl
     337              :                cprj => gpotential(kkind)%gth_potential%cprj
     338              :                lppnl = gpotential(kkind)%gth_potential%lppnl
     339              :                nppnl = gpotential(kkind)%gth_potential%nppnl
     340              :                nprj_ppnl => gpotential(kkind)%gth_potential%nprj_ppnl
     341              :                ppnl_radius = gpotential(kkind)%gth_potential%ppnl_radius
     342              :                vprj_ppnl => gpotential(kkind)%gth_potential%vprj_ppnl
     343              :                wprj_ppnl => gpotential(kkind)%gth_potential%wprj_ppnl
     344              :             ELSE IF (ASSOCIATED(spotential(kkind)%sgp_potential)) THEN
     345              :                ! SGP potential
     346              :                do_gth = .FALSE.
     347              :                nprjc = spotential(kkind)%sgp_potential%nppnl
     348              :                IF (nprjc == 0) CYCLE
     349              :                nnl = spotential(kkind)%sgp_potential%n_nonlocal
     350              :                lppnl = spotential(kkind)%sgp_potential%lmax
     351              :                a_nl => spotential(kkind)%sgp_potential%a_nonlocal
     352              :                ppnl_radius = spotential(kkind)%sgp_potential%ppnl_radius
     353              :                ALLOCATE (radp(nnl))
     354              :                radp(:) = ppnl_radius
     355              :                cprj => spotential(kkind)%sgp_potential%cprj_ppnl
     356              :                hprj => spotential(kkind)%sgp_potential%vprj_ppnl
     357              :                nppnl = SIZE(cprj, 2)
     358              :             ELSE
     359              :                CYCLE
     360              :             END IF
     361              : 
     362              :             dac = SQRT(SUM(rac*rac))
     363              :             clist => sap_int(iac)%alist(ilist)%clist(jneighbor)
     364              :             clist%catom = katom
     365              :             clist%cell = cell_c
     366              :             clist%rac = rac
     367              :             ALLOCATE (clist%acint(nsgfa, nppnl, maxder), &
     368              :                       clist%achint(nsgfa, nppnl, maxder), &
     369              :                       clist%alint(nsgfa, nppnl, 3), &
     370              :                       clist%alkint(nsgfa, nppnl, 3))
     371              :             clist%acint = 0.0_dp
     372              :             clist%achint = 0.0_dp
     373              :             clist%alint = 0.0_dp
     374              :             clist%alkint = 0.0_dp
     375              : 
     376              :             clist%nsgf_cnt = 0
     377              :             NULLIFY (clist%sgf_list)
     378              :             DO iset = 1, nseta
     379              :                ncoa = npgfa(iset)*ncoset(la_max(iset))
     380              :                sgfa = first_sgfa(1, iset)
     381              :                IF (do_gth) THEN
     382              :                   ! GTH potential
     383              :                   prjc = 1
     384              :                   work = 0.0_dp
     385              :                   DO l = 0, lppnl
     386              :                      nprjc = nprj_ppnl(l)*nco(l)
     387              :                      IF (nprjc == 0) CYCLE
     388              :                      rprjc(1) = ppnl_radius
     389              :                      IF (set_radius_a(iset) + rprjc(1) < dac) CYCLE
     390              :                      lc_max = l + 2*(nprj_ppnl(l) - 1)
     391              :                      lc_min = l
     392              :                      zetc(1) = alpha_ppnl(l)
     393              :                      ncoc = ncoset(lc_max)
     394              : 
     395              :                      ! Calculate the primitive overlap integrals
     396              :                      CALL overlap(la_max(iset), la_min(iset), npgfa(iset), rpgfa(:, iset), zeta(:, iset), &
     397              :                                   lc_max, lc_min, 1, rprjc, zetc, rac, dac, sab, nder, .TRUE., ai_work, ldai)
     398              :                      ! Transformation step projector functions (Cartesian -> spherical)
     399              :                      na = ncoa
     400              :                      nb = nprjc
     401              :                      np = ncoc
     402              :                      DO i = 1, maxder
     403              :                         first_col = (i - 1)*ldsab
     404              :                         ! CALL dgemm("N", "N", ncoa, nprjc, ncoc, 1.0_dp, sab(1, first_col + 1), SIZE(sab, 1), &
     405              :                         !            cprj(1, prjc), SIZE(cprj, 1), 0.0_dp, work(1, first_col + prjc), ldsab)
     406              :                         work(1:na, first_col + prjc:first_col + prjc + nb - 1) = &
     407              :                            MATMUL(sab(1:na, first_col + 1:first_col + np), cprj(1:np, prjc:prjc + nb - 1))
     408              :                      END DO
     409              : 
     410              :                      IF (do_soc) THEN
     411              :                         ! Calculate the primitive angular momentum integrals needed for spin-orbit coupling
     412              :                         lab = 0.0_dp
     413              :                         CALL angmom(la_max(iset), npgfa(iset), zeta(:, iset), rpgfa(:, iset), la_min(iset), &
     414              :                                     lc_max, 1, zetc, rprjc, -rac, [0._dp, 0._dp, 0._dp], lab)
     415              :                         DO i_dim = 1, 3
     416              :                            work_l(1:na, prjc:prjc + nb - 1, i_dim) = &
     417              :                               MATMUL(lab(1:na, 1:np, i_dim), cprj(1:np, prjc:prjc + nb - 1))
     418              :                         END DO
     419              :                      END IF
     420              : 
     421              :                      prjc = prjc + nprjc
     422              : 
     423              :                   END DO
     424              :                   na = nsgf_seta(iset)
     425              :                   nb = nppnl
     426              :                   np = ncoa
     427              :                   DO i = 1, maxder
     428              :                      first_col = (i - 1)*ldsab + 1
     429              :                      ! Contraction step (basis functions)
     430              :                      ! CALL dgemm("T", "N", nsgf_seta(iset), nppnl, ncoa, 1.0_dp, sphi_a(1, sgfa), SIZE(sphi_a, 1), &
     431              :                      !           work(1, first_col), ldsab, 0.0_dp, clist%acint(sgfa, 1, i), nsgfa)
     432              :                      clist%acint(sgfa:sgfa + na - 1, 1:nb, i) = &
     433              :                         MATMUL(TRANSPOSE(sphi_a(1:np, sgfa:sgfa + na - 1)), work(1:np, first_col:first_col + nb - 1))
     434              :                      ! Multiply with interaction matrix(h)
     435              :                      ! CALL dgemm("N", "N", nsgf_seta(iset), nppnl, nppnl, 1.0_dp, clist%acint(sgfa, 1, i), nsgfa, &
     436              :                      !            vprj_ppnl(1, 1), SIZE(vprj_ppnl, 1), 0.0_dp, clist%achint(sgfa, 1, i), nsgfa)
     437              :                      clist%achint(sgfa:sgfa + na - 1, 1:nb, i) = &
     438              :                         MATMUL(clist%acint(sgfa:sgfa + na - 1, 1:nb, i), vprj_ppnl(1:nb, 1:nb))
     439              :                   END DO
     440              :                   IF (do_soc) THEN
     441              :                      DO i_dim = 1, 3
     442              :                         clist%alint(sgfa:sgfa + na - 1, 1:nb, i_dim) = &
     443              :                            MATMUL(TRANSPOSE(sphi_a(1:np, sgfa:sgfa + na - 1)), work_l(1:np, 1:nb, i_dim))
     444              :                         clist%alkint(sgfa:sgfa + na - 1, 1:nb, i_dim) = &
     445              :                            MATMUL(clist%alint(sgfa:sgfa + na - 1, 1:nb, i_dim), wprj_ppnl(1:nb, 1:nb))
     446              :                      END DO
     447              :                   END IF
     448              :                ELSE
     449              :                   ! SGP potential
     450              :                   ! Calculate the primitive overlap integrals
     451              :                   CALL overlap(la_max(iset), la_min(iset), npgfa(iset), rpgfa(:, iset), zeta(:, iset), &
     452              :                                lppnl, 0, nnl, radp, a_nl, rac, dac, sab, nder, .TRUE., ai_work, ldai)
     453              :                   na = nsgf_seta(iset)
     454              :                   nb = nppnl
     455              :                   np = ncoa
     456              :                   DO i = 1, maxder
     457              :                      first_col = (i - 1)*ldsab + 1
     458              :                      ! Transformation step projector functions (cartesian->spherical)
     459              :                      ! CALL dgemm("N", "N", ncoa, nppnl, nprjc, 1.0_dp, sab(1, first_col), ldsab, &
     460              :                      !            cprj(1, 1), SIZE(cprj, 1), 0.0_dp, work(1, 1), ldsab)
     461              :                      work(1:np, 1:nb) = MATMUL(sab(1:np, first_col:first_col + nprjc - 1), cprj(1:nprjc, 1:nb))
     462              :                      ! Contraction step (basis functions)
     463              :                      ! CALL dgemm("T", "N", nsgf_seta(iset), nppnl, ncoa, 1.0_dp, sphi_a(1, sgfa), SIZE(sphi_a, 1), &
     464              :                      !            work(1, 1), ldsab, 0.0_dp, clist%acint(sgfa, 1, i), nsgfa)
     465              :                      clist%acint(sgfa:sgfa + na - 1, 1:nb, i) = &
     466              :                         MATMUL(TRANSPOSE(sphi_a(1:np, sgfa:sgfa + na - 1)), work(1:np, 1:nb))
     467              :                      ! *** Multiply with interaction matrix(h) ***
     468              :                      ncoc = sgfa + nsgf_seta(iset) - 1
     469              :                      DO j = 1, nppnl
     470              :                         clist%achint(sgfa:ncoc, j, i) = clist%acint(sgfa:ncoc, j, i)*hprj(j)
     471              :                      END DO
     472              :                   END DO
     473              :                END IF
     474              :             END DO
     475              :             clist%maxac = MAXVAL(ABS(clist%acint(:, :, 1)))
     476              :             clist%maxach = MAXVAL(ABS(clist%achint(:, :, 1)))
     477              :             IF (.NOT. do_gth) DEALLOCATE (radp)
     478              :          END DO
     479              : 
     480              :          DEALLOCATE (sab, ai_work, work)
     481              :          IF (do_soc) DEALLOCATE (lab, work_l)
     482              : !$OMP END PARALLEL
     483              : 
     484              :          ! Set up a sorting index
     485        15018 :          CALL sap_sort(sap_int)
     486              :          ! All integrals needed have been calculated and stored in sap_int
     487              :          ! We now calculate the Hamiltonian matrix elements
     488              : 
     489       240034 :          force_thread = 0.0_dp
     490        71272 :          at_thread = 0.0_dp
     491        15018 :          pv_thread = 0.0_dp
     492              : 
     493              : !$OMP PARALLEL &
     494              : !$OMP DEFAULT (NONE) &
     495              : !$OMP SHARED  (do_kp, basis_set, matrix_h, matrix_l, cell_to_index,&
     496              : !$OMP          sab_orb, matrix_p, sap_int, nkind, eps_ppnl, force, &
     497              : !$OMP          doat, do_dR, deltaR, maxder, nder, &
     498              : !$OMP          locks, virial, use_virial, calculate_forces, do_soc, natom) &
     499              : !$OMP PRIVATE (ikind, jkind, iatom, jatom, cell_b, rab, &
     500              : !$OMP          slot, iab, atom_a, f0, irow, icol, h_block, &
     501              : !$OMP          l_block_x, l_block_y, l_block_z, lock_num, &
     502              : !$OMP          r_2block, r_3block, atk, &
     503              : !$OMP          found,p_block, iac, ibc, alist_ac, alist_bc, acint, bcint, &
     504              : !$OMP          achint, bchint, alkint, blkint, &
     505              : !$OMP          na, np, nb, katom, j, fa, fb, rbc, rac, &
     506              : !$OMP          kkind, kac, kbc, i, img, hash, iatom8) &
     507        15018 : !$OMP REDUCTION (+ : at_thread, pv_thread, force_thread )
     508              : 
     509              : !$OMP SINGLE
     510              : !$       ALLOCATE (locks(nlock))
     511              : !$OMP END SINGLE
     512              : 
     513              : !$OMP DO
     514              : !$       DO lock_num = 1, nlock
     515              : !$          call omp_init_lock(locks(lock_num))
     516              : !$       END DO
     517              : !$OMP END DO
     518              : 
     519              : !$OMP DO SCHEDULE(GUIDED)
     520              :          DO slot = 1, sab_orb(1)%nl_size
     521              : 
     522              :             ikind = sab_orb(1)%nlist_task(slot)%ikind
     523              :             jkind = sab_orb(1)%nlist_task(slot)%jkind
     524              :             iatom = sab_orb(1)%nlist_task(slot)%iatom
     525              :             jatom = sab_orb(1)%nlist_task(slot)%jatom
     526              :             cell_b(:) = sab_orb(1)%nlist_task(slot)%cell(:)
     527              :             rab(1:3) = sab_orb(1)%nlist_task(slot)%r(1:3)
     528              : 
     529              :             IF (.NOT. ASSOCIATED(basis_set(ikind)%gto_basis_set)) CYCLE
     530              :             IF (.NOT. ASSOCIATED(basis_set(jkind)%gto_basis_set)) CYCLE
     531              : 
     532              :             iab = ikind + nkind*(jkind - 1)
     533              : 
     534              :             ! Use the symmetry of the first derivatives
     535              :             IF (iatom == jatom) THEN
     536              :                f0 = 1.0_dp
     537              :             ELSE
     538              :                f0 = 2.0_dp
     539              :             END IF
     540              : 
     541              :             IF (do_kp) THEN
     542              :                img = cell_to_index(cell_b(1), cell_b(2), cell_b(3))
     543              :             ELSE
     544              :                img = 1
     545              :             END IF
     546              : 
     547              :             ! Create matrix blocks for a new matrix block column
     548              :             IF (iatom <= jatom) THEN
     549              :                irow = iatom
     550              :                icol = jatom
     551              :             ELSE
     552              :                irow = jatom
     553              :                icol = iatom
     554              :             END IF
     555              :             NULLIFY (h_block)
     556              :             CALL dbcsr_get_block_p(matrix_h(1, img)%matrix, irow, icol, h_block, found)
     557              :             IF (do_soc) THEN
     558              :                NULLIFY (l_block_x, l_block_y, l_block_z)
     559              :                CALL dbcsr_get_block_p(matrix_l(1, img)%matrix, irow, icol, l_block_x, found)
     560              :                CALL dbcsr_get_block_p(matrix_l(2, img)%matrix, irow, icol, l_block_y, found)
     561              :                CALL dbcsr_get_block_p(matrix_l(3, img)%matrix, irow, icol, l_block_z, found)
     562              :             END IF
     563              : 
     564              :             IF (do_dR) THEN
     565              :                NULLIFY (r_2block, r_3block)
     566              :                CALL dbcsr_get_block_p(matrix_h(2, img)%matrix, irow, icol, r_2block, found)
     567              :                CALL dbcsr_get_block_p(matrix_h(3, img)%matrix, irow, icol, r_3block, found)
     568              :             END IF
     569              : 
     570              :             IF (calculate_forces .OR. doat) THEN
     571              :                NULLIFY (p_block)
     572              :                CALL dbcsr_get_block_p(matrix_p(1, img)%matrix, irow, icol, p_block, found)
     573              :             END IF
     574              : 
     575              :             ! loop over all kinds for projector atom
     576              :             IF (ASSOCIATED(h_block)) THEN
     577              : !$             iatom8 = INT(iatom - 1, int_8)*INT(natom, int_8) + INT(jatom, int_8)
     578              : !$             hash = INT(MOD(iatom8, INT(nlock, int_8)) + 1)
     579              : 
     580              :                DO kkind = 1, nkind
     581              :                   iac = ikind + nkind*(kkind - 1)
     582              :                   ibc = jkind + nkind*(kkind - 1)
     583              :                   IF (.NOT. ASSOCIATED(sap_int(iac)%alist)) CYCLE
     584              :                   IF (.NOT. ASSOCIATED(sap_int(ibc)%alist)) CYCLE
     585              :                   CALL get_alist(sap_int(iac), alist_ac, iatom)
     586              :                   CALL get_alist(sap_int(ibc), alist_bc, jatom)
     587              :                   IF (.NOT. ASSOCIATED(alist_ac)) CYCLE
     588              :                   IF (.NOT. ASSOCIATED(alist_bc)) CYCLE
     589              :                   DO kac = 1, alist_ac%nclist
     590              :                      DO kbc = 1, alist_bc%nclist
     591              :                         IF (alist_ac%clist(kac)%catom /= alist_bc%clist(kbc)%catom) CYCLE
     592              :                         IF (ALL(cell_b + alist_bc%clist(kbc)%cell - alist_ac%clist(kac)%cell == 0)) THEN
     593              :                            IF (alist_ac%clist(kac)%maxac*alist_bc%clist(kbc)%maxach < eps_ppnl) CYCLE
     594              :                            acint => alist_ac%clist(kac)%acint
     595              :                            bcint => alist_bc%clist(kbc)%acint
     596              :                            achint => alist_ac%clist(kac)%achint
     597              :                            bchint => alist_bc%clist(kbc)%achint
     598              :                            IF (do_soc) THEN
     599              :                               alkint => alist_ac%clist(kac)%alkint
     600              :                               blkint => alist_bc%clist(kbc)%alkint
     601              :                            END IF
     602              :                            na = SIZE(acint, 1)
     603              :                            np = SIZE(acint, 2)
     604              :                            nb = SIZE(bcint, 1)
     605              : !$                         CALL omp_set_lock(locks(hash))
     606              :                            IF (.NOT. do_dR) THEN
     607              :                               IF (iatom <= jatom) THEN
     608              :                                  h_block(1:na, 1:nb) = h_block(1:na, 1:nb) + &
     609              :                                                        MATMUL(achint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, 1)))
     610              :                               ELSE
     611              :                                  h_block(1:nb, 1:na) = h_block(1:nb, 1:na) + &
     612              :                                                        MATMUL(bchint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, 1)))
     613              :                               END IF
     614              :                            END IF
     615              :                            IF (do_soc) THEN
     616              :                               IF (iatom <= jatom) THEN
     617              :                                  l_block_x(1:na, 1:nb) = l_block_x(1:na, 1:nb) + &
     618              :                                                          MATMUL(alkint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, 1)))
     619              :                                  l_block_y(1:na, 1:nb) = l_block_y(1:na, 1:nb) + &
     620              :                                                          MATMUL(alkint(1:na, 1:np, 2), TRANSPOSE(bcint(1:nb, 1:np, 1)))
     621              :                                  l_block_z(1:na, 1:nb) = l_block_z(1:na, 1:nb) + &
     622              :                                                          MATMUL(alkint(1:na, 1:np, 3), TRANSPOSE(bcint(1:nb, 1:np, 1)))
     623              : 
     624              :                               ELSE
     625              :                                  l_block_x(1:nb, 1:na) = l_block_x(1:nb, 1:na) - &
     626              :                                                          MATMUL(blkint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, 1)))
     627              :                                  l_block_y(1:nb, 1:na) = l_block_y(1:nb, 1:na) - &
     628              :                                                          MATMUL(blkint(1:nb, 1:np, 2), TRANSPOSE(acint(1:na, 1:np, 1)))
     629              :                                  l_block_z(1:nb, 1:na) = l_block_z(1:nb, 1:na) - &
     630              :                                                          MATMUL(blkint(1:nb, 1:np, 3), TRANSPOSE(acint(1:na, 1:np, 1)))
     631              :                               END IF
     632              :                            END IF
     633              : !$                         CALL omp_unset_lock(locks(hash))
     634              :                            IF (calculate_forces) THEN
     635              :                               IF (ASSOCIATED(p_block)) THEN
     636              :                                  katom = alist_ac%clist(kac)%catom
     637              :                                  DO i = 1, 3
     638              :                                     j = i + 1
     639              :                                     IF (iatom <= jatom) THEN
     640              :                                        fa(i) = SUM(p_block(1:na, 1:nb)* &
     641              :                                                    MATMUL(acint(1:na, 1:np, j), TRANSPOSE(bchint(1:nb, 1:np, 1))))
     642              :                                        fb(i) = SUM(p_block(1:na, 1:nb)* &
     643              :                                                    MATMUL(achint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, j))))
     644              :                                     ELSE
     645              :                                        fa(i) = SUM(p_block(1:nb, 1:na)* &
     646              :                                                    MATMUL(bchint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, j))))
     647              :                                        fb(i) = SUM(p_block(1:nb, 1:na)* &
     648              :                                                    MATMUL(bcint(1:nb, 1:np, j), TRANSPOSE(achint(1:na, 1:np, 1))))
     649              :                                     END IF
     650              :                                     force_thread(i, iatom) = force_thread(i, iatom) + f0*fa(i)
     651              :                                     force_thread(i, katom) = force_thread(i, katom) - f0*fa(i)
     652              :                                     force_thread(i, jatom) = force_thread(i, jatom) + f0*fb(i)
     653              :                                     force_thread(i, katom) = force_thread(i, katom) - f0*fb(i)
     654              :                                  END DO
     655              : 
     656              :                                  IF (use_virial) THEN
     657              :                                     rac = alist_ac%clist(kac)%rac
     658              :                                     rbc = alist_bc%clist(kbc)%rac
     659              :                                     CALL virial_pair_force(pv_thread, f0, fa, rac)
     660              :                                     CALL virial_pair_force(pv_thread, f0, fb, rbc)
     661              :                                  END IF
     662              :                               END IF
     663              :                            END IF
     664              : 
     665              :                            IF (do_dR) THEN
     666              :                               i = 1; j = 2;
     667              :                               katom = alist_ac%clist(kac)%catom
     668              :                               IF (iatom <= jatom) THEN
     669              :                                  h_block(1:na, 1:nb) = h_block(1:na, 1:nb) + &
     670              :                                                        (deltaR(i, iatom) - deltaR(i, katom))* &
     671              :                                                        MATMUL(acint(1:na, 1:np, j), TRANSPOSE(bchint(1:nb, 1:np, 1)))
     672              : 
     673              :                                  h_block(1:na, 1:nb) = h_block(1:na, 1:nb) + &
     674              :                                                        (deltaR(i, jatom) - deltaR(i, katom))* &
     675              :                                                        MATMUL(achint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, j)))
     676              :                               ELSE
     677              :                                  h_block(1:nb, 1:na) = h_block(1:nb, 1:na) + &
     678              :                                                        (deltaR(i, iatom) - deltaR(i, katom))* &
     679              :                                                        MATMUL(bchint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, j)))
     680              :                                  h_block(1:nb, 1:na) = h_block(1:nb, 1:na) + &
     681              :                                                        (deltaR(i, jatom) - deltaR(i, katom))* &
     682              :                                                        MATMUL(bcint(1:nb, 1:np, j), TRANSPOSE(achint(1:na, 1:np, 1)))
     683              :                               END IF
     684              : 
     685              :                               i = 2; j = 3;
     686              :                               katom = alist_ac%clist(kac)%catom
     687              :                               IF (iatom <= jatom) THEN
     688              :                                  r_2block(1:na, 1:nb) = r_2block(1:na, 1:nb) + &
     689              :                                                         (deltaR(i, iatom) - deltaR(i, katom))* &
     690              :                                                         MATMUL(acint(1:na, 1:np, j), TRANSPOSE(bchint(1:nb, 1:np, 1)))
     691              : 
     692              :                                  r_2block(1:na, 1:nb) = r_2block(1:na, 1:nb) + &
     693              :                                                         (deltaR(i, jatom) - deltaR(i, katom))* &
     694              :                                                         MATMUL(achint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, j)))
     695              :                               ELSE
     696              :                                  r_2block(1:nb, 1:na) = r_2block(1:nb, 1:na) + &
     697              :                                                         (deltaR(i, iatom) - deltaR(i, katom))* &
     698              :                                                         MATMUL(bchint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, j)))
     699              :                                  r_2block(1:nb, 1:na) = r_2block(1:nb, 1:na) + &
     700              :                                                         (deltaR(i, jatom) - deltaR(i, katom))* &
     701              :                                                         MATMUL(bcint(1:nb, 1:np, j), TRANSPOSE(achint(1:na, 1:np, 1)))
     702              :                               END IF
     703              : 
     704              :                               i = 3; j = 4;
     705              :                               katom = alist_ac%clist(kac)%catom
     706              :                               IF (iatom <= jatom) THEN
     707              :                                  r_3block(1:na, 1:nb) = r_3block(1:na, 1:nb) + &
     708              :                                                         (deltaR(i, iatom) - deltaR(i, katom))* &
     709              :                                                         MATMUL(acint(1:na, 1:np, j), TRANSPOSE(bchint(1:nb, 1:np, 1)))
     710              : 
     711              :                                  r_3block(1:na, 1:nb) = r_3block(1:na, 1:nb) + &
     712              :                                                         (deltaR(i, jatom) - deltaR(i, katom))* &
     713              :                                                         MATMUL(achint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, j)))
     714              :                               ELSE
     715              :                                  r_3block(1:nb, 1:na) = r_3block(1:nb, 1:na) + &
     716              :                                                         (deltaR(i, iatom) - deltaR(i, katom))* &
     717              :                                                         MATMUL(bchint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, j)))
     718              :                                  r_3block(1:nb, 1:na) = r_3block(1:nb, 1:na) + &
     719              :                                                         (deltaR(i, jatom) - deltaR(i, katom))* &
     720              :                                                         MATMUL(bcint(1:nb, 1:np, j), TRANSPOSE(achint(1:na, 1:np, 1)))
     721              :                               END IF
     722              : 
     723              :                            END IF
     724              :                            IF (doat) THEN
     725              :                               IF (ASSOCIATED(p_block)) THEN
     726              :                                  katom = alist_ac%clist(kac)%catom
     727              :                                  IF (iatom <= jatom) THEN
     728              :                                     atk = SUM(p_block(1:na, 1:nb)* &
     729              :                                               MATMUL(achint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, 1))))
     730              :                                  ELSE
     731              :                                     atk = SUM(p_block(1:nb, 1:na)* &
     732              :                                               MATMUL(bchint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, 1))))
     733              :                                  END IF
     734              :                                  at_thread(katom) = at_thread(katom) + f0*atk
     735              :                               END IF
     736              :                            END IF
     737              :                            EXIT ! We have found a match and there can be only one single match
     738              :                         END IF
     739              :                      END DO
     740              :                   END DO
     741              :                END DO
     742              :             END IF
     743              :          END DO
     744              : 
     745              : !$OMP DO
     746              : !$       DO lock_num = 1, nlock
     747              : !$          call omp_destroy_lock(locks(lock_num))
     748              : !$       END DO
     749              : !$OMP END DO
     750              : 
     751              : !$OMP SINGLE
     752              : !$       DEALLOCATE (locks)
     753              : !$OMP END SINGLE NOWAIT
     754              : 
     755              : !$OMP END PARALLEL
     756              : 
     757        15018 :          CALL release_sap_int(sap_int)
     758              : 
     759        15018 :          DEALLOCATE (basis_set, gpotential, spotential)
     760        15018 :          IF (calculate_forces) THEN
     761         6431 :             CALL get_atomic_kind_set(atomic_kind_set, atom_of_kind=atom_of_kind, kind_of=kind_of)
     762              : !$OMP DO
     763              :             DO iatom = 1, natom
     764        23053 :                atom_a = atom_of_kind(iatom)
     765        23053 :                ikind = kind_of(iatom)
     766        92212 :                force(ikind)%gth_ppnl(:, atom_a) = force(ikind)%gth_ppnl(:, atom_a) + force_thread(:, iatom)
     767              :             END DO
     768              : !$OMP END DO
     769         6431 :             DEALLOCATE (atom_of_kind, kind_of)
     770              :          END IF
     771              : 
     772        15018 :          IF (calculate_forces .AND. use_virial) THEN
     773        10582 :             virial%pv_ppnl = virial%pv_ppnl + pv_thread
     774        10582 :             virial%pv_virial = virial%pv_virial + pv_thread
     775              :          END IF
     776              : 
     777        15018 :          IF (doat) THEN
     778          280 :             atcore(1:natom) = atcore(1:natom) + at_thread
     779              :          END IF
     780              : 
     781        30036 :          IF (calculate_forces .OR. doat) THEN
     782              :             ! If LSD, then recover alpha density and beta density
     783              :             ! from the total density (1) and the spin density (2)
     784         6493 :             IF (SIZE(matrix_p, 1) == 2) THEN
     785         1982 :                DO img = 1, nimages
     786              :                   CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
     787         1292 :                                  alpha_scalar=0.5_dp, beta_scalar=0.5_dp)
     788              :                   CALL dbcsr_add(matrix_p(2, img)%matrix, matrix_p(1, img)%matrix, &
     789         1982 :                                  alpha_scalar=-1.0_dp, beta_scalar=1.0_dp)
     790              :                END DO
     791              :             END IF
     792              :          END IF
     793              : 
     794              :       END IF !ppnl_present
     795              : 
     796        15018 :       CALL timestop(handle)
     797              : 
     798        30036 :    END SUBROUTINE build_core_ppnl
     799              : 
     800              : END MODULE core_ppnl
        

Generated by: LCOV version 2.0-1