LCOV - code coverage report
Current view: top level - src - qs_active_space_methods.F (source / functions) Coverage Total Hit
Test: CP2K Regtests (git:51fc4cd) Lines: 72.3 % 1602 1158
Test Date: 2026-02-04 06:28:27 Functions: 60.7 % 28 17

            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 Determine active space Hamiltonian
      10              : !> \par History
      11              : !>      04.2016 created [JGH]
      12              : !> \author JGH
      13              : ! **************************************************************************************************
      14              : MODULE qs_active_space_methods
      15              :    USE admm_types, ONLY: admm_type, &
      16              :                          get_admm_env, &
      17              :                          admm_env_release
      18              :    USE atomic_kind_types, ONLY: atomic_kind_type
      19              :    USE basis_set_types, ONLY: allocate_sto_basis_set, &
      20              :                               create_gto_from_sto_basis, &
      21              :                               deallocate_sto_basis_set, &
      22              :                               gto_basis_set_type, &
      23              :                               init_orb_basis_set, &
      24              :                               set_sto_basis_set, &
      25              :                               srules, &
      26              :                               sto_basis_set_type
      27              :    USE cell_types, ONLY: cell_type, use_perd_none, use_perd_xyz
      28              :    USE cell_methods, ONLY: init_cell, set_cell_param, write_cell_low
      29              :    USE cp_blacs_env, ONLY: cp_blacs_env_type, cp_blacs_env_create, cp_blacs_env_release, BLACS_GRID_SQUARE
      30              :    USE cp_control_types, ONLY: dft_control_type, qs_control_type
      31              :    USE cp_dbcsr_operations, ONLY: cp_dbcsr_plus_fm_fm_t, &
      32              :                                   cp_dbcsr_sm_fm_multiply, &
      33              :                                   dbcsr_allocate_matrix_set, &
      34              :                                   cp_dbcsr_m_by_n_from_template, copy_dbcsr_to_fm
      35              :    USE cp_dbcsr_output, ONLY: cp_dbcsr_write_sparse_matrix
      36              :    USE cp_files, ONLY: close_file, &
      37              :                        file_exists, &
      38              :                        open_file
      39              :    USE cp_fm_basic_linalg, ONLY: cp_fm_column_scale
      40              :    USE cp_fm_struct, ONLY: cp_fm_struct_create, &
      41              :                            cp_fm_struct_release, &
      42              :                            cp_fm_struct_type
      43              :    USE cp_fm_types, ONLY: &
      44              :       cp_fm_create, cp_fm_get_element, cp_fm_get_info, cp_fm_init_random, cp_fm_release, &
      45              :       cp_fm_set_all, cp_fm_set_element, cp_fm_to_fm, cp_fm_type, cp_fm_write_formatted
      46              :    USE cp_log_handling, ONLY: cp_get_default_logger, &
      47              :                               cp_logger_get_default_io_unit, &
      48              :                               cp_logger_type
      49              :    USE cp_output_handling, ONLY: &
      50              :       cp_p_file, cp_print_key_finished_output, cp_print_key_should_output, cp_print_key_unit_nr, &
      51              :       debug_print_level, high_print_level, low_print_level, medium_print_level, &
      52              :       silent_print_level
      53              :    USE cp_realspace_grid_cube, ONLY: cp_pw_to_cube
      54              :    USE cp_dbcsr_api, ONLY: &
      55              :       dbcsr_copy, dbcsr_csr_create, dbcsr_csr_type, dbcsr_p_type, dbcsr_type, dbcsr_release, &
      56              :       dbcsr_type_no_symmetry, dbcsr_create, dbcsr_set, dbcsr_multiply, dbcsr_iterator_next_block, &
      57              :       dbcsr_iterator_start, dbcsr_iterator_stop, dbcsr_iterator_blocks_left, &
      58              :       dbcsr_iterator_type, dbcsr_type_symmetric, dbcsr_get_occupation, dbcsr_get_info
      59              :    USE erf_complex, ONLY: erfz_fast
      60              :    USE group_dist_types, ONLY: get_group_dist, release_group_dist, group_dist_d1_type
      61              :    USE input_constants, ONLY: &
      62              :       casci_canonical, eri_method_full_gpw, eri_method_gpw_ht, eri_operator_coulomb, &
      63              :       eri_operator_erf, eri_operator_erfc, eri_operator_gaussian, eri_operator_yukawa, &
      64              :       eri_operator_trunc, eri_operator_lr_trunc, &
      65              :       manual_selection, mao_projection, no_solver, qiskit_solver, wannier_projection, &
      66              :       eri_poisson_analytic, eri_poisson_periodic, eri_poisson_mt
      67              :    USE input_section_types, ONLY: section_vals_get, section_vals_get_subs_vals, &
      68              :                                   section_vals_set_subs_vals, section_vals_type, &
      69              :                                   section_vals_val_get, &
      70              :                                   section_vals_val_set
      71              :    USE ISO_C_BINDING, ONLY: c_null_char
      72              :    USE kinds, ONLY: default_path_length, &
      73              :                     default_string_length, &
      74              :                     dp, &
      75              :                     int_8
      76              :    USE hfx_types, ONLY: hfx_create, hfx_release
      77              :    USE machine, ONLY: m_walltime, m_flush
      78              :    USE mathlib, ONLY: diamat_all
      79              :    USE mathconstants, ONLY: fourpi, twopi, pi, rootpi
      80              :    USE memory_utilities, ONLY: reallocate
      81              :    USE message_passing, ONLY: mp_comm_type, &
      82              :                               mp_para_env_type, &
      83              :                               mp_para_env_release
      84              :    USE mp2_gpw, ONLY: create_mat_munu, grep_rows_in_subgroups, build_dbcsr_from_rows
      85              :    USE mt_util, ONLY: MT0D
      86              :    USE parallel_gemm_api, ONLY: parallel_gemm
      87              :    USE particle_list_types, ONLY: particle_list_type
      88              :    USE particle_types, ONLY: particle_type
      89              :    USE periodic_table, ONLY: ptable
      90              :    USE physcon, ONLY: angstrom, bohr
      91              :    USE preconditioner_types, ONLY: preconditioner_type
      92              :    USE pw_env_methods, ONLY: pw_env_create, &
      93              :                              pw_env_rebuild
      94              :    USE pw_env_types, ONLY: pw_env_get, &
      95              :                            pw_env_release, &
      96              :                            pw_env_type
      97              :    USE pw_methods, ONLY: pw_integrate_function, &
      98              :                          pw_multiply, &
      99              :                          pw_multiply_with, &
     100              :                          pw_transfer, &
     101              :                          pw_zero, pw_integral_ab, pw_scale, &
     102              :                          pw_gauss_damp, pw_compl_gauss_damp
     103              :    USE pw_poisson_methods, ONLY: pw_poisson_rebuild, &
     104              :                                  pw_poisson_solve
     105              :    USE pw_poisson_types, ONLY: ANALYTIC0D, &
     106              :                                PERIODIC3D, &
     107              :                                greens_fn_type, &
     108              :                                pw_poisson_analytic, &
     109              :                                pw_poisson_periodic, &
     110              :                                pw_poisson_type
     111              :    USE pw_pool_types, ONLY: &
     112              :       pw_pool_type
     113              :    USE pw_types, ONLY: &
     114              :       pw_c1d_gs_type, &
     115              :       pw_r3d_rs_type
     116              :    USE qcschema, ONLY: qcschema_env_create, &
     117              :                        qcschema_env_release, &
     118              :                        qcschema_to_hdf5, &
     119              :                        qcschema_type
     120              :    USE qs_active_space_types, ONLY: active_space_type, &
     121              :                                     create_active_space_type, &
     122              :                                     csr_idx_from_combined, &
     123              :                                     csr_idx_to_combined, &
     124              :                                     eri_type, &
     125              :                                     eri_type_eri_element_func, &
     126              :                                     get_irange_csr
     127              :    USE qs_active_space_utils, ONLY: eri_to_array, &
     128              :                                     subspace_matrix_to_array
     129              :    USE qs_collocate_density, ONLY: calculate_wavefunction
     130              :    USE qs_density_matrices, ONLY: calculate_density_matrix
     131              :    USE qs_energy_types, ONLY: qs_energy_type
     132              :    USE qs_environment_types, ONLY: get_qs_env, &
     133              :                                    qs_environment_type, &
     134              :                                    set_qs_env
     135              :    USE qs_integrate_potential, ONLY: integrate_v_rspace
     136              :    USE qs_kind_types, ONLY: qs_kind_type
     137              :    USE qs_ks_methods, ONLY: qs_ks_update_qs_env, qs_ks_build_kohn_sham_matrix, &
     138              :                             evaluate_core_matrix_traces
     139              :    USE qs_ks_types, ONLY: qs_ks_did_change, &
     140              :                           qs_ks_env_type, set_ks_env
     141              :    USE qs_mo_io, ONLY: write_mo_set_to_output_unit
     142              :    USE qs_mo_methods, ONLY: calculate_subspace_eigenvalues
     143              :    USE qs_mo_types, ONLY: allocate_mo_set, &
     144              :                           get_mo_set, &
     145              :                           init_mo_set, &
     146              :                           mo_set_type
     147              :    USE qs_neighbor_list_types, ONLY: neighbor_list_set_p_type, release_neighbor_list_sets
     148              :    USE qs_ot_eigensolver, ONLY: ot_eigensolver
     149              :    USE qs_rho_methods, ONLY: qs_rho_update_rho
     150              :    USE qs_rho_types, ONLY: qs_rho_get, &
     151              :                            qs_rho_type
     152              :    USE qs_subsys_types, ONLY: qs_subsys_get, &
     153              :                               qs_subsys_type
     154              :    USE qs_scf_post_scf, ONLY: qs_scf_compute_properties
     155              :    USE scf_control_types, ONLY: scf_control_type
     156              : #ifndef __NO_SOCKETS
     157              :    USE sockets_interface, ONLY: accept_socket, &
     158              :                                 close_socket, &
     159              :                                 listen_socket, &
     160              :                                 open_bind_socket, &
     161              :                                 readbuffer, &
     162              :                                 remove_socket_file, &
     163              :                                 writebuffer
     164              : #endif
     165              :    USE task_list_methods, ONLY: generate_qs_task_list
     166              :    USE task_list_types, ONLY: allocate_task_list, &
     167              :                               deallocate_task_list, &
     168              :                               task_list_type
     169              :    USE util, ONLY: get_limit
     170              : #include "./base/base_uses.f90"
     171              : 
     172              :    IMPLICIT NONE
     173              :    PRIVATE
     174              : 
     175              :    CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_active_space_methods'
     176              : 
     177              :    PUBLIC :: active_space_main
     178              : 
     179              :    TYPE, EXTENDS(eri_type_eri_element_func) :: eri_fcidump_print
     180              :       INTEGER :: unit_nr = -1, bra_start = -1, ket_start = -1
     181              :    CONTAINS
     182              :       PROCEDURE :: func => eri_fcidump_print_func
     183              :    END TYPE eri_fcidump_print
     184              : 
     185              :    TYPE, EXTENDS(eri_type_eri_element_func) :: eri_fcidump_checksum
     186              :       INTEGER :: bra_start = 0, ket_start = 0
     187              :       REAL(KIND=dp) :: checksum = 0.0_dp
     188              :    CONTAINS
     189              :       PROCEDURE, PASS :: set => eri_fcidump_set
     190              :       PROCEDURE :: func => eri_fcidump_checksum_func
     191              :    END TYPE eri_fcidump_checksum
     192              : 
     193              : CONTAINS
     194              : 
     195              : ! **************************************************************************************************
     196              : !> \brief Sets the starting indices of the bra and ket.
     197              : !> \param this object reference
     198              : !> \param bra_start starting index of the bra
     199              : !> \param ket_start starting index of the ket
     200              : ! **************************************************************************************************
     201           90 :    SUBROUTINE eri_fcidump_set(this, bra_start, ket_start)
     202              :       CLASS(eri_fcidump_checksum) :: this
     203              :       INTEGER, INTENT(IN) :: bra_start, ket_start
     204           90 :       this%bra_start = bra_start
     205           90 :       this%ket_start = ket_start
     206           90 :    END SUBROUTINE eri_fcidump_set
     207              : 
     208              : ! **************************************************************************************************
     209              : !> \brief Main method for determining the active space Hamiltonian
     210              : !> \param qs_env ...
     211              : ! **************************************************************************************************
     212        22103 :    SUBROUTINE active_space_main(qs_env)
     213              :       TYPE(qs_environment_type), POINTER                 :: qs_env
     214              : 
     215              :       CHARACTER(len=*), PARAMETER                        :: routineN = 'active_space_main'
     216              : 
     217              :       CHARACTER(len=10)                                  :: cshell, lnam(5)
     218              :       CHARACTER(len=default_path_length)                 :: qcschema_filename
     219              :       CHARACTER(LEN=default_string_length)               :: basis_type
     220              :       INTEGER :: as_solver, eri_method, eri_operator, eri_print, group_size, handle, i, iatom, &
     221              :          ishell, isp, ispin, iw, j, jm, m, max_orb_ind, mselect, n1, n2, nao, natom, nel, &
     222              :          nelec_active, nelec_inactive, nelec_total, nmo, nmo_active, nmo_available, nmo_inactive, &
     223              :          nmo_inactive_remaining, nmo_occ, nmo_virtual, nn1, nn2, nrow_global, nspins
     224              :       INTEGER, DIMENSION(5)                              :: nshell
     225        22103 :       INTEGER, DIMENSION(:), POINTER                     :: invals
     226              :       LOGICAL                                            :: do_ddapc, do_kpoints, ex_omega, &
     227              :                                                             ex_operator, ex_perd, ex_rcut, &
     228              :                                                             explicit, stop_after_print, store_wfn
     229              :       REAL(KIND=dp) :: alpha, eri_eps_filter, eri_eps_grid, eri_eps_int, eri_gpw_cutoff, &
     230              :          eri_op_omega, eri_rcut, eri_rel_cutoff, fel, focc, maxocc, nze_percentage
     231        22103 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :)        :: eigenvalues
     232        22103 :       REAL(KIND=dp), DIMENSION(:), POINTER               :: evals_virtual
     233              :       TYPE(active_space_type), POINTER                   :: active_space_env
     234        22103 :       TYPE(atomic_kind_type), DIMENSION(:), POINTER      :: atomic_kind_set
     235              :       TYPE(cell_type), POINTER                           :: cell
     236              :       TYPE(cp_blacs_env_type), POINTER                   :: context
     237              :       TYPE(cp_fm_struct_type), POINTER                   :: fm_struct_tmp
     238              :       TYPE(cp_fm_type)                                   :: fm_dummy, mo_virtual
     239              :       TYPE(cp_fm_type), POINTER                          :: fm_target_active, fm_target_inactive, &
     240              :                                                             fmat, mo_coeff, mo_ref, mo_target
     241              :       TYPE(cp_logger_type), POINTER                      :: logger
     242              :       TYPE(dbcsr_csr_type), POINTER                      :: eri_mat
     243        44206 :       TYPE(dbcsr_p_type), DIMENSION(:), POINTER          :: ks_matrix, rho_ao, s_matrix
     244              :       TYPE(dbcsr_type), POINTER                          :: denmat
     245              :       TYPE(dft_control_type), POINTER                    :: dft_control
     246        22103 :       TYPE(mo_set_type), DIMENSION(:), POINTER           :: mos
     247              :       TYPE(mo_set_type), POINTER                         :: mo_set, mo_set_active, mo_set_inactive
     248              :       TYPE(mp_para_env_type), POINTER                    :: para_env
     249        22103 :       TYPE(particle_type), DIMENSION(:), POINTER         :: particle_set
     250              :       TYPE(preconditioner_type), POINTER                 :: local_preconditioner
     251        88412 :       TYPE(qcschema_type)                                :: qcschema_env
     252              :       TYPE(qs_energy_type), POINTER                      :: energy
     253        22103 :       TYPE(qs_kind_type), DIMENSION(:), POINTER          :: qs_kind_set
     254              :       TYPE(qs_ks_env_type), POINTER                      :: ks_env
     255              :       TYPE(qs_rho_type), POINTER                         :: rho
     256              :       TYPE(scf_control_type), POINTER                    :: scf_control
     257              :       TYPE(section_vals_type), POINTER                   :: adiabatic_rescaling, as_input, &
     258              :                                                             hfx_section, input, loc_print, &
     259              :                                                             loc_section, print_orb, xc_section
     260              : 
     261              :       !--------------------------------------------------------------------------------------------!
     262              : 
     263        22103 :       CALL get_qs_env(qs_env, input=input)
     264        22103 :       as_input => section_vals_get_subs_vals(input, "DFT%ACTIVE_SPACE")
     265        22103 :       CALL section_vals_get(as_input, explicit=explicit)
     266        22103 :       IF (.NOT. explicit) RETURN
     267           66 :       CALL timeset(routineN, handle)
     268              : 
     269           66 :       logger => cp_get_default_logger()
     270           66 :       iw = cp_logger_get_default_io_unit(logger)
     271              : 
     272           66 :       IF (iw > 0) THEN
     273              :          WRITE (iw, '(/,T2,A)') &
     274           33 :             '!-----------------------------------------------------------------------------!'
     275           33 :          WRITE (iw, '(T26,A)') "Active Space Embedding Module"
     276              :          WRITE (iw, '(T2,A)') &
     277           33 :             '!-----------------------------------------------------------------------------!'
     278              :       END IF
     279              : 
     280              :       ! k-points?
     281           66 :       CALL get_qs_env(qs_env, do_kpoints=do_kpoints, dft_control=dft_control)
     282           66 :       IF (do_kpoints) THEN
     283            0 :          CALL cp_abort(__LOCATION__, "k-points not supported in active space module")
     284              :       END IF
     285              : 
     286              :       ! adiabatic rescaling?
     287           66 :       adiabatic_rescaling => section_vals_get_subs_vals(input, "DFT%XC%ADIABATIC_RESCALING")
     288           66 :       CALL section_vals_get(adiabatic_rescaling, explicit=explicit)
     289           66 :       IF (explicit) THEN
     290            0 :          CALL cp_abort(__LOCATION__, "Adiabatic rescaling not supported in active space module")
     291              :       END IF
     292              : 
     293              :       ! Setup the possible usage of DDAPC charges
     294              :       do_ddapc = dft_control%qs_control%ddapc_restraint .OR. &
     295              :                  qs_env%cp_ddapc_ewald%do_decoupling .OR. &
     296              :                  qs_env%cp_ddapc_ewald%do_qmmm_periodic_decpl .OR. &
     297           66 :                  qs_env%cp_ddapc_ewald%do_solvation
     298              :       IF (do_ddapc) THEN
     299            0 :          CALL cp_abort(__LOCATION__, "DDAPC charges are not supported in the active space module")
     300              :       END IF
     301           66 :       IF (dft_control%do_sccs) THEN
     302            0 :          CALL cp_abort(__LOCATION__, "SCCS is not supported in the active space module")
     303              :       END IF
     304           66 :       IF (dft_control%correct_surf_dip) THEN
     305            0 :          IF (dft_control%surf_dip_correct_switch) THEN
     306            0 :             CALL cp_abort(__LOCATION__, "Surface dipole correction not supported in the AS module")
     307              :          END IF
     308              :       END IF
     309           66 :       IF (dft_control%smeagol_control%smeagol_enabled) THEN
     310            0 :          CALL cp_abort(__LOCATION__, "SMEAGOL is not supported in the active space module")
     311              :       END IF
     312           66 :       IF (dft_control%qs_control%do_kg) THEN
     313            0 :          CALL cp_abort(__LOCATION__, "KG correction not supported in the active space module")
     314              :       END IF
     315              : 
     316           66 :       NULLIFY (active_space_env)
     317           66 :       CALL create_active_space_type(active_space_env)
     318           66 :       active_space_env%energy_total = 0.0_dp
     319           66 :       active_space_env%energy_ref = 0.0_dp
     320           66 :       active_space_env%energy_inactive = 0.0_dp
     321           66 :       active_space_env%energy_active = 0.0_dp
     322              : 
     323              :       ! input options
     324              : 
     325              :       ! figure out what needs to be printed/stored
     326           66 :       IF (BTEST(cp_print_key_should_output(logger%iter_info, as_input, "FCIDUMP"), cp_p_file)) THEN
     327           66 :          active_space_env%fcidump = .TRUE.
     328              :       END IF
     329              : 
     330           66 :       CALL section_vals_val_get(as_input, "QCSCHEMA", c_val=qcschema_filename, explicit=explicit)
     331           66 :       IF (explicit) THEN
     332            0 :          active_space_env%qcschema = .TRUE.
     333            0 :          active_space_env%qcschema_filename = qcschema_filename
     334              :       END IF
     335              : 
     336           66 :       CALL section_vals_val_get(as_input, "ACTIVE_ELECTRONS", i_val=nelec_active)
     337           66 :       CALL get_qs_env(qs_env, nelectron_total=nelec_total)
     338              : 
     339           66 :       IF (nelec_active <= 0) CPABORT("Specify a positive number of active electrons.")
     340           66 :       IF (nelec_active > nelec_total) CPABORT("More active electrons than total electrons.")
     341              : 
     342           66 :       nelec_inactive = nelec_total - nelec_active
     343           66 :       IF (MOD(nelec_inactive, 2) /= 0) THEN
     344            0 :          CPABORT("The remaining number of inactive electrons has to be even.")
     345              :       END IF
     346              : 
     347           66 :       CALL section_vals_val_get(as_input, "ALPHA", r_val=alpha)
     348           66 :       IF (alpha < 0.0 .OR. alpha > 1.0) CPABORT("Specify a damping factor between 0 and 1.")
     349           66 :       active_space_env%alpha = alpha
     350              : 
     351           66 :       IF (iw > 0) THEN
     352           33 :          WRITE (iw, '(T3,A,T70,I10)') "Total number of electrons", nelec_total
     353           33 :          WRITE (iw, '(T3,A,T70,I10)') "Number of inactive electrons", nelec_inactive
     354           33 :          WRITE (iw, '(T3,A,T70,I10)') "Number of active electrons", nelec_active
     355              :       END IF
     356              : 
     357           66 :       CALL get_qs_env(qs_env, dft_control=dft_control)
     358           66 :       nspins = dft_control%nspins
     359              : 
     360           66 :       active_space_env%nelec_active = nelec_active
     361           66 :       active_space_env%nelec_inactive = nelec_inactive
     362           66 :       active_space_env%nelec_total = nelec_total
     363           66 :       active_space_env%nspins = nspins
     364           66 :       active_space_env%multiplicity = dft_control%multiplicity
     365              : 
     366              :       ! define the active/inactive space orbitals
     367           66 :       CALL section_vals_val_get(as_input, "ACTIVE_ORBITALS", explicit=explicit, i_val=nmo_active)
     368           66 :       IF (.NOT. explicit) THEN
     369            0 :          CALL cp_abort(__LOCATION__, "Number of Active Orbitals has to be specified.")
     370              :       END IF
     371           66 :       active_space_env%nmo_active = nmo_active
     372              :       ! this is safe because nelec_inactive is always even
     373           66 :       nmo_inactive = nelec_inactive/2
     374           66 :       active_space_env%nmo_inactive = nmo_inactive
     375              : 
     376           66 :       CALL section_vals_val_get(as_input, "ORBITAL_SELECTION", i_val=mselect)
     377           66 :       IF (iw > 0) THEN
     378            0 :          SELECT CASE (mselect)
     379              :          CASE DEFAULT
     380            0 :             CPABORT("Unknown orbital selection method")
     381              :          CASE (casci_canonical)
     382              :             WRITE (iw, '(/,T3,A)') &
     383           27 :                "Active space orbitals selected using energy ordered canonical orbitals"
     384              :          CASE (wannier_projection)
     385              :             WRITE (iw, '(/,T3,A)') &
     386            0 :                "Active space orbitals selected using projected Wannier orbitals"
     387              :          CASE (mao_projection)
     388              :             WRITE (iw, '(/,T3,A)') &
     389            0 :                "Active space orbitals selected using modified atomic orbitals (MAO)"
     390              :          CASE (manual_selection)
     391              :             WRITE (iw, '(/,T3,A)') &
     392           33 :                "Active space orbitals selected manually"
     393              :          END SELECT
     394              : 
     395           33 :          WRITE (iw, '(T3,A,T70,I10)') "Number of inactive orbitals", nmo_inactive
     396           33 :          WRITE (iw, '(T3,A,T70,I10)') "Number of active orbitals", nmo_active
     397              :       END IF
     398              : 
     399              :       ! get projection spaces
     400           66 :       CALL section_vals_val_get(as_input, "SUBSPACE_ATOM", i_val=iatom, explicit=explicit)
     401           66 :       IF (explicit) THEN
     402            0 :          CALL get_qs_env(qs_env, natom=natom)
     403            0 :          IF (iatom <= 0 .OR. iatom > natom) THEN
     404            0 :             IF (iw > 0) THEN
     405            0 :                WRITE (iw, '(/,T3,A,I3)') "ERROR: SUBSPACE_ATOM number is not valid", iatom
     406              :             END IF
     407            0 :             CPABORT("Select a valid SUBSPACE_ATOM")
     408              :          END IF
     409              :       END IF
     410           66 :       CALL section_vals_val_get(as_input, "SUBSPACE_SHELL", c_val=cshell, explicit=explicit)
     411           66 :       nshell = 0
     412          396 :       lnam = ""
     413           66 :       IF (explicit) THEN
     414            0 :          cshell = ADJUSTL(cshell)
     415            0 :          n1 = 1
     416            0 :          DO i = 1, 5
     417            0 :             ishell = i
     418            0 :             IF (cshell(n1:n1) == " ") THEN
     419           66 :                ishell = ishell - 1
     420              :                EXIT
     421              :             END IF
     422            0 :             READ (cshell(n1:), "(I1,A1)") nshell(i), lnam(i)
     423            0 :             n1 = n1 + 2
     424              :          END DO
     425              :       END IF
     426              : 
     427              :       ! generate orbitals
     428            0 :       SELECT CASE (mselect)
     429              :       CASE DEFAULT
     430            0 :          CPABORT("Unknown orbital selection method")
     431              :       CASE (casci_canonical)
     432           54 :          CALL get_qs_env(qs_env, mos=mos)
     433              : 
     434              :          ! total number of occupied orbitals, i.e. inactive plus active MOs
     435           54 :          nmo_occ = nmo_inactive + nmo_active
     436              : 
     437              :          ! set inactive orbital indices, these are trivially 1...nmo_inactive
     438          172 :          ALLOCATE (active_space_env%inactive_orbitals(nmo_inactive, nspins))
     439          114 :          DO ispin = 1, nspins
     440          154 :             DO i = 1, nmo_inactive
     441          100 :                active_space_env%inactive_orbitals(i, ispin) = i
     442              :             END DO
     443              :          END DO
     444              : 
     445              :          ! set active orbital indices, these are shifted by nmo_inactive
     446          216 :          ALLOCATE (active_space_env%active_orbitals(nmo_active, nspins))
     447          114 :          DO ispin = 1, nspins
     448          282 :             DO i = 1, nmo_active
     449          228 :                active_space_env%active_orbitals(i, ispin) = nmo_inactive + i
     450              :             END DO
     451              :          END DO
     452              : 
     453              :          ! allocate and initialize inactive and active mo coefficients.
     454              :          ! These are stored in a data structure for the full occupied space:
     455              :          ! for inactive mos, the active subset is set to zero, vice versa for the active mos
     456              :          ! TODO: allocate data structures only for the eaxct number MOs
     457           54 :          maxocc = 2.0_dp
     458           54 :          IF (nspins > 1) maxocc = 1.0_dp
     459          222 :          ALLOCATE (active_space_env%mos_active(nspins))
     460          168 :          ALLOCATE (active_space_env%mos_inactive(nspins))
     461          114 :          DO ispin = 1, nspins
     462           60 :             CALL get_mo_set(mos(ispin), mo_coeff=mo_ref, nao=nao)
     463           60 :             CALL cp_fm_get_info(mo_ref, context=context, para_env=para_env, nrow_global=nrow_global)
     464              :             ! the right number of active electrons per spin channel is initialized further down
     465           60 :             CALL allocate_mo_set(active_space_env%mos_active(ispin), nao, nmo_occ, 0, 0.0_dp, maxocc, 0.0_dp)
     466              :             CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=context, &
     467           60 :                                      nrow_global=nrow_global, ncol_global=nmo_occ)
     468           60 :             CALL init_mo_set(active_space_env%mos_active(ispin), fm_struct=fm_struct_tmp, name="Active Space MO")
     469           60 :             CALL cp_fm_struct_release(fm_struct_tmp)
     470           60 :             IF (nspins == 2) THEN
     471           12 :                nel = nelec_inactive/2
     472              :             ELSE
     473           48 :                nel = nelec_inactive
     474              :             END IF
     475              :             CALL allocate_mo_set(active_space_env%mos_inactive(ispin), nao, nmo_occ, nel, &
     476           60 :                                  REAL(nel, KIND=dp), maxocc, 0.0_dp)
     477              :             CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=context, &
     478           60 :                                      nrow_global=nrow_global, ncol_global=nmo_occ)
     479           60 :             CALL init_mo_set(active_space_env%mos_inactive(ispin), fm_struct=fm_struct_tmp, name="Inactive Space MO")
     480          174 :             CALL cp_fm_struct_release(fm_struct_tmp)
     481              :          END DO
     482              : 
     483              :          ! create canonical orbitals
     484           54 :          IF (dft_control%restricted) THEN
     485            0 :             CPABORT("Unclear how we define MOs in the restricted case ... stopping")
     486              :          ELSE
     487           54 :             IF (dft_control%do_admm) THEN
     488            0 :                IF (dft_control%do_admm_mo) THEN
     489            0 :                   CPABORT("ADMM currently possible only with purification none_dm")
     490              :                END IF
     491              :             END IF
     492              : 
     493          216 :             ALLOCATE (eigenvalues(nmo_occ, nspins))
     494          322 :             eigenvalues = 0.0_dp
     495           54 :             CALL get_qs_env(qs_env, matrix_ks=ks_matrix, matrix_s=s_matrix, scf_control=scf_control)
     496              : 
     497              :             ! calculate virtual MOs and copy inactive and active orbitals
     498           54 :             IF (iw > 0) THEN
     499           27 :                WRITE (iw, '(/,T3,A)') "Calculating virtual MOs..."
     500              :             END IF
     501          114 :             DO ispin = 1, nspins
     502              :                ! nmo_available is the number of MOs available from the SCF calculation:
     503              :                ! this is at least the number of occupied orbitals in the SCF, plus
     504              :                ! any number of added MOs (virtuals) requested in the SCF section
     505           60 :                CALL get_mo_set(mos(ispin), mo_coeff=mo_ref, nmo=nmo_available)
     506              : 
     507              :                ! calculate how many extra MOs we still have to compute
     508           60 :                nmo_virtual = nmo_occ - nmo_available
     509           60 :                nmo_virtual = MAX(nmo_virtual, 0)
     510              : 
     511              :                NULLIFY (evals_virtual)
     512          120 :                ALLOCATE (evals_virtual(nmo_virtual))
     513              : 
     514              :                CALL cp_fm_get_info(mo_ref, context=context, para_env=para_env, &
     515           60 :                                    nrow_global=nrow_global)
     516              : 
     517              :                CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=context, &
     518           60 :                                         nrow_global=nrow_global, ncol_global=nmo_virtual)
     519           60 :                CALL cp_fm_create(mo_virtual, fm_struct_tmp, name="virtual")
     520           60 :                CALL cp_fm_struct_release(fm_struct_tmp)
     521           60 :                CALL cp_fm_init_random(mo_virtual, nmo_virtual)
     522              : 
     523           60 :                NULLIFY (local_preconditioner)
     524              : 
     525              :                ! compute missing virtual MOs
     526              :                CALL ot_eigensolver(matrix_h=ks_matrix(ispin)%matrix, matrix_s=s_matrix(1)%matrix, &
     527              :                                    matrix_c_fm=mo_virtual, matrix_orthogonal_space_fm=mo_ref, &
     528              :                                    eps_gradient=scf_control%eps_lumos, &
     529              :                                    preconditioner=local_preconditioner, &
     530              :                                    iter_max=scf_control%max_iter_lumos, &
     531           60 :                                    size_ortho_space=nmo_available)
     532              : 
     533              :                ! get the eigenvalues
     534           60 :                CALL calculate_subspace_eigenvalues(mo_virtual, ks_matrix(ispin)%matrix, evals_virtual)
     535              : 
     536              :                ! we need to send the copy of MOs to preserve the sign
     537           60 :                CALL cp_fm_create(fm_dummy, mo_ref%matrix_struct)
     538           60 :                CALL cp_fm_to_fm(mo_ref, fm_dummy)
     539              :                CALL calculate_subspace_eigenvalues(fm_dummy, ks_matrix(ispin)%matrix, &
     540           60 :                                                    evals_arg=eigenvalues(:, ispin), do_rotation=.TRUE.)
     541              : 
     542              :                ! copy inactive orbitals
     543           60 :                mo_set => active_space_env%mos_inactive(ispin)
     544           60 :                CALL get_mo_set(mo_set, mo_coeff=mo_target)
     545          100 :                DO i = 1, SIZE(active_space_env%inactive_orbitals, 1)
     546           40 :                   m = active_space_env%inactive_orbitals(i, ispin)
     547           40 :                   CALL cp_fm_to_fm(mo_ref, mo_target, 1, m, m)
     548           40 :                   mo_set%eigenvalues(m) = eigenvalues(m, ispin)
     549          100 :                   IF (nspins > 1) THEN
     550           28 :                      mo_set%occupation_numbers(m) = 1.0
     551              :                   ELSE
     552           12 :                      mo_set%occupation_numbers(m) = 2.0
     553              :                   END IF
     554              :                END DO
     555              : 
     556              :                ! copy active orbitals
     557           60 :                mo_set => active_space_env%mos_active(ispin)
     558           60 :                CALL get_mo_set(mo_set, mo_coeff=mo_target)
     559              :                ! for mult > 1, put the polarized electrons in the alpha channel
     560           60 :                IF (nspins == 2) THEN
     561           12 :                   IF (ispin == 1) THEN
     562            6 :                      nel = (nelec_active + active_space_env%multiplicity - 1)/2
     563              :                   ELSE
     564            6 :                      nel = (nelec_active - active_space_env%multiplicity + 1)/2
     565              :                   END IF
     566              :                ELSE
     567           48 :                   nel = nelec_active
     568              :                END IF
     569           60 :                mo_set%nelectron = nel
     570           60 :                mo_set%n_el_f = REAL(nel, KIND=dp)
     571          228 :                DO i = 1, nmo_active
     572          168 :                   m = active_space_env%active_orbitals(i, ispin)
     573          168 :                   IF (m > nmo_available) THEN
     574            0 :                      CALL cp_fm_to_fm(mo_virtual, mo_target, 1, m - nmo_available, m)
     575            0 :                      eigenvalues(m, ispin) = evals_virtual(m - nmo_available)
     576            0 :                      mo_set%occupation_numbers(m) = 0.0
     577              :                   ELSE
     578          168 :                      CALL cp_fm_to_fm(mo_ref, mo_target, 1, m, m)
     579          168 :                      mo_set%occupation_numbers(m) = mos(ispin)%occupation_numbers(m)
     580              :                   END IF
     581          228 :                   mo_set%eigenvalues(m) = eigenvalues(m, ispin)
     582              :                END DO
     583              :                ! Release
     584           60 :                DEALLOCATE (evals_virtual)
     585           60 :                CALL cp_fm_release(fm_dummy)
     586          354 :                CALL cp_fm_release(mo_virtual)
     587              :             END DO
     588              : 
     589           54 :             IF (iw > 0) THEN
     590           57 :                DO ispin = 1, nspins
     591           30 :                   WRITE (iw, '(/,T3,A,I3,T66,A)') "Canonical Orbital Selection for spin", ispin, &
     592           60 :                      "[atomic units]"
     593           38 :                   DO i = 1, nmo_inactive, 4
     594            8 :                      jm = MIN(3, nmo_inactive - i)
     595           58 :                      WRITE (iw, '(T3,4(F14.6,A5))') (eigenvalues(i + j, ispin), " [I]", j=0, jm)
     596              :                   END DO
     597           61 :                   DO i = nmo_inactive + 1, nmo_inactive + nmo_active, 4
     598           31 :                      jm = MIN(3, nmo_inactive + nmo_active - i)
     599          145 :                      WRITE (iw, '(T3,4(F14.6,A5))') (eigenvalues(i + j, ispin), " [A]", j=0, jm)
     600              :                   END DO
     601           30 :                   WRITE (iw, '(/,T3,A,I3)') "Active Orbital Indices for spin", ispin
     602           88 :                   DO i = 1, SIZE(active_space_env%active_orbitals, 1), 4
     603           31 :                      jm = MIN(3, SIZE(active_space_env%active_orbitals, 1) - i)
     604          145 :                      WRITE (iw, '(T3,4(I4))') (active_space_env%active_orbitals(i + j, ispin), j=0, jm)
     605              :                   END DO
     606              :                END DO
     607              :             END IF
     608           54 :             DEALLOCATE (eigenvalues)
     609              :          END IF
     610              : 
     611              :       CASE (manual_selection)
     612              :          ! create canonical orbitals
     613           12 :          IF (dft_control%restricted) THEN
     614            0 :             CPABORT("Unclear how we define MOs in the restricted case ... stopping")
     615              :          ELSE
     616           12 :             IF (dft_control%do_admm) THEN
     617              :                ! For admm_mo, the auxiliary density is computed from the MOs, which never change
     618              :                ! in the rs-dft embedding, therefore the energy is wrong as the LR HFX never changes.
     619              :                ! For admm_dm, the auxiliary density is computed from the density matrix, which is
     620              :                ! updated at each iteration and therefore works.
     621            0 :                IF (dft_control%do_admm_mo) THEN
     622            0 :                   CPABORT("ADMM currently possible only with purification none_dm")
     623              :                END IF
     624              :             END IF
     625              : 
     626           12 :             CALL section_vals_val_get(as_input, "ACTIVE_ORBITAL_INDICES", explicit=explicit, i_vals=invals)
     627           12 :             IF (.NOT. explicit) THEN
     628              :                CALL cp_abort(__LOCATION__, "Manual orbital selection requires to explicitly "// &
     629            0 :                              "set the active orbital indices via ACTIVE_ORBITAL_INDICES")
     630              :             END IF
     631              : 
     632           12 :             IF (nspins == 1) THEN
     633            6 :                CPASSERT(SIZE(invals) == nmo_active)
     634              :             ELSE
     635            6 :                CPASSERT(SIZE(invals) == 2*nmo_active)
     636              :             END IF
     637           36 :             ALLOCATE (active_space_env%inactive_orbitals(nmo_inactive, nspins))
     638           48 :             ALLOCATE (active_space_env%active_orbitals(nmo_active, nspins))
     639              : 
     640           30 :             DO ispin = 1, nspins
     641           66 :                DO i = 1, nmo_active
     642           54 :                   active_space_env%active_orbitals(i, ispin) = invals(i + (ispin - 1)*nmo_active)
     643              :                END DO
     644              :             END DO
     645              : 
     646           12 :             CALL get_qs_env(qs_env, mos=mos)
     647              : 
     648              :             ! include MOs up to the largest index in the list
     649           48 :             max_orb_ind = MAXVAL(invals)
     650           12 :             maxocc = 2.0_dp
     651           12 :             IF (nspins > 1) maxocc = 1.0_dp
     652           54 :             ALLOCATE (active_space_env%mos_active(nspins))
     653           42 :             ALLOCATE (active_space_env%mos_inactive(nspins))
     654           30 :             DO ispin = 1, nspins
     655              :                ! init active orbitals
     656           18 :                CALL get_mo_set(mos(ispin), mo_coeff=mo_ref, nao=nao)
     657           18 :                CALL cp_fm_get_info(mo_ref, context=context, para_env=para_env, nrow_global=nrow_global)
     658           18 :                CALL allocate_mo_set(active_space_env%mos_active(ispin), nao, max_orb_ind, 0, 0.0_dp, maxocc, 0.0_dp)
     659              :                CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=context, &
     660           18 :                                         nrow_global=nrow_global, ncol_global=max_orb_ind)
     661           18 :                CALL init_mo_set(active_space_env%mos_active(ispin), fm_struct=fm_struct_tmp, name="Active Space MO")
     662           18 :                CALL cp_fm_struct_release(fm_struct_tmp)
     663              : 
     664              :                ! init inactive orbitals
     665           18 :                IF (nspins == 2) THEN
     666           12 :                   nel = nelec_inactive/2
     667              :                ELSE
     668            6 :                   nel = nelec_inactive
     669              :                END IF
     670           18 :                CALL allocate_mo_set(active_space_env%mos_inactive(ispin), nao, max_orb_ind, nel, REAL(nel, KIND=dp), maxocc, 0.0_dp)
     671              :                CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=context, &
     672           18 :                                         nrow_global=nrow_global, ncol_global=max_orb_ind)
     673           18 :                CALL init_mo_set(active_space_env%mos_inactive(ispin), fm_struct=fm_struct_tmp, name="Inactive Space MO")
     674              :                ! small hack: set the correct inactive occupations down below
     675           68 :                active_space_env%mos_inactive(ispin)%occupation_numbers = 0.0_dp
     676           48 :                CALL cp_fm_struct_release(fm_struct_tmp)
     677              :             END DO
     678              : 
     679           48 :             ALLOCATE (eigenvalues(max_orb_ind, nspins))
     680           80 :             eigenvalues = 0.0_dp
     681           12 :             CALL get_qs_env(qs_env, matrix_ks=ks_matrix, matrix_s=s_matrix, scf_control=scf_control)
     682              : 
     683              :             ! calculate virtual MOs and copy inactive and active orbitals
     684           12 :             IF (iw > 0) THEN
     685            6 :                WRITE (iw, '(/,T3,A)') "Calculating virtual MOs..."
     686              :             END IF
     687           30 :             DO ispin = 1, nspins
     688           18 :                CALL get_mo_set(mos(ispin), mo_coeff=mo_ref, nmo=nmo_available)
     689           18 :                nmo_virtual = max_orb_ind - nmo_available
     690           18 :                nmo_virtual = MAX(nmo_virtual, 0)
     691              : 
     692              :                NULLIFY (evals_virtual)
     693           36 :                ALLOCATE (evals_virtual(nmo_virtual))
     694              : 
     695              :                CALL cp_fm_get_info(mo_ref, context=context, para_env=para_env, &
     696           18 :                                    nrow_global=nrow_global)
     697              : 
     698              :                CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=context, &
     699           18 :                                         nrow_global=nrow_global, ncol_global=nmo_virtual)
     700           18 :                CALL cp_fm_create(mo_virtual, fm_struct_tmp, name="virtual")
     701           18 :                CALL cp_fm_struct_release(fm_struct_tmp)
     702           18 :                CALL cp_fm_init_random(mo_virtual, nmo_virtual)
     703              : 
     704           18 :                NULLIFY (local_preconditioner)
     705              : 
     706              :                CALL ot_eigensolver(matrix_h=ks_matrix(ispin)%matrix, matrix_s=s_matrix(1)%matrix, &
     707              :                                    matrix_c_fm=mo_virtual, matrix_orthogonal_space_fm=mo_ref, &
     708              :                                    eps_gradient=scf_control%eps_lumos, &
     709              :                                    preconditioner=local_preconditioner, &
     710              :                                    iter_max=scf_control%max_iter_lumos, &
     711           18 :                                    size_ortho_space=nmo_available)
     712              : 
     713              :                CALL calculate_subspace_eigenvalues(mo_virtual, ks_matrix(ispin)%matrix, &
     714           18 :                                                    evals_virtual)
     715              : 
     716              :                ! We need to send the copy of MOs to preserve the sign
     717           18 :                CALL cp_fm_create(fm_dummy, mo_ref%matrix_struct)
     718           18 :                CALL cp_fm_to_fm(mo_ref, fm_dummy)
     719              : 
     720              :                CALL calculate_subspace_eigenvalues(fm_dummy, ks_matrix(ispin)%matrix, &
     721           18 :                                                    evals_arg=eigenvalues(:, ispin), do_rotation=.TRUE.)
     722              : 
     723           18 :                mo_set_active => active_space_env%mos_active(ispin)
     724           18 :                CALL get_mo_set(mo_set_active, mo_coeff=fm_target_active)
     725           18 :                mo_set_inactive => active_space_env%mos_inactive(ispin)
     726           18 :                CALL get_mo_set(mo_set_inactive, mo_coeff=fm_target_inactive)
     727              : 
     728              :                ! copy orbitals
     729           18 :                nmo_inactive_remaining = nmo_inactive
     730           68 :                DO i = 1, max_orb_ind
     731              :                   ! case for i being an active orbital
     732          114 :                   IF (ANY(active_space_env%active_orbitals(:, ispin) == i)) THEN
     733           36 :                      IF (i > nmo_available) THEN
     734            0 :                         CALL cp_fm_to_fm(mo_virtual, fm_target_active, 1, i - nmo_available, i)
     735            0 :                         eigenvalues(i, ispin) = evals_virtual(i - nmo_available)
     736            0 :                         mo_set_active%occupation_numbers(i) = 0.0
     737              :                      ELSE
     738           36 :                         CALL cp_fm_to_fm(fm_dummy, fm_target_active, 1, i, i)
     739           36 :                         mo_set_active%occupation_numbers(i) = mos(ispin)%occupation_numbers(i)
     740              :                      END IF
     741           36 :                      mo_set_active%eigenvalues(i) = eigenvalues(i, ispin)
     742              :                      ! if it was not an active orbital, check whether it is an inactive orbital
     743           14 :                   ELSEIF (nmo_inactive_remaining > 0) THEN
     744            0 :                      CALL cp_fm_to_fm(fm_dummy, fm_target_inactive, 1, i, i)
     745              :                      ! store on the fly the mapping of inactive orbitals
     746            0 :                      active_space_env%inactive_orbitals(nmo_inactive - nmo_inactive_remaining + 1, ispin) = i
     747            0 :                      mo_set_inactive%eigenvalues(i) = eigenvalues(i, ispin)
     748            0 :                      mo_set_inactive%occupation_numbers(i) = mos(ispin)%occupation_numbers(i)
     749              :                      ! hack: set homo and lumo manually
     750            0 :                      IF (nmo_inactive_remaining == 1) THEN
     751            0 :                         mo_set_inactive%homo = i
     752            0 :                         mo_set_inactive%lfomo = i + 1
     753              :                      END IF
     754            0 :                      nmo_inactive_remaining = nmo_inactive_remaining - 1
     755              :                   ELSE
     756           14 :                      CYCLE
     757              :                   END IF
     758              :                END DO
     759              : 
     760              :                ! Release
     761           18 :                DEALLOCATE (evals_virtual)
     762           18 :                CALL cp_fm_release(fm_dummy)
     763          102 :                CALL cp_fm_release(mo_virtual)
     764              :             END DO
     765              : 
     766           12 :             IF (iw > 0) THEN
     767           15 :                DO ispin = 1, nspins
     768            9 :                   WRITE (iw, '(/,T3,A,I3,T66,A)') "Orbital Energies and Selection for spin", ispin, "[atomic units]"
     769              : 
     770           18 :                   DO i = 1, max_orb_ind, 4
     771            9 :                      jm = MIN(3, max_orb_ind - i)
     772            9 :                      WRITE (iw, '(T4)', advance="no")
     773           34 :                      DO j = 0, jm
     774           57 :                         IF (ANY(active_space_env%active_orbitals(:, ispin) == i + j)) THEN
     775           18 :                            WRITE (iw, '(T3,F12.6,A5)', advance="no") eigenvalues(i + j, ispin), " [A]"
     776            7 :                         ELSEIF (ANY(active_space_env%inactive_orbitals(:, ispin) == i + j)) THEN
     777            0 :                            WRITE (iw, '(T3,F12.6,A5)', advance="no") eigenvalues(i + j, ispin), " [I]"
     778              :                         ELSE
     779            7 :                            WRITE (iw, '(T3,F12.6,A5)', advance="no") eigenvalues(i + j, ispin), " [V]"
     780              :                         END IF
     781              :                      END DO
     782           18 :                      WRITE (iw, *)
     783              :                   END DO
     784            9 :                   WRITE (iw, '(/,T3,A,I3)') "Active Orbital Indices for spin", ispin
     785           24 :                   DO i = 1, SIZE(active_space_env%active_orbitals, 1), 4
     786            9 :                      jm = MIN(3, SIZE(active_space_env%active_orbitals, 1) - i)
     787           36 :                      WRITE (iw, '(T3,4(I4))') (active_space_env%active_orbitals(i + j, ispin), j=0, jm)
     788              :                   END DO
     789              :                END DO
     790              :             END IF
     791           24 :             DEALLOCATE (eigenvalues)
     792              :          END IF
     793              : 
     794              :       CASE (wannier_projection)
     795            0 :          NULLIFY (loc_section, loc_print)
     796            0 :          loc_section => section_vals_get_subs_vals(as_input, "LOCALIZE")
     797            0 :          CPASSERT(ASSOCIATED(loc_section))
     798            0 :          loc_print => section_vals_get_subs_vals(as_input, "LOCALIZE%PRINT")
     799              :          !
     800            0 :          CPABORT("not yet available")
     801              :          !
     802              :       CASE (mao_projection)
     803              :          !
     804           66 :          CPABORT("not yet available")
     805              :          !
     806              :       END SELECT
     807              : 
     808              :       ! Print orbitals on Cube files
     809           66 :       print_orb => section_vals_get_subs_vals(as_input, "PRINT_ORBITAL_CUBES")
     810           66 :       CALL section_vals_get(print_orb, explicit=explicit)
     811           66 :       CALL section_vals_val_get(print_orb, "STOP_AFTER_CUBES", l_val=stop_after_print)
     812           66 :       IF (explicit) THEN
     813              :          !
     814            2 :          CALL print_orbital_cubes(print_orb, qs_env, active_space_env%mos_active)
     815              :          !
     816            2 :          IF (stop_after_print) THEN
     817              : 
     818            0 :             IF (iw > 0) THEN
     819              :                WRITE (iw, '(/,T2,A)') &
     820            0 :                   '!----------------- Early End of Active Space Interface -----------------------!'
     821              :             END IF
     822              : 
     823            0 :             CALL timestop(handle)
     824              : 
     825            0 :             RETURN
     826              :          END IF
     827              :       END IF
     828              : 
     829              :       ! calculate inactive density matrix
     830           66 :       CALL get_qs_env(qs_env, rho=rho)
     831           66 :       CALL qs_rho_get(rho, rho_ao=rho_ao)
     832           66 :       CPASSERT(ASSOCIATED(rho_ao))
     833           66 :       CALL dbcsr_allocate_matrix_set(active_space_env%pmat_inactive, nspins)
     834          144 :       DO ispin = 1, nspins
     835           78 :          ALLOCATE (denmat)
     836           78 :          CALL dbcsr_copy(denmat, rho_ao(ispin)%matrix)
     837           78 :          mo_set => active_space_env%mos_inactive(ispin)
     838           78 :          CALL calculate_density_matrix(mo_set, denmat)
     839          144 :          active_space_env%pmat_inactive(ispin)%matrix => denmat
     840              :       END DO
     841              : 
     842              :       ! read in ERI parameters
     843           66 :       CALL section_vals_val_get(as_input, "ERI%METHOD", i_val=eri_method)
     844           66 :       active_space_env%eri%method = eri_method
     845           66 :       CALL section_vals_val_get(as_input, "ERI%OPERATOR", i_val=eri_operator, explicit=ex_operator)
     846           66 :       active_space_env%eri%operator = eri_operator
     847           66 :       CALL section_vals_val_get(as_input, "ERI%OMEGA", r_val=eri_op_omega, explicit=ex_omega)
     848           66 :       active_space_env%eri%omega = eri_op_omega
     849           66 :       CALL section_vals_val_get(as_input, "ERI%CUTOFF_RADIUS", r_val=eri_rcut, explicit=ex_rcut)
     850           66 :       active_space_env%eri%cutoff_radius = eri_rcut  ! this is already converted to bohr!
     851           66 :       CALL section_vals_val_get(as_input, "ERI%PERIODICITY", i_vals=invals, explicit=ex_perd)
     852           66 :       CALL section_vals_val_get(as_input, "ERI%EPS_INTEGRAL", r_val=eri_eps_int)
     853           66 :       active_space_env%eri%eps_integral = eri_eps_int
     854              :       ! if eri periodicity is explicitly set, we use it, otherwise we use the cell periodicity
     855           66 :       IF (ex_perd) THEN
     856           62 :          IF (SIZE(invals) == 1) THEN
     857            0 :             active_space_env%eri%periodicity(1:3) = invals(1)
     858              :          ELSE
     859          434 :             active_space_env%eri%periodicity(1:3) = invals(1:3)
     860              :          END IF
     861              :       ELSE
     862            4 :          CALL get_qs_env(qs_env, cell=cell)
     863           28 :          active_space_env%eri%periodicity(1:3) = cell%perd(1:3)
     864              :       END IF
     865           66 :       IF (iw > 0) THEN
     866           33 :          WRITE (iw, '(/,T3,A)') "Calculation of Electron Repulsion Integrals"
     867              : 
     868           27 :          SELECT CASE (eri_method)
     869              :          CASE (eri_method_full_gpw)
     870           27 :             WRITE (iw, '(T3,A,T50,A)') "Integration method", "GPW Fourier transform over MOs"
     871              :          CASE (eri_method_gpw_ht)
     872            6 :             WRITE (iw, '(T3,A,T44,A)') "Integration method", "Half transformed integrals from GPW"
     873              :          CASE DEFAULT
     874           33 :             CPABORT("Unknown ERI method")
     875              :          END SELECT
     876              : 
     877           24 :          SELECT CASE (eri_operator)
     878              :          CASE (eri_operator_coulomb)
     879           24 :             WRITE (iw, '(T3,A,T73,A)') "ERI operator", "Coulomb"
     880              : 
     881              :          CASE (eri_operator_yukawa)
     882            0 :             WRITE (iw, '(T3,A,T74,A)') "ERI operator", "Yukawa"
     883            0 :             IF (.NOT. ex_omega) CALL cp_abort(__LOCATION__, &
     884            0 :                                               "Yukawa operator requires OMEGA to be explicitly set")
     885            0 :             WRITE (iw, '(T3,A,T66,F14.3)') "ERI operator parameter OMEGA", eri_op_omega
     886              : 
     887              :          CASE (eri_operator_erf)
     888            7 :             WRITE (iw, '(T3,A,T63,A)') "ERI operator", "Longrange Coulomb"
     889            7 :             IF (.NOT. ex_omega) CALL cp_abort(__LOCATION__, &
     890            0 :                                               "Longrange operator requires OMEGA to be explicitly set")
     891            7 :             WRITE (iw, '(T3,A,T66,F14.3)') "ERI operator parameter OMEGA", eri_op_omega
     892              : 
     893              :          CASE (eri_operator_erfc)
     894            0 :             WRITE (iw, '(T3,A,T62,A)') "ERI operator", "Shortrange Coulomb"
     895            0 :             IF (.NOT. ex_omega) CALL cp_abort(__LOCATION__, &
     896            0 :                                               "Shortrange operator requires OMEGA to be explicitly set")
     897            0 :             WRITE (iw, '(T3,A,T66,F14.3)') "ERI operator parameter OMEGA", eri_op_omega
     898              : 
     899              :          CASE (eri_operator_trunc)
     900            0 :             WRITE (iw, '(T3,A,T63,A)') "ERI operator", "Truncated Coulomb"
     901            0 :             IF (.NOT. ex_rcut) CALL cp_abort(__LOCATION__, &
     902            0 :                                              "Cutoff radius not specified for trunc. Coulomb operator")
     903            0 :             WRITE (iw, '(T3,A,T66,F14.3)') "ERI operator cutoff radius (au)", eri_rcut
     904              : 
     905              :          CASE (eri_operator_lr_trunc)
     906            2 :             WRITE (iw, '(T3,A,T53,A)') "ERI operator", "Longrange truncated Coulomb"
     907            2 :             IF (.NOT. ex_rcut) CALL cp_abort(__LOCATION__, &
     908            0 :                                              "Cutoff radius not specified for trunc. longrange operator")
     909            2 :             WRITE (iw, '(T3,A,T66,F14.3)') "ERI operator cutoff radius (au)", eri_rcut
     910            2 :             IF (.NOT. ex_omega) CALL cp_abort(__LOCATION__, &
     911            0 :                                               "LR truncated operator requires OMEGA to be explicitly set")
     912            2 :             WRITE (iw, '(T3,A,T66,F14.3)') "ERI operator parameter OMEGA", eri_op_omega
     913            2 :             IF (eri_op_omega < 0.01_dp) THEN
     914            0 :                CPABORT("LR truncated operator requires OMEGA >= 0.01 to be stable")
     915              :             END IF
     916              : 
     917              :          CASE DEFAULT
     918           33 :             CPABORT("Unknown ERI operator")
     919              : 
     920              :          END SELECT
     921              : 
     922           33 :          WRITE (iw, '(T3,A,T68,E12.4)') "Accuracy of ERIs", eri_eps_int
     923           33 :          WRITE (iw, '(T3,A,T71,3I3)') "Periodicity", active_space_env%eri%periodicity(1:3)
     924              : 
     925              :          ! TODO: should be moved after ERI calculation, as it depends on screening
     926           33 :          IF (nspins < 2) THEN
     927           27 :             WRITE (iw, '(T3,A,T68,I12)') "Total Number of ERI", (nmo_active**4)/8
     928              :          ELSE
     929            6 :             WRITE (iw, '(T3,A,T68,I12)') "Total Number of ERI (aa|aa)", (nmo_active**4)/8
     930            6 :             WRITE (iw, '(T3,A,T68,I12)') "Total Number of ERI (bb|bb)", (nmo_active**4)/8
     931            6 :             WRITE (iw, '(T3,A,T68,I12)') "Total Number of ERI (aa|bb)", (nmo_active**4)/4
     932              :          END IF
     933              :       END IF
     934              : 
     935              :       ! allocate container for integrals (CSR matrix)
     936           66 :       CALL get_qs_env(qs_env, para_env=para_env)
     937           66 :       m = (nspins*(nspins + 1))/2
     938          288 :       ALLOCATE (active_space_env%eri%eri(m))
     939          156 :       DO i = 1, m
     940           90 :          CALL get_mo_set(active_space_env%mos_active(1), nmo=nmo)
     941           90 :          ALLOCATE (active_space_env%eri%eri(i)%csr_mat)
     942           90 :          eri_mat => active_space_env%eri%eri(i)%csr_mat
     943           90 :          IF (i == 1) THEN
     944           66 :             n1 = nmo
     945           66 :             n2 = nmo
     946           24 :          ELSEIF (i == 2) THEN
     947           12 :             n1 = nmo
     948           12 :             n2 = nmo
     949              :          ELSE
     950           12 :             n1 = nmo
     951           12 :             n2 = nmo
     952              :          END IF
     953           90 :          nn1 = (n1*(n1 + 1))/2
     954           90 :          nn2 = (n2*(n2 + 1))/2
     955           90 :          CALL dbcsr_csr_create(eri_mat, nn1, nn2, 0_int_8, 0, 0, para_env%get_handle())
     956          246 :          active_space_env%eri%norb = nmo
     957              :       END DO
     958              : 
     959           66 :       SELECT CASE (eri_method)
     960              :       CASE (eri_method_full_gpw, eri_method_gpw_ht)
     961           66 :          CALL section_vals_val_get(as_input, "ERI_GPW%EPS_GRID", r_val=eri_eps_grid)
     962           66 :          active_space_env%eri%eri_gpw%eps_grid = eri_eps_grid
     963           66 :          CALL section_vals_val_get(as_input, "ERI_GPW%EPS_FILTER", r_val=eri_eps_filter)
     964           66 :          active_space_env%eri%eri_gpw%eps_filter = eri_eps_filter
     965           66 :          CALL section_vals_val_get(as_input, "ERI_GPW%CUTOFF", r_val=eri_gpw_cutoff)
     966           66 :          active_space_env%eri%eri_gpw%cutoff = eri_gpw_cutoff
     967           66 :          CALL section_vals_val_get(as_input, "ERI_GPW%REL_CUTOFF", r_val=eri_rel_cutoff)
     968           66 :          active_space_env%eri%eri_gpw%rel_cutoff = eri_rel_cutoff
     969           66 :          CALL section_vals_val_get(as_input, "ERI_GPW%PRINT_LEVEL", i_val=eri_print)
     970           66 :          active_space_env%eri%eri_gpw%print_level = eri_print
     971           66 :          CALL section_vals_val_get(as_input, "ERI_GPW%STORE_WFN", l_val=store_wfn)
     972           66 :          active_space_env%eri%eri_gpw%store_wfn = store_wfn
     973           66 :          CALL section_vals_val_get(as_input, "ERI_GPW%GROUP_SIZE", i_val=group_size)
     974           66 :          active_space_env%eri%eri_gpw%group_size = group_size
     975              :          ! Always redo Poisson solver for now
     976           66 :          active_space_env%eri%eri_gpw%redo_poisson = .TRUE.
     977              :          ! active_space_env%eri%eri_gpw%redo_poisson = (ex_operator .OR. ex_perd)
     978           66 :          IF (iw > 0) THEN
     979           33 :             WRITE (iw, '(/,T2,A,T71,F10.1)') "ERI_GPW| Energy cutoff [Ry]", eri_gpw_cutoff
     980           33 :             WRITE (iw, '(T2,A,T71,F10.1)') "ERI_GPW| Relative energy cutoff [Ry]", eri_rel_cutoff
     981              :          END IF
     982              :          !
     983           66 :          CALL calculate_eri_gpw(active_space_env%mos_active, active_space_env%active_orbitals, active_space_env%eri, qs_env, iw)
     984              :          !
     985              :       CASE DEFAULT
     986           66 :          CPABORT("Unknown ERI method")
     987              :       END SELECT
     988           66 :       IF (iw > 0) THEN
     989           78 :          DO isp = 1, SIZE(active_space_env%eri%eri)
     990           45 :             eri_mat => active_space_env%eri%eri(isp)%csr_mat
     991              :             nze_percentage = 100.0_dp*(REAL(eri_mat%nze_total, KIND=dp) &
     992           45 :                                        /REAL(eri_mat%nrows_total, KIND=dp))/REAL(eri_mat%ncols_total, KIND=dp)
     993           45 :             WRITE (iw, '(/,T2,A,I2,T30,A,T68,I12)') "ERI_GPW| Spinmatrix:", isp, &
     994           90 :                "Number of  CSR non-zero elements:", eri_mat%nze_total
     995           45 :             WRITE (iw, '(T2,A,I2,T30,A,T68,F12.4)') "ERI_GPW| Spinmatrix:", isp, &
     996           90 :                "Percentage CSR non-zero elements:", nze_percentage
     997           45 :             WRITE (iw, '(T2,A,I2,T30,A,T68,I12)') "ERI_GPW| Spinmatrix:", isp, &
     998           90 :                "nrows_total", eri_mat%nrows_total
     999           45 :             WRITE (iw, '(T2,A,I2,T30,A,T68,I12)') "ERI_GPW| Spinmatrix:", isp, &
    1000           90 :                "ncols_total", eri_mat%ncols_total
    1001           45 :             WRITE (iw, '(T2,A,I2,T30,A,T68,I12)') "ERI_GPW| Spinmatrix:", isp, &
    1002          123 :                "nrows_local", eri_mat%nrows_local
    1003              :          END DO
    1004           33 :          CALL m_flush(iw)
    1005              :       END IF
    1006           66 :       CALL para_env%sync()
    1007              : 
    1008              :       ! set the reference active space density matrix
    1009           66 :       nspins = active_space_env%nspins
    1010          276 :       ALLOCATE (active_space_env%p_active(nspins))
    1011          144 :       DO isp = 1, nspins
    1012           78 :          mo_set => active_space_env%mos_active(isp)
    1013           78 :          CALL get_mo_set(mo_set, mo_coeff=mo_coeff, nmo=nmo)
    1014          144 :          CALL create_subspace_matrix(mo_coeff, active_space_env%p_active(isp), nmo)
    1015              :       END DO
    1016            0 :       SELECT CASE (mselect)
    1017              :       CASE DEFAULT
    1018            0 :          CPABORT("Unknown orbital selection method")
    1019              :       CASE (casci_canonical, manual_selection)
    1020           66 :          focc = 2.0_dp
    1021           66 :          IF (nspins == 2) focc = 1.0_dp
    1022          144 :          DO isp = 1, nspins
    1023           78 :             fmat => active_space_env%p_active(isp)
    1024           78 :             CALL cp_fm_set_all(fmat, alpha=0.0_dp)
    1025           78 :             IF (nspins == 2) THEN
    1026           24 :                IF (isp == 1) THEN
    1027           12 :                   nel = (active_space_env%nelec_active + active_space_env%multiplicity - 1)/2
    1028              :                ELSE
    1029           12 :                   nel = (active_space_env%nelec_active - active_space_env%multiplicity + 1)/2
    1030              :                END IF
    1031              :             ELSE
    1032           54 :                nel = active_space_env%nelec_active
    1033              :             END IF
    1034          348 :             DO i = 1, nmo_active
    1035          204 :                m = active_space_env%active_orbitals(i, isp)
    1036          204 :                fel = MIN(focc, REAL(nel, KIND=dp))
    1037          204 :                CALL cp_fm_set_element(fmat, m, m, fel)
    1038          204 :                nel = nel - NINT(fel)
    1039          282 :                nel = MAX(nel, 0)
    1040              :             END DO
    1041              :          END DO
    1042              :       CASE (wannier_projection)
    1043            0 :          CPABORT("NOT IMPLEMENTED")
    1044              :       CASE (mao_projection)
    1045           66 :          CPABORT("NOT IMPLEMENTED")
    1046              :       END SELECT
    1047              : 
    1048              :       ! compute alpha-beta overlap matrix in case of spin-polarized calculation
    1049           66 :       CALL calculate_spin_pol_overlap(active_space_env%mos_active, qs_env, active_space_env)
    1050              : 
    1051              :       ! figure out if we have a new xc section for the AS
    1052           66 :       xc_section => section_vals_get_subs_vals(input, "DFT%ACTIVE_SPACE%XC")
    1053           66 :       explicit = .FALSE.
    1054           66 :       IF (ASSOCIATED(xc_section)) CALL section_vals_get(xc_section, explicit=explicit)
    1055              : 
    1056              :       ! rebuild KS matrix if needed
    1057           66 :       IF (explicit) THEN
    1058              :          ! release the hfx data if it was part of the SCF functional
    1059            2 :          IF (ASSOCIATED(qs_env%x_data)) CALL hfx_release(qs_env%x_data)
    1060              :          ! also release the admm environment in case we are using admm
    1061            2 :          IF (ASSOCIATED(qs_env%admm_env)) CALL admm_env_release(qs_env%admm_env)
    1062              : 
    1063              :          CALL get_qs_env(qs_env, atomic_kind_set=atomic_kind_set, qs_kind_set=qs_kind_set, &
    1064            2 :                          particle_set=particle_set, cell=cell, ks_env=ks_env)
    1065            2 :          IF (dft_control%do_admm) THEN
    1066            0 :             basis_type = 'AUX_FIT'
    1067              :          ELSE
    1068            2 :             basis_type = 'ORB'
    1069              :          END IF
    1070            2 :          hfx_section => section_vals_get_subs_vals(xc_section, "HF")
    1071              :          CALL hfx_create(qs_env%x_data, para_env, hfx_section, atomic_kind_set, &
    1072              :                          qs_kind_set, particle_set, dft_control, cell, orb_basis=basis_type, &
    1073            2 :                          nelectron_total=nelec_total)
    1074              : 
    1075            2 :          qs_env%requires_matrix_vxc = .TRUE.  ! needs to be set only once
    1076              : 
    1077              :          ! a bit of a hack: this forces a new re-init of HFX
    1078            2 :          CALL set_ks_env(ks_env, s_mstruct_changed=.TRUE.)
    1079              :          CALL qs_ks_build_kohn_sham_matrix(qs_env, calculate_forces=.FALSE., &
    1080              :                                            just_energy=.FALSE., &
    1081            2 :                                            ext_xc_section=xc_section)
    1082              :          ! we need to reset it to false
    1083            2 :          CALL set_ks_env(ks_env, s_mstruct_changed=.FALSE.)
    1084              :       ELSE
    1085           64 :          xc_section => section_vals_get_subs_vals(input, "DFT%XC")
    1086              :       END IF
    1087              :       ! set the xc_section
    1088           66 :       active_space_env%xc_section => xc_section
    1089              : 
    1090           66 :       CALL get_qs_env(qs_env, energy=energy)
    1091              :       ! transform KS/Fock, Vxc and Hcore to AS MO basis
    1092           66 :       CALL calculate_operators(active_space_env%mos_active, qs_env, active_space_env)
    1093              :       ! set the reference energy in the active space
    1094           66 :       active_space_env%energy_ref = energy%total
    1095              :       ! calculate inactive energy and embedding potential
    1096           66 :       CALL subspace_fock_matrix(active_space_env)
    1097              : 
    1098              :       ! associate the active space environment with the qs environment
    1099           66 :       CALL set_qs_env(qs_env, active_space=active_space_env)
    1100              : 
    1101              :       ! Perform the embedding calculation only if qiskit is specified
    1102           66 :       CALL section_vals_val_get(as_input, "AS_SOLVER", i_val=as_solver)
    1103           66 :       SELECT CASE (as_solver)
    1104              :       CASE (no_solver)
    1105           66 :          IF (iw > 0) THEN
    1106           33 :             WRITE (iw, '(/,T3,A)') "No active space solver specified, skipping embedding calculation"
    1107           33 :             CALL m_flush(iw)
    1108              :          END IF
    1109           66 :          CALL para_env%sync()
    1110              :       CASE (qiskit_solver)
    1111            0 :          CALL rsdft_embedding(qs_env, active_space_env, as_input)
    1112            0 :          CALL qs_scf_compute_properties(qs_env, wf_type="MC-DFT", do_mp2=.FALSE.)
    1113              :       CASE DEFAULT
    1114           66 :          CPABORT("Unknown active space solver")
    1115              :       END SELECT
    1116              : 
    1117              :       ! Output a FCIDUMP file if requested
    1118           66 :       IF (active_space_env%fcidump) CALL fcidump(active_space_env, as_input)
    1119              : 
    1120              :       ! Output a QCSchema file if requested
    1121           66 :       IF (active_space_env%qcschema) THEN
    1122            0 :          CALL qcschema_env_create(qcschema_env, qs_env)
    1123            0 :          CALL qcschema_to_hdf5(qcschema_env, active_space_env%qcschema_filename)
    1124            0 :          CALL qcschema_env_release(qcschema_env)
    1125              :       END IF
    1126              : 
    1127           66 :       IF (iw > 0) THEN
    1128              :          WRITE (iw, '(/,T2,A)') &
    1129           33 :             '!-------------------- End of Active Space Interface --------------------------!'
    1130           33 :          CALL m_flush(iw)
    1131              :       END IF
    1132           66 :       CALL para_env%sync()
    1133              : 
    1134           66 :       CALL timestop(handle)
    1135              : 
    1136        89204 :    END SUBROUTINE active_space_main
    1137              : 
    1138              : ! **************************************************************************************************
    1139              : !> \brief computes the alpha-beta overlap within the active subspace
    1140              : !> \param mos the molecular orbital set within the active subspace
    1141              : !> \param qs_env ...
    1142              : !> \param active_space_env ...
    1143              : !> \par History
    1144              : !>      04.2016 created [JGH]
    1145              : ! **************************************************************************************************
    1146           66 :    SUBROUTINE calculate_spin_pol_overlap(mos, qs_env, active_space_env)
    1147              : 
    1148              :       TYPE(mo_set_type), DIMENSION(:), INTENT(IN)        :: mos
    1149              :       TYPE(qs_environment_type), POINTER                 :: qs_env
    1150              :       TYPE(active_space_type), POINTER                   :: active_space_env
    1151              : 
    1152              :       CHARACTER(len=*), PARAMETER :: routineN = 'calculate_spin_pol_overlap'
    1153              : 
    1154              :       INTEGER                                            :: handle, nmo, nspins
    1155              :       TYPE(cp_fm_type), POINTER                          :: mo_coeff_a, mo_coeff_b
    1156           66 :       TYPE(dbcsr_p_type), DIMENSION(:), POINTER          :: s_matrix
    1157              : 
    1158           66 :       CALL timeset(routineN, handle)
    1159              : 
    1160           66 :       nspins = active_space_env%nspins
    1161              : 
    1162              :       ! overlap in AO
    1163           66 :       IF (nspins > 1) THEN
    1164           12 :          CALL get_qs_env(qs_env, matrix_s=s_matrix)
    1165           24 :          ALLOCATE (active_space_env%sab_sub(1))
    1166              : 
    1167           12 :          CALL get_mo_set(mo_set=mos(1), mo_coeff=mo_coeff_a, nmo=nmo)
    1168           12 :          CALL get_mo_set(mo_set=mos(2), mo_coeff=mo_coeff_b, nmo=nmo)
    1169           12 :          CALL subspace_operator(mo_coeff_a, nmo, s_matrix(1)%matrix, active_space_env%sab_sub(1), mo_coeff_b)
    1170              :       END IF
    1171              : 
    1172           66 :       CALL timestop(handle)
    1173              : 
    1174           66 :    END SUBROUTINE calculate_spin_pol_overlap
    1175              : 
    1176              : ! **************************************************************************************************
    1177              : !> \brief computes the one-electron operators in the subspace of the provided orbital set
    1178              : !> \param mos the molecular orbital set within the active subspace
    1179              : !> \param qs_env ...
    1180              : !> \param active_space_env ...
    1181              : !> \par History
    1182              : !>      04.2016 created [JGH]
    1183              : ! **************************************************************************************************
    1184           66 :    SUBROUTINE calculate_operators(mos, qs_env, active_space_env)
    1185              : 
    1186              :       TYPE(mo_set_type), DIMENSION(:), INTENT(IN)        :: mos
    1187              :       TYPE(qs_environment_type), POINTER                 :: qs_env
    1188              :       TYPE(active_space_type), POINTER                   :: active_space_env
    1189              : 
    1190              :       CHARACTER(len=*), PARAMETER :: routineN = 'calculate_operators'
    1191              : 
    1192              :       INTEGER                                            :: handle, ispin, nmo, nspins
    1193              :       TYPE(cp_fm_type), POINTER                          :: mo_coeff
    1194           66 :       TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER       :: h_matrix, ks_matrix
    1195              : 
    1196           66 :       CALL timeset(routineN, handle)
    1197              : 
    1198           66 :       nspins = active_space_env%nspins
    1199              : 
    1200              :       ! Kohn-Sham / Fock operator
    1201           66 :       CALL cp_fm_release(active_space_env%ks_sub)
    1202           66 :       CALL get_qs_env(qs_env, matrix_ks_kp=ks_matrix)
    1203          276 :       ALLOCATE (active_space_env%ks_sub(nspins))
    1204          144 :       DO ispin = 1, nspins
    1205           78 :          CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
    1206          144 :          CALL subspace_operator(mo_coeff, nmo, ks_matrix(ispin, 1)%matrix, active_space_env%ks_sub(ispin))
    1207              :       END DO
    1208              : 
    1209              :       ! Core Hamiltonian
    1210           66 :       CALL cp_fm_release(active_space_env%h_sub)
    1211              : 
    1212           66 :       NULLIFY (h_matrix)
    1213           66 :       CALL get_qs_env(qs_env=qs_env, matrix_h_kp=h_matrix)
    1214          276 :       ALLOCATE (active_space_env%h_sub(nspins))
    1215          144 :       DO ispin = 1, nspins
    1216           78 :          CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
    1217          144 :          CALL subspace_operator(mo_coeff, nmo, h_matrix(1, 1)%matrix, active_space_env%h_sub(ispin))
    1218              :       END DO
    1219              : 
    1220           66 :       CALL timestop(handle)
    1221              : 
    1222           66 :    END SUBROUTINE calculate_operators
    1223              : 
    1224              : ! **************************************************************************************************
    1225              : !> \brief computes a one-electron operator in the subspace of the provided orbital set
    1226              : !> \param mo_coeff the orbital coefficient matrix
    1227              : !> \param nmo the number of subspace orbitals
    1228              : !> \param op_matrix operator matrix in AO basis
    1229              : !> \param op_sub operator in orbital basis
    1230              : !> \param mo_coeff_b the beta orbital coefficients
    1231              : !> \par History
    1232              : !>      04.2016 created [JGH]
    1233              : ! **************************************************************************************************
    1234          336 :    SUBROUTINE subspace_operator(mo_coeff, nmo, op_matrix, op_sub, mo_coeff_b)
    1235              : 
    1236              :       TYPE(cp_fm_type), INTENT(IN)                       :: mo_coeff
    1237              :       INTEGER, INTENT(IN)                                :: nmo
    1238              :       TYPE(dbcsr_type), POINTER                          :: op_matrix
    1239              :       TYPE(cp_fm_type), INTENT(INOUT)                    :: op_sub
    1240              :       TYPE(cp_fm_type), INTENT(IN), OPTIONAL             :: mo_coeff_b
    1241              : 
    1242              :       CHARACTER(len=*), PARAMETER                        :: routineN = 'subspace_operator'
    1243              : 
    1244              :       INTEGER                                            :: handle, ncol, nrow
    1245              :       TYPE(cp_fm_type)                                   :: vectors
    1246              : 
    1247          168 :       CALL timeset(routineN, handle)
    1248              : 
    1249          168 :       CALL cp_fm_get_info(matrix=mo_coeff, ncol_global=ncol, nrow_global=nrow)
    1250          168 :       CPASSERT(nmo <= ncol)
    1251              : 
    1252          168 :       IF (nmo > 0) THEN
    1253          168 :          CALL cp_fm_create(vectors, mo_coeff%matrix_struct, "vectors")
    1254          168 :          CALL create_subspace_matrix(mo_coeff, op_sub, nmo)
    1255              : 
    1256          168 :          IF (PRESENT(mo_coeff_b)) THEN
    1257              :             ! if beta orbitals are present, compute the cross alpha_beta term
    1258           12 :             CALL cp_dbcsr_sm_fm_multiply(op_matrix, mo_coeff_b, vectors, nmo)
    1259              :          ELSE
    1260              :             ! otherwise the same spin, whatever that is
    1261          156 :             CALL cp_dbcsr_sm_fm_multiply(op_matrix, mo_coeff, vectors, nmo)
    1262              :          END IF
    1263              : 
    1264          168 :          CALL parallel_gemm('T', 'N', nmo, nmo, nrow, 1.0_dp, mo_coeff, vectors, 0.0_dp, op_sub)
    1265          168 :          CALL cp_fm_release(vectors)
    1266              :       END IF
    1267              : 
    1268          168 :       CALL timestop(handle)
    1269              : 
    1270          168 :    END SUBROUTINE subspace_operator
    1271              : 
    1272              : ! **************************************************************************************************
    1273              : !> \brief creates a matrix of subspace size
    1274              : !> \param orbitals the orbital coefficient matrix
    1275              : !> \param op_sub operator in orbital basis
    1276              : !> \param n the number of orbitals
    1277              : !> \par History
    1278              : !>      04.2016 created [JGH]
    1279              : ! **************************************************************************************************
    1280          246 :    SUBROUTINE create_subspace_matrix(orbitals, op_sub, n)
    1281              : 
    1282              :       TYPE(cp_fm_type), INTENT(IN)                       :: orbitals
    1283              :       TYPE(cp_fm_type), INTENT(OUT)                      :: op_sub
    1284              :       INTEGER, INTENT(IN)                                :: n
    1285              : 
    1286              :       TYPE(cp_fm_struct_type), POINTER                   :: fm_struct
    1287              : 
    1288          246 :       IF (n > 0) THEN
    1289              : 
    1290          246 :          NULLIFY (fm_struct)
    1291              :          CALL cp_fm_struct_create(fm_struct, nrow_global=n, ncol_global=n, &
    1292              :                                   para_env=orbitals%matrix_struct%para_env, &
    1293          246 :                                   context=orbitals%matrix_struct%context)
    1294          246 :          CALL cp_fm_create(op_sub, fm_struct, name="Subspace operator")
    1295          246 :          CALL cp_fm_struct_release(fm_struct)
    1296              : 
    1297              :       END IF
    1298              : 
    1299          246 :    END SUBROUTINE create_subspace_matrix
    1300              : 
    1301              : ! **************************************************************************************************
    1302              : !> \brief computes the electron repulsion integrals using the GPW technology
    1303              : !> \param mos the molecular orbital set within the active subspace
    1304              : !> \param orbitals ...
    1305              : !> \param eri_env ...
    1306              : !> \param qs_env ...
    1307              : !> \param iw ...
    1308              : !> \par History
    1309              : !>      04.2016 created [JGH]
    1310              : ! **************************************************************************************************
    1311           66 :    SUBROUTINE calculate_eri_gpw(mos, orbitals, eri_env, qs_env, iw)
    1312              : 
    1313              :       TYPE(mo_set_type), DIMENSION(:), INTENT(IN)        :: mos
    1314              :       INTEGER, DIMENSION(:, :), POINTER                  :: orbitals
    1315              :       TYPE(eri_type)                                     :: eri_env
    1316              :       TYPE(qs_environment_type), POINTER                 :: qs_env
    1317              :       INTEGER, INTENT(IN)                                :: iw
    1318              : 
    1319              :       CHARACTER(len=*), PARAMETER                        :: routineN = 'calculate_eri_gpw'
    1320              : 
    1321              :       INTEGER :: col_local, color, handle, i1, i2, i3, i4, i_multigrid, icount2, intcount, &
    1322              :          irange(2), isp, isp1, isp2, ispin, iwa1, iwa12, iwa2, iwb1, iwb12, iwb2, iwbs, iwbt, &
    1323              :          iwfn, n_multigrid, ncol_global, ncol_local, nmm, nmo, nmo1, nmo2, nrow_global, &
    1324              :          nrow_local, nspins, number_of_subgroups, nx, row_local, stored_integrals
    1325           66 :       INTEGER, ALLOCATABLE, DIMENSION(:)                 :: eri_index
    1326           66 :       INTEGER, DIMENSION(:), POINTER                     :: col_indices, row_indices
    1327              :       LOGICAL                                            :: print1, print2, &
    1328              :                                                             skip_load_balance_distributed
    1329              :       REAL(KIND=dp)                                      :: dvol, erint, pair_int, &
    1330              :                                                             progression_factor, rc, rsize, t1, t2
    1331           66 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:)           :: eri
    1332           66 :       TYPE(atomic_kind_type), DIMENSION(:), POINTER      :: atomic_kind_set
    1333              :       TYPE(cell_type), POINTER                           :: cell
    1334              :       TYPE(cp_blacs_env_type), POINTER                   :: blacs_env, blacs_env_sub
    1335              :       TYPE(cp_fm_struct_type), POINTER                   :: fm_struct
    1336           66 :       TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:)        :: fm_matrix_pq_rnu, fm_matrix_pq_rs, &
    1337           66 :                                                             fm_mo_coeff_as
    1338              :       TYPE(cp_fm_type), POINTER                          :: mo_coeff
    1339              :       TYPE(dbcsr_p_type)                                 :: mat_munu
    1340           66 :       TYPE(dbcsr_type), ALLOCATABLE, DIMENSION(:)        :: matrix_pq_rnu, mo_coeff_as
    1341              :       TYPE(dft_control_type), POINTER                    :: dft_control
    1342              :       TYPE(mp_para_env_type), POINTER                    :: para_env
    1343              :       TYPE(neighbor_list_set_p_type), DIMENSION(:), &
    1344           66 :          POINTER                                         :: sab_orb_sub
    1345           66 :       TYPE(particle_type), DIMENSION(:), POINTER         :: particle_set
    1346              :       TYPE(pw_c1d_gs_type)                               :: pot_g, rho_g
    1347              :       TYPE(pw_env_type), POINTER                         :: pw_env_sub
    1348              :       TYPE(pw_poisson_type), POINTER                     :: poisson_env
    1349              :       TYPE(pw_pool_type), POINTER                        :: auxbas_pw_pool
    1350              :       TYPE(pw_r3d_rs_type)                               :: rho_r, wfn_r
    1351              :       TYPE(pw_r3d_rs_type), ALLOCATABLE, &
    1352           66 :          DIMENSION(:, :), TARGET                         :: wfn_a
    1353              :       TYPE(pw_r3d_rs_type), POINTER                      :: wfn1, wfn2, wfn3, wfn4
    1354              :       TYPE(qs_control_type), POINTER                     :: qs_control, qs_control_old
    1355           66 :       TYPE(qs_kind_type), DIMENSION(:), POINTER          :: qs_kind_set
    1356              :       TYPE(qs_ks_env_type), POINTER                      :: ks_env
    1357              :       TYPE(task_list_type), POINTER                      :: task_list_sub
    1358              : 
    1359           66 :       CALL timeset(routineN, handle)
    1360              : 
    1361           66 :       IF (iw > 0) t1 = m_walltime()
    1362              : 
    1363              :       ! print levels
    1364          128 :       SELECT CASE (eri_env%eri_gpw%print_level)
    1365              :       CASE (silent_print_level)
    1366           62 :          print1 = .FALSE.
    1367           62 :          print2 = .FALSE.
    1368              :       CASE (low_print_level)
    1369            2 :          print1 = .FALSE.
    1370            2 :          print2 = .FALSE.
    1371              :       CASE (medium_print_level)
    1372            2 :          print1 = .TRUE.
    1373            2 :          print2 = .FALSE.
    1374              :       CASE (high_print_level)
    1375            0 :          print1 = .TRUE.
    1376            0 :          print2 = .TRUE.
    1377              :       CASE (debug_print_level)
    1378            0 :          print1 = .TRUE.
    1379           66 :          print2 = .TRUE.
    1380              :       CASE DEFAULT
    1381              :          ! do nothing
    1382              :       END SELECT
    1383              : 
    1384              :       ! Check the input group
    1385           66 :       CALL get_qs_env(qs_env, para_env=para_env, blacs_env=blacs_env)
    1386           66 :       IF (eri_env%eri_gpw%group_size < 1) eri_env%eri_gpw%group_size = para_env%num_pe
    1387           66 :       IF (MOD(para_env%num_pe, eri_env%eri_gpw%group_size) /= 0) &
    1388            0 :          CPABORT("Group size must be a divisor of the total number of processes!")
    1389              :       ! Create a new para_env or reuse the old one
    1390           66 :       IF (eri_env%eri_gpw%group_size == para_env%num_pe) THEN
    1391           60 :          eri_env%para_env_sub => para_env
    1392           60 :          CALL eri_env%para_env_sub%retain()
    1393           60 :          blacs_env_sub => blacs_env
    1394           60 :          CALL blacs_env_sub%retain()
    1395           60 :          number_of_subgroups = 1
    1396           60 :          color = 0
    1397              :       ELSE
    1398            6 :          number_of_subgroups = para_env%num_pe/eri_env%eri_gpw%group_size
    1399            6 :          color = para_env%mepos/eri_env%eri_gpw%group_size
    1400            6 :          ALLOCATE (eri_env%para_env_sub)
    1401            6 :          CALL eri_env%para_env_sub%from_split(para_env, color)
    1402            6 :          NULLIFY (blacs_env_sub)
    1403            6 :          CALL cp_blacs_env_create(blacs_env_sub, eri_env%para_env_sub, BLACS_GRID_SQUARE, .TRUE.)
    1404              :       END IF
    1405           66 :       CALL eri_env%comm_exchange%from_split(para_env, eri_env%para_env_sub%mepos)
    1406              : 
    1407              :       ! This should be done differently! Copied from MP2 code
    1408           66 :       CALL get_qs_env(qs_env, dft_control=dft_control)
    1409          264 :       ALLOCATE (qs_control)
    1410           66 :       qs_control_old => dft_control%qs_control
    1411           66 :       qs_control = qs_control_old
    1412           66 :       dft_control%qs_control => qs_control
    1413           66 :       progression_factor = qs_control%progression_factor
    1414           66 :       n_multigrid = SIZE(qs_control%e_cutoff)
    1415           66 :       nspins = SIZE(mos)
    1416              :       ! Allocate new cutoffs (just in private qs_control, not in qs_control_old)
    1417          198 :       ALLOCATE (qs_control%e_cutoff(n_multigrid))
    1418              : 
    1419           66 :       qs_control%cutoff = eri_env%eri_gpw%cutoff*0.5_dp
    1420           66 :       qs_control%e_cutoff(1) = qs_control%cutoff
    1421          264 :       DO i_multigrid = 2, n_multigrid
    1422              :          qs_control%e_cutoff(i_multigrid) = qs_control%e_cutoff(i_multigrid - 1) &
    1423          264 :                                             /progression_factor
    1424              :       END DO
    1425           66 :       qs_control%relative_cutoff = eri_env%eri_gpw%rel_cutoff*0.5_dp
    1426              : 
    1427              :       ! For now, we will distribute neighbor lists etc. within the global communicator
    1428           66 :       CALL get_qs_env(qs_env, ks_env=ks_env)
    1429              :       CALL create_mat_munu(mat_munu, qs_env, eri_env%eri_gpw%eps_grid, blacs_env_sub, sab_orb_sub=sab_orb_sub, &
    1430           66 :                            do_alloc_blocks_from_nbl=.TRUE., dbcsr_sym_type=dbcsr_type_symmetric)
    1431           66 :       CALL dbcsr_set(mat_munu%matrix, 0.0_dp)
    1432              : 
    1433              :       ! Generate the appropriate pw_env
    1434           66 :       NULLIFY (pw_env_sub)
    1435           66 :       CALL pw_env_create(pw_env_sub)
    1436           66 :       CALL pw_env_rebuild(pw_env_sub, qs_env, external_para_env=eri_env%para_env_sub)
    1437           66 :       CALL pw_env_get(pw_env_sub, auxbas_pw_pool=auxbas_pw_pool, poisson_env=poisson_env)
    1438              : 
    1439              :       ! TODO: maybe we can let `pw_env_rebuild` do what we manually overwrite here?
    1440           66 :       IF (eri_env%eri_gpw%redo_poisson) THEN
    1441              :          ! We need to rebuild the Poisson solver on the fly
    1442          264 :          IF (SUM(eri_env%periodicity) /= 0) THEN
    1443            6 :             poisson_env%parameters%solver = pw_poisson_periodic
    1444              :          ELSE
    1445           60 :             poisson_env%parameters%solver = pw_poisson_analytic
    1446              :          END IF
    1447          264 :          poisson_env%parameters%periodic = eri_env%periodicity
    1448              : 
    1449              :          ! Rebuilds the poisson green (influence) function according
    1450              :          ! to the poisson solver and parameters set so far.
    1451              :          ! Also sets the variable poisson_env%rebuild to .FALSE.
    1452           66 :          CALL pw_poisson_rebuild(poisson_env)
    1453              : 
    1454              :          ! set the cutoff radius for the Greens function in case we use ANALYTIC Poisson solver
    1455           66 :          CALL get_qs_env(qs_env, cell=cell)
    1456           66 :          rc = cell%hmat(1, 1)
    1457          264 :          DO iwa1 = 1, 3
    1458              :             ! TODO: I think this is not the largest possible radius inscribed in the cell
    1459          264 :             rc = MIN(rc, 0.5_dp*cell%hmat(iwa1, iwa1))
    1460              :          END DO
    1461           66 :          poisson_env%green_fft%radius = rc
    1462              : 
    1463              :          ! Overwrite the Greens function with the one we want
    1464           66 :          CALL pw_eri_green_create(poisson_env%green_fft, eri_env)
    1465              : 
    1466           66 :          IF (iw > 0) THEN
    1467           33 :             CALL get_qs_env(qs_env, cell=cell)
    1468          264 :             IF (SUM(cell%perd) /= SUM(eri_env%periodicity)) THEN
    1469            0 :                IF (SUM(eri_env%periodicity) /= 0) THEN
    1470              :                   WRITE (UNIT=iw, FMT="(/,T2,A,T51,A30)") &
    1471            0 :                      "ERI_GPW| Switching Poisson solver to", "PERIODIC"
    1472              :                ELSE
    1473              :                   WRITE (UNIT=iw, FMT="(/,T2,A,T51,A30)") &
    1474            0 :                      "ERI_GPW| Switching Poisson solver to", "ANALYTIC"
    1475              :                END IF
    1476              :             END IF
    1477              :             ! print out the Greens function to check it matches the Poisson solver
    1478           36 :             SELECT CASE (poisson_env%green_fft%method)
    1479              :             CASE (PERIODIC3D)
    1480              :                WRITE (UNIT=iw, FMT="(T2,A,T51,A30)") &
    1481            3 :                   "ERI_GPW| Poisson Greens function", "PERIODIC"
    1482              :             CASE (ANALYTIC0D)
    1483              :                WRITE (UNIT=iw, FMT="(T2,A,T51,A30)") &
    1484           30 :                   "ERI_GPW| Poisson Greens function", "ANALYTIC"
    1485           30 :                WRITE (UNIT=iw, FMT="(T2,A,T71,F10.4)") "ERI_GPW| Poisson cutoff radius", &
    1486           60 :                   poisson_env%green_fft%radius*angstrom
    1487              :             CASE DEFAULT
    1488           33 :                CPABORT("Wrong Greens function setup")
    1489              :             END SELECT
    1490              :          END IF
    1491              :       END IF
    1492              : 
    1493          486 :       ALLOCATE (mo_coeff_as(nspins), fm_mo_coeff_as(nspins))
    1494          144 :       DO ispin = 1, nspins
    1495          234 :          BLOCK
    1496           78 :             REAL(KIND=dp), DIMENSION(:, :), ALLOCATABLE :: C
    1497              :             INTEGER :: nmo
    1498           78 :             TYPE(group_dist_d1_type) :: gd_array
    1499              :             TYPE(cp_fm_type), POINTER :: mo_coeff
    1500           78 :             CALL get_mo_set(mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
    1501              :             CALL grep_rows_in_subgroups(para_env, eri_env%para_env_sub, mo_coeff, gd_array, C)
    1502              : 
    1503              :             CALL build_dbcsr_from_rows(eri_env%para_env_sub, mo_coeff_as(ispin), &
    1504           78 :                                        C(:, :nmo), mat_munu%matrix, gd_array, eri_env%eri_gpw%eps_filter)
    1505           78 :             CALL release_group_dist(gd_array)
    1506          234 :             DEALLOCATE (C)
    1507              :          END BLOCK
    1508              : 
    1509           78 :          CALL dbcsr_get_info(mo_coeff_as(ispin), nfullrows_total=nrow_global, nfullcols_total=ncol_global)
    1510              : 
    1511           78 :          NULLIFY (fm_struct)
    1512              :          CALL cp_fm_struct_create(fm_struct, context=blacs_env_sub, para_env=eri_env%para_env_sub, &
    1513           78 :                                   nrow_global=nrow_global, ncol_global=ncol_global)
    1514           78 :          CALL cp_fm_create(fm_mo_coeff_as(ispin), fm_struct)
    1515           78 :          CALL cp_fm_struct_release(fm_struct)
    1516              : 
    1517          222 :          CALL copy_dbcsr_to_fm(mo_coeff_as(ispin), fm_mo_coeff_as(ispin))
    1518              :       END DO
    1519              : 
    1520           66 :       IF (eri_env%method == eri_method_gpw_ht) THEN
    1521              :          ! We need a task list
    1522           12 :          NULLIFY (task_list_sub)
    1523           12 :          skip_load_balance_distributed = dft_control%qs_control%skip_load_balance_distributed
    1524           12 :          CALL allocate_task_list(task_list_sub)
    1525              :          CALL generate_qs_task_list(ks_env, task_list_sub, basis_type="ORB", &
    1526              :                                     reorder_rs_grid_ranks=.TRUE., &
    1527              :                                     skip_load_balance_distributed=skip_load_balance_distributed, &
    1528           12 :                                     pw_env_external=pw_env_sub, sab_orb_external=sab_orb_sub)
    1529              : 
    1530              :          ! Create sparse matrices carrying the matrix products, Code borrowed from the MP2 GPW method
    1531              :          ! Create equal distributions for them (no sparsity present)
    1532              :          ! We use the routines from mp2 suggesting that one may replicate the grids later for better performance
    1533           84 :          ALLOCATE (matrix_pq_rnu(nspins), fm_matrix_pq_rnu(nspins), fm_matrix_pq_rs(nspins))
    1534           24 :          DO ispin = 1, nspins
    1535           12 :             CALL dbcsr_create(matrix_pq_rnu(ispin), template=mo_coeff_as(ispin))
    1536           12 :             CALL dbcsr_set(matrix_pq_rnu(ispin), 0.0_dp)
    1537              : 
    1538           12 :             CALL dbcsr_get_info(matrix_pq_rnu(ispin), nfullrows_total=nrow_global, nfullcols_total=ncol_global)
    1539              : 
    1540           12 :             NULLIFY (fm_struct)
    1541              :             CALL cp_fm_struct_create(fm_struct, context=blacs_env_sub, para_env=eri_env%para_env_sub, &
    1542           12 :                                      nrow_global=nrow_global, ncol_global=ncol_global)
    1543           12 :             CALL cp_fm_create(fm_matrix_pq_rnu(ispin), fm_struct)
    1544           12 :             CALL cp_fm_struct_release(fm_struct)
    1545              : 
    1546           12 :             NULLIFY (fm_struct)
    1547              :             CALL cp_fm_struct_create(fm_struct, context=blacs_env_sub, para_env=eri_env%para_env_sub, &
    1548           12 :                                      nrow_global=ncol_global, ncol_global=ncol_global)
    1549           12 :             CALL cp_fm_create(fm_matrix_pq_rs(ispin), fm_struct)
    1550           36 :             CALL cp_fm_struct_release(fm_struct)
    1551              :          END DO
    1552              : 
    1553              :          ! Copy the active space of the MOs into DBCSR matrices
    1554              :       END IF
    1555              : 
    1556           66 :       CALL auxbas_pw_pool%create_pw(wfn_r)
    1557           66 :       CALL auxbas_pw_pool%create_pw(rho_g)
    1558              :       CALL get_qs_env(qs_env, qs_kind_set=qs_kind_set, cell=cell, &
    1559           66 :                       particle_set=particle_set, atomic_kind_set=atomic_kind_set)
    1560              : 
    1561              :       ! pre-calculate wavefunctions on reals space grid
    1562           66 :       nspins = SIZE(mos)
    1563           66 :       IF (eri_env%eri_gpw%store_wfn) THEN
    1564              :          ! pre-calculate wavefunctions on reals space grid
    1565           54 :          rsize = 0.0_dp
    1566           54 :          nmo = 0
    1567          120 :          DO ispin = 1, nspins
    1568           66 :             CALL get_mo_set(mo_set=mos(ispin), nmo=nx)
    1569           66 :             nmo = MAX(nmo, nx)
    1570          318 :             rsize = REAL(SIZE(wfn_r%array), KIND=dp)*nx
    1571              :          END DO
    1572           54 :          IF (print1 .AND. iw > 0) THEN
    1573            1 :             rsize = rsize*8._dp/1000000._dp
    1574            1 :             WRITE (iw, "(T2,'ERI_GPW|',' Store active orbitals on real space grid ',T66,F12.3,' MB')") rsize
    1575              :          END IF
    1576          516 :          ALLOCATE (wfn_a(nmo, nspins))
    1577          120 :          DO ispin = 1, nspins
    1578           66 :             CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
    1579          300 :             DO i1 = 1, SIZE(orbitals, 1)
    1580          180 :                iwfn = orbitals(i1, ispin)
    1581          180 :                CALL auxbas_pw_pool%create_pw(wfn_a(iwfn, ispin))
    1582              :                CALL calculate_wavefunction(mo_coeff, iwfn, wfn_a(iwfn, ispin), rho_g, atomic_kind_set, &
    1583          180 :                                            qs_kind_set, cell, dft_control, particle_set, pw_env_sub)
    1584          246 :                IF (print2 .AND. iw > 0) THEN
    1585            0 :                   WRITE (iw, "(T2,'ERI_GPW|',' Orbital stored ',I4,'  Spin ',i1)") iwfn, ispin
    1586              :                END IF
    1587              :             END DO
    1588              :          END DO
    1589              :       ELSE
    1590              :          ! Even if we do not store all WFNs, we still need containers for the functions to store
    1591           12 :          ALLOCATE (wfn1, wfn2)
    1592           12 :          CALL auxbas_pw_pool%create_pw(wfn1)
    1593           12 :          CALL auxbas_pw_pool%create_pw(wfn2)
    1594           12 :          IF (eri_env%method /= eri_method_gpw_ht) THEN
    1595            6 :             ALLOCATE (wfn3, wfn4)
    1596            6 :             CALL auxbas_pw_pool%create_pw(wfn3)
    1597            6 :             CALL auxbas_pw_pool%create_pw(wfn4)
    1598              :          END IF
    1599              :       END IF
    1600              : 
    1601              :       ! get some of the grids ready
    1602           66 :       CALL auxbas_pw_pool%create_pw(rho_r)
    1603           66 :       CALL auxbas_pw_pool%create_pw(pot_g)
    1604              : 
    1605              :       ! run the FFT once, to set up buffers and to take into account the memory
    1606           66 :       CALL pw_zero(rho_r)
    1607           66 :       CALL pw_transfer(rho_r, rho_g)
    1608           66 :       dvol = rho_r%pw_grid%dvol
    1609              : 
    1610           66 :       IF (iw > 0) THEN
    1611           33 :          CALL m_flush(iw)
    1612              :       END IF
    1613              :       ! calculate the integrals
    1614           66 :       stored_integrals = 0
    1615          144 :       DO isp1 = 1, nspins
    1616           78 :          CALL get_mo_set(mo_set=mos(isp1), nmo=nmo1)
    1617           78 :          nmm = (nmo1*(nmo1 + 1))/2
    1618           78 :          irange = get_irange_csr(nmm, eri_env%comm_exchange)
    1619          426 :          DO i1 = 1, SIZE(orbitals, 1)
    1620          204 :             iwa1 = orbitals(i1, isp1)
    1621          204 :             IF (eri_env%eri_gpw%store_wfn) THEN
    1622          180 :                wfn1 => wfn_a(iwa1, isp1)
    1623              :             ELSE
    1624              :                CALL calculate_wavefunction(fm_mo_coeff_as(isp1), iwa1, wfn1, rho_g, atomic_kind_set, &
    1625           24 :                                            qs_kind_set, cell, dft_control, particle_set, pw_env_sub)
    1626              :             END IF
    1627          682 :             DO i2 = i1, SIZE(orbitals, 1)
    1628          400 :                iwa2 = orbitals(i2, isp1)
    1629          400 :                iwa12 = csr_idx_to_combined(iwa1, iwa2, nmo1)
    1630              :                ! Skip calculation directly if the pair is not part of our subgroup
    1631          400 :                IF (iwa12 < irange(1) .OR. iwa12 > irange(2)) CYCLE
    1632          391 :                iwa12 = iwa12 - irange(1) + 1
    1633          391 :                IF (eri_env%eri_gpw%store_wfn) THEN
    1634          361 :                   wfn2 => wfn_a(iwa2, isp1)
    1635              :                ELSE
    1636              :                   CALL calculate_wavefunction(fm_mo_coeff_as(isp1), iwa2, wfn2, rho_g, atomic_kind_set, &
    1637           30 :                                               qs_kind_set, cell, dft_control, particle_set, pw_env_sub)
    1638              :                END IF
    1639              :                ! calculate charge distribution and potential
    1640          391 :                CALL pw_zero(rho_r)
    1641          391 :                CALL pw_multiply(rho_r, wfn1, wfn2)
    1642          391 :                CALL pw_transfer(rho_r, rho_g)
    1643          391 :                CALL pw_poisson_solve(poisson_env, rho_g, pair_int, pot_g)
    1644              : 
    1645              :                ! screening using pair_int
    1646          391 :                IF (pair_int < eri_env%eps_integral) CYCLE
    1647          391 :                CALL pw_transfer(pot_g, rho_r)
    1648              :                !
    1649          986 :                IF (eri_env%method == eri_method_gpw_ht) THEN
    1650           30 :                   CALL pw_scale(rho_r, dvol)
    1651           60 :                   DO isp2 = isp1, nspins
    1652           30 :                      CALL get_mo_set(mo_set=mos(isp2), nmo=nmo2)
    1653           30 :                      nx = (nmo2*(nmo2 + 1))/2
    1654          150 :                      ALLOCATE (eri(nx), eri_index(nx))
    1655           30 :                      CALL dbcsr_set(mat_munu%matrix, 0.0_dp)
    1656              :                      CALL integrate_v_rspace(rho_r, hmat=mat_munu, qs_env=qs_env, &
    1657              :                                              calculate_forces=.FALSE., compute_tau=.FALSE., gapw=.FALSE., &
    1658           30 :                                              pw_env_external=pw_env_sub, task_list_external=task_list_sub)
    1659              : 
    1660              :                      CALL dbcsr_multiply("N", "N", 1.0_dp, mat_munu%matrix, mo_coeff_as(isp2), &
    1661           30 :                                          0.0_dp, matrix_pq_rnu(isp2), filter_eps=eri_env%eri_gpw%eps_filter)
    1662           30 :                      CALL copy_dbcsr_to_fm(matrix_pq_rnu(isp2), fm_matrix_pq_rnu(isp2))
    1663              : 
    1664           30 :                      CALL cp_fm_get_info(fm_matrix_pq_rnu(isp2), ncol_global=ncol_global, nrow_global=nrow_global)
    1665              : 
    1666              :                      CALL parallel_gemm("T", "N", ncol_global, ncol_global, nrow_global, 0.5_dp, &
    1667              :                                         fm_matrix_pq_rnu(isp2), fm_mo_coeff_as(isp2), &
    1668           30 :                                         0.0_dp, fm_matrix_pq_rs(isp2))
    1669              :                      CALL parallel_gemm("T", "N", ncol_global, ncol_global, nrow_global, 0.5_dp, &
    1670              :                                         fm_mo_coeff_as(isp2), fm_matrix_pq_rnu(isp2), &
    1671           30 :                                         1.0_dp, fm_matrix_pq_rs(isp2))
    1672              : 
    1673              :                      CALL cp_fm_get_info(fm_matrix_pq_rs(isp2), ncol_local=ncol_local, nrow_local=nrow_local, &
    1674           30 :                                          col_indices=col_indices, row_indices=row_indices)
    1675              : 
    1676           30 :                      icount2 = 0
    1677           90 :                      DO col_local = 1, ncol_local
    1678           60 :                         iwb2 = orbitals(col_indices(col_local), isp2)
    1679          162 :                         DO row_local = 1, nrow_local
    1680           72 :                            iwb1 = orbitals(row_indices(row_local), isp2)
    1681              : 
    1682          132 :                            IF (iwb1 <= iwb2) THEN
    1683           54 :                               iwb12 = csr_idx_to_combined(iwb1, iwb2, nmo2)
    1684           54 :                               erint = fm_matrix_pq_rs(isp2)%local_data(row_local, col_local)
    1685           54 :                               IF (ABS(erint) > eri_env%eps_integral .AND. (irange(1) - 1 + iwa12 <= iwb12 .OR. isp1 /= isp2)) THEN
    1686           30 :                                  icount2 = icount2 + 1
    1687           30 :                                  eri(icount2) = erint
    1688           30 :                                  eri_index(icount2) = iwb12
    1689              :                               END IF
    1690              :                            END IF
    1691              :                         END DO
    1692              :                      END DO
    1693           30 :                      stored_integrals = stored_integrals + icount2
    1694              :                      !
    1695           30 :                      isp = (isp1 - 1)*isp2 + (isp2 - isp1 + 1)
    1696           30 :                      CALL update_csr_matrix(eri_env%eri(isp)%csr_mat, icount2, eri, eri_index, iwa12)
    1697              :                      !
    1698          150 :                      DEALLOCATE (eri, eri_index)
    1699              :                   END DO
    1700          361 :                ELSEIF (eri_env%method == eri_method_full_gpw) THEN
    1701          776 :                   DO isp2 = isp1, nspins
    1702          415 :                      CALL get_mo_set(mo_set=mos(isp2), nmo=nmo2)
    1703          415 :                      nx = (nmo2*(nmo2 + 1))/2
    1704         2075 :                      ALLOCATE (eri(nx), eri_index(nx))
    1705          415 :                      icount2 = 0
    1706          415 :                      iwbs = 1
    1707          415 :                      IF (isp1 == isp2) iwbs = i1
    1708          415 :                      isp = (isp1 - 1)*isp2 + (isp2 - isp1 + 1)
    1709         1486 :                      DO i3 = iwbs, SIZE(orbitals, 1)
    1710         1071 :                         iwb1 = orbitals(i3, isp2)
    1711         1071 :                         IF (eri_env%eri_gpw%store_wfn) THEN
    1712         1046 :                            wfn3 => wfn_a(iwb1, isp2)
    1713              :                         ELSE
    1714              :                            CALL calculate_wavefunction(fm_mo_coeff_as(isp1), iwb1, wfn3, rho_g, atomic_kind_set, &
    1715           25 :                                                        qs_kind_set, cell, dft_control, particle_set, pw_env_sub)
    1716              :                         END IF
    1717         1071 :                         CALL pw_zero(wfn_r)
    1718         1071 :                         CALL pw_multiply(wfn_r, rho_r, wfn3)
    1719         1071 :                         iwbt = i3
    1720         1071 :                         IF (isp1 == isp2 .AND. i1 == i3) iwbt = i2
    1721         3386 :                         DO i4 = iwbt, SIZE(orbitals, 1)
    1722         1900 :                            iwb2 = orbitals(i4, isp2)
    1723         1900 :                            IF (eri_env%eri_gpw%store_wfn) THEN
    1724         1870 :                               wfn4 => wfn_a(iwb2, isp2)
    1725              :                            ELSE
    1726              :                               CALL calculate_wavefunction(fm_mo_coeff_as(isp1), iwb2, wfn4, rho_g, atomic_kind_set, &
    1727           30 :                                                           qs_kind_set, cell, dft_control, particle_set, pw_env_sub)
    1728              :                            END IF
    1729              :                            ! We reduce the amount of communication by collecting the local sums first and sum globally later
    1730         1900 :                            erint = pw_integral_ab(wfn_r, wfn4, local_only=.TRUE.)
    1731         1900 :                            icount2 = icount2 + 1
    1732         1900 :                            eri(icount2) = erint
    1733         2971 :                            eri_index(icount2) = csr_idx_to_combined(iwb1, iwb2, nmo2)
    1734              :                         END DO
    1735              :                      END DO
    1736              :                      ! Now, we sum the integrals globally
    1737          415 :                      CALL eri_env%para_env_sub%sum(eri)
    1738              :                      ! and we reorder the integrals to prevent storing too small integrals
    1739          415 :                      intcount = 0
    1740          415 :                      icount2 = 0
    1741              :                      iwbs = 1
    1742              :                      IF (isp1 == isp2) iwbs = i1
    1743          415 :                      isp = (isp1 - 1)*isp2 + (isp2 - isp1 + 1)
    1744         1486 :                      DO i3 = iwbs, SIZE(orbitals, 1)
    1745         1071 :                         iwb1 = orbitals(i3, isp2)
    1746         1071 :                         iwbt = i3
    1747         1071 :                         IF (isp1 == isp2 .AND. i1 == i3) iwbt = i2
    1748         3386 :                         DO i4 = iwbt, SIZE(orbitals, 1)
    1749         1900 :                            iwb2 = orbitals(i4, isp2)
    1750         1900 :                            intcount = intcount + 1
    1751         1900 :                            erint = eri(intcount)
    1752         2971 :                            IF (ABS(erint) > eri_env%eps_integral) THEN
    1753         1526 :                               IF (MOD(intcount, eri_env%para_env_sub%num_pe) == eri_env%para_env_sub%mepos) THEN
    1754          765 :                                  icount2 = icount2 + 1
    1755          765 :                                  eri(icount2) = erint
    1756          765 :                                  eri_index(icount2) = eri_index(intcount)
    1757              :                               END IF
    1758              :                            END IF
    1759              :                         END DO
    1760              :                      END DO
    1761          415 :                      stored_integrals = stored_integrals + icount2
    1762              :                      !
    1763          415 :                      CALL update_csr_matrix(eri_env%eri(isp)%csr_mat, icount2, eri, eri_index, iwa12)
    1764              :                      !
    1765         1191 :                      DEALLOCATE (eri, eri_index)
    1766              :                   END DO
    1767              :                ELSE
    1768            0 :                   CPABORT("Unknown option")
    1769              :                END IF
    1770              :             END DO
    1771              :          END DO
    1772              :       END DO
    1773              : 
    1774           66 :       IF (print1 .AND. iw > 0) THEN
    1775            1 :          WRITE (iw, "(T2,'ERI_GPW|',' Number of Integrals stored locally',T71,I10)") stored_integrals
    1776              :       END IF
    1777              : 
    1778           66 :       IF (eri_env%eri_gpw%store_wfn) THEN
    1779          120 :          DO ispin = 1, nspins
    1780          300 :             DO i1 = 1, SIZE(orbitals, 1)
    1781          180 :                iwfn = orbitals(i1, ispin)
    1782          246 :                CALL wfn_a(iwfn, ispin)%release()
    1783              :             END DO
    1784              :          END DO
    1785           54 :          DEALLOCATE (wfn_a)
    1786              :       ELSE
    1787           12 :          CALL wfn1%release()
    1788           12 :          CALL wfn2%release()
    1789           12 :          DEALLOCATE (wfn1, wfn2)
    1790           12 :          IF (eri_env%method /= eri_method_gpw_ht) THEN
    1791            6 :             CALL wfn3%release()
    1792            6 :             CALL wfn4%release()
    1793            6 :             DEALLOCATE (wfn3, wfn4)
    1794              :          END IF
    1795              :       END IF
    1796           66 :       CALL auxbas_pw_pool%give_back_pw(wfn_r)
    1797           66 :       CALL auxbas_pw_pool%give_back_pw(rho_g)
    1798           66 :       CALL auxbas_pw_pool%give_back_pw(rho_r)
    1799           66 :       CALL auxbas_pw_pool%give_back_pw(pot_g)
    1800              : 
    1801           66 :       IF (eri_env%method == eri_method_gpw_ht) THEN
    1802           24 :          DO ispin = 1, nspins
    1803           12 :             CALL dbcsr_release(mo_coeff_as(ispin))
    1804           12 :             CALL dbcsr_release(matrix_pq_rnu(ispin))
    1805           12 :             CALL cp_fm_release(fm_matrix_pq_rnu(ispin))
    1806           24 :             CALL cp_fm_release(fm_matrix_pq_rs(ispin))
    1807              :          END DO
    1808           12 :          DEALLOCATE (matrix_pq_rnu, fm_matrix_pq_rnu, fm_matrix_pq_rs)
    1809           12 :          CALL deallocate_task_list(task_list_sub)
    1810              :       END IF
    1811          144 :       DO ispin = 1, nspins
    1812           78 :          CALL dbcsr_release(mo_coeff_as(ispin))
    1813          144 :          CALL cp_fm_release(fm_mo_coeff_as(ispin))
    1814              :       END DO
    1815           66 :       DEALLOCATE (mo_coeff_as, fm_mo_coeff_as)
    1816           66 :       CALL release_neighbor_list_sets(sab_orb_sub)
    1817           66 :       CALL cp_blacs_env_release(blacs_env_sub)
    1818           66 :       CALL dbcsr_release(mat_munu%matrix)
    1819           66 :       DEALLOCATE (mat_munu%matrix)
    1820           66 :       CALL pw_env_release(pw_env_sub)
    1821              :       ! Return to the old qs_control
    1822           66 :       dft_control%qs_control => qs_control_old
    1823           66 :       DEALLOCATE (qs_control%e_cutoff)
    1824           66 :       DEALLOCATE (qs_control)
    1825              : 
    1826              :       ! print out progress
    1827           66 :       IF (iw > 0) THEN
    1828           33 :          t2 = m_walltime()
    1829           33 :          WRITE (iw, '(/,T2,A,T66,F14.2)') "ERI_GPW| ERI calculation took (sec)", t2 - t1
    1830           33 :          CALL m_flush(iw)
    1831              :       END IF
    1832              : 
    1833           66 :       CALL timestop(handle)
    1834              : 
    1835          132 :    END SUBROUTINE calculate_eri_gpw
    1836              : 
    1837              : ! **************************************************************************************************
    1838              : !> \brief Sets the Green's function for the ERI calculation. Here we deal with the G=0 case!
    1839              : !> \param green ...
    1840              : !> \param eri_env ...
    1841              : !> \par History
    1842              : !>      04.2016 created [JGH]
    1843              : !>      08.2025 added support for the LR truncation [SB]
    1844              : ! **************************************************************************************************
    1845           66 :    SUBROUTINE pw_eri_green_create(green, eri_env)
    1846              : 
    1847              :       TYPE(greens_fn_type), INTENT(INOUT)                :: green
    1848              :       TYPE(eri_type)                                     :: eri_env
    1849              : 
    1850              :       COMPLEX(KIND=dp)                                   :: erf_fac_p, z_p
    1851              :       INTEGER                                            :: ig
    1852              :       REAL(KIND=dp)                                      :: cossin_fac, ea, erfcos_fac, exp_prefac, &
    1853              :                                                             g, G0, g2, g3d, ga, Ginf, omega, &
    1854              :                                                             omega2, Rc, Rc2
    1855              : 
    1856              :       ! initialize influence function
    1857              :       ASSOCIATE (gf => green%influence_fn, grid => green%influence_fn%pw_grid)
    1858           72 :          SELECT CASE (green%method)
    1859              :          CASE (PERIODIC3D)
    1860              : 
    1861           68 :             SELECT CASE (eri_env%operator)
    1862              :             CASE (eri_operator_coulomb)
    1863       262145 :                DO ig = grid%first_gne0, grid%ngpts_cut_local
    1864       262143 :                   g2 = grid%gsq(ig)
    1865       262145 :                   gf%array(ig) = fourpi/g2
    1866              :                END DO
    1867            2 :                IF (grid%have_g0) gf%array(1) = 0.0_dp
    1868              : 
    1869              :             CASE (eri_operator_yukawa)
    1870            0 :                CALL cp_warn(__LOCATION__, "Yukawa operator has not been tested")
    1871            0 :                omega2 = eri_env%omega**2
    1872            0 :                DO ig = grid%first_gne0, grid%ngpts_cut_local
    1873            0 :                   g2 = grid%gsq(ig)
    1874            0 :                   gf%array(ig) = fourpi/(omega2 + g2)
    1875              :                END DO
    1876            0 :                IF (grid%have_g0) gf%array(1) = fourpi/omega2
    1877              : 
    1878              :             CASE (eri_operator_erf)
    1879            0 :                omega2 = eri_env%omega**2
    1880            0 :                DO ig = grid%first_gne0, grid%ngpts_cut_local
    1881            0 :                   g2 = grid%gsq(ig)
    1882            0 :                   gf%array(ig) = fourpi/g2*EXP(-0.25_dp*g2/omega2)
    1883              :                END DO
    1884            0 :                IF (grid%have_g0) gf%array(1) = 0.0_dp
    1885              : 
    1886              :             CASE (eri_operator_erfc)
    1887            0 :                omega2 = eri_env%omega**2
    1888            0 :                DO ig = grid%first_gne0, grid%ngpts_cut_local
    1889            0 :                   g2 = grid%gsq(ig)
    1890            0 :                   gf%array(ig) = fourpi/g2*(1.0_dp - EXP(-0.25_dp*g2/omega2))
    1891              :                END DO
    1892            0 :                IF (grid%have_g0) gf%array(1) = pi/omega2
    1893              : 
    1894              :             CASE (eri_operator_trunc)
    1895            0 :                Rc = eri_env%cutoff_radius
    1896            0 :                DO ig = grid%first_gne0, grid%ngpts_cut_local
    1897            0 :                   g2 = grid%gsq(ig)
    1898            0 :                   g = SQRT(g2)
    1899              :                   ! Taylor expansion around zero
    1900            0 :                   IF (g*Rc >= 0.005_dp) THEN
    1901            0 :                      gf%array(ig) = fourpi/g2*(1.0_dp - COS(g*Rc))
    1902              :                   ELSE
    1903            0 :                      gf%array(ig) = fourpi/g2*(g*Rc)**2/2.0_dp*(1.0_dp - (g*Rc)**2/12.0_dp)
    1904              :                   END IF
    1905              :                END DO
    1906            0 :                IF (grid%have_g0) gf%array(1) = twopi*Rc**2
    1907              : 
    1908              :             CASE (eri_operator_lr_trunc)
    1909            4 :                omega = eri_env%omega
    1910            4 :                omega2 = omega**2
    1911            4 :                Rc = eri_env%cutoff_radius
    1912            4 :                Rc2 = Rc**2
    1913            4 :                G0 = 0.001_dp ! threshold for the G=0 case
    1914            4 :                Ginf = 20.0_dp  ! threshold for the Taylor exapnsion arounf G=∞
    1915       843752 :                DO ig = grid%first_gne0, grid%ngpts_cut_local
    1916       843748 :                   g2 = grid%gsq(ig)
    1917       843748 :                   g = SQRT(g2)
    1918       843752 :                   IF (g <= 2.0_dp*G0) THEN
    1919              :                      gf%array(ig) = -pi/omega2*erf(omega*Rc) &
    1920              :                                     + twopi*Rc2*erf(omega*Rc) &
    1921            0 :                                     + 2*rootpi*Rc*EXP(-omega2*Rc2)/omega
    1922       843748 :                   ELSE IF (g >= 2.0_dp*Ginf*omega) THEN
    1923              :                      ! exponential prefactor
    1924         1488 :                      exp_prefac = EXP(-omega2*Rc2)/(rootpi*(omega2*Rc2 + 0.25_dp*g2/omega2))
    1925              :                      ! cos sin factor
    1926         1488 :                      cossin_fac = omega*Rc*COS(g*Rc) - 0.5_dp*g/omega*SIN(g*Rc)
    1927              :                      ! real erf term with cosine
    1928         1488 :                      erfcos_fac = ERF(omega*Rc)*COS(g*Rc)
    1929              :                      ! Combine terms
    1930         1488 :                      gf%array(ig) = fourpi/g2*(-exp_prefac*cossin_fac - erfcos_fac)
    1931              :                   ELSE
    1932              :                      ! exponential prefactor
    1933       842260 :                      exp_prefac = twopi/g2*EXP(-0.25_dp*g2/omega2)
    1934              :                      ! Compute complex arguments for erf
    1935       842260 :                      z_p = CMPLX(omega*Rc, 0.5_dp*g/omega, kind=dp)
    1936              :                      ! Evaluate complex error functions
    1937       842260 :                      erf_fac_p = 2.0_dp*REAL(erfz_fast(z_p))
    1938              :                      ! Real erf term with cosine
    1939       842260 :                      erfcos_fac = fourpi/g2*ERF(omega*Rc)*COS(g*Rc)
    1940              :                      ! Combine terms
    1941       842260 :                      gf%array(ig) = exp_prefac*erf_fac_p - erfcos_fac
    1942              :                   END IF
    1943              :                END DO
    1944            4 :                IF (grid%have_g0) THEN
    1945              :                   gf%array(1) = -pi/omega2*ERF(omega*Rc) &
    1946              :                                 + twopi*Rc2*ERF(omega*Rc) &
    1947            2 :                                 + 2*rootpi*Rc*EXP(-omega2*Rc2)/omega
    1948              :                END IF
    1949              : 
    1950              :             CASE DEFAULT
    1951            6 :                CPABORT("Please specify a valid operator for the periodic Poisson solver")
    1952              :             END SELECT
    1953              : 
    1954              :             ! The analytic Poisson solver simply limits the domain of integration
    1955              :             ! of the Fourier transform to a sphere of radius Rc, rather than integrating
    1956              :             ! over all space (-∞,∞)
    1957              :          CASE (ANALYTIC0D)
    1958              : 
    1959          106 :             SELECT CASE (eri_env%operator)
    1960              :                ! This is identical to the truncated Coulomb operator integrated
    1961              :                ! over all space, when the truncation radius is equal to the radius of
    1962              :                ! the Poisson solver
    1963              :             CASE (eri_operator_coulomb, eri_operator_trunc)
    1964           46 :                IF (eri_env%operator == eri_operator_coulomb) THEN
    1965           46 :                   Rc = green%radius
    1966              :                ELSE
    1967            0 :                   Rc = eri_env%cutoff_radius
    1968              :                END IF
    1969     14545757 :                DO ig = grid%first_gne0, grid%ngpts_cut_local
    1970     14545711 :                   g2 = grid%gsq(ig)
    1971     14545711 :                   g = SQRT(g2)
    1972              :                   ! Taylor expansion around zero
    1973     14545757 :                   IF (g*Rc >= 0.005_dp) THEN
    1974     14545711 :                      gf%array(ig) = fourpi/g2*(1.0_dp - COS(g*Rc))
    1975              :                   ELSE
    1976            0 :                      gf%array(ig) = fourpi/g2*(g*Rc)**2/2.0_dp*(1.0_dp - (g*Rc)**2/12.0_dp)
    1977              :                   END IF
    1978              :                END DO
    1979           46 :                IF (grid%have_g0) gf%array(1) = twopi*Rc**2
    1980              : 
    1981              :                ! Not tested
    1982              :             CASE (eri_operator_yukawa)
    1983            0 :                CALL cp_warn(__LOCATION__, "Yukawa operator has not been tested")
    1984            0 :                Rc = green%radius
    1985            0 :                omega = eri_env%omega
    1986            0 :                ea = EXP(-omega*Rc)
    1987            0 :                DO ig = grid%first_gne0, grid%ngpts_cut_local
    1988            0 :                   g2 = grid%gsq(ig)
    1989            0 :                   g = SQRT(g2)
    1990            0 :                   g3d = fourpi/(omega**2 + g2)
    1991            0 :                   gf%array(ig) = g3d*(1.0_dp - ea*(COS(g*Rc) + omega/g*SIN(g*Rc)))
    1992              :                END DO
    1993            0 :                IF (grid%have_g0) gf%array(1) = fourpi/(omega**2)*(1.0_dp - ea*(1.0_dp + omega*Rc))
    1994              : 
    1995              :                ! Long-range Coulomb
    1996              :                ! TODO: this should be equivalent to LR truncated Coulomb from above!
    1997              :             CASE (eri_operator_erf, eri_operator_lr_trunc)
    1998           14 :                IF (eri_env%operator == eri_operator_erf) THEN
    1999           14 :                   Rc = green%radius
    2000              :                ELSE
    2001            0 :                   Rc = eri_env%cutoff_radius
    2002              :                END IF
    2003           14 :                omega2 = eri_env%omega**2
    2004      1512007 :                DO ig = grid%first_gne0, grid%ngpts_cut_local
    2005      1511993 :                   g2 = grid%gsq(ig)
    2006      1511993 :                   g = SQRT(g2)
    2007      1511993 :                   ga = -0.25_dp*g2/omega2
    2008      1512007 :                   gf%array(ig) = fourpi/g2*EXP(ga)*(1.0_dp - COS(g*Rc))
    2009              :                END DO
    2010           14 :                IF (grid%have_g0) gf%array(1) = twopi*Rc**2
    2011              : 
    2012              :                ! Short-range Coulomb
    2013              :                ! TODO: this should actually be properly derived and see whether it is correct
    2014              :             CASE (eri_operator_erfc)
    2015              :                CALL cp_warn(__LOCATION__, &
    2016            0 :                             "Short-range Coulomb operator may be incorrect with ANALYTIC0D Poisson solver")
    2017            0 :                Rc = green%radius
    2018            0 :                omega2 = eri_env%omega**2
    2019            0 :                DO ig = grid%first_gne0, grid%ngpts_cut_local
    2020            0 :                   g2 = grid%gsq(ig)
    2021            0 :                   g = SQRT(g2)
    2022            0 :                   ga = -0.25_dp*g2/omega2
    2023            0 :                   gf%array(ig) = fourpi/g2*(1.0_dp - EXP(ga))*(1.0_dp - COS(g*Rc))
    2024              :                END DO
    2025            0 :                IF (grid%have_g0) gf%array(1) = pi/omega2
    2026              : 
    2027              :             CASE DEFAULT
    2028           60 :                CPABORT("Unsupported operator")
    2029              :             END SELECT
    2030              : 
    2031              :          CASE DEFAULT
    2032           66 :             CPABORT("Unsupported Poisson solver")
    2033              :          END SELECT
    2034              :       END ASSOCIATE
    2035              : 
    2036           66 :    END SUBROUTINE pw_eri_green_create
    2037              : 
    2038              : ! **************************************************************************************************
    2039              : !> \brief Adds data for a new row to the csr matrix
    2040              : !> \param csr_mat ...
    2041              : !> \param nnz ...
    2042              : !> \param rdat ...
    2043              : !> \param rind ...
    2044              : !> \param irow ...
    2045              : !> \par History
    2046              : !>      04.2016 created [JGH]
    2047              : ! **************************************************************************************************
    2048          445 :    SUBROUTINE update_csr_matrix(csr_mat, nnz, rdat, rind, irow)
    2049              : 
    2050              :       TYPE(dbcsr_csr_type), INTENT(INOUT)                :: csr_mat
    2051              :       INTEGER, INTENT(IN)                                :: nnz
    2052              :       REAL(KIND=dp), DIMENSION(:), INTENT(IN)            :: rdat
    2053              :       INTEGER, DIMENSION(:), INTENT(IN)                  :: rind
    2054              :       INTEGER, INTENT(IN)                                :: irow
    2055              : 
    2056              :       INTEGER                                            :: k, nrow, nze, nze_new
    2057              : 
    2058          445 :       IF (irow /= 0) THEN
    2059          445 :          nze = csr_mat%nze_local
    2060          445 :          nze_new = nze + nnz
    2061              :          ! values
    2062          445 :          CALL reallocate(csr_mat%nzval_local%r_dp, 1, nze_new)
    2063         1240 :          csr_mat%nzval_local%r_dp(nze + 1:nze_new) = rdat(1:nnz)
    2064              :          ! col indices
    2065          445 :          CALL reallocate(csr_mat%colind_local, 1, nze_new)
    2066         1240 :          csr_mat%colind_local(nze + 1:nze_new) = rind(1:nnz)
    2067              :          ! rows
    2068          445 :          nrow = csr_mat%nrows_local
    2069          445 :          CALL reallocate(csr_mat%rowptr_local, 1, irow + 1)
    2070         1198 :          csr_mat%rowptr_local(nrow + 1:irow) = nze + 1
    2071          445 :          csr_mat%rowptr_local(irow + 1) = nze_new + 1
    2072              :          ! nzerow
    2073          445 :          CALL reallocate(csr_mat%nzerow_local, 1, irow)
    2074         1198 :          DO k = nrow + 1, irow
    2075         1198 :             csr_mat%nzerow_local(k) = csr_mat%rowptr_local(k + 1) - csr_mat%rowptr_local(k)
    2076              :          END DO
    2077          445 :          csr_mat%nrows_local = irow
    2078          445 :          csr_mat%nze_local = csr_mat%nze_local + nnz
    2079              :       END IF
    2080          445 :       csr_mat%nze_total = csr_mat%nze_total + nnz
    2081          445 :       csr_mat%has_indices = .TRUE.
    2082              : 
    2083          445 :    END SUBROUTINE update_csr_matrix
    2084              : 
    2085              : ! **************************************************************************************************
    2086              : !> \brief Computes and prints the active orbitals on Cube Files
    2087              : !> \param input ...
    2088              : !> \param qs_env the qs_env in which the qs_env lives
    2089              : !> \param mos ...
    2090              : ! **************************************************************************************************
    2091            2 :    SUBROUTINE print_orbital_cubes(input, qs_env, mos)
    2092              :       TYPE(section_vals_type), POINTER                   :: input
    2093              :       TYPE(qs_environment_type), POINTER                 :: qs_env
    2094              :       TYPE(mo_set_type), DIMENSION(:), INTENT(IN)        :: mos
    2095              : 
    2096              :       CHARACTER(LEN=default_path_length)                 :: filebody, filename, title
    2097              :       INTEGER                                            :: i, imo, isp, nmo, str(3), unit_nr
    2098            2 :       INTEGER, DIMENSION(:), POINTER                     :: alist, blist, istride
    2099              :       LOGICAL                                            :: do_mo, explicit_a, explicit_b
    2100            2 :       TYPE(atomic_kind_type), DIMENSION(:), POINTER      :: atomic_kind_set
    2101              :       TYPE(cell_type), POINTER                           :: cell
    2102              :       TYPE(cp_fm_type), POINTER                          :: mo_coeff
    2103              :       TYPE(dft_control_type), POINTER                    :: dft_control
    2104              :       TYPE(mp_para_env_type), POINTER                    :: para_env
    2105              :       TYPE(particle_list_type), POINTER                  :: particles
    2106            2 :       TYPE(particle_type), DIMENSION(:), POINTER         :: particle_set
    2107              :       TYPE(pw_c1d_gs_type)                               :: wf_g
    2108              :       TYPE(pw_env_type), POINTER                         :: pw_env
    2109              :       TYPE(pw_pool_type), POINTER                        :: auxbas_pw_pool
    2110              :       TYPE(pw_r3d_rs_type)                               :: wf_r
    2111            2 :       TYPE(qs_kind_type), DIMENSION(:), POINTER          :: qs_kind_set
    2112              :       TYPE(qs_subsys_type), POINTER                      :: subsys
    2113              :       TYPE(section_vals_type), POINTER                   :: dft_section, scf_input
    2114              : 
    2115            2 :       CALL section_vals_val_get(input, "FILENAME", c_val=filebody)
    2116            2 :       CALL section_vals_val_get(input, "STRIDE", i_vals=istride)
    2117            2 :       IF (SIZE(istride) == 1) THEN
    2118            8 :          str(1:3) = istride(1)
    2119            0 :       ELSEIF (SIZE(istride) == 3) THEN
    2120            0 :          str(1:3) = istride(1:3)
    2121              :       ELSE
    2122            0 :          CPABORT("STRIDE arguments inconsistent")
    2123              :       END IF
    2124            2 :       CALL section_vals_val_get(input, "ALIST", i_vals=alist, explicit=explicit_a)
    2125            2 :       CALL section_vals_val_get(input, "BLIST", i_vals=blist, explicit=explicit_b)
    2126              : 
    2127              :       CALL get_qs_env(qs_env=qs_env, &
    2128              :                       dft_control=dft_control, &
    2129              :                       para_env=para_env, &
    2130              :                       subsys=subsys, &
    2131              :                       atomic_kind_set=atomic_kind_set, &
    2132              :                       qs_kind_set=qs_kind_set, &
    2133              :                       cell=cell, &
    2134              :                       particle_set=particle_set, &
    2135              :                       pw_env=pw_env, &
    2136            2 :                       input=scf_input)
    2137              : 
    2138            2 :       CALL qs_subsys_get(subsys, particles=particles)
    2139              :       !
    2140            2 :       CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
    2141            2 :       CALL auxbas_pw_pool%create_pw(wf_r)
    2142            2 :       CALL auxbas_pw_pool%create_pw(wf_g)
    2143              :       !
    2144            2 :       dft_section => section_vals_get_subs_vals(scf_input, "DFT")
    2145              :       !
    2146            4 :       DO isp = 1, SIZE(mos)
    2147            2 :          CALL get_mo_set(mo_set=mos(isp), mo_coeff=mo_coeff, nmo=nmo)
    2148              : 
    2149            2 :          IF (SIZE(mos) > 1) THEN
    2150            0 :             SELECT CASE (isp)
    2151              :             CASE (1)
    2152              :                CALL write_mo_set_to_output_unit(mos(isp), qs_kind_set, particle_set, &
    2153            0 :                                                 dft_section, 4, 0, final_mos=.TRUE., spin="ALPHA")
    2154              :             CASE (2)
    2155              :                CALL write_mo_set_to_output_unit(mos(isp), qs_kind_set, particle_set, &
    2156            0 :                                                 dft_section, 4, 0, final_mos=.TRUE., spin="BETA")
    2157              :             CASE DEFAULT
    2158            0 :                CPABORT("Invalid spin")
    2159              :             END SELECT
    2160              :          ELSE
    2161              :             CALL write_mo_set_to_output_unit(mos(isp), qs_kind_set, particle_set, &
    2162            2 :                                              dft_section, 4, 0, final_mos=.TRUE.)
    2163              :          END IF
    2164              : 
    2165           22 :          DO imo = 1, nmo
    2166           16 :             IF (isp == 1 .AND. explicit_a) THEN
    2167           16 :                IF (alist(1) == -1) THEN
    2168              :                   do_mo = .TRUE.
    2169              :                ELSE
    2170           16 :                   do_mo = .FALSE.
    2171           64 :                   DO i = 1, SIZE(alist)
    2172           64 :                      IF (imo == alist(i)) do_mo = .TRUE.
    2173              :                   END DO
    2174              :                END IF
    2175            0 :             ELSE IF (isp == 2 .AND. explicit_b) THEN
    2176            0 :                IF (blist(1) == -1) THEN
    2177              :                   do_mo = .TRUE.
    2178              :                ELSE
    2179            0 :                   do_mo = .FALSE.
    2180            0 :                   DO i = 1, SIZE(blist)
    2181            0 :                      IF (imo == blist(i)) do_mo = .TRUE.
    2182              :                   END DO
    2183              :                END IF
    2184              :             ELSE
    2185              :                do_mo = .TRUE.
    2186              :             END IF
    2187           16 :             IF (.NOT. do_mo) CYCLE
    2188              :             CALL calculate_wavefunction(mo_coeff, imo, wf_r, wf_g, atomic_kind_set, &
    2189            6 :                                         qs_kind_set, cell, dft_control, particle_set, pw_env)
    2190            6 :             IF (para_env%is_source()) THEN
    2191            3 :                WRITE (filename, '(A,A1,I4.4,A1,I1.1,A)') TRIM(filebody), "_", imo, "_", isp, ".cube"
    2192            3 :                CALL open_file(filename, unit_number=unit_nr, file_status="UNKNOWN", file_action="WRITE")
    2193            3 :                WRITE (title, *) "Active Orbital ", imo, " spin ", isp
    2194              :             ELSE
    2195            3 :                unit_nr = -1
    2196              :             END IF
    2197            6 :             CALL cp_pw_to_cube(wf_r, unit_nr, title, particles=particles, stride=istride)
    2198            8 :             IF (para_env%is_source()) THEN
    2199           13 :                CALL close_file(unit_nr)
    2200              :             END IF
    2201              :          END DO
    2202              :       END DO
    2203              : 
    2204            2 :       CALL auxbas_pw_pool%give_back_pw(wf_r)
    2205            2 :       CALL auxbas_pw_pool%give_back_pw(wf_g)
    2206              : 
    2207            2 :    END SUBROUTINE print_orbital_cubes
    2208              : 
    2209              : ! **************************************************************************************************
    2210              : !> \brief Writes a FCIDUMP file
    2211              : !> \param active_space_env ...
    2212              : !> \param as_input ...
    2213              : !> \par History
    2214              : !>      04.2016 created [JGH]
    2215              : ! **************************************************************************************************
    2216           66 :    SUBROUTINE fcidump(active_space_env, as_input)
    2217              : 
    2218              :       TYPE(active_space_type), POINTER                   :: active_space_env
    2219              :       TYPE(section_vals_type), POINTER                   :: as_input
    2220              : 
    2221              :       INTEGER                                            :: i, i1, i2, i3, i4, isym, iw, m1, m2, &
    2222              :                                                             nmo, norb, nspins
    2223              :       REAL(KIND=dp)                                      :: checksum, esub
    2224           66 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :)        :: fmat
    2225              :       TYPE(cp_logger_type), POINTER                      :: logger
    2226              :       TYPE(eri_fcidump_checksum)                         :: eri_checksum
    2227              : 
    2228           66 :       checksum = 0.0_dp
    2229              : 
    2230          132 :       logger => cp_get_default_logger()
    2231              :       iw = cp_print_key_unit_nr(logger, as_input, "FCIDUMP", &
    2232           66 :                                 extension=".fcidump", file_status="REPLACE", file_action="WRITE", file_form="FORMATTED")
    2233              :       !
    2234           66 :       nspins = active_space_env%nspins
    2235           66 :       norb = SIZE(active_space_env%active_orbitals, 1)
    2236           66 :       IF (nspins == 1) THEN
    2237              :          ASSOCIATE (ms2 => active_space_env%multiplicity, &
    2238              :                     nelec => active_space_env%nelec_active)
    2239              : 
    2240           54 :             IF (iw > 0) THEN
    2241           27 :                WRITE (iw, "(A,A,I4,A,I4,A,I2,A)") "&FCI", " NORB=", norb, ",NELEC=", nelec, ",MS2=", ms2, ","
    2242           27 :                isym = 1
    2243           99 :                WRITE (iw, "(A,1000(I1,','))") "  ORBSYM=", (isym, i=1, norb)
    2244           27 :                isym = 0
    2245           27 :                WRITE (iw, "(A,I1,A)") "  ISYM=", isym, ","
    2246           27 :                WRITE (iw, "(A)") " /"
    2247              :             END IF
    2248              :             !
    2249              :             ! Print integrals: ERI
    2250              :             CALL active_space_env%eri%eri_foreach(1, active_space_env%active_orbitals, &
    2251           54 :                                                   eri_fcidump_print(iw, 1, 1), 1, 1)
    2252           54 :             CALL eri_checksum%set(1, 1)
    2253           54 :             CALL active_space_env%eri%eri_foreach(1, active_space_env%active_orbitals, eri_checksum, 1, 1)
    2254              : 
    2255              :             ! Print integrals: Fij
    2256              :             ! replicate Fock matrix
    2257           54 :             nmo = active_space_env%eri%norb
    2258          216 :             ALLOCATE (fmat(nmo, nmo))
    2259           54 :             CALL replicate_and_symmetrize_matrix(nmo, active_space_env%fock_sub(1), fmat)
    2260           54 :             IF (iw > 0) THEN
    2261           27 :                i3 = 0; i4 = 0
    2262           99 :                DO m1 = 1, SIZE(active_space_env%active_orbitals, 1)
    2263           72 :                   i1 = active_space_env%active_orbitals(m1, 1)
    2264          245 :                   DO m2 = m1, SIZE(active_space_env%active_orbitals, 1)
    2265          146 :                      i2 = active_space_env%active_orbitals(m2, 1)
    2266          146 :                      checksum = checksum + ABS(fmat(i1, i2))
    2267          218 :                      WRITE (iw, "(ES23.16,4I4)") fmat(i1, i2), m1, m2, i3, i4
    2268              :                   END DO
    2269              :                END DO
    2270              :             END IF
    2271           54 :             DEALLOCATE (fmat)
    2272              :             ! Print energy
    2273           54 :             esub = active_space_env%energy_inactive
    2274           54 :             i1 = 0; i2 = 0; i3 = 0; i4 = 0
    2275           54 :             checksum = checksum + ABS(esub)
    2276          108 :             IF (iw > 0) WRITE (iw, "(ES23.16,4I4)") esub, i1, i2, i3, i4
    2277              :          END ASSOCIATE
    2278              : 
    2279              :       ELSE
    2280              :          ASSOCIATE (ms2 => active_space_env%multiplicity, &
    2281              :                     nelec => active_space_env%nelec_active)
    2282              : 
    2283           12 :             IF (iw > 0) THEN
    2284            6 :                WRITE (iw, "(A,A,I4,A,I4,A,I2,A)") "&FCI", " NORB=", norb, ",NELEC=", nelec, ",MS2=", ms2, ","
    2285            6 :                isym = 1
    2286           21 :                WRITE (iw, "(A,1000(I1,','))") "  ORBSYM=", (isym, i=1, norb)
    2287            6 :                isym = 0
    2288            6 :                WRITE (iw, "(A,I1,A)") "  ISYM=", isym, ","
    2289            6 :                WRITE (iw, "(A,I1,A)") "  UHF=", 1, ","
    2290            6 :                WRITE (iw, "(A)") " /"
    2291              :             END IF
    2292              :             !
    2293              :             ! Print integrals: ERI
    2294              :             ! alpha-alpha
    2295              :             CALL active_space_env%eri%eri_foreach(1, active_space_env%active_orbitals, &
    2296           12 :                                                   eri_fcidump_print(iw, 1, 1), 1, 1)
    2297           12 :             CALL eri_checksum%set(1, 1)
    2298           12 :             CALL active_space_env%eri%eri_foreach(1, active_space_env%active_orbitals, eri_checksum, 1, 1)
    2299              :             ! alpha-beta
    2300              :             CALL active_space_env%eri%eri_foreach(2, active_space_env%active_orbitals, &
    2301           12 :                                                   eri_fcidump_print(iw, 1, norb + 1), 1, 2)
    2302           12 :             CALL eri_checksum%set(1, norb + 1)
    2303           12 :             CALL active_space_env%eri%eri_foreach(2, active_space_env%active_orbitals, eri_checksum, 1, 2)
    2304              :             ! beta-beta
    2305              :             CALL active_space_env%eri%eri_foreach(3, active_space_env%active_orbitals, &
    2306           12 :                                                   eri_fcidump_print(iw, norb + 1, norb + 1), 2, 2)
    2307           12 :             CALL eri_checksum%set(norb + 1, norb + 1)
    2308           12 :             CALL active_space_env%eri%eri_foreach(3, active_space_env%active_orbitals, eri_checksum, 2, 2)
    2309              :             ! Print integrals: Fij
    2310              :             ! alpha
    2311           12 :             nmo = active_space_env%eri%norb
    2312           48 :             ALLOCATE (fmat(nmo, nmo))
    2313           12 :             CALL replicate_and_symmetrize_matrix(nmo, active_space_env%fock_sub(1), fmat)
    2314           12 :             IF (iw > 0) THEN
    2315            6 :                i3 = 0; i4 = 0
    2316           21 :                DO m1 = 1, norb
    2317           15 :                   i1 = active_space_env%active_orbitals(m1, 1)
    2318           48 :                   DO m2 = m1, norb
    2319           27 :                      i2 = active_space_env%active_orbitals(m2, 1)
    2320           27 :                      checksum = checksum + ABS(fmat(i1, i2))
    2321           42 :                      WRITE (iw, "(ES23.16,4I4)") fmat(i1, i2), m1, m2, i3, i4
    2322              :                   END DO
    2323              :                END DO
    2324              :             END IF
    2325           12 :             DEALLOCATE (fmat)
    2326              :             ! beta
    2327           48 :             ALLOCATE (fmat(nmo, nmo))
    2328           12 :             CALL replicate_and_symmetrize_matrix(nmo, active_space_env%fock_sub(2), fmat)
    2329           12 :             IF (iw > 0) THEN
    2330            6 :                i3 = 0; i4 = 0
    2331           21 :                DO m1 = 1, SIZE(active_space_env%active_orbitals, 1)
    2332           15 :                   i1 = active_space_env%active_orbitals(m1, 2)
    2333           48 :                   DO m2 = m1, SIZE(active_space_env%active_orbitals, 1)
    2334           27 :                      i2 = active_space_env%active_orbitals(m2, 2)
    2335           27 :                      checksum = checksum + ABS(fmat(i1, i2))
    2336           42 :                      WRITE (iw, "(ES23.16,4I4)") fmat(i1, i2), m1 + norb, m2 + norb, i3, i4
    2337              :                   END DO
    2338              :                END DO
    2339              :             END IF
    2340           12 :             DEALLOCATE (fmat)
    2341              :             ! Print energy
    2342           12 :             esub = active_space_env%energy_inactive
    2343           12 :             i1 = 0; i2 = 0; i3 = 0; i4 = 0
    2344           12 :             checksum = checksum + ABS(esub)
    2345           24 :             IF (iw > 0) WRITE (iw, "(ES23.16,4I4)") esub, i1, i2, i3, i4
    2346              :          END ASSOCIATE
    2347              :       END IF
    2348              :       !
    2349           66 :       CALL cp_print_key_finished_output(iw, logger, as_input, "FCIDUMP")
    2350              : 
    2351              :       !>>
    2352           66 :       iw = cp_logger_get_default_io_unit(logger)
    2353           66 :       IF (iw > 0) WRITE (iw, '(T4,A,T66,F12.8)') "FCIDUMP| Checksum:", eri_checksum%checksum + checksum
    2354              :       !<<
    2355              : 
    2356          132 :    END SUBROUTINE fcidump
    2357              : 
    2358              : ! **************************************************************************************************
    2359              : !> \brief replicate and symmetrize a matrix
    2360              : !> \param norb the number of orbitals
    2361              : !> \param distributed_matrix ...
    2362              : !> \param replicated_matrix ...
    2363              : ! **************************************************************************************************
    2364          234 :    SUBROUTINE replicate_and_symmetrize_matrix(norb, distributed_matrix, replicated_matrix)
    2365              :       INTEGER, INTENT(IN)                                :: norb
    2366              :       TYPE(cp_fm_type), INTENT(IN)                       :: distributed_matrix
    2367              :       REAL(dp), DIMENSION(:, :), INTENT(INOUT)           :: replicated_matrix
    2368              : 
    2369              :       INTEGER                                            :: i1, i2
    2370              :       REAL(dp)                                           :: mval
    2371              : 
    2372         4086 :       replicated_matrix(:, :) = 0.0_dp
    2373         1008 :       DO i1 = 1, norb
    2374         2934 :          DO i2 = i1, norb
    2375         1926 :             CALL cp_fm_get_element(distributed_matrix, i1, i2, mval)
    2376         1926 :             replicated_matrix(i1, i2) = mval
    2377         2700 :             replicated_matrix(i2, i1) = mval
    2378              :          END DO
    2379              :       END DO
    2380          234 :    END SUBROUTINE replicate_and_symmetrize_matrix
    2381              : 
    2382              : ! **************************************************************************************************
    2383              : !> \brief Calculates active space Fock matrix and inactive energy
    2384              : !> \param active_space_env ...
    2385              : !> \par History
    2386              : !>      06.2016 created [JGH]
    2387              : ! **************************************************************************************************
    2388           66 :    SUBROUTINE subspace_fock_matrix(active_space_env)
    2389              : 
    2390              :       TYPE(active_space_type), POINTER                   :: active_space_env
    2391              : 
    2392              :       INTEGER                                            :: i1, i2, is, norb, nspins
    2393              :       REAL(KIND=dp)                                      :: eeri, eref, esub, mval
    2394           66 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :)        :: ks_a_mat, ks_a_ref, ks_b_mat, ks_b_ref, &
    2395           66 :                                                             ks_mat, ks_ref, p_a_mat, p_b_mat, p_mat
    2396              :       TYPE(cp_fm_type), POINTER                          :: matrix, mo_coef
    2397              :       TYPE(dbcsr_csr_type), POINTER                      :: eri, eri_aa, eri_ab, eri_bb
    2398              : 
    2399           66 :       eref = active_space_env%energy_ref
    2400           66 :       nspins = active_space_env%nspins
    2401              : 
    2402           66 :       IF (nspins == 1) THEN
    2403           54 :          CALL get_mo_set(active_space_env%mos_active(1), nmo=norb, mo_coeff=mo_coef)
    2404              :          !
    2405              :          ! Loop over ERI, calculate subspace HF energy and Fock matrix
    2406              :          !
    2407              :          ! replicate KS, Core, and P matrices
    2408          540 :          ALLOCATE (ks_mat(norb, norb), ks_ref(norb, norb), p_mat(norb, norb))
    2409          810 :          ks_ref = 0.0_dp
    2410              : 
    2411              :          ! ks_mat contains the KS/Fock matrix (of full density) projected onto the AS MO subspace (f_ref in eq. 19)
    2412           54 :          CALL replicate_and_symmetrize_matrix(norb, active_space_env%ks_sub(1), ks_mat)
    2413           54 :          CALL replicate_and_symmetrize_matrix(norb, active_space_env%p_active(1), p_mat)
    2414              : 
    2415              :          ! compute ks_ref = V_H[rho^A] + V_HFX[rho^A]
    2416           54 :          eri => active_space_env%eri%eri(1)%csr_mat
    2417              :          CALL build_subspace_fock_matrix(active_space_env%active_orbitals, eri, p_mat, ks_ref, &
    2418           54 :                                          active_space_env%eri%comm_exchange)
    2419              : 
    2420              :          ! compute eeri = E_H[rho^A] + E_HFX[rho^A] as
    2421              :          ! eeri = 1/2 * (SUM_pq (V_H[rho^A] + V_HFX[rho^A])_pq * D^A_pq)
    2422          810 :          eeri = 0.5_dp*SUM(ks_ref*p_mat)
    2423              : 
    2424              :          ! now calculate the inactive energy acoording to eq. 19, that is
    2425              :          ! esub = E^I = E_ref - f_ref .* D^A + E_H[rho^A] + E_HFX[rho^A]
    2426              :          ! where f^ref = ks_mat, which is the KS/Fock matrix in MO basis, transformed previously
    2427              :          ! and is equal to ks_mat = h^0 + V_core + V_H[rho] + V_HFX[rho]
    2428          810 :          esub = eref - SUM(ks_mat(1:norb, 1:norb)*p_mat(1:norb, 1:norb)) + eeri
    2429              : 
    2430              :          ! reuse ks_mat to store f^I = f^ref - (V_H[rho^A] + V_HFX[rho^A]) according to eq. 20
    2431          810 :          ks_mat(1:norb, 1:norb) = ks_mat(1:norb, 1:norb) - ks_ref(1:norb, 1:norb)
    2432              :          ! this is now the embedding potential for the AS calculation!
    2433              : 
    2434           54 :          active_space_env%energy_inactive = esub
    2435              : 
    2436           54 :          CALL cp_fm_release(active_space_env%fock_sub)
    2437          216 :          ALLOCATE (active_space_env%fock_sub(nspins))
    2438          108 :          DO is = 1, nspins
    2439           54 :             matrix => active_space_env%ks_sub(is)
    2440              :             CALL cp_fm_create(active_space_env%fock_sub(is), matrix%matrix_struct, &
    2441          108 :                               name="Active Fock operator")
    2442              :          END DO
    2443           54 :          matrix => active_space_env%fock_sub(1)
    2444          270 :          DO i1 = 1, norb
    2445          810 :             DO i2 = 1, norb
    2446          594 :                mval = ks_mat(i1, i2)
    2447          756 :                CALL cp_fm_set_element(matrix, i1, i2, mval)
    2448              :             END DO
    2449              :          END DO
    2450              :       ELSE
    2451              : 
    2452           12 :          CALL get_mo_set(active_space_env%mos_active(1), nmo=norb)
    2453              :          !
    2454              :          ! Loop over ERI, calculate subspace HF energy and Fock matrix
    2455              :          !
    2456              :          ! replicate KS, Core, and P matrices
    2457              :          ALLOCATE (ks_a_mat(norb, norb), ks_b_mat(norb, norb), &
    2458              :               &    ks_a_ref(norb, norb), ks_b_ref(norb, norb), &
    2459          228 :               &     p_a_mat(norb, norb), p_b_mat(norb, norb))
    2460          552 :          ks_a_ref(:, :) = 0.0_dp; ks_b_ref(:, :) = 0.0_dp
    2461              : 
    2462           12 :          CALL replicate_and_symmetrize_matrix(norb, active_space_env%p_active(1), p_a_mat)
    2463           12 :          CALL replicate_and_symmetrize_matrix(norb, active_space_env%p_active(2), p_b_mat)
    2464           12 :          CALL replicate_and_symmetrize_matrix(norb, active_space_env%ks_sub(1), ks_a_mat)
    2465           12 :          CALL replicate_and_symmetrize_matrix(norb, active_space_env%ks_sub(2), ks_b_mat)
    2466              :          !
    2467              :          !
    2468           12 :          eri_aa => active_space_env%eri%eri(1)%csr_mat
    2469           12 :          eri_ab => active_space_env%eri%eri(2)%csr_mat
    2470           12 :          eri_bb => active_space_env%eri%eri(3)%csr_mat
    2471              :          CALL build_subspace_spin_fock_matrix(active_space_env%active_orbitals, eri_aa, eri_ab, p_a_mat, p_b_mat, ks_a_ref, &
    2472           12 :                                               tr_mixed_eri=.FALSE., comm_exchange=active_space_env%eri%comm_exchange)
    2473              :          CALL build_subspace_spin_fock_matrix(active_space_env%active_orbitals, eri_bb, eri_ab, p_b_mat, p_a_mat, ks_b_ref, &
    2474           12 :                                               tr_mixed_eri=.TRUE., comm_exchange=active_space_env%eri%comm_exchange)
    2475              :          !
    2476              :          ! calculate energy
    2477           12 :          eeri = 0.0_dp
    2478          552 :          eeri = 0.5_dp*(SUM(ks_a_ref*p_a_mat) + SUM(ks_b_ref*p_b_mat))
    2479          552 :          esub = eref - SUM(ks_a_mat*p_a_mat) - SUM(ks_b_mat*p_b_mat) + eeri
    2480          276 :          ks_a_mat(:, :) = ks_a_mat(:, :) - ks_a_ref(:, :)
    2481          276 :          ks_b_mat(:, :) = ks_b_mat(:, :) - ks_b_ref(:, :)
    2482              :          !
    2483           12 :          active_space_env%energy_inactive = esub
    2484              :          !
    2485           12 :          CALL cp_fm_release(active_space_env%fock_sub)
    2486           60 :          ALLOCATE (active_space_env%fock_sub(nspins))
    2487           36 :          DO is = 1, nspins
    2488           24 :             matrix => active_space_env%ks_sub(is)
    2489              :             CALL cp_fm_create(active_space_env%fock_sub(is), matrix%matrix_struct, &
    2490           36 :                               name="Active Fock operator")
    2491              :          END DO
    2492              : 
    2493           12 :          matrix => active_space_env%fock_sub(1)
    2494           60 :          DO i1 = 1, norb
    2495          276 :             DO i2 = 1, norb
    2496          216 :                mval = ks_a_mat(i1, i2)
    2497          264 :                CALL cp_fm_set_element(matrix, i1, i2, mval)
    2498              :             END DO
    2499              :          END DO
    2500           12 :          matrix => active_space_env%fock_sub(2)
    2501           72 :          DO i1 = 1, norb
    2502          276 :             DO i2 = 1, norb
    2503          216 :                mval = ks_b_mat(i1, i2)
    2504          264 :                CALL cp_fm_set_element(matrix, i1, i2, mval)
    2505              :             END DO
    2506              :          END DO
    2507              : 
    2508              :       END IF
    2509              : 
    2510           66 :    END SUBROUTINE subspace_fock_matrix
    2511              : 
    2512              : ! **************************************************************************************************
    2513              : !> \brief build subspace fockian
    2514              : !> \param active_orbitals the active orbital indices
    2515              : !> \param eri two electon integrals in MO
    2516              : !> \param p_mat density matrix
    2517              : !> \param ks_ref fockian matrix
    2518              : !> \param comm_exchange ...
    2519              : ! **************************************************************************************************
    2520           54 :    SUBROUTINE build_subspace_fock_matrix(active_orbitals, eri, p_mat, ks_ref, comm_exchange)
    2521              :       INTEGER, DIMENSION(:, :), INTENT(IN)               :: active_orbitals
    2522              :       TYPE(dbcsr_csr_type), INTENT(IN)                   :: eri
    2523              :       REAL(dp), DIMENSION(:, :), INTENT(IN)              :: p_mat
    2524              :       REAL(dp), DIMENSION(:, :), INTENT(INOUT)           :: ks_ref
    2525              :       TYPE(mp_comm_type), INTENT(IN)                     :: comm_exchange
    2526              : 
    2527              :       CHARACTER(LEN=*), PARAMETER :: routineN = 'build_subspace_fock_matrix'
    2528              : 
    2529              :       INTEGER                                            :: handle, i1, i12, i12l, i2, i3, i34, &
    2530              :                                                             i34l, i4, irptr, m1, m2, nindex, &
    2531              :                                                             nmo_total, norb
    2532              :       INTEGER, DIMENSION(2)                              :: irange
    2533              :       REAL(dp)                                           :: erint
    2534              :       TYPE(mp_comm_type)                                 :: mp_group
    2535              : 
    2536           54 :       CALL timeset(routineN, handle)
    2537              : 
    2538              :       ! Nothing to do
    2539           54 :       norb = SIZE(active_orbitals, 1)
    2540           54 :       nmo_total = SIZE(p_mat, 1)
    2541           54 :       nindex = (nmo_total*(nmo_total + 1))/2
    2542           54 :       CALL mp_group%set_handle(eri%mp_group%get_handle())
    2543           54 :       irange = get_irange_csr(nindex, comm_exchange)
    2544          198 :       DO m1 = 1, norb
    2545          144 :          i1 = active_orbitals(m1, 1)
    2546          490 :          DO m2 = m1, norb
    2547          292 :             i2 = active_orbitals(m2, 1)
    2548          292 :             i12 = csr_idx_to_combined(i1, i2, nmo_total)
    2549          436 :             IF (i12 >= irange(1) .AND. i12 <= irange(2)) THEN
    2550          283 :                i12l = i12 - irange(1) + 1
    2551          283 :                irptr = eri%rowptr_local(i12l) - 1
    2552          806 :                DO i34l = 1, eri%nzerow_local(i12l)
    2553          523 :                   i34 = eri%colind_local(irptr + i34l)
    2554          523 :                   CALL csr_idx_from_combined(i34, nmo_total, i3, i4)
    2555          523 :                   erint = eri%nzval_local%r_dp(irptr + i34l)
    2556              :                   ! Coulomb
    2557          523 :                   ks_ref(i1, i2) = ks_ref(i1, i2) + erint*p_mat(i3, i4)
    2558          523 :                   IF (i3 /= i4) THEN
    2559          278 :                      ks_ref(i1, i2) = ks_ref(i1, i2) + erint*p_mat(i3, i4)
    2560              :                   END IF
    2561          523 :                   IF (i12 /= i34) THEN
    2562          377 :                      ks_ref(i3, i4) = ks_ref(i3, i4) + erint*p_mat(i1, i2)
    2563          377 :                      IF (i1 /= i2) THEN
    2564          244 :                         ks_ref(i3, i4) = ks_ref(i3, i4) + erint*p_mat(i1, i2)
    2565              :                      END IF
    2566              :                   END IF
    2567              :                   ! Exchange
    2568          523 :                   erint = -0.5_dp*erint
    2569          523 :                   ks_ref(i1, i3) = ks_ref(i1, i3) + erint*p_mat(i2, i4)
    2570          523 :                   IF (i1 /= i2) THEN
    2571          318 :                      ks_ref(i2, i3) = ks_ref(i2, i3) + erint*p_mat(i1, i4)
    2572              :                   END IF
    2573          523 :                   IF (i3 /= i4) THEN
    2574          278 :                      ks_ref(i1, i4) = ks_ref(i1, i4) + erint*p_mat(i2, i3)
    2575              :                   END IF
    2576         1329 :                   IF (i1 /= i2 .AND. i3 /= i4) THEN
    2577          219 :                      ks_ref(i2, i4) = ks_ref(i2, i4) + erint*p_mat(i1, i3)
    2578              :                   END IF
    2579              :                END DO
    2580              :             END IF
    2581              :          END DO
    2582              :       END DO
    2583              :       !
    2584          198 :       DO m1 = 1, norb
    2585          144 :          i1 = active_orbitals(m1, 1)
    2586          490 :          DO m2 = m1, norb
    2587          292 :             i2 = active_orbitals(m2, 1)
    2588          436 :             ks_ref(i2, i1) = ks_ref(i1, i2)
    2589              :          END DO
    2590              :       END DO
    2591         1566 :       CALL mp_group%sum(ks_ref)
    2592              : 
    2593           54 :       CALL timestop(handle)
    2594              : 
    2595           54 :    END SUBROUTINE build_subspace_fock_matrix
    2596              : 
    2597              : ! **************************************************************************************************
    2598              : !> \brief build subspace fockian for unrestricted spins
    2599              : !> \param active_orbitals the active orbital indices
    2600              : !> \param eri_aa two electon integrals in MO with parallel spins
    2601              : !> \param eri_ab two electon integrals in MO with anti-parallel spins
    2602              : !> \param p_a_mat density matrix for up-spin
    2603              : !> \param p_b_mat density matrix for down-spin
    2604              : !> \param ks_a_ref fockian matrix for up-spin
    2605              : !> \param tr_mixed_eri boolean to indicate Coulomb interaction alignment
    2606              : !> \param comm_exchange ...
    2607              : ! **************************************************************************************************
    2608           24 :    SUBROUTINE build_subspace_spin_fock_matrix(active_orbitals, eri_aa, eri_ab, p_a_mat, p_b_mat, ks_a_ref, tr_mixed_eri, &
    2609              :                                               comm_exchange)
    2610              :       INTEGER, DIMENSION(:, :), INTENT(IN)               :: active_orbitals
    2611              :       TYPE(dbcsr_csr_type), INTENT(IN)                   :: eri_aa, eri_ab
    2612              :       REAL(dp), DIMENSION(:, :), INTENT(IN)              :: p_a_mat, p_b_mat
    2613              :       REAL(dp), DIMENSION(:, :), INTENT(INOUT)           :: ks_a_ref
    2614              :       LOGICAL, INTENT(IN)                                :: tr_mixed_eri
    2615              :       TYPE(mp_comm_type), INTENT(IN)                     :: comm_exchange
    2616              : 
    2617              :       CHARACTER(LEN=*), PARAMETER :: routineN = 'build_subspace_spin_fock_matrix'
    2618              : 
    2619              :       INTEGER                                            :: handle, i1, i12, i12l, i2, i3, i34, &
    2620              :                                                             i34l, i4, irptr, m1, m2, nindex, &
    2621              :                                                             nmo_total, norb, spin1, spin2
    2622              :       INTEGER, DIMENSION(2)                              :: irange
    2623              :       REAL(dp)                                           :: erint
    2624              :       TYPE(mp_comm_type)                                 :: mp_group
    2625              : 
    2626           24 :       CALL timeset(routineN, handle)
    2627              : 
    2628           24 :       norb = SIZE(active_orbitals, 1)
    2629           24 :       nmo_total = SIZE(p_a_mat, 1)
    2630           24 :       nindex = (nmo_total*(nmo_total + 1))/2
    2631           24 :       irange = get_irange_csr(nindex, comm_exchange)
    2632           24 :       IF (tr_mixed_eri) THEN
    2633              :          spin1 = 2
    2634           24 :          spin2 = 1
    2635              :       ELSE
    2636           12 :          spin1 = 1
    2637           12 :          spin2 = 2
    2638              :       END IF
    2639           84 :       DO m1 = 1, norb
    2640           60 :          i1 = active_orbitals(m1, spin1)
    2641          192 :          DO m2 = m1, norb
    2642          108 :             i2 = active_orbitals(m2, spin1)
    2643          108 :             i12 = csr_idx_to_combined(i1, i2, nmo_total)
    2644          168 :             IF (i12 >= irange(1) .AND. i12 <= irange(2)) THEN
    2645          108 :                i12l = i12 - irange(1) + 1
    2646          108 :                irptr = eri_aa%rowptr_local(i12l) - 1
    2647          259 :                DO i34l = 1, eri_aa%nzerow_local(i12l)
    2648          151 :                   i34 = eri_aa%colind_local(irptr + i34l)
    2649          151 :                   CALL csr_idx_from_combined(i34, nmo_total, i3, i4)
    2650          151 :                   erint = eri_aa%nzval_local%r_dp(irptr + i34l)
    2651              :                   ! Coulomb
    2652              :                   !F_ij += (ij|kl)*d_kl
    2653          151 :                   ks_a_ref(i1, i2) = ks_a_ref(i1, i2) + erint*p_a_mat(i3, i4)
    2654          151 :                   IF (i12 /= i34) THEN
    2655              :                      !F_kl += (ij|kl)*d_ij
    2656           97 :                      ks_a_ref(i3, i4) = ks_a_ref(i3, i4) + erint*p_a_mat(i1, i2)
    2657              :                   END IF
    2658              :                   ! Exchange
    2659          151 :                   erint = -1.0_dp*erint
    2660              :                   !F_ik -= (ij|kl)*d_jl
    2661          151 :                   ks_a_ref(i1, i3) = ks_a_ref(i1, i3) + erint*p_a_mat(i2, i4)
    2662          151 :                   IF (i1 /= i2) THEN
    2663              :                      !F_jk -= (ij|kl)*d_il
    2664           71 :                      ks_a_ref(i2, i3) = ks_a_ref(i2, i3) + erint*p_a_mat(i1, i4)
    2665              :                   END IF
    2666          151 :                   IF (i3 /= i4) THEN
    2667              :                      !F_il -= (ij|kl)*d_jk
    2668           68 :                      ks_a_ref(i1, i4) = ks_a_ref(i1, i4) + erint*p_a_mat(i2, i3)
    2669              :                   END IF
    2670          410 :                   IF (i1 /= i2 .AND. i3 /= i4) THEN
    2671              :                      !F_jl -= (ij|kl)*d_ik
    2672           42 :                      ks_a_ref(i2, i4) = ks_a_ref(i2, i4) + erint*p_a_mat(i1, i3)
    2673              :                   END IF
    2674              :                END DO
    2675              :             END IF
    2676              :          END DO
    2677              :       END DO
    2678              :       !
    2679              : 
    2680           24 :       irange = get_irange_csr(nindex, comm_exchange)
    2681           84 :       DO m1 = 1, norb
    2682           60 :          i1 = active_orbitals(m1, 1)
    2683          192 :          DO m2 = m1, norb
    2684          108 :             i2 = active_orbitals(m2, 1)
    2685          108 :             i12 = csr_idx_to_combined(i1, i2, nmo_total)
    2686          168 :             IF (i12 >= irange(1) .AND. i12 <= irange(2)) THEN
    2687          108 :                i12l = i12 - irange(1) + 1
    2688          108 :                irptr = eri_ab%rowptr_local(i12l) - 1
    2689          350 :                DO i34l = 1, eri_ab%nzerow_local(i12l)
    2690          242 :                   i34 = eri_ab%colind_local(irptr + i34l)
    2691          242 :                   CALL csr_idx_from_combined(i34, nmo_total, i3, i4)
    2692          242 :                   erint = eri_ab%nzval_local%r_dp(irptr + i34l)
    2693              :                   ! Coulomb
    2694          350 :                   IF (tr_mixed_eri) THEN
    2695              :                      !F_kl += (kl beta|ij alpha )*d_alpha_ij
    2696          121 :                      ks_a_ref(i3, i4) = ks_a_ref(i3, i4) + erint*p_b_mat(i1, i2)
    2697              :                   ELSE
    2698              :                      !F_ij += (ij alpha|kl beta )*d_beta_kl
    2699          121 :                      ks_a_ref(i1, i2) = ks_a_ref(i1, i2) + erint*p_b_mat(i3, i4)
    2700              :                   END IF
    2701              :                END DO
    2702              :             END IF
    2703              :          END DO
    2704              :       END DO
    2705              :       !
    2706           84 :       DO m1 = 1, norb
    2707           60 :          i1 = active_orbitals(m1, spin1)
    2708          192 :          DO m2 = m1, norb
    2709          108 :             i2 = active_orbitals(m2, spin1)
    2710          168 :             ks_a_ref(i2, i1) = ks_a_ref(i1, i2)
    2711              :          END DO
    2712              :       END DO
    2713           24 :       CALL mp_group%set_handle(eri_aa%mp_group%get_handle())
    2714         1080 :       CALL mp_group%sum(ks_a_ref)
    2715              : 
    2716           24 :       CALL timestop(handle)
    2717              : 
    2718           24 :    END SUBROUTINE build_subspace_spin_fock_matrix
    2719              : 
    2720              : ! **************************************************************************************************
    2721              : !> \brief Creates a local basis
    2722              : !> \param pro_basis_set ...
    2723              : !> \param zval ...
    2724              : !> \param ishell ...
    2725              : !> \param nshell ...
    2726              : !> \param lnam ...
    2727              : !> \par History
    2728              : !>      05.2016 created [JGH]
    2729              : ! **************************************************************************************************
    2730            0 :    SUBROUTINE create_pro_basis(pro_basis_set, zval, ishell, nshell, lnam)
    2731              :       TYPE(gto_basis_set_type), POINTER                  :: pro_basis_set
    2732              :       INTEGER, INTENT(IN)                                :: zval, ishell
    2733              :       INTEGER, DIMENSION(:), INTENT(IN)                  :: nshell
    2734              :       CHARACTER(len=*), DIMENSION(:), INTENT(IN)         :: lnam
    2735              : 
    2736            0 :       CHARACTER(len=6), DIMENSION(:), POINTER            :: sym
    2737              :       INTEGER                                            :: i, l, nj
    2738              :       INTEGER, DIMENSION(4, 7)                           :: ne
    2739            0 :       INTEGER, DIMENSION(:), POINTER                     :: lq, nq
    2740            0 :       REAL(KIND=dp), DIMENSION(:), POINTER               :: zet
    2741              :       TYPE(sto_basis_set_type), POINTER                  :: sto_basis_set
    2742              : 
    2743            0 :       CPASSERT(.NOT. ASSOCIATED(pro_basis_set))
    2744            0 :       NULLIFY (sto_basis_set)
    2745              : 
    2746              :       ! electronic configuration
    2747            0 :       ne = 0
    2748            0 :       DO l = 1, 4 !lq(1)+1
    2749            0 :          nj = 2*(l - 1) + 1
    2750            0 :          DO i = l, 7 ! nq(1)
    2751            0 :             ne(l, i) = ptable(zval)%e_conv(l - 1) - 2*nj*(i - l)
    2752            0 :             ne(l, i) = MAX(ne(l, i), 0)
    2753            0 :             ne(l, i) = MIN(ne(l, i), 2*nj)
    2754              :          END DO
    2755              :       END DO
    2756            0 :       ALLOCATE (nq(ishell), lq(ishell), zet(ishell), sym(ishell))
    2757            0 :       DO i = 1, ishell
    2758            0 :          nq(i) = nshell(i)
    2759            0 :          SELECT CASE (lnam(i))
    2760              :          CASE ('S', 's')
    2761            0 :             lq(i) = 0
    2762              :          CASE ('P', 'p')
    2763            0 :             lq(i) = 1
    2764              :          CASE ('D', 'd')
    2765            0 :             lq(i) = 2
    2766              :          CASE ('F', 'f')
    2767            0 :             lq(i) = 3
    2768              :          CASE DEFAULT
    2769            0 :             CPABORT("Wrong l QN")
    2770              :          END SELECT
    2771            0 :          sym(i) = lnam(i)
    2772            0 :          zet(i) = srules(zval, ne, nq(1), lq(1))
    2773              :       END DO
    2774            0 :       CALL allocate_sto_basis_set(sto_basis_set)
    2775            0 :       CALL set_sto_basis_set(sto_basis_set, nshell=1, nq=nq, lq=lq, zet=zet, symbol=sym)
    2776            0 :       CALL create_gto_from_sto_basis(sto_basis_set, pro_basis_set, 6)
    2777            0 :       pro_basis_set%norm_type = 2
    2778            0 :       CALL init_orb_basis_set(pro_basis_set)
    2779            0 :       CALL deallocate_sto_basis_set(sto_basis_set)
    2780              : 
    2781            0 :    END SUBROUTINE create_pro_basis
    2782              : 
    2783              : ! **************************************************************************************************
    2784              : !> \brief Update the density matrix in AO basis with the active density contribution
    2785              : !> \param active_space_env the active space environment
    2786              : !> \param rho_ao the density matrix in AO basis
    2787              : ! **************************************************************************************************
    2788            0 :    SUBROUTINE update_density_ao(active_space_env, rho_ao)
    2789              :       TYPE(active_space_type), POINTER                   :: active_space_env
    2790              :       TYPE(dbcsr_p_type), DIMENSION(:), POINTER          :: rho_ao
    2791              : 
    2792              :       INTEGER                                            :: ispin, nao, nmo, nspins
    2793              :       TYPE(cp_fm_type)                                   :: R, U
    2794              :       TYPE(cp_fm_type), POINTER                          :: C_active, p_active_mo
    2795              :       TYPE(dbcsr_type), POINTER                          :: p_inactive_ao
    2796            0 :       TYPE(mo_set_type), DIMENSION(:), POINTER           :: mos_active
    2797              : 
    2798              :       ! Transform the AS density matrix P_MO to the atomic orbital basis,
    2799              :       ! this is simply C * P_MO * C^T
    2800            0 :       nspins = active_space_env%nspins
    2801            0 :       mos_active => active_space_env%mos_active
    2802            0 :       DO ispin = 1, nspins
    2803              :          ! size of p_inactive_ao is (nao x nao)
    2804            0 :          p_inactive_ao => active_space_env%pmat_inactive(ispin)%matrix
    2805              : 
    2806              :          ! copy p_inactive_ao to rho_ao
    2807            0 :          CALL dbcsr_copy(rho_ao(ispin)%matrix, p_inactive_ao)
    2808              : 
    2809              :          ! size of p_active_mo is (nmo x nmo)
    2810            0 :          p_active_mo => active_space_env%p_active(ispin)
    2811              : 
    2812              :          ! calculate R = p_mo
    2813            0 :          CALL cp_fm_create(R, p_active_mo%matrix_struct)
    2814            0 :          CALL cp_fm_to_fm(p_active_mo, R)
    2815              : 
    2816              :          ! calculate U = C * p_mo
    2817            0 :          CALL get_mo_set(mos_active(ispin), mo_coeff=C_active, nao=nao, nmo=nmo)
    2818            0 :          CALL cp_fm_create(U, C_active%matrix_struct)
    2819            0 :          CALL parallel_gemm("N", "N", nao, nmo, nmo, 1.0_dp, C_active, R, 0.0_dp, U)
    2820              : 
    2821              :          CALL cp_dbcsr_plus_fm_fm_t(sparse_matrix=rho_ao(ispin)%matrix, &
    2822            0 :                                     matrix_v=U, matrix_g=C_active, ncol=nmo, alpha=1.0_dp)
    2823              : 
    2824            0 :          CALL cp_fm_release(R)
    2825            0 :          CALL cp_fm_release(U)
    2826              :       END DO
    2827              : 
    2828            0 :    END SUBROUTINE update_density_ao
    2829              : 
    2830              : ! **************************************************************************************************
    2831              : !> \brief Print each value on the master node
    2832              : !> \param this object reference
    2833              : !> \param i i-index
    2834              : !> \param j j-index
    2835              : !> \param k k-index
    2836              : !> \param l l-index
    2837              : !> \param val value of the integral at (i,j,k.l)
    2838              : !> \return always true to dump all integrals
    2839              : ! **************************************************************************************************
    2840         1590 :    LOGICAL FUNCTION eri_fcidump_print_func(this, i, j, k, l, val) RESULT(cont)
    2841              :       CLASS(eri_fcidump_print), INTENT(inout) :: this
    2842              :       INTEGER, INTENT(in)                     :: i, j, k, l
    2843              :       REAL(KIND=dp), INTENT(in)               :: val
    2844              : 
    2845              :       ! write to the actual file only on the master
    2846         1590 :       IF (this%unit_nr > 0) THEN
    2847          795 :          WRITE (this%unit_nr, "(ES23.16,4I4)") val, i + this%bra_start - 1, j + this%bra_start - 1, &
    2848         1590 :               &                                     k + this%ket_start - 1, l + this%ket_start - 1
    2849              :       END IF
    2850              : 
    2851         1590 :       cont = .TRUE.
    2852         1590 :    END FUNCTION eri_fcidump_print_func
    2853              : 
    2854              : ! **************************************************************************************************
    2855              : !> \brief checksum each value on the master node
    2856              : !> \param this object reference
    2857              : !> \param i i-index
    2858              : !> \param j j-index
    2859              : !> \param k k-index
    2860              : !> \param l l-index
    2861              : !> \param val value of the integral at (i,j,k.l)
    2862              : !> \return always true to dump all integrals
    2863              : ! **************************************************************************************************
    2864         1590 :    LOGICAL FUNCTION eri_fcidump_checksum_func(this, i, j, k, l, val) RESULT(cont)
    2865              :       CLASS(eri_fcidump_checksum), INTENT(inout) :: this
    2866              :       INTEGER, INTENT(in)                     :: i, j, k, l
    2867              :       REAL(KIND=dp), INTENT(in)               :: val
    2868              :       MARK_USED(i)
    2869              :       MARK_USED(j)
    2870              :       MARK_USED(k)
    2871              :       MARK_USED(l)
    2872              : 
    2873         1590 :       this%checksum = this%checksum + ABS(val)
    2874              : 
    2875         1590 :       cont = .TRUE.
    2876         1590 :    END FUNCTION eri_fcidump_checksum_func
    2877              : 
    2878              : ! **************************************************************************************************
    2879              : !> \brief Update active space density matrix from a fortran array
    2880              : !> \param p_act_mo density matrix in active space MO basis
    2881              : !> \param active_space_env active space environment
    2882              : !> \param ispin spin index
    2883              : !> \author Vladimir Rybkin
    2884              : ! **************************************************************************************************
    2885            0 :    SUBROUTINE update_active_density(p_act_mo, active_space_env, ispin)
    2886              :       REAL(KIND=dp), DIMENSION(:)                        :: p_act_mo
    2887              :       TYPE(active_space_type), POINTER                   :: active_space_env
    2888              :       INTEGER                                            :: ispin
    2889              : 
    2890              :       INTEGER                                            :: i1, i2, m1, m2, nmo_active
    2891              :       REAL(KIND=dp)                                      :: alpha, pij_new, pij_old
    2892              :       TYPE(cp_fm_type), POINTER                          :: p_active
    2893              : 
    2894            0 :       p_active => active_space_env%p_active(ispin)
    2895            0 :       nmo_active = active_space_env%nmo_active
    2896            0 :       alpha = active_space_env%alpha
    2897              : 
    2898            0 :       DO i1 = 1, nmo_active
    2899            0 :          m1 = active_space_env%active_orbitals(i1, ispin)
    2900            0 :          DO i2 = 1, nmo_active
    2901            0 :             m2 = active_space_env%active_orbitals(i2, ispin)
    2902            0 :             CALL cp_fm_get_element(p_active, m1, m2, pij_old)
    2903            0 :             pij_new = p_act_mo(i2 + (i1 - 1)*nmo_active)
    2904            0 :             pij_new = alpha*pij_new + (1.0_dp - alpha)*pij_old
    2905            0 :             CALL cp_fm_set_element(p_active, m1, m2, pij_new)
    2906              :          END DO
    2907              :       END DO
    2908              : 
    2909            0 :    END SUBROUTINE update_active_density
    2910              : 
    2911              : ! **************************************************************************************************
    2912              : !> \brief Compute and print the AS rdm and the natural orbitals occupation numbers
    2913              : !> \param active_space_env active space environment
    2914              : !> \param iw output unit
    2915              : !> \author Stefano Battaglia
    2916              : ! **************************************************************************************************
    2917            0 :    SUBROUTINE print_pmat_noon(active_space_env, iw)
    2918              :       TYPE(active_space_type), POINTER                   :: active_space_env
    2919              :       INTEGER                                            :: iw
    2920              : 
    2921              :       INTEGER                                            :: i1, i2, ii, ispin, jm, m1, m2, &
    2922              :                                                             nmo_active, nspins
    2923            0 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :)        :: noon, pmat
    2924              :       TYPE(cp_fm_type), POINTER                          :: p_active
    2925              : 
    2926            0 :       nspins = active_space_env%nspins
    2927            0 :       nmo_active = active_space_env%nmo_active
    2928              : 
    2929            0 :       ALLOCATE (noon(nmo_active, nspins))
    2930            0 :       ALLOCATE (pmat(nmo_active, nmo_active))
    2931              : 
    2932            0 :       DO ispin = 1, nspins
    2933            0 :          p_active => active_space_env%p_active(ispin)
    2934            0 :          noon(:, ispin) = 0.0_dp
    2935            0 :          pmat = 0.0_dp
    2936              : 
    2937            0 :          DO i1 = 1, nmo_active
    2938            0 :             m1 = active_space_env%active_orbitals(i1, ispin)
    2939            0 :             DO i2 = 1, nmo_active
    2940            0 :                m2 = active_space_env%active_orbitals(i2, ispin)
    2941            0 :                CALL cp_fm_get_element(p_active, m1, m2, pmat(i1, i2))
    2942              :             END DO
    2943              :          END DO
    2944              : 
    2945            0 :          IF (iw > 0) THEN
    2946            0 :             WRITE (iw, '(/,T3,A,I2,A)') "Active space density matrix for spin ", ispin
    2947            0 :             DO i1 = 1, nmo_active
    2948            0 :                DO ii = 1, nmo_active, 8
    2949            0 :                   jm = MIN(7, nmo_active - ii)
    2950            0 :                   WRITE (iw, '(T3,6(F9.4))') (pmat(i1, ii + i2), i2=0, jm)
    2951              :                END DO
    2952              :             END DO
    2953              :          END IF
    2954              : 
    2955              :          ! diagonalize the density matrix
    2956            0 :          CALL diamat_all(pmat, noon(:, ispin))
    2957              : 
    2958            0 :          IF (iw > 0) THEN
    2959            0 :             WRITE (iw, '(/,T3,A,I2,A)') "Natural orbitals occupation numbers for spin ", ispin
    2960            0 :             DO i1 = 1, nmo_active, 8
    2961            0 :                jm = MIN(7, nmo_active - i1)
    2962              :                ! noons are stored in ascending order, so reverse-print them
    2963            0 :                WRITE (iw, '(T3,6(F9.4))') (noon(nmo_active - i1 - i2 + 1, ispin), i2=0, jm)
    2964              :             END DO
    2965              :          END IF
    2966              : 
    2967              :       END DO
    2968              : 
    2969            0 :       DEALLOCATE (noon)
    2970            0 :       DEALLOCATE (pmat)
    2971              : 
    2972            0 :    END SUBROUTINE print_pmat_noon
    2973              : 
    2974              : ! **************************************************************************************************
    2975              : !> \brief ...
    2976              : !> \param qs_env ...
    2977              : !> \param active_space_env ...
    2978              : !> \param as_input ...
    2979              : ! **************************************************************************************************
    2980            0 :    SUBROUTINE rsdft_embedding(qs_env, active_space_env, as_input)
    2981              :       TYPE(qs_environment_type), POINTER                 :: qs_env
    2982              :       TYPE(active_space_type), POINTER                   :: active_space_env
    2983              :       TYPE(section_vals_type), POINTER                   :: as_input
    2984              : 
    2985              :       CHARACTER(len=*), PARAMETER                        :: routineN = 'rsdft_embedding'
    2986              :       INTEGER                                            :: handle
    2987              : 
    2988              : #ifdef __NO_SOCKETS
    2989              :       CALL timeset(routineN, handle)
    2990              :       CPABORT("CP2K was compiled with the __NO_SOCKETS option!")
    2991              :       MARK_USED(qs_env)
    2992              :       MARK_USED(active_space_env)
    2993              :       MARK_USED(as_input)
    2994              : #else
    2995              : 
    2996              :       INTEGER                                            :: iw, client_fd, socket_fd, iter, max_iter
    2997              :       LOGICAL                                            :: converged, do_scf_embedding, ionode
    2998              :       REAL(KIND=dp)                                      :: delta_E, energy_corr, energy_new, &
    2999              :                                                             energy_old, energy_scf, eps_iter, t1, t2
    3000              :       TYPE(cp_logger_type), POINTER                      :: logger
    3001            0 :       TYPE(dbcsr_p_type), DIMENSION(:), POINTER          :: rho_ao
    3002            0 :       TYPE(mo_set_type), DIMENSION(:), POINTER           :: mos_active
    3003              :       TYPE(mp_para_env_type), POINTER                    :: para_env
    3004              :       TYPE(qs_energy_type), POINTER                      :: energy
    3005              :       TYPE(qs_ks_env_type), POINTER                      :: ks_env
    3006              :       TYPE(qs_rho_type), POINTER                         :: rho
    3007              : 
    3008            0 :       CALL timeset(routineN, handle)
    3009              : 
    3010            0 :       t1 = m_walltime()
    3011              : 
    3012            0 :       logger => cp_get_default_logger()
    3013            0 :       iw = cp_logger_get_default_io_unit(logger)
    3014              : 
    3015            0 :       CALL get_qs_env(qs_env, para_env=para_env)
    3016            0 :       ionode = para_env%is_source()
    3017              : 
    3018              :       ! get info from the input
    3019            0 :       CALL section_vals_val_get(as_input, "SCF_EMBEDDING", l_val=do_scf_embedding)
    3020            0 :       active_space_env%do_scf_embedding = do_scf_embedding
    3021            0 :       CALL section_vals_val_get(as_input, "MAX_ITER", i_val=max_iter)
    3022            0 :       IF (max_iter < 0) CPABORT("Specify a non-negative number of max iterations.")
    3023            0 :       CALL section_vals_val_get(as_input, "EPS_ITER", r_val=eps_iter)
    3024            0 :       IF (eps_iter < 0.0) CPABORT("Specify a non-negative convergence threshold.")
    3025              : 
    3026              :       ! create the socket and wait for the client to connect
    3027            0 :       CALL initialize_socket(socket_fd, client_fd, as_input, ionode)
    3028            0 :       CALL para_env%sync()
    3029              : 
    3030              :       ! send two-electron integrals to the client
    3031            0 :       CALL send_eri_to_client(client_fd, active_space_env, para_env)
    3032              : 
    3033              :       ! get pointer to density in ao basis
    3034            0 :       CALL get_qs_env(qs_env, rho=rho, energy=energy, ks_env=ks_env)
    3035            0 :       CALL qs_rho_get(rho, rho_ao=rho_ao)
    3036              : 
    3037            0 :       IF (iw > 0) THEN
    3038              :          WRITE (UNIT=iw, FMT="(/,T2,A,/)") &
    3039            0 :             "RANGE-SEPARATED DFT EMBEDDING SELF-CONSISTENT OPTIMIZATION"
    3040              : 
    3041            0 :          WRITE (iw, '(T3,A,T68,I12)') "Max. iterations", max_iter
    3042            0 :          WRITE (iw, '(T3,A,T68,E12.4)') "Conv. threshold", eps_iter
    3043            0 :          WRITE (iw, '(T3,A,T66,F14.2)') "Density damping", active_space_env%alpha
    3044              : 
    3045              :          WRITE (UNIT=iw, FMT="(/,T3,A,T11,A,T21,A,T34,A,T55,A,T75,A,/,T3,A)") &
    3046            0 :             "Iter", "Update", "Time", "Corr. energy", "Total energy", "Change", REPEAT("-", 78)
    3047              :       END IF
    3048              :       ! CALL cp_add_iter_level(logger%iter_info, "QS_SCF")
    3049              : 
    3050            0 :       iter = 0
    3051            0 :       converged = .FALSE.
    3052              :       ! store the scf energy
    3053            0 :       energy_scf = active_space_env%energy_ref
    3054            0 :       energy_new = energy_scf
    3055            0 :       mos_active => active_space_env%mos_active
    3056              :       ! CALL set_qs_env(qs_env, active_space=active_space_env)
    3057              : 
    3058              :       ! start the self-consistent embedding loop
    3059            0 :       DO WHILE (iter < max_iter)
    3060            0 :          iter = iter + 1
    3061              : 
    3062              :          ! send V_emb and E_ina to the active space solver and update
    3063              :          ! the active space environment with the new active energy and density
    3064            0 :          CALL send_fock_to_client(client_fd, active_space_env, para_env)
    3065              : 
    3066              :          ! update energies
    3067            0 :          energy_old = energy_new
    3068            0 :          energy_new = active_space_env%energy_total
    3069            0 :          energy_corr = energy_new - energy_scf
    3070            0 :          delta_E = energy_new - energy_old
    3071              : 
    3072              :          ! get timer
    3073            0 :          t2 = t1
    3074            0 :          t1 = m_walltime()
    3075              :          ! print out progress
    3076            0 :          IF ((iw > 0)) THEN
    3077              :             WRITE (UNIT=iw, &
    3078              :                    FMT="(T3,I4,T11,A,T19,F6.1,T28,F18.10,T49,F18.10,T70,ES11.2)") &
    3079            0 :                iter, 'P_Mix', t1 - t2, energy_corr, energy_new, delta_E
    3080            0 :             CALL m_flush(iw)
    3081              :          END IF
    3082              : 
    3083              :          ! update total density in AO basis with the AS contribution
    3084            0 :          CALL update_density_ao(active_space_env, rho_ao) ! rho_ao is updated
    3085              : 
    3086              :          ! calculate F_ks in AO basis (which contains Vxc) with the new density
    3087            0 :          CALL qs_rho_update_rho(rho, qs_env=qs_env) ! updates rho_r and rho_g using rho_ao
    3088            0 :          CALL qs_ks_did_change(qs_env%ks_env, rho_changed=.TRUE.) ! set flags about the change
    3089              :          ! Re-evaluate the traces between the density matrix and the core Hamiltonians
    3090            0 :          CALL evaluate_core_matrix_traces(qs_env)
    3091              :          ! the ks matrix will be rebuilt so this is fine now
    3092              :          ! CALL set_ks_env(qs_env%ks_env, potential_changed=.FALSE.)
    3093              :          CALL qs_ks_build_kohn_sham_matrix(qs_env, calculate_forces=.FALSE., &
    3094              :                                            just_energy=.FALSE., &
    3095            0 :                                            ext_xc_section=active_space_env%xc_section)
    3096              : 
    3097              :          ! update the reference energy
    3098            0 :          active_space_env%energy_ref = energy%total
    3099              : 
    3100              :          ! transform KS/Fock, Vxc and Hcore from AO to MO basis
    3101            0 :          CALL calculate_operators(mos_active, qs_env, active_space_env)
    3102              : 
    3103              :          ! calculate the new inactive energy and embedding potential
    3104            0 :          CALL subspace_fock_matrix(active_space_env)
    3105              : 
    3106              :          ! check if it is a one-shot correction
    3107            0 :          IF (.NOT. active_space_env%do_scf_embedding) THEN
    3108            0 :             IF (iw > 0) THEN
    3109              :                WRITE (UNIT=iw, FMT="(/,T3,A,I5,A)") &
    3110            0 :                   "*** one-shot embedding correction finished ***"
    3111              :             END IF
    3112              :             converged = .TRUE.
    3113              :             EXIT
    3114              :             ! check for convergence
    3115            0 :          ELSEIF (ABS(delta_E) <= eps_iter) THEN
    3116            0 :             IF (iw > 0) THEN
    3117              :                WRITE (UNIT=iw, FMT="(/,T3,A,I5,A)") &
    3118            0 :                   "*** rs-DFT embedding run converged in ", iter, " iteration(s) ***"
    3119              :             END IF
    3120              :             converged = .TRUE.
    3121              :             EXIT
    3122              :          END IF
    3123              :       END DO
    3124              : 
    3125              :       IF (.NOT. converged) THEN
    3126            0 :          IF (iw > 0) THEN
    3127              :             WRITE (UNIT=iw, FMT="(/,T3,A,I5,A)") &
    3128            0 :                "*** rs-DFT embedding did not converged after ", iter, " iteration(s) ***"
    3129              :          END IF
    3130              :       END IF
    3131              : 
    3132              :       ! update qs total energy to the final rs-DFT energy
    3133            0 :       energy%total = active_space_env%energy_total
    3134              : 
    3135              :       ! print final energy contributions
    3136            0 :       IF (iw > 0) THEN
    3137              :          WRITE (UNIT=iw, FMT="(/,T3,A)") &
    3138            0 :             "Final energy contributions:"
    3139              :          WRITE (UNIT=iw, FMT="(T6,A,T56,F20.10)") &
    3140            0 :             "Inactive energy:", active_space_env%energy_inactive
    3141              :          WRITE (UNIT=iw, FMT="(T6,A,T56,F20.10)") &
    3142            0 :             "Active energy:", active_space_env%energy_active
    3143              :          WRITE (UNIT=iw, FMT="(T6,A,T56,F20.10)") &
    3144            0 :             "Correlation energy:", energy_corr
    3145              :          WRITE (UNIT=iw, FMT="(T6,A,T56,F20.10)") &
    3146            0 :             "Total rs-DFT energy:", active_space_env%energy_total
    3147              :       END IF
    3148              : 
    3149              :       ! print the AS rdm and the natural orbital occupation numbers
    3150            0 :       CALL print_pmat_noon(active_space_env, iw)
    3151              : 
    3152            0 :       CALL finalize_socket(socket_fd, client_fd, as_input, ionode)
    3153            0 :       CALL para_env%sync()
    3154              : #endif
    3155              : 
    3156            0 :       CALL timestop(handle)
    3157              : 
    3158            0 :    END SUBROUTINE rsdft_embedding
    3159              : 
    3160              : #ifndef __NO_SOCKETS
    3161              : ! **************************************************************************************************
    3162              : !> \brief Creates the socket, spawns the client and connects to it
    3163              : !> \param socket_fd the socket file descriptor
    3164              : !> \param client_fd the client file descriptor
    3165              : !> \param as_input active space inpute section
    3166              : !> \param ionode logical flag indicating if the process is the master
    3167              : ! **************************************************************************************************
    3168            0 :    SUBROUTINE initialize_socket(socket_fd, client_fd, as_input, ionode)
    3169              :       INTEGER, INTENT(OUT)                               :: socket_fd, client_fd
    3170              :       TYPE(section_vals_type), INTENT(IN), POINTER       :: as_input
    3171              :       LOGICAL, INTENT(IN)                                :: ionode
    3172              : 
    3173              :       CHARACTER(len=*), PARAMETER                        :: routineN = 'initialize_socket'
    3174              :       INTEGER, PARAMETER                                 :: backlog = 10
    3175              : 
    3176              :       CHARACTER(len=default_path_length)                 :: hostname
    3177              :       INTEGER                                            :: handle, iw, port, protocol
    3178              :       LOGICAL                                            :: inet
    3179              :       TYPE(cp_logger_type), POINTER                      :: logger
    3180              : 
    3181            0 :       CALL timeset(routineN, handle)
    3182              : 
    3183            0 :       logger => cp_get_default_logger()
    3184            0 :       iw = cp_logger_get_default_io_unit(logger)
    3185              : 
    3186              :       ! protocol == 0 for UNIX, protocol > 0 for INET
    3187            0 :       CALL section_vals_val_get(as_input, "SOCKET%INET", l_val=inet)
    3188            0 :       IF (inet) THEN
    3189            0 :          protocol = 1
    3190              :       ELSE
    3191            0 :          protocol = 0
    3192              :       END IF
    3193            0 :       CALL section_vals_val_get(as_input, "SOCKET%HOST", c_val=hostname)
    3194            0 :       CALL section_vals_val_get(as_input, "SOCKET%PORT", i_val=port)
    3195              : 
    3196            0 :       IF (ionode) THEN
    3197            0 :          CALL open_bind_socket(socket_fd, protocol, port, TRIM(hostname)//C_NULL_CHAR)
    3198            0 :          WRITE (iw, '(/,T2,A,A)') "@SERVER: Created socket with address ", TRIM(hostname)
    3199            0 :          CALL listen_socket(socket_fd, backlog)
    3200              : 
    3201              :          ! wait until a connetion request arrives
    3202            0 :          WRITE (iw, '(T2,A)') "@SERVER: Waiting for requests..."
    3203            0 :          CALL accept_socket(socket_fd, client_fd)
    3204            0 :          WRITE (iw, '(T2,A,I2)') "@SERVER: Accepted socket with fd ", client_fd
    3205              :       END IF
    3206              : 
    3207            0 :       CALL timestop(handle)
    3208              : 
    3209            0 :    END SUBROUTINE initialize_socket
    3210              : 
    3211              : ! **************************************************************************************************
    3212              : !> \brief Closes the connection to the socket and deletes the file
    3213              : !> \param socket_fd the socket file descriptor
    3214              : !> \param client_fd the client file descriptor
    3215              : !> \param as_input active space inpute section
    3216              : !> \param ionode logical flag indicating if the process is the master
    3217              : ! **************************************************************************************************
    3218            0 :    SUBROUTINE finalize_socket(socket_fd, client_fd, as_input, ionode)
    3219              :       INTEGER, INTENT(IN)                                :: socket_fd, client_fd
    3220              :       TYPE(section_vals_type), INTENT(IN), POINTER       :: as_input
    3221              :       LOGICAL, INTENT(IN)                                :: ionode
    3222              : 
    3223              :       CHARACTER(len=*), PARAMETER                        :: routineN = 'finalize_socket'
    3224              :       INTEGER, PARAMETER                                 :: header_len = 12
    3225              : 
    3226              :       CHARACTER(len=default_path_length)                 :: hostname
    3227              :       INTEGER                                            :: handle
    3228              : 
    3229            0 :       CALL timeset(routineN, handle)
    3230              : 
    3231            0 :       CALL section_vals_val_get(as_input, "SOCKET%HOST", c_val=hostname)
    3232              : 
    3233            0 :       IF (ionode) THEN
    3234              :          ! signal the client to quit
    3235            0 :          CALL writebuffer(client_fd, "QUIT        ", header_len)
    3236              :          ! close the connection
    3237            0 :          CALL close_socket(client_fd)
    3238            0 :          CALL close_socket(socket_fd)
    3239              : 
    3240              :          ! delete the socket file
    3241            0 :          IF (file_exists(TRIM(hostname))) THEN
    3242            0 :             CALL remove_socket_file(TRIM(hostname)//C_NULL_CHAR)
    3243              :          END IF
    3244              :       END IF
    3245              : 
    3246            0 :       CALL timestop(handle)
    3247              : 
    3248            0 :    END SUBROUTINE finalize_socket
    3249              : 
    3250              : ! **************************************************************************************************
    3251              : !> \brief Sends the two-electron integrals to the client vie the socket
    3252              : !> \param client_fd the client file descriptor
    3253              : !> \param active_space_env active space environment
    3254              : !> \param para_env parallel environment
    3255              : ! **************************************************************************************************
    3256            0 :    SUBROUTINE send_eri_to_client(client_fd, active_space_env, para_env)
    3257              :       INTEGER, INTENT(IN)                                :: client_fd
    3258              :       TYPE(active_space_type), INTENT(IN), POINTER       :: active_space_env
    3259              :       TYPE(mp_para_env_type), INTENT(IN), POINTER        :: para_env
    3260              : 
    3261              :       CHARACTER(len=*), PARAMETER :: routineN = 'send_eri_to_client'
    3262              :       INTEGER, PARAMETER                                 :: header_len = 12
    3263              : 
    3264              :       CHARACTER(len=default_string_length)               :: header
    3265              :       INTEGER                                            :: handle, iw
    3266              :       LOGICAL                                            :: ionode
    3267            0 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:)           :: eri_aa, eri_ab, eri_bb, s_ab
    3268              :       TYPE(cp_logger_type), POINTER                      :: logger
    3269              : 
    3270            0 :       CALL timeset(routineN, handle)
    3271              : 
    3272            0 :       logger => cp_get_default_logger()
    3273            0 :       iw = cp_logger_get_default_io_unit(logger)
    3274            0 :       ionode = para_env%is_source()
    3275              : 
    3276            0 :       ALLOCATE (eri_aa(active_space_env%nmo_active**4))
    3277            0 :       CALL eri_to_array(active_space_env%eri, eri_aa, active_space_env%active_orbitals, 1, 1)
    3278            0 :       IF (active_space_env%nspins == 2) THEN
    3279            0 :          ALLOCATE (eri_ab(active_space_env%nmo_active**4))
    3280            0 :          CALL eri_to_array(active_space_env%eri, eri_ab, active_space_env%active_orbitals, 1, 2)
    3281            0 :          ALLOCATE (eri_bb(active_space_env%nmo_active**4))
    3282            0 :          CALL eri_to_array(active_space_env%eri, eri_bb, active_space_env%active_orbitals, 2, 2)
    3283              :          ! get the overlap_ab matrix into Fortran array
    3284            0 :          ALLOCATE (s_ab(active_space_env%nmo_active**2))
    3285              :          ASSOCIATE (act_indices_a => active_space_env%active_orbitals(:, 1), &
    3286              :                     act_indices_b => active_space_env%active_orbitals(:, 2))
    3287            0 :             CALL subspace_matrix_to_array(active_space_env%sab_sub(1), s_ab, act_indices_a, act_indices_b)
    3288              :          END ASSOCIATE
    3289              :       END IF
    3290              : 
    3291              :       ! ask the status of the client
    3292            0 :       IF (ionode) CALL writebuffer(client_fd, "STATUS      ", header_len)
    3293              :       DO
    3294            0 :          header = ""
    3295            0 :          CALL para_env%sync()
    3296            0 :          IF (ionode) THEN
    3297              :             ! IF (iw > 0) WRITE(iw, *) "@SERVER: Waiting for messages..."
    3298            0 :             CALL readbuffer(client_fd, header, header_len)
    3299              :          END IF
    3300            0 :          CALL para_env%bcast(header, para_env%source)
    3301              : 
    3302              :          ! IF (iw > 0) WRITE(iw, *) "@SERVER: Message from client: ", TRIM(header)
    3303              : 
    3304            0 :          IF (TRIM(header) == "READY") THEN
    3305              :             ! if the client is ready, send the data
    3306            0 :             CALL para_env%sync()
    3307            0 :             IF (ionode) THEN
    3308            0 :                CALL writebuffer(client_fd, "TWOBODY     ", header_len)
    3309            0 :                CALL writebuffer(client_fd, active_space_env%nspins)
    3310            0 :                CALL writebuffer(client_fd, active_space_env%nmo_active)
    3311            0 :                CALL writebuffer(client_fd, active_space_env%nelec_active)
    3312            0 :                CALL writebuffer(client_fd, active_space_env%multiplicity)
    3313              :                ! send the alpha component
    3314            0 :                CALL writebuffer(client_fd, eri_aa, SIZE(eri_aa))
    3315              :                ! send the beta part for unrestricted calculations
    3316            0 :                IF (active_space_env%nspins == 2) THEN
    3317            0 :                   CALL writebuffer(client_fd, eri_ab, SIZE(eri_ab))
    3318            0 :                   CALL writebuffer(client_fd, eri_bb, SIZE(eri_bb))
    3319            0 :                   CALL writebuffer(client_fd, s_ab, SIZE(s_ab))
    3320              :                END IF
    3321              :             END IF
    3322            0 :          ELSE IF (TRIM(header) == "RECEIVED") THEN
    3323              :             EXIT
    3324              :          END IF
    3325              :       END DO
    3326              : 
    3327            0 :       DEALLOCATE (eri_aa)
    3328            0 :       IF (active_space_env%nspins == 2) THEN
    3329            0 :          DEALLOCATE (eri_ab)
    3330            0 :          DEALLOCATE (eri_bb)
    3331            0 :          DEALLOCATE (s_ab)
    3332              :       END IF
    3333              : 
    3334            0 :       CALL para_env%sync()
    3335              : 
    3336            0 :       CALL timestop(handle)
    3337              : 
    3338            0 :    END SUBROUTINE send_eri_to_client
    3339              : 
    3340              : ! **************************************************************************************************
    3341              : !> \brief Sends the one-electron embedding potential and the inactive energy to the client
    3342              : !> \param client_fd the client file descriptor
    3343              : !> \param active_space_env active space environment
    3344              : !> \param para_env parallel environment
    3345              : ! **************************************************************************************************
    3346            0 :    SUBROUTINE send_fock_to_client(client_fd, active_space_env, para_env)
    3347              :       INTEGER, INTENT(IN)                                :: client_fd
    3348              :       TYPE(active_space_type), INTENT(IN), POINTER       :: active_space_env
    3349              :       TYPE(mp_para_env_type), INTENT(IN), POINTER        :: para_env
    3350              : 
    3351              :       CHARACTER(len=*), PARAMETER :: routineN = 'send_fock_to_client'
    3352              :       INTEGER, PARAMETER                                 :: header_len = 12
    3353              : 
    3354              :       CHARACTER(len=default_string_length)               :: header
    3355              :       INTEGER                                            :: handle, iw
    3356              :       LOGICAL                                            :: debug, ionode
    3357            0 :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:)           :: fock_a, fock_b, p_act_mo_a, p_act_mo_b
    3358              :       TYPE(cp_logger_type), POINTER                      :: logger
    3359              : 
    3360            0 :       CALL timeset(routineN, handle)
    3361              : 
    3362              :       ! Set to .TRUE. to activate debug output
    3363            0 :       debug = .FALSE.
    3364              : 
    3365            0 :       logger => cp_get_default_logger()
    3366            0 :       iw = cp_logger_get_default_io_unit(logger)
    3367            0 :       ionode = para_env%is_source()
    3368              : 
    3369            0 :       ALLOCATE (p_act_mo_a(active_space_env%nmo_active**2))
    3370            0 :       ALLOCATE (fock_a(active_space_env%nmo_active**2))
    3371            0 :       IF (active_space_env%nspins == 2) THEN
    3372            0 :          ALLOCATE (p_act_mo_b(active_space_env%nmo_active**2))
    3373            0 :          ALLOCATE (fock_b(active_space_env%nmo_active**2))
    3374              :       END IF
    3375              : 
    3376              :       ! get the fock matrix into Fortran arrays
    3377              :       ASSOCIATE (act_indices => active_space_env%active_orbitals(:, 1))
    3378            0 :          CALL subspace_matrix_to_array(active_space_env%fock_sub(1), fock_a, act_indices, act_indices)
    3379              :       END ASSOCIATE
    3380              : 
    3381            0 :       IF (active_space_env%nspins == 2) THEN
    3382              :          ASSOCIATE (act_indices => active_space_env%active_orbitals(:, 2))
    3383            0 :             CALL subspace_matrix_to_array(active_space_env%fock_sub(2), fock_b, act_indices, act_indices)
    3384              :          END ASSOCIATE
    3385              :       END IF
    3386              : 
    3387              :       ! ask the status of the client
    3388            0 :       IF (ionode) CALL writebuffer(client_fd, "STATUS      ", header_len)
    3389              :       DO
    3390            0 :          header = ""
    3391              : 
    3392            0 :          CALL para_env%sync()
    3393            0 :          IF (ionode) THEN
    3394              :             IF (debug .AND. iw > 0) WRITE (iw, *) "@SERVER: Waiting for messages..."
    3395            0 :             CALL readbuffer(client_fd, header, header_len)
    3396              :          END IF
    3397            0 :          CALL para_env%bcast(header, para_env%source)
    3398              : 
    3399              :          IF (debug .AND. iw > 0) WRITE (iw, *) "@SERVER: Message from client: ", TRIM(header)
    3400              : 
    3401            0 :          IF (TRIM(header) == "READY") THEN
    3402              :             ! if the client is ready, send the data
    3403            0 :             CALL para_env%sync()
    3404            0 :             IF (ionode) THEN
    3405            0 :                CALL writebuffer(client_fd, "ONEBODY     ", header_len)
    3406            0 :                CALL writebuffer(client_fd, active_space_env%energy_inactive)
    3407              :                ! send the alpha component
    3408            0 :                CALL writebuffer(client_fd, fock_a, SIZE(fock_a))
    3409              :                ! send the beta part for unrestricted calculations
    3410            0 :                IF (active_space_env%nspins == 2) THEN
    3411            0 :                   CALL writebuffer(client_fd, fock_b, SIZE(fock_b))
    3412              :                END IF
    3413              :             END IF
    3414              : 
    3415            0 :          ELSE IF (TRIM(header) == "HAVEDATA") THEN
    3416              :             ! qiskit has data to transfer, let them know we want it and wait for it
    3417            0 :             CALL para_env%sync()
    3418            0 :             IF (ionode) THEN
    3419              :                IF (debug .AND. iw > 0) WRITE (iw, *) "@SERVER: Qiskit has data to transfer"
    3420            0 :                CALL writebuffer(client_fd, "GETDENSITY  ", header_len)
    3421              : 
    3422              :                ! read the active energy and density
    3423            0 :                CALL readbuffer(client_fd, active_space_env%energy_active)
    3424            0 :                CALL readbuffer(client_fd, p_act_mo_a, SIZE(p_act_mo_a))
    3425            0 :                IF (active_space_env%nspins == 2) THEN
    3426            0 :                   CALL readbuffer(client_fd, p_act_mo_b, SIZE(p_act_mo_b))
    3427              :                END IF
    3428              :             END IF
    3429              : 
    3430              :             ! broadcast the data to all processors
    3431            0 :             CALL para_env%bcast(active_space_env%energy_active, para_env%source)
    3432            0 :             CALL para_env%bcast(p_act_mo_a, para_env%source)
    3433            0 :             IF (active_space_env%nspins == 2) THEN
    3434            0 :                CALL para_env%bcast(p_act_mo_b, para_env%source)
    3435              :             END IF
    3436              : 
    3437              :             ! update total and reference energies in active space enviornment
    3438            0 :             active_space_env%energy_total = active_space_env%energy_inactive + active_space_env%energy_active
    3439              : 
    3440              :             ! update the active density matrix in the active space environment
    3441            0 :             CALL update_active_density(p_act_mo_a, active_space_env, 1)
    3442            0 :             IF (active_space_env%nspins == 2) THEN
    3443            0 :                CALL update_active_density(p_act_mo_b, active_space_env, 2)
    3444              :             END IF
    3445              : 
    3446              :             ! the non-iterative part is done, we can continue
    3447              :             EXIT
    3448              :          END IF
    3449              : 
    3450              :       END DO
    3451              : 
    3452            0 :       DEALLOCATE (p_act_mo_a)
    3453            0 :       DEALLOCATE (fock_a)
    3454            0 :       IF (active_space_env%nspins == 2) THEN
    3455            0 :          DEALLOCATE (p_act_mo_b)
    3456            0 :          DEALLOCATE (fock_b)
    3457              :       END IF
    3458              : 
    3459            0 :       CALL para_env%sync()
    3460              : 
    3461            0 :       CALL timestop(handle)
    3462              : 
    3463            0 :    END SUBROUTINE send_fock_to_client
    3464              : #endif
    3465              : 
    3466            0 : END MODULE qs_active_space_methods
        

Generated by: LCOV version 2.0-1