LCOV - code coverage report
Current view: top level - src - iao_analysis.F (source / functions) Coverage Total Hit
Test: CP2K Regtests (git:1155b05) Lines: 98.2 % 876 860
Test Date: 2026-03-21 06:31:29 Functions: 100.0 % 14 14

            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              : ! **************************************************************************************************
       9              : !> \brief Calculate intrinsic atomic orbitals and analyze wavefunctions
      10              : !> \par History
      11              : !>      03.2023 created [JGH]
      12              : !> \author JGH
      13              : ! **************************************************************************************************
      14              : MODULE iao_analysis
      15              :    USE ai_contraction,                  ONLY: block_add,&
      16              :                                               contraction
      17              :    USE ai_overlap,                      ONLY: overlap_ab
      18              :    USE atomic_charges,                  ONLY: print_atomic_charges
      19              :    USE atomic_kind_types,               ONLY: atomic_kind_type,&
      20              :                                               get_atomic_kind
      21              :    USE auto_basis,                      ONLY: create_oce_basis
      22              :    USE basis_set_types,                 ONLY: deallocate_gto_basis_set,&
      23              :                                               get_gto_basis_set,&
      24              :                                               gto_basis_set_p_type,&
      25              :                                               gto_basis_set_type,&
      26              :                                               init_orb_basis_set
      27              :    USE bibliography,                    ONLY: Knizia2013,&
      28              :                                               cite_reference
      29              :    USE cell_types,                      ONLY: cell_type,&
      30              :                                               pbc
      31              :    USE cp_array_utils,                  ONLY: cp_2d_r_p_type
      32              :    USE cp_control_types,                ONLY: dft_control_type
      33              :    USE cp_dbcsr_api,                    ONLY: &
      34              :         dbcsr_copy, dbcsr_create, dbcsr_distribution_type, dbcsr_get_block_p, dbcsr_get_info, &
      35              :         dbcsr_iterator_blocks_left, dbcsr_iterator_next_block, dbcsr_iterator_start, &
      36              :         dbcsr_iterator_stop, dbcsr_iterator_type, dbcsr_multiply, dbcsr_p_type, dbcsr_release, &
      37              :         dbcsr_replicate_all, dbcsr_set, dbcsr_type, dbcsr_type_no_symmetry
      38              :    USE cp_dbcsr_contrib,                ONLY: dbcsr_reserve_diag_blocks,&
      39              :                                               dbcsr_trace
      40              :    USE cp_dbcsr_operations,             ONLY: copy_dbcsr_to_fm,&
      41              :                                               copy_fm_to_dbcsr,&
      42              :                                               cp_dbcsr_plus_fm_fm_t,&
      43              :                                               cp_dbcsr_sm_fm_multiply,&
      44              :                                               dbcsr_deallocate_matrix_set
      45              :    USE cp_fm_basic_linalg,              ONLY: cp_fm_column_scale,&
      46              :                                               cp_fm_geadd,&
      47              :                                               cp_fm_invert,&
      48              :                                               cp_fm_rot_cols,&
      49              :                                               cp_fm_rot_rows
      50              :    USE cp_fm_struct,                    ONLY: cp_fm_struct_create,&
      51              :                                               cp_fm_struct_release,&
      52              :                                               cp_fm_struct_type
      53              :    USE cp_fm_types,                     ONLY: cp_fm_create,&
      54              :                                               cp_fm_get_diag,&
      55              :                                               cp_fm_get_element,&
      56              :                                               cp_fm_get_info,&
      57              :                                               cp_fm_release,&
      58              :                                               cp_fm_set_all,&
      59              :                                               cp_fm_to_fm,&
      60              :                                               cp_fm_type
      61              :    USE cp_log_handling,                 ONLY: cp_get_default_logger,&
      62              :                                               cp_logger_type
      63              :    USE cp_output_handling,              ONLY: cp_print_key_finished_output,&
      64              :                                               cp_print_key_unit_nr
      65              :    USE cp_realspace_grid_cube,          ONLY: cp_pw_to_cube
      66              :    USE iao_types,                       ONLY: iao_env_type
      67              :    USE input_constants,                 ONLY: do_iaoloc_enone,&
      68              :                                               do_iaoloc_pm2
      69              :    USE input_section_types,             ONLY: section_get_ivals,&
      70              :                                               section_vals_get,&
      71              :                                               section_vals_type,&
      72              :                                               section_vals_val_get
      73              :    USE kinds,                           ONLY: default_path_length,&
      74              :                                               default_string_length,&
      75              :                                               dp
      76              :    USE machine,                         ONLY: m_walltime
      77              :    USE mathconstants,                   ONLY: twopi
      78              :    USE mathlib,                         ONLY: invmat_symm,&
      79              :                                               jacobi
      80              :    USE message_passing,                 ONLY: mp_comm_type
      81              :    USE min_basis_set,                   ONLY: create_minbas_set
      82              :    USE molden_utils,                    ONLY: write_mos_molden
      83              :    USE orbital_pointers,                ONLY: ncoset
      84              :    USE parallel_gemm_api,               ONLY: parallel_gemm
      85              :    USE particle_list_types,             ONLY: particle_list_type
      86              :    USE particle_methods,                ONLY: get_particle_set
      87              :    USE particle_types,                  ONLY: particle_type
      88              :    USE pw_env_types,                    ONLY: pw_env_get,&
      89              :                                               pw_env_type
      90              :    USE pw_pool_types,                   ONLY: pw_pool_type
      91              :    USE pw_types,                        ONLY: pw_c1d_gs_type,&
      92              :                                               pw_r3d_rs_type
      93              :    USE qs_collocate_density,            ONLY: calculate_wavefunction
      94              :    USE qs_environment_types,            ONLY: get_qs_env,&
      95              :                                               qs_environment_type
      96              :    USE qs_interactions,                 ONLY: init_interaction_radii_orb_basis
      97              :    USE qs_kind_types,                   ONLY: get_qs_kind,&
      98              :                                               qs_kind_type
      99              :    USE qs_ks_types,                     ONLY: get_ks_env,&
     100              :                                               qs_ks_env_type
     101              :    USE qs_loc_utils,                    ONLY: compute_berry_operator
     102              :    USE qs_localization_methods,         ONLY: initialize_weights,&
     103              :                                               rotate_orbitals,&
     104              :                                               scdm_qrfact
     105              :    USE qs_mo_methods,                   ONLY: make_basis_lowdin
     106              :    USE qs_mo_types,                     ONLY: allocate_mo_set,&
     107              :                                               deallocate_mo_set,&
     108              :                                               get_mo_set,&
     109              :                                               init_mo_set,&
     110              :                                               mo_set_type,&
     111              :                                               set_mo_set
     112              :    USE qs_moments,                      ONLY: build_local_moment_matrix
     113              :    USE qs_neighbor_list_types,          ONLY: neighbor_list_set_p_type,&
     114              :                                               release_neighbor_list_sets
     115              :    USE qs_neighbor_lists,               ONLY: setup_neighbor_list
     116              :    USE qs_overlap,                      ONLY: build_overlap_matrix_simple
     117              :    USE qs_subsys_types,                 ONLY: qs_subsys_get,&
     118              :                                               qs_subsys_type
     119              : #include "./base/base_uses.f90"
     120              : 
     121              :    IMPLICIT NONE
     122              :    PRIVATE
     123              : 
     124              :    CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'iao_analysis'
     125              : 
     126              :    PUBLIC ::  iao_wfn_analysis, iao_charges, iao_calculate_dmat, print_center_spread
     127              : 
     128              :    INTERFACE iao_calculate_dmat
     129              :       MODULE PROCEDURE iao_calculate_dmat_diag, & ! diagonal occupation matrix
     130              :          iao_calculate_dmat_full    ! full occupation matrix
     131              :    END INTERFACE
     132              : 
     133              : ! **************************************************************************************************
     134              : 
     135              : CONTAINS
     136              : 
     137              : ! **************************************************************************************************
     138              : !> \brief ...
     139              : !> \param qs_env ...
     140              : !> \param iao_env ...
     141              : !> \param unit_nr ...
     142              : !> \param c_iao_coef ...
     143              : !> \param mos ...
     144              : !> \param bond_centers ...
     145              : ! **************************************************************************************************
     146           34 :    SUBROUTINE iao_wfn_analysis(qs_env, iao_env, unit_nr, c_iao_coef, mos, bond_centers)
     147              :       TYPE(qs_environment_type), POINTER                 :: qs_env
     148              :       TYPE(iao_env_type), INTENT(IN)                     :: iao_env
     149              :       INTEGER, INTENT(IN)                                :: unit_nr
     150              :       TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:), &
     151              :          OPTIONAL                                        :: c_iao_coef
     152              :       TYPE(mo_set_type), DIMENSION(:), OPTIONAL, POINTER :: mos
     153              :       REAL(KIND=dp), DIMENSION(:, :, :), OPTIONAL        :: bond_centers
     154              : 
     155              :       CHARACTER(len=*), PARAMETER                        :: routineN = 'iao_wfn_analysis'
     156              : 
     157              :       CHARACTER(LEN=2)                                   :: element_symbol
     158              :       CHARACTER(LEN=default_string_length)               :: bname
     159              :       INTEGER                                            :: handle, i, iatom, ikind, isgf, ispin, &
     160              :                                                             nao, natom, nimages, nkind, no, norb, &
     161              :                                                             nref, ns, nsgf, nspin, nvec, nx, order
     162           34 :       INTEGER, ALLOCATABLE, DIMENSION(:)                 :: first_sgf, last_sgf, nsgfat
     163              :       INTEGER, DIMENSION(2)                              :: nocc
     164              :       LOGICAL                                            :: cubes_iao, cubes_ibo, molden_iao, &
     165              :                                                             molden_ibo, uniform_occupation
     166              :       REAL(KIND=dp)                                      :: fin, fout, t1, t2, trace, zval
     167           34 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:)           :: fdiag
     168           34 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :)        :: mcharge
     169           34 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :)     :: moments
     170           34 :       REAL(KIND=dp), DIMENSION(:), POINTER               :: occupation_numbers
     171              :       TYPE(cell_type), POINTER                           :: cell
     172           34 :       TYPE(cp_2d_r_p_type), ALLOCATABLE, DIMENSION(:)    :: smat_kind
     173           34 :       TYPE(cp_2d_r_p_type), ALLOCATABLE, DIMENSION(:, :) :: oce_atom
     174              :       TYPE(cp_fm_struct_type), POINTER                   :: fm_struct
     175              :       TYPE(cp_fm_type)                                   :: p_orb_ref, p_ref_orb, smo
     176           34 :       TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:)        :: c_loc_orb, ciao, cloc, cvec, iao_coef
     177           34 :       TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :)     :: zij_atom
     178              :       TYPE(cp_fm_type), POINTER                          :: mo_coeff
     179           34 :       TYPE(dbcsr_p_type), DIMENSION(:), POINTER          :: sref, sro, sxo
     180           34 :       TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER       :: matrix_s
     181              :       TYPE(dbcsr_type)                                   :: dmat
     182              :       TYPE(dbcsr_type), POINTER                          :: smat
     183              :       TYPE(dft_control_type), POINTER                    :: dft_control
     184           34 :       TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER  :: oce_basis_set_list, orb_basis_set_list, &
     185           34 :                                                             ref_basis_set_list
     186              :       TYPE(gto_basis_set_type), POINTER                  :: ocebasis, orbbasis, refbasis
     187           34 :       TYPE(mo_set_type), ALLOCATABLE, DIMENSION(:)       :: mo_iao, mo_loc
     188           34 :       TYPE(mo_set_type), DIMENSION(:), POINTER           :: my_mos
     189              :       TYPE(neighbor_list_set_p_type), DIMENSION(:), &
     190           34 :          POINTER                                         :: sro_list, srr_list, sxo_list
     191           34 :       TYPE(particle_type), DIMENSION(:), POINTER         :: particle_set
     192           34 :       TYPE(qs_kind_type), DIMENSION(:), POINTER          :: qs_kind_set
     193              :       TYPE(qs_kind_type), POINTER                        :: qs_kind
     194              :       TYPE(qs_ks_env_type), POINTER                      :: ks_env
     195              :       TYPE(section_vals_type), POINTER                   :: iao_cubes_section, iao_molden_section, &
     196              :                                                             ibo_cc_section, ibo_cubes_section, &
     197              :                                                             ibo_molden_section
     198              : 
     199              :       ! only do IAO analysis if explicitly requested
     200            0 :       CPASSERT(iao_env%do_iao)
     201              : 
     202              :       ! k-points?
     203           34 :       CALL get_qs_env(qs_env, dft_control=dft_control)
     204           34 :       nspin = dft_control%nspins
     205           34 :       nimages = dft_control%nimages
     206           34 :       IF (nimages > 1) THEN
     207            0 :          IF (unit_nr > 0) THEN
     208              :             WRITE (UNIT=unit_nr, FMT="(T2,A)") &
     209            0 :                "K-Points: Intrinsic Atomic Orbitals Analysis not available."
     210              :          END IF
     211              :       END IF
     212            0 :       IF (nimages > 1) RETURN
     213              : 
     214           34 :       CALL timeset(routineN, handle)
     215              : 
     216           34 :       IF (unit_nr > 0) THEN
     217           17 :          WRITE (unit_nr, '(/,T2,A)') '!-----------------------------------------------------------------------------!'
     218           17 :          WRITE (UNIT=unit_nr, FMT="(T24,A)") "INTRINSIC ATOMIC ORBITALS ANALYSIS"
     219           17 :          WRITE (UNIT=unit_nr, FMT="(T13,A)") "G. Knizia, J. Chem. Theory Comput. 9, 4834-4843 (2013)"
     220           17 :          WRITE (unit_nr, '(T2,A)') '!-----------------------------------------------------------------------------!'
     221              :       END IF
     222           34 :       CALL cite_reference(Knizia2013)
     223              : 
     224              :       ! input sections
     225           34 :       iao_molden_section => iao_env%iao_molden_section
     226           34 :       iao_cubes_section => iao_env%iao_cubes_section
     227           34 :       ibo_molden_section => iao_env%ibo_molden_section
     228           34 :       ibo_cubes_section => iao_env%ibo_cubes_section
     229           34 :       ibo_cc_section => iao_env%ibo_cc_section
     230              : 
     231              :       !
     232           34 :       molden_iao = .FALSE.
     233           34 :       IF (ASSOCIATED(iao_molden_section)) THEN
     234            4 :          CALL section_vals_get(iao_molden_section, explicit=molden_iao)
     235              :       END IF
     236           34 :       cubes_iao = .FALSE.
     237           34 :       IF (ASSOCIATED(iao_cubes_section)) THEN
     238            4 :          CALL section_vals_get(iao_cubes_section, explicit=cubes_iao)
     239              :       END IF
     240           34 :       molden_ibo = .FALSE.
     241           34 :       IF (ASSOCIATED(ibo_molden_section)) THEN
     242            4 :          CALL section_vals_get(ibo_molden_section, explicit=molden_ibo)
     243              :       END IF
     244           34 :       cubes_ibo = .FALSE.
     245           34 :       IF (ASSOCIATED(ibo_cubes_section)) THEN
     246            4 :          CALL section_vals_get(ibo_cubes_section, explicit=cubes_ibo)
     247              :       END IF
     248              : 
     249              :       ! check for or generate reference basis
     250           34 :       CALL create_minbas_set(qs_env, unit_nr)
     251              : 
     252              :       ! overlap matrices
     253           34 :       CALL get_qs_env(qs_env, matrix_s_kp=matrix_s)
     254           34 :       smat => matrix_s(1, 1)%matrix
     255              :       !
     256           34 :       CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set)
     257           34 :       nkind = SIZE(qs_kind_set)
     258          276 :       ALLOCATE (ref_basis_set_list(nkind), orb_basis_set_list(nkind))
     259          104 :       DO ikind = 1, nkind
     260           70 :          qs_kind => qs_kind_set(ikind)
     261           70 :          NULLIFY (ref_basis_set_list(ikind)%gto_basis_set)
     262           70 :          NULLIFY (orb_basis_set_list(ikind)%gto_basis_set)
     263           70 :          NULLIFY (refbasis, orbbasis)
     264           70 :          CALL get_qs_kind(qs_kind=qs_kind, basis_set=orbbasis, basis_type="ORB")
     265           70 :          IF (ASSOCIATED(orbbasis)) orb_basis_set_list(ikind)%gto_basis_set => orbbasis
     266           70 :          CALL get_qs_kind(qs_kind=qs_kind, basis_set=refbasis, basis_type="MIN")
     267          104 :          IF (ASSOCIATED(refbasis)) ref_basis_set_list(ikind)%gto_basis_set => refbasis
     268              :       END DO
     269              : 
     270              :       ! neighbor lists
     271           34 :       NULLIFY (srr_list, sro_list)
     272           34 :       CALL setup_neighbor_list(srr_list, ref_basis_set_list, qs_env=qs_env)
     273           34 :       CALL setup_neighbor_list(sro_list, ref_basis_set_list, orb_basis_set_list, qs_env=qs_env)
     274              : 
     275              :       ! Srr and Sro overlap matrices
     276           34 :       NULLIFY (sref, sro)
     277           34 :       CALL get_qs_env(qs_env, ks_env=ks_env)
     278              :       CALL build_overlap_matrix_simple(ks_env, sref, &
     279           34 :                                        ref_basis_set_list, ref_basis_set_list, srr_list)
     280              :       CALL build_overlap_matrix_simple(ks_env, sro, &
     281           34 :                                        ref_basis_set_list, orb_basis_set_list, sro_list)
     282              :       !
     283           34 :       IF (PRESENT(mos)) THEN
     284           30 :          my_mos => mos
     285              :       ELSE
     286            4 :          CALL get_qs_env(qs_env, mos=my_mos)
     287              :       END IF
     288           34 :       CALL get_mo_set(my_mos(1), mo_coeff=mo_coeff)
     289           34 :       CALL dbcsr_get_info(sro(1)%matrix, nfullrows_total=nref, nfullcols_total=nao)
     290              : 
     291              :       ! Projectors
     292           34 :       NULLIFY (fm_struct)
     293              :       CALL cp_fm_struct_create(fm_struct, context=mo_coeff%matrix_struct%context, nrow_global=nref, &
     294           34 :                                ncol_global=nao, para_env=mo_coeff%matrix_struct%para_env)
     295           34 :       CALL cp_fm_create(p_ref_orb, fm_struct)
     296           34 :       CALL cp_fm_struct_release(fm_struct)
     297           34 :       NULLIFY (fm_struct)
     298              :       CALL cp_fm_struct_create(fm_struct, context=mo_coeff%matrix_struct%context, nrow_global=nao, &
     299           34 :                                ncol_global=nref, para_env=mo_coeff%matrix_struct%para_env)
     300           34 :       CALL cp_fm_create(p_orb_ref, fm_struct)
     301           34 :       CALL cp_fm_struct_release(fm_struct)
     302              :       !
     303           34 :       CALL iao_projectors(smat, sref(1)%matrix, sro(1)%matrix, p_orb_ref, p_ref_orb, iao_env%eps_svd)
     304              : 
     305              :       ! occupied orbitals, orbitals considered for IAO generation
     306           34 :       nocc(1:2) = 0
     307           70 :       DO ispin = 1, nspin
     308              :          CALL get_mo_set(my_mos(ispin), mo_coeff=mo_coeff, nmo=norb, uniform_occupation=uniform_occupation, &
     309           36 :                          occupation_numbers=occupation_numbers)
     310           36 :          IF (uniform_occupation) THEN
     311           34 :             nvec = norb
     312              :          ELSE
     313            2 :             nvec = 0
     314           34 :             DO i = 1, norb
     315           34 :                IF (occupation_numbers(i) > iao_env%eps_occ) THEN
     316           32 :                   nvec = i
     317              :                ELSE
     318              :                   EXIT
     319              :                END IF
     320              :             END DO
     321              :          END IF
     322          106 :          nocc(ispin) = nvec
     323              :       END DO
     324              :       ! generate IAOs
     325          208 :       ALLOCATE (iao_coef(nspin), cvec(nspin))
     326           70 :       DO ispin = 1, nspin
     327           36 :          CALL get_mo_set(my_mos(ispin), mo_coeff=mo_coeff)
     328           36 :          nvec = nocc(ispin)
     329           70 :          IF (nvec > 0) THEN
     330           36 :             NULLIFY (fm_struct)
     331           36 :             CALL cp_fm_struct_create(fm_struct, ncol_global=nvec, template_fmstruct=mo_coeff%matrix_struct)
     332           36 :             CALL cp_fm_create(cvec(ispin), fm_struct)
     333           36 :             CALL cp_fm_struct_release(fm_struct)
     334           36 :             CALL cp_fm_to_fm(mo_coeff, cvec(ispin), nvec)
     335              :             !
     336           36 :             NULLIFY (fm_struct)
     337           36 :             CALL cp_fm_struct_create(fm_struct, ncol_global=nref, template_fmstruct=mo_coeff%matrix_struct)
     338           36 :             CALL cp_fm_create(iao_coef(ispin), fm_struct)
     339           36 :             CALL cp_fm_struct_release(fm_struct)
     340              :             !
     341           36 :             CALL intrinsic_ao_calc(smat, p_orb_ref, p_ref_orb, cvec(ispin), iao_coef(ispin))
     342              :          END IF
     343              :       END DO
     344              :       !
     345           34 :       IF (iao_env%do_charges) THEN
     346              :          ! MOs in IAO basis
     347          104 :          ALLOCATE (ciao(nspin))
     348           70 :          DO ispin = 1, nspin
     349           36 :             CALL get_mo_set(my_mos(ispin), mo_coeff=mo_coeff, nao=nao, nmo=norb)
     350           36 :             CALL cp_fm_create(smo, mo_coeff%matrix_struct)
     351           36 :             NULLIFY (fm_struct)
     352           36 :             CALL cp_fm_struct_create(fm_struct, nrow_global=nref, template_fmstruct=mo_coeff%matrix_struct)
     353           36 :             CALL cp_fm_create(ciao(ispin), fm_struct)
     354           36 :             CALL cp_fm_struct_release(fm_struct)
     355              :             ! CIAO = A(T)*S*C
     356           36 :             CALL cp_dbcsr_sm_fm_multiply(smat, mo_coeff, smo, ncol=norb)
     357           36 :             CALL parallel_gemm('T', 'N', nref, norb, nao, 1.0_dp, iao_coef(ispin), smo, 0.0_dp, ciao(ispin))
     358          106 :             CALL cp_fm_release(smo)
     359              :          END DO
     360              :          !
     361              :          ! population analysis
     362           34 :          IF (unit_nr > 0) THEN
     363           17 :             WRITE (unit_nr, '(/,T2,A)') 'Intrinsic AO Population Analysis '
     364              :          END IF
     365           34 :          CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set, natom=natom)
     366           34 :          CALL get_ks_env(ks_env=ks_env, particle_set=particle_set)
     367          136 :          ALLOCATE (mcharge(natom, nspin))
     368           34 :          CALL dbcsr_get_info(sref(1)%matrix, nfullrows_total=nref)
     369           70 :          DO ispin = 1, nspin
     370              :             CALL get_mo_set(my_mos(ispin), &
     371              :                             uniform_occupation=uniform_occupation, &
     372           36 :                             occupation_numbers=occupation_numbers)
     373              :             ! diagonal block matrix
     374           36 :             CALL dbcsr_create(dmat, template=sref(1)%matrix)
     375           36 :             CALL dbcsr_reserve_diag_blocks(dmat)
     376              :             !
     377           36 :             CALL iao_calculate_dmat(ciao(ispin), dmat, occupation_numbers, uniform_occupation)
     378           36 :             CALL dbcsr_trace(dmat, trace)
     379           36 :             IF (unit_nr > 0) THEN
     380           18 :                WRITE (unit_nr, '(T2,A,I2,T66,F15.4)') 'Number of Electrons: Trace(Piao) Spin ', ispin, trace
     381              :             END IF
     382           36 :             CALL iao_charges(dmat, mcharge(:, ispin))
     383          142 :             CALL dbcsr_release(dmat)
     384              :          END DO
     385              :          CALL print_atomic_charges(particle_set, qs_kind_set, unit_nr, "Intrinsic Atomic Orbital Charges", &
     386           34 :                                    electronic_charges=mcharge)
     387           34 :          DEALLOCATE (mcharge)
     388           68 :          CALL cp_fm_release(ciao)
     389              :       END IF
     390              :       !
     391              :       ! Deallocate the neighbor list structure
     392           34 :       CALL release_neighbor_list_sets(srr_list)
     393           34 :       CALL release_neighbor_list_sets(sro_list)
     394           34 :       IF (ASSOCIATED(sref)) CALL dbcsr_deallocate_matrix_set(sref)
     395           34 :       IF (ASSOCIATED(sro)) CALL dbcsr_deallocate_matrix_set(sro)
     396           34 :       CALL cp_fm_release(p_ref_orb)
     397           34 :       CALL cp_fm_release(p_orb_ref)
     398           34 :       CALL cp_fm_release(cvec)
     399           34 :       DEALLOCATE (orb_basis_set_list)
     400              :       !
     401           34 :       IF (iao_env%do_oce) THEN
     402              :          ! generate OCE basis
     403            4 :          IF (unit_nr > 0) THEN
     404            2 :             WRITE (unit_nr, '(T2,A)') "IAO One-Center Expansion: OCE Basis Set"
     405              :          END IF
     406            4 :          CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set)
     407            4 :          nkind = SIZE(qs_kind_set)
     408           40 :          ALLOCATE (oce_basis_set_list(nkind), orb_basis_set_list(nkind))
     409           16 :          DO ikind = 1, nkind
     410           12 :             qs_kind => qs_kind_set(ikind)
     411           12 :             NULLIFY (orb_basis_set_list(ikind)%gto_basis_set)
     412           12 :             NULLIFY (orbbasis)
     413           12 :             CALL get_qs_kind(qs_kind=qs_kind, basis_set=orbbasis, basis_type="ORB")
     414           16 :             IF (ASSOCIATED(orbbasis)) orb_basis_set_list(ikind)%gto_basis_set => orbbasis
     415              :          END DO
     416           16 :          DO ikind = 1, nkind
     417           12 :             orbbasis => orb_basis_set_list(ikind)%gto_basis_set
     418           12 :             CPASSERT(ASSOCIATED(orbbasis))
     419           12 :             NULLIFY (ocebasis)
     420           12 :             CALL create_oce_basis(ocebasis, orbbasis, iao_env%lmax_oce, iao_env%nbas_oce)
     421           12 :             CALL init_orb_basis_set(ocebasis)
     422           12 :             CALL init_interaction_radii_orb_basis(ocebasis, dft_control%qs_control%eps_pgf_orb)
     423           12 :             oce_basis_set_list(ikind)%gto_basis_set => ocebasis
     424           16 :             IF (unit_nr > 0) THEN
     425            6 :                qs_kind => qs_kind_set(ikind)
     426            6 :                CALL get_qs_kind(qs_kind, zeff=zval, element_symbol=element_symbol)
     427            6 :                CALL get_gto_basis_set(ocebasis, name=bname, nsgf=nsgf)
     428            6 :                WRITE (unit_nr, '(T2,A,A,T14,A,I4,T40,A,A30)') "Kind: ", element_symbol, "NBasFun: ", nsgf, &
     429           12 :                   "OCE Basis: ", ADJUSTL(TRIM(bname))
     430              :             END IF
     431              :          END DO
     432            4 :          IF (unit_nr > 0) WRITE (unit_nr, *)
     433              :          ! OCE basis overlap
     434           24 :          ALLOCATE (smat_kind(nkind))
     435           16 :          DO ikind = 1, nkind
     436           12 :             ocebasis => oce_basis_set_list(ikind)%gto_basis_set
     437           12 :             CALL get_gto_basis_set(gto_basis_set=ocebasis, nsgf=nsgf)
     438           52 :             ALLOCATE (smat_kind(ikind)%array(nsgf, nsgf))
     439              :          END DO
     440            4 :          CALL oce_overlap_matrix(smat_kind, oce_basis_set_list)
     441              :          ! overlap between IAO OCE basis and orbital basis
     442            4 :          NULLIFY (sxo, sxo_list)
     443            4 :          CALL setup_neighbor_list(sxo_list, oce_basis_set_list, orb_basis_set_list, qs_env=qs_env)
     444            4 :          CALL get_qs_env(qs_env, ks_env=ks_env)
     445              :          CALL build_overlap_matrix_simple(ks_env, sxo, &
     446            4 :                                           oce_basis_set_list, orb_basis_set_list, sxo_list)
     447              :          !
     448              :          ! One Center Expansion of IAOs
     449            4 :          CALL get_qs_env(qs_env=qs_env, natom=natom)
     450           36 :          ALLOCATE (oce_atom(natom, nspin))
     451            4 :          CALL iao_oce_expansion(qs_env, oce_basis_set_list, smat_kind, sxo, iao_coef, oce_atom, unit_nr)
     452              :          !
     453           16 :          DO ikind = 1, nkind
     454           12 :             ocebasis => oce_basis_set_list(ikind)%gto_basis_set
     455           16 :             CALL deallocate_gto_basis_set(ocebasis)
     456              :          END DO
     457            4 :          DEALLOCATE (oce_basis_set_list, orb_basis_set_list)
     458              :          !
     459            4 :          CALL release_neighbor_list_sets(sxo_list)
     460            4 :          IF (ASSOCIATED(sxo)) CALL dbcsr_deallocate_matrix_set(sxo)
     461           16 :          DO ikind = 1, nkind
     462           16 :             DEALLOCATE (smat_kind(ikind)%array)
     463              :          END DO
     464            4 :          DEALLOCATE (smat_kind)
     465            8 :          DO ispin = 1, nspin
     466           24 :             DO iatom = 1, natom
     467           20 :                DEALLOCATE (oce_atom(iatom, ispin)%array)
     468              :             END DO
     469              :          END DO
     470            4 :          DEALLOCATE (oce_atom)
     471              :       END IF
     472              :       !
     473           34 :       IF (molden_iao) THEN
     474              :          ! Print basis functions: molden file
     475            4 :          CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set)
     476            4 :          CALL get_qs_env(qs_env, cell=cell)
     477            4 :          CALL get_ks_env(ks_env=ks_env, particle_set=particle_set)
     478            4 :          IF (unit_nr > 0) THEN
     479            2 :             WRITE (unit_nr, '(T2,A)') '  Write IAO in MOLDEN Format'
     480              :          END IF
     481           16 :          ALLOCATE (mo_iao(nspin))
     482            8 :          DO ispin = 1, nspin
     483            4 :             CALL get_mo_set(my_mos(ispin), nao=nao)
     484            4 :             CALL allocate_mo_set(mo_iao(ispin), nao, nref, nref, 0.0_dp, 1.0_dp, 0.0_dp)
     485            4 :             CALL init_mo_set(mo_iao(ispin), fm_ref=iao_coef(ispin), name="iao_set")
     486            4 :             CALL cp_fm_to_fm(iao_coef(ispin), mo_iao(ispin)%mo_coeff)
     487            4 :             CALL set_mo_set(mo_iao(ispin), homo=nref)
     488           56 :             mo_iao(ispin)%occupation_numbers = 1.0_dp
     489              :          END DO
     490              :          ! Print IAO basis functions: molden format
     491            4 :          CALL write_mos_molden(mo_iao, qs_kind_set, particle_set, iao_molden_section, cell=cell)
     492            8 :          DO ispin = 1, nspin
     493            8 :             CALL deallocate_mo_set(mo_iao(ispin))
     494              :          END DO
     495            4 :          DEALLOCATE (mo_iao)
     496              :       END IF
     497           34 :       IF (cubes_iao) THEN
     498            4 :          IF (unit_nr > 0) THEN
     499            2 :             WRITE (unit_nr, '(T2,A)') '  Write IAO as CUBE Files'
     500              :          END IF
     501              :          ! Print basis functions: cube file
     502            4 :          CALL print_iao_cubes(qs_env, iao_cubes_section, iao_coef, ref_basis_set_list)
     503              :       END IF
     504              :       !
     505              :       ! Intrinsic bond orbitals
     506           34 :       IF (iao_env%do_bondorbitals) THEN
     507           32 :          IF (unit_nr > 0) THEN
     508           16 :             WRITE (unit_nr, '(/,T2,A)') 'Intrinsic Bond Orbital Generation'
     509              :          END IF
     510              :          ! MOs in IAO basis -> ciao
     511           98 :          ALLOCATE (ciao(nspin))
     512           66 :          DO ispin = 1, nspin
     513           34 :             CALL get_mo_set(my_mos(ispin), mo_coeff=mo_coeff, nao=nao, nmo=norb)
     514           34 :             CALL cp_fm_create(smo, mo_coeff%matrix_struct)
     515           34 :             NULLIFY (fm_struct)
     516           34 :             CALL cp_fm_struct_create(fm_struct, nrow_global=nref, template_fmstruct=mo_coeff%matrix_struct)
     517           34 :             CALL cp_fm_create(ciao(ispin), fm_struct)
     518           34 :             CALL cp_fm_struct_release(fm_struct)
     519              :             ! CIAO = A(T)*S*C
     520           34 :             CALL cp_dbcsr_sm_fm_multiply(smat, mo_coeff, smo, ncol=norb)
     521           34 :             CALL parallel_gemm('T', 'N', nref, norb, nao, 1.0_dp, iao_coef(ispin), smo, 0.0_dp, ciao(ispin))
     522          100 :             CALL cp_fm_release(smo)
     523              :          END DO
     524              :          !
     525              :          ! localize occupied orbitals nocc(ispin), see IAO generation
     526              :          !
     527          164 :          ALLOCATE (cloc(nspin), c_loc_orb(nspin))
     528           66 :          DO ispin = 1, nspin
     529           34 :             NULLIFY (fm_struct)
     530              :             CALL cp_fm_struct_create(fm_struct, ncol_global=nocc(ispin), &
     531           34 :                                      template_fmstruct=ciao(ispin)%matrix_struct)
     532           34 :             CALL cp_fm_create(cloc(ispin), fm_struct)
     533           34 :             CALL cp_fm_struct_release(fm_struct)
     534           66 :             CALL cp_fm_to_fm(ciao(ispin), cloc(ispin), ncol=nocc(ispin))
     535              :          END DO
     536           32 :          CALL cp_fm_release(ciao)
     537           32 :          CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set, natom=natom)
     538          160 :          ALLOCATE (first_sgf(natom), last_sgf(natom), nsgfat(natom))
     539           32 :          CALL get_ks_env(ks_env=ks_env, particle_set=particle_set)
     540              :          CALL get_particle_set(particle_set, qs_kind_set, first_sgf=first_sgf, last_sgf=last_sgf, &
     541           32 :                                nsgf=nsgfat, basis=ref_basis_set_list)
     542              : 
     543           32 :          IF (iao_env%loc_operator /= do_iaoloc_pm2) THEN
     544            0 :             CPABORT("IAO localization operator NYA")
     545              :          END IF
     546           32 :          IF (iao_env%eloc_function /= do_iaoloc_enone) THEN
     547            0 :             CPABORT("IAO energy weight localization NYA")
     548              :          END IF
     549              : 
     550           66 :          DO ispin = 1, nspin
     551           34 :             nvec = nocc(ispin)
     552           66 :             IF (nvec > 0) THEN
     553          326 :                ALLOCATE (zij_atom(natom, 1), fdiag(nvec))
     554           34 :                NULLIFY (fm_struct)
     555              :                CALL cp_fm_struct_create(fm_struct, nrow_global=nvec, ncol_global=nvec, &
     556           34 :                                         template_fmstruct=cloc(ispin)%matrix_struct)
     557          156 :                DO iatom = 1, natom
     558          122 :                   CALL cp_fm_create(zij_atom(iatom, 1), fm_struct)
     559          122 :                   isgf = first_sgf(iatom)
     560              :                   CALL parallel_gemm('T', 'N', nvec, nvec, nsgfat(iatom), 1.0_dp, cloc(ispin), cloc(ispin), &
     561              :                                      0.0_dp, zij_atom(iatom, 1), &
     562          156 :                                      a_first_col=1, a_first_row=isgf, b_first_col=1, b_first_row=isgf)
     563              :                END DO
     564           34 :                CALL cp_fm_struct_release(fm_struct)
     565              :                !
     566           34 :                t1 = m_walltime()
     567           34 :                order = 4
     568           34 :                fin = 0.0_dp
     569          156 :                DO iatom = 1, natom
     570          122 :                   CALL cp_fm_get_diag(zij_atom(iatom, 1), fdiag)
     571          932 :                   fin = fin + SUM(fdiag**order)
     572              :                END DO
     573           34 :                fin = fin**(1._dp/order)
     574              :                ! QR localization
     575           34 :                CALL scdm_qrfact(cloc(ispin))
     576              :                !
     577          156 :                DO iatom = 1, natom
     578          122 :                   isgf = first_sgf(iatom)
     579              :                   CALL parallel_gemm('T', 'N', nvec, nvec, nsgfat(iatom), 1.0_dp, cloc(ispin), cloc(ispin), &
     580              :                                      0.0_dp, zij_atom(iatom, 1), &
     581          156 :                                      a_first_col=1, a_first_row=isgf, b_first_col=1, b_first_row=isgf)
     582              :                END DO
     583           34 :                fout = 0.0_dp
     584          156 :                DO iatom = 1, natom
     585          122 :                   CALL cp_fm_get_diag(zij_atom(iatom, 1), fdiag)
     586          932 :                   fout = fout + SUM(fdiag**order)
     587              :                END DO
     588           34 :                fout = fout**(1._dp/order)
     589           34 :                DEALLOCATE (fdiag)
     590           34 :                t2 = m_walltime()
     591           34 :                IF (unit_nr > 0) THEN
     592              :                   WRITE (unit_nr, '(T2,A,F14.8,A,F14.8,A,F8.2)') &
     593           17 :                      'SCDM pre-localization: fin=', fin, "  fout=", fout, "    Time=", t2 - t1
     594              :                END IF
     595              :                !
     596           34 :                CALL pm_localization(zij_atom, cloc(ispin), order, 1.E-12_dp, unit_nr)
     597              :                !
     598           34 :                CALL cp_fm_release(zij_atom)
     599           34 :                CALL get_mo_set(my_mos(ispin), nao=nao)
     600           34 :                NULLIFY (fm_struct)
     601              :                CALL cp_fm_struct_create(fm_struct, nrow_global=nao, &
     602           34 :                                         template_fmstruct=cloc(ispin)%matrix_struct)
     603           34 :                CALL cp_fm_create(c_loc_orb(ispin), fm_struct)
     604           34 :                CALL cp_fm_struct_release(fm_struct)
     605              :                CALL parallel_gemm('N', 'N', nao, nvec, nref, 1.0_dp, iao_coef(ispin), cloc(ispin), &
     606           68 :                                   0.0_dp, c_loc_orb(ispin))
     607              :             END IF
     608              :          END DO
     609           32 :          DEALLOCATE (first_sgf, last_sgf, nsgfat)
     610           32 :          CALL cp_fm_release(cloc)
     611              :          !
     612           32 :          IF (iao_env%do_center) THEN
     613           32 :             IF (unit_nr > 0) THEN
     614           16 :                WRITE (unit_nr, '(T2,A)') '  Calculate Localized Orbital Centers and Spread'
     615           16 :                IF (iao_env%pos_periodic) THEN
     616           13 :                   WRITE (unit_nr, '(T2,A)') '    Use Berry Phase Position Operator'
     617              :                ELSE
     618            3 :                   WRITE (unit_nr, '(T2,A)') '    Use Local Position Operator'
     619              :                END IF
     620              :             END IF
     621           96 :             nvec = MAXVAL(nocc)
     622              :             ! x  y  z  m2(1)  m2(2)
     623          128 :             ALLOCATE (moments(5, nvec, nspin))
     624         1158 :             moments = 0.0_dp
     625           32 :             IF (iao_env%pos_periodic) THEN
     626           26 :                CALL center_spread_berry(qs_env, c_loc_orb, moments)
     627              :             ELSE
     628            6 :                CALL center_spread_loc(qs_env, c_loc_orb, moments)
     629              :             END IF
     630              :             !
     631           32 :             IF (ASSOCIATED(ibo_cc_section)) THEN
     632            4 :                CALL print_center_spread(moments, nocc, ibo_cc_section)
     633              :             END IF
     634              :             !
     635           32 :             IF (PRESENT(bond_centers)) THEN
     636           28 :                nx = SIZE(bond_centers, 1)
     637           28 :                no = SIZE(bond_centers, 2)
     638           28 :                ns = SIZE(bond_centers, 3)
     639           28 :                CPASSERT(no >= nvec)
     640           28 :                CPASSERT(ns == nspin)
     641           28 :                CPASSERT(nx >= 3)
     642           28 :                nx = MIN(nx, 5)
     643          982 :                bond_centers(1:nx, 1:nvec, 1:ns) = moments(1:nx, 1:nvec, 1:ns)
     644              :             END IF
     645           32 :             DEALLOCATE (moments)
     646              :          END IF
     647              :          !
     648           32 :          IF (molden_ibo) THEN
     649            4 :             CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set)
     650            4 :             CALL get_ks_env(ks_env=ks_env, particle_set=particle_set)
     651            4 :             IF (unit_nr > 0) THEN
     652            2 :                WRITE (unit_nr, '(T2,A)') '  Write IBO in MOLDEN Format'
     653              :             END IF
     654           16 :             ALLOCATE (mo_loc(nspin))
     655            8 :             DO ispin = 1, nspin
     656            4 :                CALL get_mo_set(my_mos(ispin), nao=nao)
     657            4 :                nvec = nocc(ispin)
     658            4 :                CALL allocate_mo_set(mo_loc(ispin), nao, nvec, nvec, 0.0_dp, 1.0_dp, 0.0_dp)
     659            4 :                CALL init_mo_set(mo_loc(ispin), fm_ref=c_loc_orb(ispin), name="ibo_orb")
     660            4 :                CALL cp_fm_to_fm(c_loc_orb(ispin), mo_loc(ispin)%mo_coeff)
     661            4 :                CALL set_mo_set(mo_loc(ispin), homo=nvec)
     662           40 :                mo_loc(ispin)%occupation_numbers = 1.0_dp
     663              :             END DO
     664              :             ! Print IBO functions: molden format
     665            4 :             CALL write_mos_molden(mo_loc, qs_kind_set, particle_set, ibo_molden_section, cell=cell)
     666            8 :             DO ispin = 1, nspin
     667            8 :                CALL deallocate_mo_set(mo_loc(ispin))
     668              :             END DO
     669            4 :             DEALLOCATE (mo_loc)
     670              :          END IF
     671              :          !
     672           32 :          IF (cubes_ibo) THEN
     673            4 :             IF (unit_nr > 0) THEN
     674            2 :                WRITE (unit_nr, '(T2,A)') '  Write IBO on CUBE files'
     675              :             END IF
     676              :             ! Print localized orbital functions: cube file
     677            4 :             CALL print_ibo_cubes(qs_env, ibo_cubes_section, c_loc_orb)
     678              :          END IF
     679              :          !
     680           32 :          IF (PRESENT(mos) .AND. ALLOCATED(c_loc_orb)) THEN
     681              :             ! c_loc_orb -> mos
     682           58 :             DO ispin = 1, nspin
     683           30 :                CALL get_mo_set(mos(ispin), mo_coeff=mo_coeff)
     684           30 :                nvec = nocc(ispin)
     685           58 :                CALL cp_fm_to_fm(c_loc_orb(ispin), mo_coeff, ncol=nvec)
     686              :             END DO
     687              :          END IF
     688           32 :          CALL cp_fm_release(c_loc_orb)
     689              :       END IF
     690           34 :       DEALLOCATE (ref_basis_set_list)
     691              : 
     692           34 :       IF (PRESENT(c_iao_coef)) THEN
     693           30 :          CALL cp_fm_release(c_iao_coef)
     694           92 :          ALLOCATE (c_iao_coef(nspin))
     695           62 :          DO ispin = 1, nspin
     696           32 :             CALL cp_fm_create(c_iao_coef(ispin), iao_coef(ispin)%matrix_struct)
     697           62 :             CALL cp_fm_to_fm(iao_coef(ispin), c_iao_coef(ispin))
     698              :          END DO
     699              :       END IF
     700           34 :       CALL cp_fm_release(iao_coef)
     701              : 
     702           34 :       IF (unit_nr > 0) THEN
     703              :          WRITE (unit_nr, '(T2,A)') &
     704           17 :             '!----------------------------END OF IAO ANALYSIS------------------------------!'
     705              :       END IF
     706              : 
     707           34 :       CALL timestop(handle)
     708              : 
     709          204 :    END SUBROUTINE iao_wfn_analysis
     710              : 
     711              : ! **************************************************************************************************
     712              : !> \brief Computes projector matrices for ref to orb basis and reverse
     713              : !> \param smat ...
     714              : !> \param sref ...
     715              : !> \param s_r_o ...
     716              : !> \param p_o_r ...
     717              : !> \param p_r_o ...
     718              : !> \param eps_svd ...
     719              : ! **************************************************************************************************
     720          238 :    SUBROUTINE iao_projectors(smat, sref, s_r_o, p_o_r, p_r_o, eps_svd)
     721              :       TYPE(dbcsr_type), INTENT(INOUT)                    :: smat, sref, s_r_o
     722              :       TYPE(cp_fm_type), INTENT(INOUT)                    :: p_o_r, p_r_o
     723              :       REAL(KIND=dp), INTENT(IN)                          :: eps_svd
     724              : 
     725              :       CHARACTER(len=*), PARAMETER                        :: routineN = 'iao_projectors'
     726              : 
     727              :       INTEGER                                            :: handle, norb, nref
     728              :       TYPE(cp_fm_struct_type), POINTER                   :: fm_struct
     729              :       TYPE(cp_fm_type)                                   :: fm_inv, fm_s, fm_sro
     730              : 
     731           34 :       CALL timeset(routineN, handle)
     732              : 
     733           34 :       CALL cp_fm_get_info(p_r_o, nrow_global=nref, ncol_global=norb)
     734           34 :       CALL cp_fm_create(fm_sro, p_r_o%matrix_struct)
     735           34 :       CALL copy_dbcsr_to_fm(s_r_o, fm_sro)
     736              : 
     737           34 :       NULLIFY (fm_struct)
     738              :       CALL cp_fm_struct_create(fm_struct, nrow_global=norb, ncol_global=norb, &
     739           34 :                                template_fmstruct=p_r_o%matrix_struct)
     740           34 :       CALL cp_fm_create(fm_s, fm_struct, name="smat")
     741           34 :       CALL cp_fm_create(fm_inv, fm_struct, name="sinv")
     742           34 :       CALL cp_fm_struct_release(fm_struct)
     743           34 :       CALL copy_dbcsr_to_fm(smat, fm_s)
     744           34 :       CALL cp_fm_invert(fm_s, fm_inv, eps_svd=eps_svd)
     745           34 :       CALL parallel_gemm('N', 'T', norb, nref, norb, 1.0_dp, fm_inv, fm_sro, 0.0_dp, p_o_r)
     746           34 :       CALL cp_fm_release(fm_s)
     747           34 :       CALL cp_fm_release(fm_inv)
     748              : 
     749           34 :       NULLIFY (fm_struct)
     750              :       CALL cp_fm_struct_create(fm_struct, nrow_global=nref, ncol_global=nref, &
     751           34 :                                template_fmstruct=p_r_o%matrix_struct)
     752           34 :       CALL cp_fm_create(fm_s, fm_struct, name="sref")
     753           34 :       CALL cp_fm_create(fm_inv, fm_struct, name="sinv")
     754           34 :       CALL cp_fm_struct_release(fm_struct)
     755           34 :       CALL copy_dbcsr_to_fm(sref, fm_s)
     756           34 :       CALL cp_fm_invert(fm_s, fm_inv, eps_svd=eps_svd)
     757           34 :       CALL parallel_gemm('N', 'N', nref, norb, nref, 1.0_dp, fm_inv, fm_sro, 0.0_dp, p_r_o)
     758           34 :       CALL cp_fm_release(fm_s)
     759           34 :       CALL cp_fm_release(fm_inv)
     760              : 
     761           34 :       CALL cp_fm_release(fm_sro)
     762              : 
     763           34 :       CALL timestop(handle)
     764              : 
     765           34 :    END SUBROUTINE iao_projectors
     766              : 
     767              : ! **************************************************************************************************
     768              : !> \brief Computes intrinsic orbitals for a given MO vector set
     769              : !> \param smat ...
     770              : !> \param p_o_r ...
     771              : !> \param p_r_o ...
     772              : !> \param cvec ...
     773              : !> \param avec ...
     774              : ! **************************************************************************************************
     775          324 :    SUBROUTINE intrinsic_ao_calc(smat, p_o_r, p_r_o, cvec, avec)
     776              :       TYPE(dbcsr_type), INTENT(INOUT)                    :: smat
     777              :       TYPE(cp_fm_type), INTENT(INOUT)                    :: p_o_r, p_r_o, cvec, avec
     778              : 
     779              :       CHARACTER(len=*), PARAMETER                        :: routineN = 'intrinsic_ao_calc'
     780              : 
     781              :       INTEGER                                            :: handle, nao, norb, nref
     782              :       TYPE(cp_fm_struct_type), POINTER                   :: fm_struct
     783              :       TYPE(cp_fm_type)                                   :: ctvec, pc, sc, sct, vec1
     784              : 
     785           36 :       CALL timeset(routineN, handle)
     786              : 
     787              :       ! number of orbitals
     788           36 :       CALL cp_fm_get_info(cvec, ncol_global=norb)
     789              :       ! basis set sizes
     790           36 :       CALL cp_fm_get_info(p_r_o, nrow_global=nref, ncol_global=nao)
     791              :       ! temp array
     792           36 :       NULLIFY (fm_struct)
     793              :       CALL cp_fm_struct_create(fm_struct, nrow_global=nref, ncol_global=norb, &
     794           36 :                                template_fmstruct=cvec%matrix_struct)
     795           36 :       CALL cp_fm_create(pc, fm_struct)
     796           36 :       CALL cp_fm_struct_release(fm_struct)
     797              :       ! CT = orth(Por.Pro.C)
     798           36 :       CALL parallel_gemm('N', 'N', nref, norb, nao, 1.0_dp, p_r_o, cvec, 0.0_dp, pc)
     799           36 :       CALL cp_fm_create(ctvec, cvec%matrix_struct)
     800           36 :       CALL parallel_gemm('N', 'N', nao, norb, nref, 1.0_dp, p_o_r, pc, 0.0_dp, ctvec)
     801           36 :       CALL cp_fm_release(pc)
     802           36 :       CALL make_basis_lowdin(ctvec, norb, smat)
     803              :       ! S*C and S*CT
     804           36 :       CALL cp_fm_create(sc, cvec%matrix_struct)
     805           36 :       CALL cp_dbcsr_sm_fm_multiply(smat, cvec, sc, ncol=norb)
     806           36 :       CALL cp_fm_create(sct, cvec%matrix_struct)
     807           36 :       CALL cp_dbcsr_sm_fm_multiply(smat, ctvec, sct, ncol=norb)
     808              :       CALL cp_fm_struct_create(fm_struct, nrow_global=norb, ncol_global=nref, &
     809           36 :                                template_fmstruct=cvec%matrix_struct)
     810           36 :       CALL cp_fm_create(pc, fm_struct)
     811           36 :       CALL cp_fm_struct_release(fm_struct)
     812              :       ! V1 = (CT*SCT(T))Por
     813           36 :       CALL parallel_gemm('T', 'N', norb, nref, nao, 1.0_dp, sct, p_o_r, 0.0_dp, pc)
     814           36 :       NULLIFY (fm_struct)
     815              :       CALL cp_fm_struct_create(fm_struct, nrow_global=nao, ncol_global=nref, &
     816           36 :                                template_fmstruct=cvec%matrix_struct)
     817           36 :       CALL cp_fm_create(vec1, fm_struct)
     818           36 :       CALL cp_fm_struct_release(fm_struct)
     819           36 :       CALL parallel_gemm('N', 'N', nao, nref, norb, 1.0_dp, ctvec, pc, 0.0_dp, vec1)
     820              :       ! A = C*SC(T)*V1
     821           36 :       CALL parallel_gemm('T', 'N', norb, nref, nao, 1.0_dp, sc, vec1, 0.0_dp, pc)
     822           36 :       CALL parallel_gemm('N', 'N', nao, nref, norb, 1.0_dp, cvec, pc, 0.0_dp, avec)
     823              :       ! V1 = Por - V1
     824           36 :       CALL cp_fm_geadd(1.0_dp, 'N', p_o_r, -1.0_dp, vec1)
     825              :       ! A = A + V1
     826           36 :       CALL cp_fm_geadd(1.0_dp, 'N', vec1, 1.0_dp, avec)
     827              :       ! A = A - C*SC(T)*V1
     828           36 :       CALL parallel_gemm('T', 'N', norb, nref, nao, 1.0_dp, sc, vec1, 0.0_dp, pc)
     829           36 :       CALL parallel_gemm('N', 'N', nao, nref, norb, -1.0_dp, cvec, pc, 1.0_dp, avec)
     830              :       ! A = orth(A)
     831           36 :       CALL make_basis_lowdin(avec, nref, smat)
     832              : 
     833              :       ! clean up
     834           36 :       CALL cp_fm_release(pc)
     835           36 :       CALL cp_fm_release(vec1)
     836           36 :       CALL cp_fm_release(sc)
     837           36 :       CALL cp_fm_release(sct)
     838           36 :       CALL cp_fm_release(ctvec)
     839              : 
     840           36 :       CALL timestop(handle)
     841              : 
     842           36 :    END SUBROUTINE intrinsic_ao_calc
     843              : 
     844              : ! **************************************************************************************************
     845              : !> \brief Calculate the density matrix from fm vectors using occupation numbers
     846              : !> \param cvec ...
     847              : !> \param density_matrix ...
     848              : !> \param occupation ...
     849              : !> \param uniform_occupation ...
     850              : ! **************************************************************************************************
     851          132 :    SUBROUTINE iao_calculate_dmat_diag(cvec, density_matrix, occupation, uniform_occupation)
     852              : 
     853              :       TYPE(cp_fm_type), INTENT(IN)                       :: cvec
     854              :       TYPE(dbcsr_type)                                   :: density_matrix
     855              :       REAL(KIND=dp), DIMENSION(:), INTENT(IN)            :: occupation
     856              :       LOGICAL, INTENT(IN)                                :: uniform_occupation
     857              : 
     858              :       CHARACTER(len=*), PARAMETER :: routineN = 'iao_calculate_dmat_diag'
     859              : 
     860              :       INTEGER                                            :: handle, ncol
     861              :       REAL(KIND=dp)                                      :: alpha
     862              :       TYPE(cp_fm_type)                                   :: fm_tmp
     863              : 
     864          132 :       CALL timeset(routineN, handle)
     865              : 
     866          132 :       CALL dbcsr_set(density_matrix, 0.0_dp)
     867              : 
     868          132 :       CALL cp_fm_get_info(cvec, ncol_global=ncol)
     869          132 :       IF (.NOT. uniform_occupation) THEN
     870           98 :          CALL cp_fm_create(fm_tmp, cvec%matrix_struct)
     871           98 :          CALL cp_fm_to_fm(cvec, fm_tmp)
     872           98 :          CALL cp_fm_column_scale(fm_tmp, occupation(1:ncol))
     873           98 :          alpha = 1.0_dp
     874              :          CALL cp_dbcsr_plus_fm_fm_t(sparse_matrix=density_matrix, &
     875              :                                     matrix_v=cvec, matrix_g=fm_tmp, &
     876           98 :                                     ncol=ncol, alpha=alpha)
     877           98 :          CALL cp_fm_release(fm_tmp)
     878              :       ELSE
     879           34 :          alpha = occupation(1)
     880              :          CALL cp_dbcsr_plus_fm_fm_t(sparse_matrix=density_matrix, &
     881           34 :                                     matrix_v=cvec, ncol=ncol, alpha=alpha)
     882              :       END IF
     883              : 
     884          132 :       CALL timestop(handle)
     885              : 
     886          132 :    END SUBROUTINE iao_calculate_dmat_diag
     887              : 
     888              : ! **************************************************************************************************
     889              : !> \brief Calculate the density matrix from fm vectors using an occupation matrix
     890              : !> \param cvec ...
     891              : !> \param density_matrix ...
     892              : !> \param weight ...
     893              : !> \param occmat ...
     894              : ! **************************************************************************************************
     895           16 :    SUBROUTINE iao_calculate_dmat_full(cvec, density_matrix, weight, occmat)
     896              : 
     897              :       TYPE(cp_fm_type), INTENT(IN)                       :: cvec
     898              :       TYPE(dbcsr_type)                                   :: density_matrix
     899              :       REAL(KIND=dp), DIMENSION(:), INTENT(IN)            :: weight
     900              :       TYPE(cp_fm_type), INTENT(IN)                       :: occmat
     901              : 
     902              :       CHARACTER(len=*), PARAMETER :: routineN = 'iao_calculate_dmat_full'
     903              : 
     904              :       INTEGER                                            :: handle, ic, jc, ncol
     905              :       REAL(KIND=dp)                                      :: alpha
     906              :       TYPE(cp_fm_struct_type), POINTER                   :: fm_struct
     907              :       TYPE(cp_fm_type)                                   :: fm1, fm2
     908              : 
     909           16 :       CALL timeset(routineN, handle)
     910              : 
     911           16 :       CALL dbcsr_set(density_matrix, 0.0_dp)
     912              : 
     913              :       CALL cp_fm_struct_create(fm_struct, ncol_global=1, &
     914           16 :                                template_fmstruct=cvec%matrix_struct)
     915           16 :       CALL cp_fm_create(fm1, fm_struct)
     916           16 :       CALL cp_fm_create(fm2, fm_struct)
     917           16 :       CALL cp_fm_struct_release(fm_struct)
     918              : 
     919           16 :       CALL cp_fm_get_info(cvec, ncol_global=ncol)
     920          432 :       DO ic = 1, ncol
     921          416 :          CALL cp_fm_to_fm(cvec, fm1, ncol=1, source_start=ic, target_start=1)
     922        11248 :          DO jc = 1, ncol
     923        10816 :             CALL cp_fm_to_fm(cvec, fm2, ncol=1, source_start=jc, target_start=1)
     924        10816 :             CALL cp_fm_get_element(occmat, ic, jc, alpha)
     925        10816 :             alpha = weight(ic)*alpha
     926              :             CALL cp_dbcsr_plus_fm_fm_t(density_matrix, fm1, fm2, ncol=1, &
     927        11232 :                                        alpha=alpha, symmetry_mode=1)
     928              :          END DO
     929              :       END DO
     930           16 :       CALL cp_fm_release(fm1)
     931           16 :       CALL cp_fm_release(fm2)
     932              : 
     933           16 :       CALL timestop(handle)
     934              : 
     935           16 :    END SUBROUTINE iao_calculate_dmat_full
     936              : 
     937              : ! **************************************************************************************************
     938              : !> \brief compute the atomic charges (orthogonal basis)
     939              : !> \param p_matrix ...
     940              : !> \param charges ...
     941              : ! **************************************************************************************************
     942          208 :    SUBROUTINE iao_charges(p_matrix, charges)
     943              :       TYPE(dbcsr_type)                                   :: p_matrix
     944              :       REAL(KIND=dp), DIMENSION(:)                        :: charges
     945              : 
     946              :       INTEGER                                            :: i, iblock_col, iblock_row
     947              :       REAL(kind=dp)                                      :: trace
     948          208 :       REAL(KIND=dp), DIMENSION(:, :), POINTER            :: p_block
     949              :       TYPE(dbcsr_iterator_type)                          :: iter
     950              :       TYPE(mp_comm_type)                                 :: group
     951              : 
     952         1120 :       charges = 0.0_dp
     953              : 
     954          208 :       CALL dbcsr_iterator_start(iter, p_matrix)
     955          664 :       DO WHILE (dbcsr_iterator_blocks_left(iter))
     956          456 :          NULLIFY (p_block)
     957          456 :          CALL dbcsr_iterator_next_block(iter, iblock_row, iblock_col, p_block)
     958          456 :          CPASSERT(iblock_row == iblock_col)
     959          456 :          trace = 0.0_dp
     960         1758 :          DO i = 1, SIZE(p_block, 1)
     961         1758 :             trace = trace + p_block(i, i)
     962              :          END DO
     963          456 :          charges(iblock_row) = trace
     964              :       END DO
     965          208 :       CALL dbcsr_iterator_stop(iter)
     966              : 
     967          208 :       CALL dbcsr_get_info(p_matrix, group=group)
     968         2032 :       CALL group%sum(charges)
     969              : 
     970          208 :    END SUBROUTINE iao_charges
     971              : 
     972              : ! **************************************************************************************************
     973              : !> \brief Computes and prints the Cube Files for the intrinsic atomic orbitals
     974              : !> \param qs_env ...
     975              : !> \param print_section ...
     976              : !> \param iao_coef ...
     977              : !> \param basis_set_list ...
     978              : ! **************************************************************************************************
     979            4 :    SUBROUTINE print_iao_cubes(qs_env, print_section, iao_coef, basis_set_list)
     980              :       TYPE(qs_environment_type), POINTER                 :: qs_env
     981              :       TYPE(section_vals_type), POINTER                   :: print_section
     982              :       TYPE(cp_fm_type), DIMENSION(:), INTENT(IN)         :: iao_coef
     983              :       TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER  :: basis_set_list
     984              : 
     985              :       CHARACTER(LEN=default_path_length)                 :: filename, title
     986              :       INTEGER                                            :: i, i_rep, ispin, ivec, iw, j, n_rep, &
     987              :                                                             nat, natom, norb, nspins, nstart
     988            4 :       INTEGER, DIMENSION(:), POINTER                     :: atom_list, blk_sizes, first_bas, stride
     989              :       LOGICAL                                            :: mpi_io
     990            4 :       TYPE(atomic_kind_type), DIMENSION(:), POINTER      :: atomic_kind_set
     991              :       TYPE(cell_type), POINTER                           :: cell
     992              :       TYPE(cp_logger_type), POINTER                      :: logger
     993              :       TYPE(dft_control_type), POINTER                    :: dft_control
     994              :       TYPE(particle_list_type), POINTER                  :: particles
     995            4 :       TYPE(particle_type), DIMENSION(:), POINTER         :: particle_set
     996              :       TYPE(pw_c1d_gs_type)                               :: wf_g
     997              :       TYPE(pw_env_type), POINTER                         :: pw_env
     998              :       TYPE(pw_pool_type), POINTER                        :: auxbas_pw_pool
     999              :       TYPE(pw_r3d_rs_type)                               :: wf_r
    1000            4 :       TYPE(qs_kind_type), DIMENSION(:), POINTER          :: qs_kind_set
    1001              :       TYPE(qs_subsys_type), POINTER                      :: subsys
    1002              : 
    1003            8 :       logger => cp_get_default_logger()
    1004            4 :       stride => section_get_ivals(print_section, "STRIDE")
    1005            4 :       CALL section_vals_val_get(print_section, "ATOM_LIST", n_rep_val=n_rep)
    1006              : 
    1007              :       CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set, qs_kind_set=qs_kind_set, &
    1008            4 :                       subsys=subsys, cell=cell, particle_set=particle_set, pw_env=pw_env, dft_control=dft_control)
    1009            4 :       CALL qs_subsys_get(subsys, particles=particles)
    1010              : 
    1011            4 :       CALL get_qs_env(qs_env=qs_env, natom=natom)
    1012           20 :       ALLOCATE (blk_sizes(natom), first_bas(0:natom))
    1013            4 :       CALL get_particle_set(particle_set, qs_kind_set, nsgf=blk_sizes, basis=basis_set_list)
    1014            4 :       first_bas(0) = 0
    1015           20 :       DO i = 1, natom
    1016           20 :          first_bas(i) = first_bas(i - 1) + blk_sizes(i)
    1017              :       END DO
    1018              : 
    1019            4 :       CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
    1020            4 :       CALL auxbas_pw_pool%create_pw(wf_r)
    1021            4 :       CALL auxbas_pw_pool%create_pw(wf_g)
    1022              : 
    1023            4 :       nspins = SIZE(iao_coef)
    1024            4 :       nstart = MIN(1, n_rep)
    1025              : 
    1026            8 :       DO ispin = 1, nspins
    1027           12 :          DO i_rep = nstart, n_rep
    1028            4 :             CALL cp_fm_get_info(iao_coef(ispin), ncol_global=norb)
    1029            4 :             IF (i_rep == 0) THEN
    1030            0 :                nat = natom
    1031              :             ELSE
    1032            4 :                CALL section_vals_val_get(print_section, "ATOM_LIST", i_rep_val=i_rep, i_vals=atom_list)
    1033            4 :                nat = SIZE(atom_list)
    1034              :             END IF
    1035           20 :             DO i = 1, nat
    1036            8 :                IF (i_rep == 0) THEN
    1037            0 :                   j = i
    1038              :                ELSE
    1039            8 :                   j = atom_list(i)
    1040              :                END IF
    1041            8 :                CPASSERT(j >= 1 .AND. j <= natom)
    1042           34 :                DO ivec = first_bas(j - 1) + 1, first_bas(j)
    1043           22 :                   WRITE (filename, '(a4,I5.5,a1,I1.1)') "IAO_", ivec, "_", ispin
    1044           22 :                   WRITE (title, *) "Intrinsic Atomic Orbitals ", ivec, " atom ", j, " spin ", ispin
    1045           22 :                   mpi_io = .TRUE.
    1046              :                   iw = cp_print_key_unit_nr(logger, print_section, "", extension=".cube", &
    1047              :                                             middle_name=TRIM(filename), file_position="REWIND", log_filename=.FALSE., &
    1048           22 :                                             mpi_io=mpi_io)
    1049              :                   CALL calculate_wavefunction(iao_coef(ispin), ivec, wf_r, wf_g, atomic_kind_set, qs_kind_set, &
    1050           22 :                                               cell, dft_control, particle_set, pw_env)
    1051           22 :                   CALL cp_pw_to_cube(wf_r, iw, title, particles=particles, stride=stride, mpi_io=mpi_io)
    1052           30 :                   CALL cp_print_key_finished_output(iw, logger, print_section, "", mpi_io=mpi_io)
    1053              :                END DO
    1054              :             END DO
    1055              :          END DO
    1056              :       END DO
    1057            4 :       DEALLOCATE (blk_sizes, first_bas)
    1058            4 :       CALL auxbas_pw_pool%give_back_pw(wf_r)
    1059            4 :       CALL auxbas_pw_pool%give_back_pw(wf_g)
    1060              : 
    1061            8 :    END SUBROUTINE print_iao_cubes
    1062              : 
    1063              : ! **************************************************************************************************
    1064              : !> \brief Computes and prints the Cube Files for the intrinsic bond orbitals
    1065              : !> \param qs_env ...
    1066              : !> \param print_section ...
    1067              : !> \param ibo_coef ...
    1068              : ! **************************************************************************************************
    1069            4 :    SUBROUTINE print_ibo_cubes(qs_env, print_section, ibo_coef)
    1070              :       TYPE(qs_environment_type), POINTER                 :: qs_env
    1071              :       TYPE(section_vals_type), POINTER                   :: print_section
    1072              :       TYPE(cp_fm_type), DIMENSION(:), INTENT(IN)         :: ibo_coef
    1073              : 
    1074              :       CHARACTER(LEN=default_path_length)                 :: filename, title
    1075              :       INTEGER                                            :: i, i_rep, ispin, iw, j, n_rep, norb, &
    1076              :                                                             nspins, nstart, nstate
    1077            4 :       INTEGER, DIMENSION(:), POINTER                     :: state_list, stride
    1078              :       LOGICAL                                            :: mpi_io
    1079            4 :       TYPE(atomic_kind_type), DIMENSION(:), POINTER      :: atomic_kind_set
    1080              :       TYPE(cell_type), POINTER                           :: cell
    1081              :       TYPE(cp_logger_type), POINTER                      :: logger
    1082              :       TYPE(dft_control_type), POINTER                    :: dft_control
    1083              :       TYPE(particle_list_type), POINTER                  :: particles
    1084            4 :       TYPE(particle_type), DIMENSION(:), POINTER         :: particle_set
    1085              :       TYPE(pw_c1d_gs_type)                               :: wf_g
    1086              :       TYPE(pw_env_type), POINTER                         :: pw_env
    1087              :       TYPE(pw_pool_type), POINTER                        :: auxbas_pw_pool
    1088              :       TYPE(pw_r3d_rs_type)                               :: wf_r
    1089            4 :       TYPE(qs_kind_type), DIMENSION(:), POINTER          :: qs_kind_set
    1090              :       TYPE(qs_subsys_type), POINTER                      :: subsys
    1091              : 
    1092            0 :       CPASSERT(ASSOCIATED(print_section))
    1093              : 
    1094            4 :       logger => cp_get_default_logger()
    1095            4 :       stride => section_get_ivals(print_section, "STRIDE")
    1096            4 :       CALL section_vals_val_get(print_section, "STATE_LIST", n_rep_val=n_rep)
    1097              : 
    1098              :       CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set, qs_kind_set=qs_kind_set, &
    1099            4 :                       subsys=subsys, cell=cell, particle_set=particle_set, pw_env=pw_env, dft_control=dft_control)
    1100            4 :       CALL qs_subsys_get(subsys, particles=particles)
    1101              : 
    1102            4 :       CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
    1103            4 :       CALL auxbas_pw_pool%create_pw(wf_r)
    1104            4 :       CALL auxbas_pw_pool%create_pw(wf_g)
    1105              : 
    1106            4 :       nspins = SIZE(ibo_coef)
    1107            4 :       nstart = MIN(1, n_rep)
    1108              : 
    1109            8 :       DO ispin = 1, nspins
    1110           12 :          DO i_rep = nstart, n_rep
    1111            4 :             CALL cp_fm_get_info(ibo_coef(ispin), ncol_global=norb)
    1112            4 :             IF (i_rep == 0) THEN
    1113            0 :                nstate = norb
    1114              :             ELSE
    1115            4 :                CALL section_vals_val_get(print_section, "STATE_LIST", i_rep_val=i_rep, i_vals=state_list)
    1116            4 :                nstate = SIZE(state_list)
    1117              :             END IF
    1118           16 :             DO i = 1, nstate
    1119            4 :                IF (i_rep == 0) THEN
    1120            0 :                   j = i
    1121              :                ELSE
    1122            4 :                   j = state_list(i)
    1123              :                END IF
    1124            4 :                CPASSERT(j >= 1 .AND. j <= norb)
    1125            4 :                WRITE (filename, '(a4,I5.5,a1,I1.1)') "IBO_", j, "_", ispin
    1126            4 :                WRITE (title, *) "Intrinsic Atomic Orbitals ", j, " spin ", ispin
    1127            4 :                mpi_io = .TRUE.
    1128              :                iw = cp_print_key_unit_nr(logger, print_section, "", extension=".cube", &
    1129              :                                          middle_name=TRIM(filename), file_position="REWIND", log_filename=.FALSE., &
    1130            4 :                                          mpi_io=mpi_io)
    1131              :                CALL calculate_wavefunction(ibo_coef(ispin), j, wf_r, wf_g, atomic_kind_set, qs_kind_set, &
    1132            4 :                                            cell, dft_control, particle_set, pw_env)
    1133            4 :                CALL cp_pw_to_cube(wf_r, iw, title, particles=particles, stride=stride, mpi_io=mpi_io)
    1134            8 :                CALL cp_print_key_finished_output(iw, logger, print_section, "", mpi_io=mpi_io)
    1135              :             END DO
    1136              :          END DO
    1137              :       END DO
    1138              : 
    1139            4 :       CALL auxbas_pw_pool%give_back_pw(wf_r)
    1140            4 :       CALL auxbas_pw_pool%give_back_pw(wf_g)
    1141              : 
    1142            4 :    END SUBROUTINE print_ibo_cubes
    1143              : 
    1144              : ! **************************************************************************************************
    1145              : !> \brief prints charge center and spreads for all orbitals
    1146              : !> \param moments ...
    1147              : !> \param nocc ...
    1148              : !> \param print_section ...
    1149              : !> \param iounit ...
    1150              : ! **************************************************************************************************
    1151           32 :    SUBROUTINE print_center_spread(moments, nocc, print_section, iounit)
    1152              :       REAL(KIND=dp), DIMENSION(:, :, :), INTENT(IN)      :: moments
    1153              :       INTEGER, DIMENSION(:), INTENT(IN)                  :: nocc
    1154              :       TYPE(section_vals_type), OPTIONAL, POINTER         :: print_section
    1155              :       INTEGER, INTENT(IN), OPTIONAL                      :: iounit
    1156              : 
    1157              :       CHARACTER(LEN=default_path_length)                 :: filename
    1158              :       INTEGER                                            :: is, ispin, iw, nspin
    1159              :       TYPE(cp_logger_type), POINTER                      :: logger
    1160              : 
    1161           32 :       logger => cp_get_default_logger()
    1162           32 :       nspin = SIZE(moments, 3)
    1163              : 
    1164           32 :       IF (PRESENT(print_section)) THEN
    1165            4 :          WRITE (filename, '(A18,I1.1)') "IBO_CENTERS_SPREAD"
    1166              :          iw = cp_print_key_unit_nr(logger, print_section, "", extension=".csp", &
    1167            4 :                                    middle_name=TRIM(filename), file_position="REWIND", log_filename=.FALSE.)
    1168           28 :       ELSEIF (PRESENT(iounit)) THEN
    1169           28 :          iw = iounit
    1170              :       ELSE
    1171            0 :          iw = -1
    1172              :       END IF
    1173           32 :       IF (iw > 0) THEN
    1174           33 :          DO ispin = 1, nspin
    1175           17 :             WRITE (iw, "(T2,A,i1)") "Intrinsic Bond Orbitals: Centers and Spread for Spin ", ispin
    1176           17 :             WRITE (iw, "(A7,T30,A6,T68,A7)") "State", "Center", "Spreads"
    1177          133 :             DO is = 1, nocc(ispin)
    1178          117 :                WRITE (iw, "(i7,3F15.8,8X,2F10.5)") is, moments(:, is, ispin)
    1179              :             END DO
    1180              :          END DO
    1181              :       END IF
    1182           32 :       IF (PRESENT(print_section)) THEN
    1183            4 :          CALL cp_print_key_finished_output(iw, logger, print_section, "")
    1184              :       END IF
    1185              : 
    1186           32 :    END SUBROUTINE print_center_spread
    1187              : 
    1188              : ! **************************************************************************************************
    1189              : !> \brief ...
    1190              : !> \param qs_env ...
    1191              : !> \param c_loc_orb ...
    1192              : !> \param moments ...
    1193              : ! **************************************************************************************************
    1194            6 :    SUBROUTINE center_spread_loc(qs_env, c_loc_orb, moments)
    1195              :       TYPE(qs_environment_type), POINTER                 :: qs_env
    1196              :       TYPE(cp_fm_type), DIMENSION(:), INTENT(IN)         :: c_loc_orb
    1197              :       REAL(KIND=dp), DIMENSION(:, :, :), INTENT(INOUT)   :: moments
    1198              : 
    1199              :       CHARACTER(len=*), PARAMETER                        :: routineN = 'center_spread_loc'
    1200              :       INTEGER, DIMENSION(6), PARAMETER                   :: list = [1, 2, 3, 4, 7, 9]
    1201              : 
    1202              :       INTEGER                                            :: handle, i, iop, iorb, ispin, nao, norb, &
    1203              :                                                             nspin
    1204              :       REAL(KIND=dp)                                      :: xmii
    1205              :       REAL(KIND=dp), DIMENSION(3)                        :: rpoint
    1206              :       TYPE(cell_type), POINTER                           :: cell
    1207              :       TYPE(cp_fm_struct_type), POINTER                   :: fm_struct
    1208              :       TYPE(cp_fm_type)                                   :: ccmat, ocvec
    1209            6 :       TYPE(dbcsr_p_type), DIMENSION(:), POINTER          :: dipmat
    1210            6 :       TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER       :: matrix_s
    1211              :       TYPE(dbcsr_type), POINTER                          :: omat, smat
    1212              : 
    1213            6 :       CALL timeset(routineN, handle)
    1214              : 
    1215            6 :       CALL get_qs_env(qs_env=qs_env, matrix_s_kp=matrix_s)
    1216            6 :       smat => matrix_s(1, 1)%matrix
    1217            6 :       rpoint = 0.0_dp
    1218            6 :       nspin = SIZE(c_loc_orb)
    1219          228 :       moments = 0.0_dp
    1220              : 
    1221           60 :       ALLOCATE (dipmat(9))
    1222           60 :       DO i = 1, 9
    1223           54 :          ALLOCATE (dipmat(i)%matrix)
    1224           54 :          CALL dbcsr_copy(dipmat(i)%matrix, smat)
    1225           60 :          CALL dbcsr_set(dipmat(i)%matrix, 0.0_dp)
    1226              :       END DO
    1227              : 
    1228            6 :       CALL build_local_moment_matrix(qs_env, dipmat, 2, rpoint)
    1229              : 
    1230           12 :       DO ispin = 1, nspin
    1231            6 :          CALL cp_fm_get_info(c_loc_orb(ispin), nrow_global=nao, ncol_global=norb)
    1232            6 :          CALL cp_fm_create(ocvec, c_loc_orb(ispin)%matrix_struct)
    1233            6 :          NULLIFY (fm_struct)
    1234              :          CALL cp_fm_struct_create(fm_struct, nrow_global=norb, ncol_global=norb, &
    1235            6 :                                   template_fmstruct=c_loc_orb(ispin)%matrix_struct)
    1236            6 :          CALL cp_fm_create(ccmat, fm_struct)
    1237            6 :          CALL cp_fm_struct_release(fm_struct)
    1238           42 :          DO i = 1, 6
    1239           36 :             iop = list(i)
    1240           36 :             omat => dipmat(iop)%matrix
    1241           36 :             CALL cp_dbcsr_sm_fm_multiply(omat, c_loc_orb(ispin), ocvec, ncol=norb)
    1242              :             CALL parallel_gemm('T', 'N', norb, norb, nao, 1.0_dp, c_loc_orb(ispin), ocvec, &
    1243           36 :                                0.0_dp, ccmat)
    1244          258 :             DO iorb = 1, norb
    1245          216 :                CALL cp_fm_get_element(ccmat, iorb, iorb, xmii)
    1246          252 :                IF (iop <= 3) THEN
    1247          108 :                   moments(iop, iorb, ispin) = moments(iop, iorb, ispin) + xmii
    1248          108 :                   moments(4, iorb, ispin) = moments(4, iorb, ispin) - xmii**2
    1249              :                ELSE
    1250          108 :                   moments(4, iorb, ispin) = moments(4, iorb, ispin) + xmii
    1251              :                END IF
    1252              :             END DO
    1253              :          END DO
    1254            6 :          CALL cp_fm_release(ocvec)
    1255           24 :          CALL cp_fm_release(ccmat)
    1256              :       END DO
    1257              : 
    1258           60 :       DO i = 1, 9
    1259           54 :          CALL dbcsr_release(dipmat(i)%matrix)
    1260           60 :          DEALLOCATE (dipmat(i)%matrix)
    1261              :       END DO
    1262            6 :       DEALLOCATE (dipmat)
    1263              : 
    1264            6 :       CALL get_qs_env(qs_env=qs_env, cell=cell)
    1265           12 :       DO ispin = 1, nspin
    1266           48 :          DO iorb = 1, norb
    1267          150 :             moments(1:3, iorb, ispin) = pbc(moments(1:3, iorb, ispin), cell)
    1268              :          END DO
    1269              :       END DO
    1270              : 
    1271            6 :       CALL timestop(handle)
    1272              : 
    1273            6 :    END SUBROUTINE center_spread_loc
    1274              : 
    1275              : ! **************************************************************************************************
    1276              : !> \brief ...
    1277              : !> \param qs_env ...
    1278              : !> \param c_loc_orb ...
    1279              : !> \param moments ...
    1280              : ! **************************************************************************************************
    1281           26 :    SUBROUTINE center_spread_berry(qs_env, c_loc_orb, moments)
    1282              :       TYPE(qs_environment_type), POINTER                 :: qs_env
    1283              :       TYPE(cp_fm_type), DIMENSION(:), INTENT(IN)         :: c_loc_orb
    1284              :       REAL(KIND=dp), DIMENSION(:, :, :), INTENT(INOUT)   :: moments
    1285              : 
    1286              :       CHARACTER(len=*), PARAMETER :: routineN = 'center_spread_berry'
    1287              : 
    1288              :       COMPLEX(KIND=dp)                                   :: z
    1289              :       INTEGER                                            :: dim_op, handle, i, idir, ispin, istate, &
    1290              :                                                             j, jdir, nao, norb, nspin
    1291              :       REAL(dp), DIMENSION(3)                             :: c, cpbc
    1292              :       REAL(dp), DIMENSION(6)                             :: weights
    1293              :       REAL(KIND=dp)                                      :: imagpart, realpart, spread_i, spread_ii
    1294              :       TYPE(cell_type), POINTER                           :: cell
    1295              :       TYPE(cp_fm_struct_type), POINTER                   :: fm_struct
    1296              :       TYPE(cp_fm_type)                                   :: opvec
    1297           26 :       TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :)     :: zij
    1298           26 :       TYPE(dbcsr_p_type), DIMENSION(:), POINTER          :: matrix_s
    1299           26 :       TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER       :: op_sm_set
    1300              : 
    1301           26 :       CALL timeset(routineN, handle)
    1302              : 
    1303           26 :       CALL get_qs_env(qs_env=qs_env, matrix_s=matrix_s, cell=cell)
    1304              : 
    1305           26 :       IF (cell%orthorhombic) THEN
    1306           26 :          dim_op = 3
    1307              :       ELSE
    1308            0 :          dim_op = 6
    1309              :       END IF
    1310          312 :       ALLOCATE (op_sm_set(2, dim_op))
    1311          104 :       DO i = 1, dim_op
    1312          260 :          DO j = 1, 2
    1313          156 :             NULLIFY (op_sm_set(j, i)%matrix)
    1314          156 :             ALLOCATE (op_sm_set(j, i)%matrix)
    1315          156 :             CALL dbcsr_copy(op_sm_set(j, i)%matrix, matrix_s(1)%matrix)
    1316          234 :             CALL dbcsr_set(op_sm_set(j, i)%matrix, 0.0_dp)
    1317              :          END DO
    1318              :       END DO
    1319           26 :       CALL initialize_weights(cell, weights)
    1320              : 
    1321           26 :       CALL compute_berry_operator(qs_env, cell, op_sm_set, dim_op)
    1322              : 
    1323           26 :       nspin = SIZE(c_loc_orb, 1)
    1324           54 :       DO ispin = 1, nspin
    1325           28 :          CALL cp_fm_get_info(c_loc_orb(ispin), nrow_global=nao, ncol_global=norb)
    1326           28 :          CALL cp_fm_create(opvec, c_loc_orb(ispin)%matrix_struct)
    1327           28 :          NULLIFY (fm_struct)
    1328              :          CALL cp_fm_struct_create(fm_struct, nrow_global=norb, ncol_global=norb, &
    1329           28 :                                   template_fmstruct=c_loc_orb(ispin)%matrix_struct)
    1330          336 :          ALLOCATE (zij(2, dim_op))
    1331              : 
    1332              :          ! Compute zij here
    1333          112 :          DO i = 1, dim_op
    1334          280 :             DO j = 1, 2
    1335          168 :                CALL cp_fm_create(zij(j, i), fm_struct)
    1336          168 :                CALL cp_fm_set_all(zij(j, i), 0.0_dp)
    1337          168 :                CALL cp_dbcsr_sm_fm_multiply(op_sm_set(j, i)%matrix, c_loc_orb(ispin), opvec, ncol=norb)
    1338          252 :                CALL parallel_gemm("T", "N", norb, norb, nao, 1.0_dp, c_loc_orb(ispin), opvec, 0.0_dp, zij(j, i))
    1339              :             END DO
    1340              :          END DO
    1341              : 
    1342           28 :          CALL cp_fm_struct_release(fm_struct)
    1343           28 :          CALL cp_fm_release(opvec)
    1344              : 
    1345          172 :          DO istate = 1, norb
    1346          144 :             c = 0.0_dp
    1347          144 :             spread_i = 0.0_dp
    1348          144 :             spread_ii = 0.0_dp
    1349          576 :             DO jdir = 1, dim_op
    1350          432 :                CALL cp_fm_get_element(zij(1, jdir), istate, istate, realpart)
    1351          432 :                CALL cp_fm_get_element(zij(2, jdir), istate, istate, imagpart)
    1352          432 :                z = CMPLX(realpart, imagpart, dp)
    1353              :                spread_i = spread_i - weights(jdir)* &
    1354          432 :                           LOG(realpart*realpart + imagpart*imagpart)/twopi/twopi
    1355              :                spread_ii = spread_ii + weights(jdir)* &
    1356          432 :                            (1.0_dp - (realpart*realpart + imagpart*imagpart))/twopi/twopi
    1357          576 :                IF (jdir < 4) THEN
    1358         1728 :                   DO idir = 1, 3
    1359              :                      c(idir) = c(idir) + &
    1360         1728 :                                (cell%hmat(idir, jdir)/twopi)*AIMAG(LOG(z))
    1361              :                   END DO
    1362              :                END IF
    1363              :             END DO
    1364          144 :             cpbc = pbc(c, cell)
    1365          576 :             moments(1:3, istate, ispin) = cpbc(1:3)
    1366          144 :             moments(4, istate, ispin) = spread_i
    1367          172 :             moments(5, istate, ispin) = spread_ii
    1368              :          END DO
    1369              : 
    1370           82 :          CALL cp_fm_release(zij)
    1371              : 
    1372              :       END DO
    1373              : 
    1374          104 :       DO i = 1, dim_op
    1375          260 :          DO j = 1, 2
    1376          156 :             CALL dbcsr_release(op_sm_set(j, i)%matrix)
    1377          234 :             DEALLOCATE (op_sm_set(j, i)%matrix)
    1378              :          END DO
    1379              :       END DO
    1380           26 :       DEALLOCATE (op_sm_set)
    1381              : 
    1382           26 :       CALL timestop(handle)
    1383              : 
    1384           52 :    END SUBROUTINE center_spread_berry
    1385              : 
    1386              : ! **************************************************************************************************
    1387              : !> \brief ...
    1388              : !> \param qs_env ...
    1389              : !> \param oce_basis_set_list ...
    1390              : !> \param smat_kind ...
    1391              : !> \param sxo ...
    1392              : !> \param iao_coef ...
    1393              : !> \param oce_atom ...
    1394              : !> \param unit_nr ...
    1395              : ! **************************************************************************************************
    1396            4 :    SUBROUTINE iao_oce_expansion(qs_env, oce_basis_set_list, smat_kind, sxo, iao_coef, oce_atom, unit_nr)
    1397              :       TYPE(qs_environment_type), POINTER                 :: qs_env
    1398              :       TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER  :: oce_basis_set_list
    1399              :       TYPE(cp_2d_r_p_type), DIMENSION(:)                 :: smat_kind
    1400              :       TYPE(dbcsr_p_type), DIMENSION(:)                   :: sxo
    1401              :       TYPE(cp_fm_type), DIMENSION(:)                     :: iao_coef
    1402              :       TYPE(cp_2d_r_p_type), DIMENSION(:, :)              :: oce_atom
    1403              :       INTEGER, INTENT(IN)                                :: unit_nr
    1404              : 
    1405              :       CHARACTER(len=*), PARAMETER                        :: routineN = 'iao_oce_expansion'
    1406              : 
    1407              :       INTEGER                                            :: handle, i, i1, i2, iao, iatom, ikind, &
    1408              :                                                             iset, ishell, ispin, l, m, maxl, n, &
    1409              :                                                             natom, nkind, noce, ns, nset, nsgf, &
    1410              :                                                             nspin
    1411            4 :       INTEGER, DIMENSION(:), POINTER                     :: iao_blk_sizes, nshell, oce_blk_sizes, &
    1412            4 :                                                             orb_blk_sizes
    1413            4 :       INTEGER, DIMENSION(:, :), POINTER                  :: first_sgf, last_sgf, lval
    1414              :       LOGICAL                                            :: found
    1415            4 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:)           :: ev, vector
    1416            4 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :)        :: amat, bmat, filter, oce_comp, prol, vec
    1417            4 :       REAL(KIND=dp), DIMENSION(:, :), POINTER            :: ablock
    1418            4 :       TYPE(cp_2d_r_p_type), ALLOCATABLE, DIMENSION(:)    :: sinv_kind
    1419              :       TYPE(dbcsr_distribution_type)                      :: dbcsr_dist
    1420              :       TYPE(dbcsr_type)                                   :: iao_vec, sx_vec
    1421              :       TYPE(gto_basis_set_type), POINTER                  :: oce_basis
    1422              :       TYPE(mp_comm_type)                                 :: group
    1423            4 :       TYPE(particle_type), DIMENSION(:), POINTER         :: particle_set
    1424            4 :       TYPE(qs_kind_type), DIMENSION(:), POINTER          :: qs_kind_set
    1425              :       TYPE(qs_ks_env_type), POINTER                      :: ks_env
    1426              : 
    1427            4 :       CALL timeset(routineN, handle)
    1428              : 
    1429              :       ! basis sets: block sizes
    1430              :       CALL dbcsr_get_info(sxo(1)%matrix, row_blk_size=oce_blk_sizes, col_blk_size=orb_blk_sizes, &
    1431            4 :                           distribution=dbcsr_dist)
    1432            4 :       CALL get_qs_env(qs_env=qs_env, ks_env=ks_env, natom=natom)
    1433            4 :       CALL get_ks_env(ks_env=ks_env, particle_set=particle_set, qs_kind_set=qs_kind_set)
    1434           12 :       ALLOCATE (iao_blk_sizes(natom))
    1435           20 :       DO iatom = 1, natom
    1436           16 :          CALL get_atomic_kind(particle_set(iatom)%atomic_kind, kind_number=ikind)
    1437           16 :          CALL get_qs_kind(qs_kind_set(ikind), nsgf=ns, basis_type="MIN")
    1438           20 :          iao_blk_sizes(iatom) = ns
    1439              :       END DO
    1440              : 
    1441              :       CALL dbcsr_create(iao_vec, name="IAO_VEC", dist=dbcsr_dist, &
    1442              :                         matrix_type=dbcsr_type_no_symmetry, row_blk_size=orb_blk_sizes, &
    1443            4 :                         col_blk_size=iao_blk_sizes)
    1444              :       CALL dbcsr_create(sx_vec, name="SX_VEC", dist=dbcsr_dist, &
    1445              :                         matrix_type=dbcsr_type_no_symmetry, row_blk_size=oce_blk_sizes, &
    1446            4 :                         col_blk_size=iao_blk_sizes)
    1447            4 :       CALL dbcsr_reserve_diag_blocks(matrix=sx_vec)
    1448              : 
    1449            4 :       nkind = SIZE(smat_kind)
    1450           24 :       ALLOCATE (sinv_kind(nkind))
    1451           16 :       DO ikind = 1, nkind
    1452           12 :          noce = SIZE(smat_kind(ikind)%array, 1)
    1453           48 :          ALLOCATE (sinv_kind(ikind)%array(noce, noce))
    1454       125116 :          sinv_kind(ikind)%array = smat_kind(ikind)%array
    1455           16 :          CALL invmat_symm(sinv_kind(ikind)%array)
    1456              :       END DO
    1457            4 :       CALL dbcsr_get_info(iao_vec, group=group)
    1458              : 
    1459            4 :       nspin = SIZE(iao_coef, 1)
    1460           16 :       ALLOCATE (oce_comp(natom, nspin))
    1461           24 :       oce_comp = 0.0_dp
    1462            8 :       DO ispin = 1, nspin
    1463            4 :          CALL copy_fm_to_dbcsr(iao_coef(ispin), iao_vec)
    1464              :          CALL dbcsr_multiply("N", "N", 1.0_dp, sxo(1)%matrix, iao_vec, 0.0_dp, sx_vec, &
    1465            4 :                              retain_sparsity=.TRUE.)
    1466           20 :          DO iatom = 1, natom
    1467           16 :             CALL get_atomic_kind(particle_set(iatom)%atomic_kind, kind_number=ikind)
    1468              :             CALL dbcsr_get_block_p(matrix=sx_vec, &
    1469           16 :                                    row=iatom, col=iatom, BLOCK=ablock, found=found)
    1470           20 :             IF (found) THEN
    1471            8 :                n = SIZE(ablock, 1)
    1472            8 :                m = SIZE(ablock, 2)
    1473       213906 :                ablock = MATMUL(sinv_kind(ikind)%array, ablock)
    1474           88 :                ALLOCATE (amat(n, m), bmat(m, m), ev(m), vec(m, m))
    1475        25462 :                amat(1:n, 1:m) = MATMUL(smat_kind(ikind)%array(1:n, 1:n), ablock(1:n, 1:m))
    1476         9964 :                bmat(1:m, 1:m) = MATMUL(TRANSPOSE(ablock(1:n, 1:m)), amat(1:n, 1:m))
    1477            8 :                CALL jacobi(bmat, ev, vec)
    1478           30 :                oce_comp(iatom, ispin) = SUM(ev(1:m))/REAL(m, KIND=dp)
    1479           30 :                DO i = 1, m
    1480           22 :                   ev(i) = 1._dp/SQRT(ev(i))
    1481          116 :                   bmat(1:m, i) = vec(1:m, i)*ev(i)
    1482              :                END DO
    1483          714 :                bmat(1:m, 1:m) = MATMUL(bmat(1:m, 1:m), TRANSPOSE(vec(1:m, 1:m)))
    1484        14548 :                ablock(1:n, 1:m) = MATMUL(ablock(1:n, 1:m), bmat(1:m, 1:m))
    1485            8 :                DEALLOCATE (amat, bmat, ev, vec)
    1486              :             END IF
    1487              :          END DO
    1488            4 :          CALL dbcsr_replicate_all(sx_vec)
    1489           24 :          DO iatom = 1, natom
    1490           16 :             CALL get_atomic_kind(particle_set(iatom)%atomic_kind, kind_number=ikind)
    1491              :             CALL dbcsr_get_block_p(matrix=sx_vec, &
    1492           16 :                                    row=iatom, col=iatom, BLOCK=ablock, found=found)
    1493           16 :             CPASSERT(found)
    1494           16 :             n = SIZE(ablock, 1)
    1495           16 :             m = SIZE(ablock, 2)
    1496           64 :             ALLOCATE (oce_atom(iatom, ispin)%array(n, m))
    1497         4728 :             oce_atom(iatom, ispin)%array = ablock
    1498              :          END DO
    1499              :       END DO
    1500              : 
    1501            4 :       CALL group%sum(oce_comp)
    1502            4 :       IF (unit_nr > 0) THEN
    1503           10 :          DO iatom = 1, natom
    1504            8 :             WRITE (unit_nr, "(T4,A,I6,T30,A,2F12.4)") "Atom:", iatom, " Completeness: ", oce_comp(iatom, 1:nspin)
    1505            8 :             CALL get_atomic_kind(particle_set(iatom)%atomic_kind, kind_number=ikind)
    1506            8 :             oce_basis => oce_basis_set_list(ikind)%gto_basis_set
    1507              :             CALL get_gto_basis_set(oce_basis, nset=nset, nshell=nshell, nsgf=nsgf, maxl=maxl, &
    1508            8 :                                    l=lval, first_sgf=first_sgf, last_sgf=last_sgf)
    1509           72 :             ALLOCATE (filter(nsgf, 0:maxl), vector(nsgf), prol(0:maxl, nspin))
    1510         4496 :             filter = 0.0_dp
    1511           70 :             DO iset = 1, nset
    1512          238 :                DO ishell = 1, nshell(iset)
    1513          168 :                   l = lval(ishell, iset)
    1514          168 :                   i1 = first_sgf(ishell, iset)
    1515          168 :                   i2 = last_sgf(ishell, iset)
    1516          970 :                   filter(i1:i2, l) = 1.0_dp
    1517              :                END DO
    1518              :             END DO
    1519              :             !
    1520            8 :             n = SIZE(oce_atom(iatom, 1)%array, 1)
    1521            8 :             m = SIZE(oce_atom(iatom, 1)%array, 2)
    1522            8 :             CPASSERT(n == nsgf)
    1523           30 :             DO iao = 1, m
    1524          176 :                prol = 0.0_dp
    1525           44 :                DO ispin = 1, nspin
    1526          176 :                   DO l = 0, maxl
    1527        14076 :                      vector(1:n) = oce_atom(iatom, ispin)%array(1:n, iao)*filter(1:n, l)
    1528      1611250 :                      prol(l, ispin) = SUM(MATMUL(smat_kind(ikind)%array(1:n, 1:n), vector(1:n))*vector(1:n))
    1529              :                   END DO
    1530              :                END DO
    1531          294 :                WRITE (unit_nr, "(T4,I3,T15,A,T39,(6F7.4))") iao, " l-contributions:", (SUM(prol(l, :)), l=0, maxl)
    1532              :             END DO
    1533           18 :             DEALLOCATE (filter, vector, prol)
    1534              :          END DO
    1535            2 :          WRITE (unit_nr, *)
    1536              :       END IF
    1537              : 
    1538            4 :       DEALLOCATE (oce_comp)
    1539           16 :       DO ikind = 1, nkind
    1540           16 :          DEALLOCATE (sinv_kind(ikind)%array)
    1541              :       END DO
    1542            4 :       DEALLOCATE (sinv_kind)
    1543            4 :       DEALLOCATE (iao_blk_sizes)
    1544            4 :       CALL dbcsr_release(iao_vec)
    1545            4 :       CALL dbcsr_release(sx_vec)
    1546              : 
    1547            4 :       CALL timestop(handle)
    1548              : 
    1549           12 :    END SUBROUTINE iao_oce_expansion
    1550              : 
    1551              : ! **************************************************************************************************
    1552              : !> \brief ...
    1553              : !> \param smat_kind ...
    1554              : !> \param basis_set_list ...
    1555              : ! **************************************************************************************************
    1556            4 :    SUBROUTINE oce_overlap_matrix(smat_kind, basis_set_list)
    1557              :       TYPE(cp_2d_r_p_type), DIMENSION(:)                 :: smat_kind
    1558              :       TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER  :: basis_set_list
    1559              : 
    1560              :       CHARACTER(len=*), PARAMETER :: routineN = 'oce_overlap_matrix'
    1561              : 
    1562              :       INTEGER                                            :: handle, ikind, iset, jset, ldsab, m1, &
    1563              :                                                             m2, n1, n2, ncoa, ncob, nkind, nseta, &
    1564              :                                                             sgfa, sgfb
    1565            4 :       INTEGER, DIMENSION(:), POINTER                     :: la_max, la_min, npgfa, nsgfa
    1566            4 :       INTEGER, DIMENSION(:, :), POINTER                  :: first_sgfa
    1567            4 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :)        :: oint, owork
    1568              :       REAL(KIND=dp), DIMENSION(3)                        :: rab
    1569            4 :       REAL(KIND=dp), DIMENSION(:, :), POINTER            :: rpgfa, scon_a, smat, zeta
    1570              :       TYPE(gto_basis_set_type), POINTER                  :: basis_set
    1571              : 
    1572            4 :       CALL timeset(routineN, handle)
    1573              : 
    1574            4 :       rab(1:3) = 0.0_dp
    1575              : 
    1576            4 :       nkind = SIZE(smat_kind)
    1577            4 :       ldsab = 0
    1578           16 :       DO ikind = 1, nkind
    1579           12 :          basis_set => basis_set_list(ikind)%gto_basis_set
    1580           12 :          CPASSERT(ASSOCIATED(basis_set))
    1581           12 :          CALL get_gto_basis_set(gto_basis_set=basis_set, maxco=m1, nsgf=m2)
    1582           16 :          ldsab = MAX(m1, m2, ldsab)
    1583              :       END DO
    1584              : 
    1585           24 :       ALLOCATE (oint(ldsab, ldsab), owork(ldsab, ldsab))
    1586              : 
    1587           16 :       DO ikind = 1, nkind
    1588           12 :          basis_set => basis_set_list(ikind)%gto_basis_set
    1589           12 :          CPASSERT(ASSOCIATED(basis_set))
    1590           12 :          smat => smat_kind(ikind)%array
    1591       125116 :          smat = 0.0_dp
    1592              :          ! basis ikind
    1593           12 :          first_sgfa => basis_set%first_sgf
    1594           12 :          la_max => basis_set%lmax
    1595           12 :          la_min => basis_set%lmin
    1596           12 :          npgfa => basis_set%npgf
    1597           12 :          nseta = basis_set%nset
    1598           12 :          nsgfa => basis_set%nsgf_set
    1599           12 :          rpgfa => basis_set%pgf_radius
    1600           12 :          scon_a => basis_set%scon
    1601           12 :          zeta => basis_set%zet
    1602              : 
    1603          118 :          DO iset = 1, nseta
    1604              : 
    1605          102 :             ncoa = npgfa(iset)*ncoset(la_max(iset))
    1606          102 :             n1 = npgfa(iset)*(ncoset(la_max(iset)) - ncoset(la_min(iset) - 1))
    1607          102 :             sgfa = first_sgfa(1, iset)
    1608              : 
    1609         1236 :             DO jset = 1, nseta
    1610              : 
    1611         1122 :                ncob = npgfa(jset)*ncoset(la_max(jset))
    1612         1122 :                n2 = npgfa(jset)*(ncoset(la_max(jset)) - ncoset(la_min(jset) - 1))
    1613         1122 :                sgfb = first_sgfa(1, jset)
    1614              : 
    1615              :                ! calculate integrals and derivatives
    1616              :                CALL overlap_ab(la_max(iset), la_min(iset), npgfa(iset), rpgfa(:, iset), zeta(:, iset), &
    1617              :                                la_max(jset), la_min(jset), npgfa(jset), rpgfa(:, jset), zeta(:, jset), &
    1618         1122 :                                rab, sab=oint(:, :))
    1619              :                ! Contraction
    1620              :                CALL contraction(oint(:, :), owork, ca=scon_a(:, sgfa:), na=n1, ma=nsgfa(iset), &
    1621         1122 :                                 cb=scon_a(:, sgfb:), nb=n2, mb=nsgfa(jset), fscale=1.0_dp, trans=.FALSE.)
    1622              :                CALL block_add("IN", owork, nsgfa(iset), nsgfa(jset), smat, &
    1623         1224 :                               sgfa, sgfb, trans=.FALSE.)
    1624              : 
    1625              :             END DO
    1626              :          END DO
    1627              : 
    1628              :       END DO
    1629            4 :       DEALLOCATE (oint, owork)
    1630              : 
    1631            4 :       CALL timestop(handle)
    1632              : 
    1633            4 :    END SUBROUTINE oce_overlap_matrix
    1634              : 
    1635              : ! **************************************************************************************************
    1636              : !> \brief 2 by 2 rotation optimization for the Pipek-Mezey operator
    1637              : !> \param zij_fm_set ...
    1638              : !> \param vectors ...
    1639              : !> \param order ...
    1640              : !> \param accuracy ...
    1641              : !> \param unit_nr ...
    1642              : !> \par History
    1643              : !>       05-2005 created
    1644              : !>       10-2023 adapted from jacobi_rotation_pipek [JGH]
    1645              : !> \author MI
    1646              : ! **************************************************************************************************
    1647           34 :    SUBROUTINE pm_localization(zij_fm_set, vectors, order, accuracy, unit_nr)
    1648              : 
    1649              :       TYPE(cp_fm_type), DIMENSION(:, :), INTENT(IN)      :: zij_fm_set
    1650              :       TYPE(cp_fm_type), INTENT(IN)                       :: vectors
    1651              :       INTEGER, INTENT(IN)                                :: order
    1652              :       REAL(KIND=dp), INTENT(IN)                          :: accuracy
    1653              :       INTEGER, INTENT(IN)                                :: unit_nr
    1654              : 
    1655              :       INTEGER, PARAMETER                                 :: max_sweeps = 250
    1656              : 
    1657              :       INTEGER                                            :: iatom, istate, jstate, natom, nstate, &
    1658              :                                                             sweeps
    1659              :       REAL(KIND=dp)                                      :: aij, bij, ct, ftarget, mii, mij, mjj, &
    1660              :                                                             st, t1, t2, theta, tolerance, tt
    1661              :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:)           :: fdiag
    1662              :       TYPE(cp_fm_type)                                   :: rmat
    1663              : 
    1664           34 :       CALL cp_fm_create(rmat, zij_fm_set(1, 1)%matrix_struct)
    1665           34 :       CALL cp_fm_set_all(rmat, 0.0_dp, 1.0_dp)
    1666              : 
    1667           34 :       CALL cp_fm_get_info(rmat, nrow_global=nstate)
    1668          102 :       ALLOCATE (fdiag(nstate))
    1669           34 :       tolerance = 1.0e10_dp
    1670           34 :       sweeps = 0
    1671           34 :       natom = SIZE(zij_fm_set, 1)
    1672           34 :       IF (unit_nr > 0) THEN
    1673              :          WRITE (unit_nr, '(A,T30,A,T45,A,T63,A,T77,A)') &
    1674           17 :             " Jacobi Localization  ", "Sweep", "Functional", "Gradient", "Time"
    1675              :       END IF
    1676              :       ! do jacobi sweeps until converged
    1677          214 :       DO sweeps = 1, max_sweeps
    1678          214 :          t1 = m_walltime()
    1679          214 :          tolerance = 0.0_dp
    1680         1464 :          DO istate = 1, nstate
    1681         5858 :             DO jstate = istate + 1, nstate
    1682              :                aij = 0.0_dp
    1683              :                bij = 0.0_dp
    1684        31548 :                DO iatom = 1, natom
    1685        27154 :                   CALL cp_fm_get_element(zij_fm_set(iatom, 1), istate, istate, mii)
    1686        27154 :                   CALL cp_fm_get_element(zij_fm_set(iatom, 1), istate, jstate, mij)
    1687        27154 :                   CALL cp_fm_get_element(zij_fm_set(iatom, 1), jstate, jstate, mjj)
    1688        31548 :                   IF (order == 2) THEN
    1689            0 :                      aij = aij + 4._dp*mij**2 - (mii - mjj)**2
    1690            0 :                      bij = bij + 4._dp*mij*(mii - mjj)
    1691        27154 :                   ELSEIF (order == 4) THEN
    1692              :                      aij = aij - mii**4 - mjj**4 + 6._dp*(mii**2 + mjj**2)*mij**2 + &
    1693        27154 :                            mii**3*mjj + mii*mjj**3
    1694        27154 :                      bij = bij + 4._dp*mij*(mii**3 - mjj**3)
    1695              :                   ELSE
    1696            0 :                      CPABORT("PM order")
    1697              :                   END IF
    1698              :                END DO
    1699         4394 :                IF ((aij**2 + bij**2) < 1.E-10_dp) CYCLE
    1700         4240 :                tolerance = tolerance + bij**2
    1701         4240 :                theta = 0.25_dp*ATAN2(bij, -aij)
    1702         4240 :                ct = COS(theta)
    1703         4240 :                st = SIN(theta)
    1704              : 
    1705         4240 :                CALL cp_fm_rot_cols(rmat, istate, jstate, ct, st)
    1706              : 
    1707        32098 :                DO iatom = 1, natom
    1708        26608 :                   CALL cp_fm_rot_cols(zij_fm_set(iatom, 1), istate, jstate, ct, st)
    1709        31002 :                   CALL cp_fm_rot_rows(zij_fm_set(iatom, 1), istate, jstate, ct, st)
    1710              :                END DO
    1711              :             END DO
    1712              :          END DO
    1713          214 :          ftarget = 0.0_dp
    1714         1032 :          DO iatom = 1, natom
    1715          818 :             CALL cp_fm_get_diag(zij_fm_set(iatom, 1), fdiag)
    1716         6966 :             ftarget = ftarget + SUM(fdiag**order)
    1717              :          END DO
    1718          214 :          ftarget = ftarget**(1._dp/order)
    1719          214 :          tolerance = SQRT(tolerance)
    1720          214 :          t2 = m_walltime()
    1721          214 :          tt = t2 - t1
    1722          214 :          IF (unit_nr > 0) THEN
    1723          107 :             WRITE (unit_nr, '(T31,I4,T39,F16.8,T55,G16.8,T71,F10.2)') sweeps, ftarget, tolerance, tt
    1724              :          END IF
    1725          214 :          IF (tolerance < accuracy) EXIT
    1726              :       END DO
    1727           34 :       DEALLOCATE (fdiag)
    1728              : 
    1729           34 :       CALL rotate_orbitals(rmat, vectors)
    1730           34 :       CALL cp_fm_release(rmat)
    1731              : 
    1732           68 :    END SUBROUTINE pm_localization
    1733              : ! **************************************************************************************************
    1734              : 
    1735          140 : END MODULE iao_analysis
        

Generated by: LCOV version 2.0-1