(git:374b731)
Loading...
Searching...
No Matches
qs_tddfpt2_densities.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
9 USE admm_types, ONLY: admm_type,&
14 USE cp_fm_types, ONLY: cp_fm_get_info,&
16 USE dbcsr_api, ONLY: dbcsr_p_type,&
17 dbcsr_scale
18 USE kinds, ONLY: default_string_length,&
19 dp
21 USE pw_env_types, ONLY: pw_env_get
23 USE pw_types, ONLY: pw_c1d_gs_type,&
33 USE qs_rho_methods, ONLY: qs_rho_copy,&
35 USE qs_rho_types, ONLY: qs_rho_get,&
39#include "./base/base_uses.f90"
40
41 IMPLICIT NONE
42
43 PRIVATE
44
45 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_tddfpt2_densities'
46
47 INTEGER, PARAMETER, PRIVATE :: maxspins = 2
48
50
51! **************************************************************************************************
52
53CONTAINS
54
55! **************************************************************************************************
56!> \brief Compute the ground-state charge density expressed in primary basis set.
57!> \param rho_orb_struct ground-state density in primary basis set
58!> \param rho_xc_struct ground-state density in primary basis set for GAPW_XC
59!> \param is_rks_triplets indicates that the triplet excited states calculation using
60!> spin-unpolarised molecular orbitals has been requested
61!> \param qs_env Quickstep environment
62!> \param sub_env parallel (sub)group environment
63!> \param wfm_rho_orb work dense matrix with shape [nao x nao] distributed among
64!> processors of the given parallel group (modified on exit)
65!> \par History
66!> * 06.2018 created by splitting the subroutine tddfpt_apply_admm_correction() in two
67!> subroutines tddfpt_construct_ground_state_orb_density() and
68!> tddfpt_construct_aux_fit_density [Sergey Chulkov]
69! **************************************************************************************************
70 SUBROUTINE tddfpt_construct_ground_state_orb_density(rho_orb_struct, rho_xc_struct, is_rks_triplets, &
71 qs_env, sub_env, wfm_rho_orb)
72 TYPE(qs_rho_type), POINTER :: rho_orb_struct, rho_xc_struct
73 LOGICAL, INTENT(in) :: is_rks_triplets
74 TYPE(qs_environment_type), POINTER :: qs_env
75 TYPE(tddfpt_subgroup_env_type), INTENT(in) :: sub_env
76 TYPE(cp_fm_type), INTENT(IN) :: wfm_rho_orb
77
78 CHARACTER(LEN=*), PARAMETER :: routinen = 'tddfpt_construct_ground_state_orb_density'
79
80 INTEGER :: handle, ispin, nao, nspins
81 INTEGER, DIMENSION(maxspins) :: nmo_occ
82 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rho_ij_ao
83 TYPE(dft_control_type), POINTER :: dft_control
84 TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
85
86 CALL timeset(routinen, handle)
87
88 nspins = SIZE(sub_env%mos_occ)
89 DO ispin = 1, nspins
90 CALL cp_fm_get_info(sub_env%mos_occ(ispin), nrow_global=nao, ncol_global=nmo_occ(ispin))
91 END DO
92
93 CALL qs_rho_get(rho_orb_struct, rho_ao=rho_ij_ao)
94 DO ispin = 1, nspins
95 CALL parallel_gemm('N', 'T', nao, nao, nmo_occ(ispin), 1.0_dp, &
96 sub_env%mos_occ(ispin), sub_env%mos_occ(ispin), &
97 0.0_dp, wfm_rho_orb)
98
99 CALL copy_fm_to_dbcsr(wfm_rho_orb, rho_ij_ao(ispin)%matrix, keep_sparsity=.true.)
100 END DO
101
102 ! take into account that all MOs are doubly occupied in spin-restricted case
103 IF (nspins == 1 .AND. (.NOT. is_rks_triplets)) &
104 CALL dbcsr_scale(rho_ij_ao(1)%matrix, 2.0_dp)
105
106 CALL get_qs_env(qs_env, dft_control=dft_control)
107
108 IF (dft_control%qs_control%gapw) THEN
109 CALL qs_rho_update_rho(rho_orb_struct, qs_env, &
110 local_rho_set=sub_env%local_rho_set, &
111 pw_env_external=sub_env%pw_env, &
112 task_list_external=sub_env%task_list_orb_soft, &
113 para_env_external=sub_env%para_env)
114 CALL prepare_gapw_den(qs_env, local_rho_set=sub_env%local_rho_set)
115 ELSEIF (dft_control%qs_control%gapw_xc) THEN
116 CALL qs_rho_update_rho(rho_orb_struct, qs_env, &
117 rho_xc_external=rho_xc_struct, &
118 local_rho_set=sub_env%local_rho_set, &
119 pw_env_external=sub_env%pw_env, &
120 task_list_external=sub_env%task_list_orb, &
121 task_list_external_soft=sub_env%task_list_orb_soft, &
122 para_env_external=sub_env%para_env)
123 CALL pw_env_get(sub_env%pw_env, auxbas_pw_pool=auxbas_pw_pool)
124 CALL qs_rho_copy(rho_xc_struct, rho_orb_struct, auxbas_pw_pool, nspins)
125 CALL prepare_gapw_den(qs_env, local_rho_set=sub_env%local_rho_set, do_rho0=.false.)
126 ELSE
127 CALL qs_rho_update_rho(rho_orb_struct, qs_env, &
128 pw_env_external=sub_env%pw_env, &
129 task_list_external=sub_env%task_list_orb, &
130 para_env_external=sub_env%para_env)
131 END IF
132
133 CALL timestop(handle)
134
136
137! **************************************************************************************************
138!> \brief Project a charge density expressed in primary basis set into the auxiliary basis set.
139!> \param rho_orb_struct response density in primary basis set
140!> \param rho_aux_fit_struct response density in auxiliary basis set (modified on exit)
141!> \param local_rho_set GAPW density of auxiliary basis set density
142!> \param qs_env Quickstep environment
143!> \param sub_env parallel (sub)group environment
144!> \param wfm_rho_orb work dense matrix with shape [nao x nao] distributed among
145!> processors of the given parallel group (modified on exit)
146!> \param wfm_rho_aux_fit work dense matrix with shape [nao_aux x nao_aux] distributed among
147!> processors of the given parallel group (modified on exit)
148!> \param wfm_aux_orb work dense matrix with shape [nao_aux x nao] distributed among
149!> processors of the given parallel group (modified on exit)
150!> \par History
151!> * 03.2017 the subroutine tddfpt_apply_admm_correction() was originally created by splitting
152!> the subroutine tddfpt_apply_hfx() in two parts [Sergey Chulkov]
153!> * 06.2018 created by splitting the subroutine tddfpt_apply_admm_correction() in two subroutines
154!> tddfpt_construct_ground_state_orb_density() and tddfpt_construct_aux_fit_density()
155!> in order to avoid code duplication [Sergey Chulkov]
156! **************************************************************************************************
157 SUBROUTINE tddfpt_construct_aux_fit_density(rho_orb_struct, rho_aux_fit_struct, local_rho_set, &
158 qs_env, sub_env, &
159 wfm_rho_orb, wfm_rho_aux_fit, wfm_aux_orb)
160 TYPE(qs_rho_type), POINTER :: rho_orb_struct, rho_aux_fit_struct
161 TYPE(local_rho_type), POINTER :: local_rho_set
162 TYPE(qs_environment_type), POINTER :: qs_env
163 TYPE(tddfpt_subgroup_env_type), INTENT(in) :: sub_env
164 TYPE(cp_fm_type), INTENT(IN) :: wfm_rho_orb, wfm_rho_aux_fit, wfm_aux_orb
165
166 CHARACTER(LEN=*), PARAMETER :: routinen = 'tddfpt_construct_aux_fit_density'
167
168 CHARACTER(LEN=default_string_length) :: basis_type
169 INTEGER :: handle, ispin, nao, nao_aux, nspins
170 REAL(kind=dp), DIMENSION(:), POINTER :: tot_rho_aux_fit_r
171 TYPE(admm_type), POINTER :: admm_env
172 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rho_ao_aux_fit, rho_ao_orb
173 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
174 POINTER :: sab_aux_fit
175 TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER :: rho_aux_fit_g
176 TYPE(pw_r3d_rs_type), DIMENSION(:), POINTER :: rho_aux_fit_r
177 TYPE(qs_ks_env_type), POINTER :: ks_env
178 TYPE(task_list_type), POINTER :: task_list
179
180 CALL timeset(routinen, handle)
181
182 cpassert(ASSOCIATED(sub_env%admm_A))
183
184 CALL get_qs_env(qs_env, ks_env=ks_env, admm_env=admm_env)
185 CALL qs_rho_get(rho_orb_struct, rho_ao=rho_ao_orb)
186 CALL qs_rho_get(rho_aux_fit_struct, rho_ao=rho_ao_aux_fit, rho_g=rho_aux_fit_g, &
187 rho_r=rho_aux_fit_r, tot_rho_r=tot_rho_aux_fit_r)
188
189 nspins = SIZE(rho_ao_orb)
190
191 IF (admm_env%do_gapw) THEN
192 basis_type = "AUX_FIT_SOFT"
193 task_list => sub_env%task_list_aux_fit_soft
194 ELSE
195 basis_type = "AUX_FIT"
196 task_list => sub_env%task_list_aux_fit
197 END IF
198
199 CALL cp_fm_get_info(sub_env%admm_A, nrow_global=nao_aux, ncol_global=nao)
200 DO ispin = 1, nspins
201 ! TO DO: consider sub_env%admm_A to be a DBCSR matrix
202 CALL copy_dbcsr_to_fm(rho_ao_orb(ispin)%matrix, wfm_rho_orb)
203 CALL parallel_gemm('N', 'N', nao_aux, nao, nao, 1.0_dp, sub_env%admm_A, &
204 wfm_rho_orb, 0.0_dp, wfm_aux_orb)
205 CALL parallel_gemm('N', 'T', nao_aux, nao_aux, nao, 1.0_dp, sub_env%admm_A, wfm_aux_orb, &
206 0.0_dp, wfm_rho_aux_fit)
207 CALL copy_fm_to_dbcsr(wfm_rho_aux_fit, rho_ao_aux_fit(ispin)%matrix, keep_sparsity=.true.)
208
209 CALL calculate_rho_elec(matrix_p=rho_ao_aux_fit(ispin)%matrix, &
210 rho=rho_aux_fit_r(ispin), rho_gspace=rho_aux_fit_g(ispin), &
211 total_rho=tot_rho_aux_fit_r(ispin), ks_env=ks_env, &
212 soft_valid=.false., basis_type=basis_type, &
213 pw_env_external=sub_env%pw_env, task_list_external=task_list)
214 END DO
215 IF (admm_env%do_gapw) THEN
216 CALL get_admm_env(qs_env%admm_env, sab_aux_fit=sab_aux_fit)
217 CALL calculate_rho_atom_coeff(qs_env, rho_ao_aux_fit, &
218 rho_atom_set=local_rho_set%rho_atom_set, &
219 qs_kind_set=admm_env%admm_gapw_env%admm_kind_set, &
220 oce=admm_env%admm_gapw_env%oce, sab=sab_aux_fit, para_env=sub_env%para_env)
221 CALL prepare_gapw_den(qs_env, local_rho_set=local_rho_set, &
222 do_rho0=.false., kind_set_external=admm_env%admm_gapw_env%admm_kind_set)
223 END IF
224
225 CALL timestop(handle)
226
228
229END MODULE qs_tddfpt2_densities
Types and set/get functions for auxiliary density matrix methods.
Definition admm_types.F:15
subroutine, public get_admm_env(admm_env, mo_derivs_aux_fit, mos_aux_fit, sab_aux_fit, sab_aux_fit_asymm, sab_aux_fit_vs_orb, matrix_s_aux_fit, matrix_s_aux_fit_kp, matrix_s_aux_fit_vs_orb, matrix_s_aux_fit_vs_orb_kp, task_list_aux_fit, matrix_ks_aux_fit, matrix_ks_aux_fit_kp, matrix_ks_aux_fit_im, matrix_ks_aux_fit_dft, matrix_ks_aux_fit_hfx, matrix_ks_aux_fit_dft_kp, matrix_ks_aux_fit_hfx_kp, rho_aux_fit, rho_aux_fit_buffer, admm_dm)
Get routine for the ADMM env.
Definition admm_types.F:593
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
DBCSR operations in CP2K.
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
subroutine, public copy_fm_to_dbcsr(fm, matrix, keep_sparsity)
Copy a BLACS matrix to a dbcsr matrix.
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
subroutine, public cp_fm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
basic linear algebra operations for full matrixes
container for various plainwaves related things
subroutine, public pw_env_get(pw_env, pw_pools, cube_info, gridlevel_info, auxbas_pw_pool, auxbas_grid, auxbas_rs_desc, auxbas_rs_grid, rs_descs, rs_grids, xc_pw_pool, vdw_pw_pool, poisson_env, interp_section)
returns the various attributes of the pw env
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Calculate the plane wave density by collocating the primitive Gaussian functions (pgf).
subroutine, public calculate_rho_elec(matrix_p, matrix_p_kp, rho, rho_gspace, total_rho, ks_env, soft_valid, compute_tau, compute_grad, basis_type, der_type, idir, task_list_external, pw_env_external)
computes the density corresponding to a given density matrix on the grid
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, rhs)
Get the QUICKSTEP environment.
subroutine, public prepare_gapw_den(qs_env, local_rho_set, do_rho0, kind_set_external)
...
Define the neighbor list data types and the corresponding functionality.
subroutine, public calculate_rho_atom_coeff(qs_env, rho_ao, rho_atom_set, qs_kind_set, oce, sab, para_env)
...
methods of the rho structure (defined in qs_rho_types)
subroutine, public qs_rho_update_rho(rho_struct, qs_env, rho_xc_external, local_rho_set, task_list_external, task_list_external_soft, pw_env_external, para_env_external)
updates rho_r and rho_g to the rhorho_ao. if use_kinetic_energy_density also computes tau_r and tau_g...
subroutine, public qs_rho_copy(rho_input, rho_output, auxbas_pw_pool, mspin)
Allocate a density structure and fill it with data from an input structure SIZE(rho_input) == mspin =...
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public qs_rho_get(rho_struct, rho_ao, rho_ao_im, rho_ao_kp, rho_ao_im_kp, rho_r, drho_r, rho_g, drho_g, tau_r, tau_g, rho_r_valid, drho_r_valid, rho_g_valid, drho_g_valid, tau_r_valid, tau_g_valid, tot_rho_r, tot_rho_g, rho_r_sccs, soft_valid, complex_rho_ao)
returns info about the density described by this object. If some representation is not available an e...
subroutine, public tddfpt_construct_aux_fit_density(rho_orb_struct, rho_aux_fit_struct, local_rho_set, qs_env, sub_env, wfm_rho_orb, wfm_rho_aux_fit, wfm_aux_orb)
Project a charge density expressed in primary basis set into the auxiliary basis set.
subroutine, public tddfpt_construct_ground_state_orb_density(rho_orb_struct, rho_xc_struct, is_rks_triplets, qs_env, sub_env, wfm_rho_orb)
Compute the ground-state charge density expressed in primary basis set.
types for task lists
stores some data used in wavefunction fitting
Definition admm_types.F:120
represent a full matrix
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...
keeps the density in various representations, keeping track of which ones are valid.