Skip to content

Commit

Permalink
Replace arrays of cp_fm_p_type with arrays of cp_fm_type (#2513)
Browse files Browse the repository at this point in the history
This replacement improves the readibility and a layer of pointers is lost. In a few cases, the arrays can be allocated on the stack instead of the heap. In one case, I kept the cp_fm_p_type because pointers to mos are totally sufficient in that case and a copy of data was not intended.

* Remove cp_fm_p_type in negf_methods

* Replace cp_fm_p_type in kpoint_types

* Refactor cp_fm_p_type in qs_ks_utils

* Refactor cp_fm_p_type in qs_efield_berry

* qs_ks_apply_restraints

* mixed_cdft_types

* Remove a pointer attribute

* negf_env_types

* mixed_cdft_methods

* Make pretty

* et_coupling, mixed_cdft_utils

* zij_fm_set

* qs_scf_env_types

* qs_active_space_types

* qs_active_space_types

* qs_linres_op

* qs_tddfpt2_types

* qs_dcdr

* qs_wannier90

* Refactor psi0_order

* Refactor qs_loc_types

* moloc_coeff

* matrix_exp

* Fix comment

* fc_psi0

* psi1_fc

* fm_work in qs_linres_current_utils

* dipole_op_mos_occ

* some work matrices

* p_op_l1

* Continue with Sb

* Continue with b

* Continue with R,X

* mo_coeff_array

* cpmos, h1_psi0, p_psi0

* psi1

* work_arrays rr_p*, r_p*

* et_coupling_proj

* mat_h

* mat_w

* invS

* More qs_tddfpt2

* qs_moments

* More work matrices

* dBerry_psi0

* work matrices qs_linres_methods

* work matrices qs_linres_methods

* et_mo_coeff

* Remove optional argument (never used)

* Last work matrices

* S_krylov, S_evects

* t_env%evects

* Aop_evects

* Aop_ritz

* Remove never used ec_env%cpmos

* mos_occ

* qs_linres_current

* efg_psi0, pso_psi0, dso_psi0

* psi1_efg, psi1_pso, psi1_dso

* qs_dcdr

* qs_diis

* mos_occ

* qs_cdft_types
  • Loading branch information
fstein93 committed Jan 18, 2023
1 parent baad364 commit 60432cd
Show file tree
Hide file tree
Showing 90 changed files with 2,370 additions and 2,900 deletions.
16 changes: 4 additions & 12 deletions src/cp_control_types.F
Original file line number Diff line number Diff line change
Expand Up @@ -10,8 +10,8 @@
!> settings for the DFT-based calculations.
! **************************************************************************************************
MODULE cp_control_types
USE cp_fm_types, ONLY: cp_fm_p_type,&
cp_fm_release
USE cp_fm_types, ONLY: cp_fm_release,&
cp_fm_type
USE input_constants, ONLY: do_full_density
USE kinds, ONLY: default_path_length,&
default_string_length,&
Expand Down Expand Up @@ -391,7 +391,7 @@ MODULE cp_control_types
! TDDFPT operator without the perturbation kernel.
! **************************************************************************************************
TYPE tddfpt_control_type
TYPE(cp_fm_p_type), DIMENSION(:), &
TYPE(cp_fm_type), DIMENSION(:), &
POINTER :: lumos
REAL(KIND=dp) :: tolerance
INTEGER :: n_ev
Expand Down Expand Up @@ -869,19 +869,11 @@ SUBROUTINE tddfpt_control_release(tddfpt_control)
TYPE(tddfpt_control_type), POINTER :: tddfpt_control

INTEGER :: ispin
LOGICAL :: dummy

IF (ASSOCIATED(tddfpt_control)) THEN
IF (ASSOCIATED(tddfpt_control%lumos)) THEN
DO ispin = 1, SIZE(tddfpt_control%lumos)
CALL cp_fm_release(tddfpt_control%lumos(ispin)%matrix)
DEALLOCATE (tddfpt_control%lumos(ispin)%matrix)
NULLIFY (tddfpt_control%lumos(ispin)%matrix)
!MK the following line just avoids a crash of TDDFT runs using
!MK the sdbg version compiled with the NAG compiler when
!MK tddfpt_control%lumos is deallocated. This is most likely a
!MK compiler bug and thus the line might become obsolete
dummy = ASSOCIATED(tddfpt_control%lumos(ispin)%matrix)
CALL cp_fm_release(tddfpt_control%lumos(ispin))
END DO
DEALLOCATE (tddfpt_control%lumos)
END IF
Expand Down
4 changes: 2 additions & 2 deletions src/dft_plus_u.F
Original file line number Diff line number Diff line change
Expand Up @@ -351,8 +351,8 @@ SUBROUTINE lowdin(qs_env, matrix_h, matrix_w, should_output, output_unit, &
fm_work2_local_alloc = .FALSE.

IF (ASSOCIATED(scf_env%scf_work1)) THEN
IF (ASSOCIATED(scf_env%scf_work1(1)%matrix)) THEN
fm_work1 => scf_env%scf_work1(1)%matrix
IF (ASSOCIATED(scf_env%scf_work1(1)%matrix_struct)) THEN
fm_work1 => scf_env%scf_work1(1)
END IF
END IF

Expand Down
12 changes: 0 additions & 12 deletions src/ec_env_types.F
Original file line number Diff line number Diff line change
Expand Up @@ -13,8 +13,6 @@
! **************************************************************************************************
MODULE ec_env_types
USE cp_dbcsr_operations, ONLY: dbcsr_deallocate_matrix_set
USE cp_fm_types, ONLY: cp_fm_p_type,&
cp_fm_release
USE dbcsr_api, ONLY: dbcsr_p_type
USE dm_ls_scf_types, ONLY: ls_scf_env_type,&
ls_scf_release
Expand Down Expand Up @@ -95,7 +93,6 @@ MODULE ec_env_types
TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: mao_coef
! CP equations
TYPE(qs_p_env_type), POINTER :: p_env
TYPE(cp_fm_p_type), DIMENSION(:), POINTER :: cpmos
TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_hz, matrix_z, matrix_wz, z_admm
! Harris (rhoout), and response density (rhoz) on grid
TYPE(pw_type), DIMENSION(:), POINTER :: rhoout_r, rhoz_r
Expand Down Expand Up @@ -146,15 +143,6 @@ SUBROUTINE ec_env_release(ec_env)
IF (ASSOCIATED(ec_env%dispersion_env)) THEN
CALL qs_dispersion_release(ec_env%dispersion_env)
END IF
! CP env
IF (ASSOCIATED(ec_env%cpmos)) THEN
DO iab = 1, SIZE(ec_env%cpmos)
CALL cp_fm_release(ec_env%cpmos(iab)%matrix)
DEALLOCATE (ec_env%cpmos(iab)%matrix)
END DO
DEALLOCATE (ec_env%cpmos)
NULLIFY (ec_env%cpmos)
END IF

IF (ASSOCIATED(ec_env%matrix_z)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_z)
IF (ASSOCIATED(ec_env%matrix_hz)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_hz)
Expand Down
1 change: 0 additions & 1 deletion src/ec_environment.F
Original file line number Diff line number Diff line change
Expand Up @@ -141,7 +141,6 @@ SUBROUTINE init_ec_env(qs_env, ec_env, dft_section, ec_section)
NULLIFY (ec_env%force)
NULLIFY (ec_env%dispersion_env)
NULLIFY (ec_env%xc_section)
NULLIFY (ec_env%cpmos)
NULLIFY (ec_env%matrix_z)
NULLIFY (ec_env%matrix_hz)
NULLIFY (ec_env%matrix_wz)
Expand Down
29 changes: 11 additions & 18 deletions src/et_coupling.F
Original file line number Diff line number Diff line change
Expand Up @@ -21,7 +21,6 @@ MODULE et_coupling
USE cp_fm_struct, ONLY: cp_fm_struct_type
USE cp_fm_types, ONLY: cp_fm_create,&
cp_fm_get_info,&
cp_fm_p_type,&
cp_fm_release,&
cp_fm_type
USE cp_log_handling, ONLY: cp_get_default_logger,&
Expand Down Expand Up @@ -74,18 +73,18 @@ SUBROUTINE calc_et_coupling(qs_env)
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: a, b, S_det
REAL(KIND=dp), DIMENSION(2) :: eigenv
REAL(KIND=dp), DIMENSION(2, 2) :: S_mat, tmp_mat, U, W_mat
TYPE(cp_fm_p_type), DIMENSION(:), POINTER :: rest_MO
TYPE(cp_fm_pool_p_type), DIMENSION(:), POINTER :: mo_mo_fm_pools
TYPE(cp_fm_struct_type), POINTER :: mo_mo_fmstruct
TYPE(cp_fm_type) :: inverse_mat, SMO, Tinverse, tmp2
TYPE(cp_fm_type), DIMENSION(2) :: rest_MO
TYPE(cp_logger_type), POINTER :: logger
TYPE(cp_para_env_type), POINTER :: para_env
TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s
TYPE(dft_control_type), POINTER :: dft_control
TYPE(qs_energy_type), POINTER :: energy
TYPE(section_vals_type), POINTER :: et_coupling_section

NULLIFY (rest_MO, mo_mo_fmstruct, rest_MO, energy, matrix_s, dft_control, para_env)
NULLIFY (mo_mo_fmstruct, energy, matrix_s, dft_control, para_env)

CALL timeset(routineN, handle)

Expand All @@ -99,7 +98,6 @@ SUBROUTINE calc_et_coupling(qs_env)
extension=".log")

is_spin_constraint = .FALSE.
ALLOCATE (rest_MO(2))
ALLOCATE (a(dft_control%nspins))
ALLOCATE (b(dft_control%nspins))
ALLOCATE (S_det(dft_control%nspins))
Expand All @@ -126,16 +124,15 @@ SUBROUTINE calc_et_coupling(qs_env)
matrix_struct=mo_mo_fmstruct, &
name="ET_SMO"//TRIM(ADJUSTL(cp_to_string(1)))//"MATRIX")
DO j = 1, 2
ALLOCATE (rest_MO(j)%matrix)
CALL cp_fm_create(matrix=rest_MO(j)%matrix, &
CALL cp_fm_create(matrix=rest_MO(j), &
matrix_struct=mo_mo_fmstruct, &
name="ET_rest_MO"//TRIM(ADJUSTL(cp_to_string(j)))//"MATRIX")
END DO

! calculate MO-overlap

CALL get_qs_env(qs_env, matrix_s=matrix_s)
CALL cp_dbcsr_sm_fm_multiply(matrix_s(1)%matrix, qs_env%et_coupling%et_mo_coeff(i)%matrix, &
CALL cp_dbcsr_sm_fm_multiply(matrix_s(1)%matrix, qs_env%et_coupling%et_mo_coeff(i), &
tmp2, nmo, 1.0_dp, 0.0_dp)
CALL parallel_gemm('T', 'N', nmo, nmo, nao, 1.0_dp, &
qs_env%mos(i)%mo_coeff, &
Expand All @@ -144,12 +141,12 @@ SUBROUTINE calc_et_coupling(qs_env)
! calculate the MO-representation of the restraint matrix A

CALL cp_dbcsr_sm_fm_multiply(qs_env%et_coupling%rest_mat(1)%matrix, &
qs_env%et_coupling%et_mo_coeff(i)%matrix, &
qs_env%et_coupling%et_mo_coeff(i), &
tmp2, nmo, 1.0_dp, 0.0_dp)

CALL parallel_gemm('T', 'N', nmo, nmo, nao, 1.0_dp, &
qs_env%mos(i)%mo_coeff, &
tmp2, 0.0_dp, rest_MO(1)%matrix)
tmp2, 0.0_dp, rest_MO(1))

! calculate the MO-representation of the restraint matrix D

Expand All @@ -158,8 +155,8 @@ SUBROUTINE calc_et_coupling(qs_env)
tmp2, nmo, 1.0_dp, 0.0_dp)

CALL parallel_gemm('T', 'N', nmo, nmo, nao, 1.0_dp, &
qs_env%et_coupling%et_mo_coeff(i)%matrix, &
tmp2, 0.0_dp, rest_MO(2)%matrix)
qs_env%et_coupling%et_mo_coeff(i), &
tmp2, 0.0_dp, rest_MO(2))

CALL cp_fm_invert(SMO, inverse_mat, S_det(i))

Expand All @@ -169,15 +166,15 @@ SUBROUTINE calc_et_coupling(qs_env)

DO j = 1, ncol_local
DO k = 1, nrow_local
b(i) = b(i) + rest_MO(2)%matrix%local_data(k, j)*inverse_mat%local_data(k, j)
b(i) = b(i) + rest_MO(2)%local_data(k, j)*inverse_mat%local_data(k, j)
END DO
END DO

CALL cp_fm_transpose(inverse_mat, Tinverse)
a(i) = 0.0_dp
DO j = 1, ncol_local
DO k = 1, nrow_local
a(i) = a(i) + rest_MO(1)%matrix%local_data(k, j)*Tinverse%local_data(k, j)
a(i) = a(i) + rest_MO(1)%local_data(k, j)*Tinverse%local_data(k, j)
END DO
END DO
IF (is_spin_constraint) THEN
Expand All @@ -190,16 +187,12 @@ SUBROUTINE calc_et_coupling(qs_env)

CALL cp_fm_release(tmp2)
DO j = 1, 2
IF (ASSOCIATED(rest_MO(j)%matrix)) THEN
CALL cp_fm_release(rest_MO(j)%matrix)
DEALLOCATE (rest_MO(j)%matrix)
END IF
CALL cp_fm_release(rest_MO(j))
END DO
CALL cp_fm_release(SMO)
CALL cp_fm_release(Tinverse)
CALL cp_fm_release(inverse_mat)
END DO
DEALLOCATE (rest_MO)

! solve eigenstates for the projector matrix

Expand Down

0 comments on commit 60432cd

Please sign in to comment.