(git:dbde302)
Loading...
Searching...
No Matches
qs_dcdr_ao.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Calculate the derivatives of the MO coefficients wrt nuclear coordinates
10!> \author Sandra Luber, Edward Ditler
11! **************************************************************************************************
12
14
18 USE cp_dbcsr_api, ONLY: dbcsr_copy,&
21 dbcsr_set,&
25 USE cp_fm_types, ONLY: cp_fm_create,&
32 USE kinds, ONLY: default_string_length,&
33 dp
34 USE orbital_pointers, ONLY: ncoset
36 USE pw_env_types, ONLY: pw_env_get,&
38 USE pw_methods, ONLY: pw_axpy,&
39 pw_copy,&
40 pw_scale,&
45 USE pw_pool_types, ONLY: pw_pool_p_type,&
47 USE pw_types, ONLY: pw_c1d_gs_type,&
57 USE qs_integrate_potential, ONLY: integrate_v_dbasis,&
58 integrate_v_rspace
60 USE qs_ks_types, ONLY: get_ks_env,&
72 USE qs_rho_types, ONLY: qs_rho_create,&
76 USE qs_vxc, ONLY: qs_vxc_create
77 USE xc, ONLY: xc_calc_2nd_deriv,&
83
84!$ USE OMP_LIB, ONLY: omp_get_max_threads, omp_get_thread_num, omp_get_num_threads
85!$ USE OMP_LIB, ONLY: omp_lock_kind, &
86!$ omp_init_lock, omp_set_lock, &
87!$ omp_unset_lock, omp_destroy_lock
88
89#include "./base/base_uses.f90"
90
91 IMPLICIT NONE
92
93 PRIVATE
96 PUBLIC :: hr_mult_by_delta_1d
97
98 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_dcdr_ao'
99
100CONTAINS
101
102! **************************************************************************************************
103!> \brief Build the perturbed density matrix correction depending on the overlap derivative
104!> \param qs_env ...
105!> \param dcdr_env ...
106!> \param overlap1 Overlap derivative in AO basis
107!> \author Edward Ditler
108! **************************************************************************************************
109 SUBROUTINE apply_op_constant_term(qs_env, dcdr_env, overlap1)
110 TYPE(qs_environment_type), POINTER :: qs_env
111 TYPE(dcdr_env_type) :: dcdr_env
112 TYPE(dbcsr_p_type), OPTIONAL :: overlap1
113
114 CHARACTER(len=*), PARAMETER :: routinen = 'apply_op_constant_term'
115
116 INTEGER :: handle, ispin
117 REAL(kind=dp) :: energy_hartree
118 TYPE(cp_fm_type) :: rho_ao_fm, rho_ao_s1, rho_ao_s1_rho_ao, &
119 s1_ao
120 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rho1_ao, rho_ao
121 TYPE(pw_c1d_gs_type) :: rho1_tot_gspace, v_hartree_gspace
122 TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER :: rho1_g, rho1_g_pw
123 TYPE(pw_env_type), POINTER :: pw_env
124 TYPE(pw_poisson_type), POINTER :: poisson_env
125 TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
126 TYPE(pw_r3d_rs_type) :: v_hartree_rspace
127 TYPE(pw_r3d_rs_type), DIMENSION(:), POINTER :: rho1_r, rho_r, tau1_r, v_rspace_new, &
128 v_xc, v_xc_tau
129 TYPE(pw_r3d_rs_type), POINTER :: weights
130 TYPE(qs_rho_type), POINTER :: perturbed_density, rho
131 TYPE(section_vals_type), POINTER :: input, xc_section
132 TYPE(xc_derivative_set_type) :: deriv_set
133 TYPE(xc_rho_set_type) :: rho_set
134
135 ! Build the perturbed density matrix correction depending on the overlap derivative
136 ! P1 = C0 C1 + C1 C0
137 ! - C0_(mu j) S1_(jk) C0_(k nu)
138 ! This routine is adapted from apply_op_2_dft. There, build_dm_response builds
139 ! C0 * dCR + dCR * C0.
140 ! build_dm_response is computing $-1 * (C^0 C^1 + C^1 C^0)$ and later on in the
141 ! integration the factor 2 is applied to account for the occupancy.
142 ! The sign is negative because the kernel is on the RHS of the Sternheimer equation.
143 !
144 ! The correction factor in this routine needs to have
145 ! the opposite sign mathematically as (C0 C1 + C1 C0)
146 ! so the same sign in the code because of the $-1$ in dCR
147 ! so the opposite sign in the code because we are on the LHS of the Sternheimer equation.
148 !
149 ! This term must not go into the kernel applied by the linear response solver, because
150 ! for the (P)CG algorithm, all constant terms have to be on one side of the equations
151 ! and all solution dependent terms must be on the other side.
152
153 CALL timeset(routinen, handle)
154
155 NULLIFY (auxbas_pw_pool, pw_env, rho1_r, rho1_g_pw, &
156 v_xc, poisson_env, input, rho, rho1_g, v_xc_tau)
157
158 CALL cp_fm_create(rho_ao_fm, dcdr_env%aoao_fm_struct)
159 CALL cp_fm_create(rho_ao_s1, dcdr_env%aoao_fm_struct)
160 CALL cp_fm_create(rho_ao_s1_rho_ao, dcdr_env%aoao_fm_struct)
161 CALL cp_fm_create(s1_ao, dcdr_env%aoao_fm_struct)
162
163 IF (PRESENT(overlap1)) THEN
164 CALL copy_dbcsr_to_fm(overlap1%matrix, s1_ao)
165 ELSE
166 CALL copy_dbcsr_to_fm(dcdr_env%matrix_s1(dcdr_env%beta + 1)%matrix, s1_ao)
167 END IF
168
169 DO ispin = 1, dcdr_env%nspins
170 CALL dbcsr_set(dcdr_env%perturbed_dm_correction(ispin)%matrix, 0._dp)
171 CALL dbcsr_set(dcdr_env%matrix_apply_op_constant(ispin)%matrix, 0.0_dp)
172
173 CALL parallel_gemm('N', 'T', dcdr_env%nao, dcdr_env%nao, dcdr_env%nmo(ispin), &
174 1.0_dp, dcdr_env%mo_coeff(ispin), dcdr_env%mo_coeff(ispin), &
175 0.0_dp, rho_ao_fm)
176
177 CALL parallel_gemm('N', 'N', dcdr_env%nao, dcdr_env%nao, dcdr_env%nao, &
178 1.0_dp, rho_ao_fm, s1_ao, &
179 0.0_dp, rho_ao_s1)
180
181 CALL parallel_gemm('N', 'N', dcdr_env%nao, dcdr_env%nao, dcdr_env%nao, &
182 -1._dp, rho_ao_s1, rho_ao_fm, & ! this is the sign mentioned above.
183 0.0_dp, rho_ao_s1_rho_ao)
184
185 CALL copy_fm_to_dbcsr(rho_ao_s1_rho_ao, dcdr_env%perturbed_dm_correction(ispin)%matrix)
186 END DO
187
188 CALL cp_fm_release(rho_ao_fm)
189 CALL cp_fm_release(rho_ao_s1)
190 CALL cp_fm_release(rho_ao_s1_rho_ao)
191 CALL cp_fm_release(s1_ao)
192 ! Done building the density matrix correction
193
194 ! Build the density struct from the environment
195 NULLIFY (perturbed_density)
196 ALLOCATE (perturbed_density)
197 CALL qs_rho_create(perturbed_density)
198 CALL qs_rho_rebuild(perturbed_density, qs_env=qs_env)
199
200 ! ... set the density matrix to be the perturbed density matrix
201 CALL qs_rho_get(perturbed_density, rho_ao=rho1_ao)
202 DO ispin = 1, dcdr_env%nspins
203 CALL dbcsr_copy(rho1_ao(ispin)%matrix, dcdr_env%perturbed_dm_correction(ispin)%matrix)
204 END DO
205
206 ! ... updates rho_r and rho_g to the rho%rho_ao.
207 CALL qs_rho_update_rho(rho_struct=perturbed_density, &
208 qs_env=qs_env)
209
210 ! Also update the qs_env%rho
211 CALL get_qs_env(qs_env, rho=rho)
212 CALL qs_rho_update_rho(rho, qs_env=qs_env)
213 CALL qs_rho_get(rho, rho_ao=rho_ao, rho_r=rho_r)
214
215 energy_hartree = 0.0_dp
216
217 CALL get_qs_env(qs_env=qs_env, &
218 pw_env=pw_env, &
219 xcint_weights=weights, &
220 input=input)
221
222 ! Create the temporary grids
223 CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool, &
224 poisson_env=poisson_env)
225
226 ! Allocate deriv_set and rho_set
227 xc_section => section_vals_get_subs_vals(input, "DFT%XC")
228
229 CALL xc_prep_2nd_deriv(deriv_set, rho_set, &
230 rho_r, auxbas_pw_pool, weights, &
231 xc_section=xc_section)
232
233 ! Done with deriv_set and rho_set
234
235 ALLOCATE (v_rspace_new(dcdr_env%nspins))
236 CALL auxbas_pw_pool%create_pw(v_hartree_gspace)
237 CALL auxbas_pw_pool%create_pw(v_hartree_rspace)
238
239 ! Calculate the Hartree potential on the total density
240 CALL auxbas_pw_pool%create_pw(rho1_tot_gspace)
241
242 CALL qs_rho_get(perturbed_density, rho_g=rho1_g, rho_r=rho1_r, tau_r=tau1_r)
243 CALL pw_copy(rho1_g(1), rho1_tot_gspace)
244 DO ispin = 2, dcdr_env%nspins
245 CALL pw_axpy(rho1_g(ispin), rho1_tot_gspace)
246 END DO
247
248 CALL pw_poisson_solve(poisson_env, rho1_tot_gspace, &
249 energy_hartree, &
250 v_hartree_gspace)
251 CALL pw_transfer(v_hartree_gspace, v_hartree_rspace)
252
253 CALL auxbas_pw_pool%give_back_pw(rho1_tot_gspace)
254
255 ! Calculate the second derivative of the exchange-correlation potential
256 CALL xc_calc_2nd_deriv(v_xc, v_xc_tau, deriv_set, rho_set, &
257 rho1_r, rho1_g_pw, tau1_r, auxbas_pw_pool, &
258 weights, xc_section, gapw=.false.)
259
260 DO ispin = 1, dcdr_env%nspins
261 v_rspace_new(ispin) = v_xc(ispin)
262 END DO
263 DEALLOCATE (v_xc)
264
265 ! Done calculating the potentials
266
267 !-------------------------------!
268 ! Add both hartree and xc terms !
269 !-------------------------------!
270 CALL pw_scale(v_hartree_rspace, v_hartree_rspace%pw_grid%dvol)
271 DO ispin = 1, dcdr_env%nspins
272 CALL pw_scale(v_rspace_new(ispin), v_rspace_new(ispin)%pw_grid%dvol)
273 END DO
274
275 DO ispin = 1, dcdr_env%nspins
276 CALL dbcsr_set(dcdr_env%matrix_apply_op_constant(ispin)%matrix, 0.0_dp)
277 CALL pw_axpy(v_hartree_rspace, v_rspace_new(ispin))
278 IF (dcdr_env%nspins == 1) THEN
279 CALL pw_scale(v_rspace_new(1), 2.0_dp)
280 END IF
281
282 CALL integrate_v_rspace(v_rspace=v_rspace_new(ispin), &
283 hmat=dcdr_env%matrix_apply_op_constant(ispin), &
284 qs_env=qs_env, &
285 calculate_forces=.false.)
286 END DO
287
288 CALL auxbas_pw_pool%give_back_pw(v_hartree_gspace)
289 CALL auxbas_pw_pool%give_back_pw(v_hartree_rspace)
290 DO ispin = 1, dcdr_env%nspins
291 CALL auxbas_pw_pool%give_back_pw(v_rspace_new(ispin))
292 END DO
293 DEALLOCATE (v_rspace_new)
294
295 IF (ASSOCIATED(v_xc_tau)) THEN
296 CALL pw_scale(v_xc_tau(1), 2._dp*v_xc_tau(1)%pw_grid%dvol)
297 CALL integrate_v_rspace(v_rspace=v_xc_tau(1), &
298 hmat=dcdr_env%matrix_apply_op_constant(1), &
299 qs_env=qs_env, &
300 compute_tau=.true., &
301 calculate_forces=.false.)
302
303 CALL auxbas_pw_pool%give_back_pw(v_xc_tau(1))
304 DEALLOCATE (v_xc_tau)
305 END IF
306
307 CALL qs_rho_release(perturbed_density)
308 DEALLOCATE (perturbed_density)
309 CALL xc_rho_set_release(rho_set, auxbas_pw_pool)
310 CALL xc_dset_release(deriv_set)
311
312 CALL timestop(handle)
313
314 END SUBROUTINE apply_op_constant_term
315
316! **************************************************************************************************
317!> \brief Calculate the derivative of the Hartree term due to the core charge density
318!> \param qs_env ...
319!> \param dcdr_env ...
320!> \author Edward Ditler
321! **************************************************************************************************
322 SUBROUTINE d_core_charge_density_dr(qs_env, dcdr_env)
323 ! drho_core contribution
324 ! sum over all directions
325 ! output in ao x ao
326 TYPE(qs_environment_type), POINTER :: qs_env
327 TYPE(dcdr_env_type) :: dcdr_env
328
329 CHARACTER(len=*), PARAMETER :: routinen = 'd_core_charge_density_dR'
330
331 INTEGER :: beta, handle
332 TYPE(cp_logger_type), POINTER :: logger
333 TYPE(dft_control_type), POINTER :: dft_control
334 TYPE(pw_c1d_gs_type) :: drho_g, v_hartree_gspace
335 TYPE(pw_env_type), POINTER :: pw_env
336 TYPE(pw_poisson_type), POINTER :: poisson_env
337 TYPE(pw_pool_p_type), DIMENSION(:), POINTER :: pw_pools
338 TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
339 TYPE(pw_r3d_rs_type) :: v_hartree_rspace
340 TYPE(qs_rho_type), POINTER :: rho
341
342 CALL timeset(routinen, handle)
343
344 logger => cp_get_default_logger()
345
346 NULLIFY (pw_env, auxbas_pw_pool, pw_pools, poisson_env, dft_control, &
347 rho)
348
349 CALL get_qs_env(qs_env=qs_env, pw_env=pw_env, rho=rho, &
350 dft_control=dft_control)
351
352 CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool, poisson_env=poisson_env, &
353 pw_pools=pw_pools)
354
355 ! Create the Hartree potential grids in real and reciprocal space.
356 CALL auxbas_pw_pool%create_pw(v_hartree_gspace)
357 CALL auxbas_pw_pool%create_pw(v_hartree_rspace)
358 ! Create the grid for the derivative of the core potential
359 CALL auxbas_pw_pool%create_pw(drho_g)
360
361 DO beta = 1, 3
362 CALL pw_zero(v_hartree_gspace)
363 CALL pw_zero(v_hartree_rspace)
364 CALL pw_zero(drho_g)
365
366 ! Calculate the Hartree potential on the perturbed density and Poisson solve it
367 CALL calculate_drho_core(drho_core=drho_g, qs_env=qs_env, &
368 beta=beta, lambda=dcdr_env%lambda)
369 CALL pw_poisson_solve(poisson_env, drho_g, &
370 vhartree=v_hartree_gspace)
371 CALL pw_transfer(v_hartree_gspace, v_hartree_rspace)
372 CALL pw_scale(v_hartree_rspace, v_hartree_rspace%pw_grid%dvol)
373
374 ! Calculate the integrals
375 CALL integrate_v_rspace(v_rspace=v_hartree_rspace, &
376 hmat=dcdr_env%matrix_core_charge_1(beta), &
377 qs_env=qs_env, &
378 calculate_forces=.false.)
379 END DO
380
381 CALL auxbas_pw_pool%give_back_pw(drho_g)
382 CALL auxbas_pw_pool%give_back_pw(v_hartree_rspace)
383 CALL auxbas_pw_pool%give_back_pw(v_hartree_gspace)
384
385 CALL timestop(handle)
386 END SUBROUTINE d_core_charge_density_dr
387
388! **************************************************************************************************
389!> \brief Core Hamiltonian contributions to the operator (the pseudopotentials)
390!> \param qs_env ...
391!> \param dcdr_env ..
392!> \author Edward Ditler
393! **************************************************************************************************
394 SUBROUTINE core_dr(qs_env, dcdr_env)
395 TYPE(qs_environment_type), POINTER :: qs_env
396 TYPE(dcdr_env_type) :: dcdr_env
397
398 CHARACTER(LEN=*), PARAMETER :: routinen = 'core_dR'
399
400 CHARACTER(LEN=default_string_length) :: my_basis_type
401 INTEGER :: handle, nder
402 LOGICAL :: calculate_forces
403 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rho_ao
404 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_h, matrix_p_pass
405 TYPE(qs_ks_env_type), POINTER :: ks_env
406 TYPE(qs_rho_type), POINTER :: rho
407
408 CALL timeset(routinen, handle)
409
410 CALL get_qs_env(qs_env=qs_env, ks_env=ks_env)
411 CALL get_ks_env(ks_env=ks_env, rho=rho)
412 CALL qs_rho_get(rho, rho_ao=rho_ao)
413
414 nder = 1
415 calculate_forces = .false.
416
417 my_basis_type = "ORB"
418
419 NULLIFY (matrix_h)
420 matrix_p_pass(1:1, 1:1) => rho_ao(1:1)
421 CALL core_matrices(qs_env, matrix_h, matrix_p_pass, calculate_forces, nder, &
422 dcdr_env=dcdr_env)
423
424 CALL timestop(handle)
425
426 END SUBROUTINE core_dr
427
428! **************************************************************************************************
429!> \brief The derivatives of the basis functions going into the HXC potential wrt nuclear positions
430!> \param qs_env ...
431!> \param dcdr_env ...
432!> \author Edward Ditler
433! **************************************************************************************************
434 SUBROUTINE d_vhxc_dr(qs_env, dcdr_env)
435 TYPE(qs_environment_type), POINTER :: qs_env
436 TYPE(dcdr_env_type) :: dcdr_env
437
438 CHARACTER(len=*), PARAMETER :: routinen = 'd_vhxc_dR'
439
440 INTEGER :: handle, idir, ispin
441 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rho_ao
442 TYPE(pw_c1d_gs_type) :: drho_g_total, v_hartree_gspace
443 TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER :: drho_g
444 TYPE(pw_env_type), POINTER :: pw_env
445 TYPE(pw_poisson_type), POINTER :: poisson_env
446 TYPE(pw_pool_p_type), DIMENSION(:), POINTER :: pw_pools
447 TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
448 TYPE(pw_r3d_rs_type) :: drho_r_total, v_hartree_rspace
449 TYPE(pw_r3d_rs_type), DIMENSION(:), POINTER :: drho_r, dtau_r, rho_r, v_xc, v_xc_tau
450 TYPE(pw_r3d_rs_type), POINTER :: weights
451 TYPE(qs_rho_type), POINTER :: rho
452 TYPE(section_vals_type), POINTER :: input, xc_section
453 TYPE(xc_derivative_set_type) :: my_deriv_set
454 TYPE(xc_rho_set_type) :: my_rho_set
455
456 CALL timeset(routinen, handle)
457
458 CALL get_qs_env(qs_env=qs_env, &
459 pw_env=pw_env, &
460 input=input, &
461 xcint_weights=weights, &
462 rho=rho)
463 CALL qs_rho_get(rho, rho_ao=rho_ao, rho_r=rho_r)
464
465 xc_section => section_vals_get_subs_vals(input, "DFT%XC")
466
467 ! get the tmp grids
468 ALLOCATE (drho_r(dcdr_env%nspins))
469 ALLOCATE (drho_g(dcdr_env%nspins))
470
471 CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool, &
472 pw_pools=pw_pools, poisson_env=poisson_env)
473 CALL auxbas_pw_pool%create_pw(v_hartree_gspace)
474 CALL auxbas_pw_pool%create_pw(v_hartree_rspace)
475
476 DO ispin = 1, dcdr_env%nspins
477 CALL auxbas_pw_pool%create_pw(drho_r(ispin))
478 CALL auxbas_pw_pool%create_pw(drho_g(ispin))
479 END DO
480 CALL auxbas_pw_pool%create_pw(drho_g_total)
481 CALL auxbas_pw_pool%create_pw(drho_r_total)
482
483 DO idir = 1, 3
484 CALL pw_zero(v_hartree_gspace)
485 CALL pw_zero(v_hartree_rspace)
486 CALL pw_zero(drho_g_total)
487 CALL pw_zero(drho_r_total)
488
489 DO ispin = 1, dcdr_env%nspins
490 CALL pw_zero(drho_r(ispin))
491 CALL pw_zero(drho_g(ispin))
492
493 ! Get the density
494 CALL calculate_drho_elec_dr(matrix_p=rho_ao(ispin)%matrix, &
495 drho=drho_r(ispin), &
496 drho_gspace=drho_g(ispin), &
497 qs_env=qs_env, &
498 beta=idir, lambda=dcdr_env%lambda)
499
500 CALL pw_axpy(drho_g(ispin), drho_g_total)
501 CALL pw_axpy(drho_r(ispin), drho_r_total)
502 END DO
503 ! Get the Hartree potential corresponding to the perturbed density
504 CALL pw_poisson_solve(poisson_env, drho_g_total, &
505 vhartree=v_hartree_gspace)
506 CALL pw_transfer(v_hartree_gspace, v_hartree_rspace)
507
508 ! Get the XC potential corresponding to the perturbed density
509 CALL xc_prep_2nd_deriv(my_deriv_set, my_rho_set, &
510 rho_r, auxbas_pw_pool, weights, &
511 xc_section=xc_section)
512
513 NULLIFY (dtau_r)
514 CALL xc_calc_2nd_deriv(v_xc, v_xc_tau, my_deriv_set, my_rho_set, &
515 drho_r, drho_g, dtau_r, auxbas_pw_pool, &
516 weights, xc_section, gapw=.false.)
517 IF (ASSOCIATED(v_xc_tau)) cpabort("Meta functionals are not supported!")
518
519 CALL xc_dset_release(my_deriv_set)
520 CALL xc_rho_set_release(my_rho_set)
521
522 !-------------------------------!
523 ! Add both hartree and xc terms !
524 !-------------------------------!
525 DO ispin = 1, dcdr_env%nspins
526 ! Can the dvol be different?
527 CALL pw_scale(v_xc(ispin), v_xc(ispin)%pw_grid%dvol)
528 CALL pw_axpy(v_hartree_rspace, v_xc(ispin), v_hartree_rspace%pw_grid%dvol)
529
530 CALL integrate_v_rspace(v_rspace=v_xc(ispin), &
531 hmat=dcdr_env%matrix_d_vhxc_dR(idir, ispin), &
532 qs_env=qs_env, &
533 calculate_forces=.false.)
534
535 ! v_xc gets allocated again in xc_calc_2nd_deriv
536 CALL auxbas_pw_pool%give_back_pw(v_xc(ispin))
537 END DO ! ispin
538 DEALLOCATE (v_xc)
539 END DO ! idir
540
541 CALL auxbas_pw_pool%give_back_pw(v_hartree_gspace)
542 CALL auxbas_pw_pool%give_back_pw(v_hartree_rspace)
543 CALL auxbas_pw_pool%give_back_pw(drho_g_total)
544 CALL auxbas_pw_pool%give_back_pw(drho_r_total)
545
546 DO ispin = 1, dcdr_env%nspins
547 CALL auxbas_pw_pool%give_back_pw(drho_g(ispin))
548 CALL auxbas_pw_pool%give_back_pw(drho_r(ispin))
549 END DO
550
551 DEALLOCATE (drho_g)
552 DEALLOCATE (drho_r)
553
554 CALL timestop(handle)
555
556 END SUBROUTINE d_vhxc_dr
557
558! **************************************************************************************************
559!> \brief The derivatives of the basis functions over which the HXC potential is integrated,
560!> so < da/dR | Vhxc | b >
561!> \param qs_env ...
562!> \param dcdr_env ...
563!> \author Edward Ditler
564! **************************************************************************************************
565 SUBROUTINE vhxc_r_perturbed_basis_functions(qs_env, dcdr_env)
566 TYPE(qs_environment_type), POINTER :: qs_env
567 TYPE(dcdr_env_type) :: dcdr_env
568
569 CHARACTER(LEN=*), PARAMETER :: routinen = 'vhxc_R_perturbed_basis_functions'
570
571 INTEGER :: handle, ispin
572 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_vhxc_dbasis
573 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_p
574 TYPE(pw_env_type), POINTER :: pw_env
575 TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
576 TYPE(pw_r3d_rs_type), DIMENSION(:), POINTER :: v_hxc_r, v_tau_rspace
577 TYPE(pw_r3d_rs_type), POINTER :: v_hartree_r
578 TYPE(qs_energy_type), POINTER :: energy
579 TYPE(qs_ks_env_type), POINTER :: ks_env
580 TYPE(qs_rho_type), POINTER :: rho_struct
581 TYPE(section_vals_type), POINTER :: input, xc_section
582
583 CALL timeset(routinen, handle)
584
585 NULLIFY (rho_struct, energy, input, ks_env, pw_env, matrix_p)
586 CALL get_qs_env(qs_env, &
587 rho=rho_struct, &
588 energy=energy, &
589 input=input, &
590 ks_env=ks_env, &
591 pw_env=pw_env, &
592 v_hartree_rspace=v_hartree_r)
593 CALL qs_rho_get(rho_struct, rho_ao_kp=matrix_p)
594 xc_section => section_vals_get_subs_vals(input, "DFT%XC")
595
596 NULLIFY (auxbas_pw_pool)
597 CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
598
599 ! *** calculate the xc potential on the pw density ***
600 ! *** associates v_hxc_r if the xc potential needs to be computed.
601 ! If we do wavefunction fitting, we need the vxc_potential in the auxiliary basis set
602 NULLIFY (v_hxc_r, v_tau_rspace)
603 CALL qs_vxc_create(ks_env=ks_env, rho_struct=rho_struct, xc_section=xc_section, &
604 vxc_rho=v_hxc_r, vxc_tau=v_tau_rspace, exc=energy%exc)
605
606 DO ispin = 1, dcdr_env%nspins
607 CALL pw_scale(v_hxc_r(ispin), v_hxc_r(ispin)%pw_grid%dvol)
608
609 ! sum up potentials and integrate
610 CALL pw_axpy(v_hartree_r, v_hxc_r(ispin), 1._dp)
611
612 matrix_vhxc_dbasis => dcdr_env%matrix_vhxc_perturbed_basis(ispin, :)
613 CALL integrate_v_dbasis(v_rspace=v_hxc_r(ispin), &
614 matrix_p=matrix_p(ispin, 1)%matrix, &
615 matrix_vhxc_dbasis=matrix_vhxc_dbasis, &
616 qs_env=qs_env, &
617 lambda=dcdr_env%lambda)
618
619 CALL auxbas_pw_pool%give_back_pw(v_hxc_r(ispin))
620 END DO
621
622 DEALLOCATE (v_hxc_r)
623
624 CALL timestop(handle)
626
627! **************************************************************************************************
628!> \brief Enforce that one of the basis functions in < a | O | b > is centered on atom lambda.
629!> \param matrix ...
630!> \param qs_kind_set ...
631!> \param basis_type ...
632!> \param sab_nl ...
633!> \param lambda Atom index
634!> \param direction_Or True: < a | O | b==lambda >, False: < a==lambda | O | b >
635! **************************************************************************************************
636 SUBROUTINE hr_mult_by_delta_1d(matrix, qs_kind_set, basis_type, sab_nl, lambda, direction_Or)
637
638 TYPE(dbcsr_type), POINTER :: matrix
639 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
640 CHARACTER(LEN=*), INTENT(IN) :: basis_type
641 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
642 POINTER :: sab_nl
643 INTEGER, INTENT(IN) :: lambda
644 LOGICAL, INTENT(IN) :: direction_or
645
646 CHARACTER(len=*), PARAMETER :: routinen = 'hr_mult_by_delta_1d'
647
648 INTEGER :: handle, iatom, icol, ikind, irow, jatom, &
649 jkind, ldsab, mepos, nkind, nseta, &
650 nsetb, nthread
651 INTEGER, DIMENSION(3) :: cell
652 INTEGER, DIMENSION(:), POINTER :: la_max, la_min, lb_max, lb_min, npgfa, &
653 npgfb, nsgfa, nsgfb
654 INTEGER, DIMENSION(:, :), POINTER :: first_sgfa, first_sgfb
655 LOGICAL :: do_symmetric, found
656 REAL(kind=dp), DIMENSION(3) :: rab
657 REAL(kind=dp), DIMENSION(:), POINTER :: set_radius_a, set_radius_b
658 REAL(kind=dp), DIMENSION(:, :), POINTER :: k_block, rpgfa, rpgfb, scon_a, scon_b, &
659 zeta, zetb
660 TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER :: basis_set_list
661 TYPE(gto_basis_set_type), POINTER :: basis_set_a, basis_set_b
663 DIMENSION(:), POINTER :: nl_iterator
664
665 CALL timeset(routinen, handle)
666
667 nkind = SIZE(qs_kind_set)
668
669 ! check for symmetry
670 cpassert(SIZE(sab_nl) > 0)
671 CALL get_neighbor_list_set_p(neighbor_list_sets=sab_nl, symmetric=do_symmetric)
672
673 ! prepare basis set
674 ALLOCATE (basis_set_list(nkind))
675 CALL basis_set_list_setup(basis_set_list, basis_type, qs_kind_set)
676
677 ! *** Allocate work storage ***
678 ldsab = get_memory_usage(qs_kind_set, basis_type)
679
680 nthread = 1
681!$ nthread = omp_get_max_threads()
682 ! Iterate of neighbor list
683 CALL neighbor_list_iterator_create(nl_iterator, sab_nl, nthread=nthread)
684
685!$OMP PARALLEL DEFAULT(NONE) &
686!$OMP SHARED (nthread,ldsab,nl_iterator, do_symmetric) &
687!$OMP SHARED (ncoset,matrix,basis_set_list) &
688!$OMP SHARED (direction_or, lambda) &
689!$OMP PRIVATE (k_block,mepos,ikind,jkind,iatom,jatom,rab,cell) &
690!$OMP PRIVATE (basis_set_a,basis_set_b) &
691!$OMP PRIVATE (first_sgfa, la_max, la_min, npgfa, nsgfa, nseta, rpgfa, set_radius_a) &
692!$OMP PRIVATE (zeta, first_sgfb, lb_max, lb_min, npgfb, nsetb, rpgfb, set_radius_b, nsgfb) &
693!$OMP PRIVATE (zetb, scon_a, scon_b, irow, icol, found)
694
695 mepos = 0
696!$ mepos = omp_get_thread_num()
697
698 DO WHILE (neighbor_list_iterate(nl_iterator, mepos=mepos) == 0)
699 CALL get_iterator_info(nl_iterator, mepos=mepos, ikind=ikind, jkind=jkind, &
700 iatom=iatom, jatom=jatom, r=rab, cell=cell)
701 basis_set_a => basis_set_list(ikind)%gto_basis_set
702 IF (.NOT. ASSOCIATED(basis_set_a)) cycle
703 basis_set_b => basis_set_list(jkind)%gto_basis_set
704 IF (.NOT. ASSOCIATED(basis_set_b)) cycle
705 ! basis ikind
706 first_sgfa => basis_set_a%first_sgf
707 la_max => basis_set_a%lmax
708 la_min => basis_set_a%lmin
709 npgfa => basis_set_a%npgf
710 nseta = basis_set_a%nset
711 nsgfa => basis_set_a%nsgf_set
712 rpgfa => basis_set_a%pgf_radius
713 set_radius_a => basis_set_a%set_radius
714 scon_a => basis_set_a%scon
715 zeta => basis_set_a%zet
716 ! basis jkind
717 first_sgfb => basis_set_b%first_sgf
718 lb_max => basis_set_b%lmax
719 lb_min => basis_set_b%lmin
720 npgfb => basis_set_b%npgf
721 nsetb = basis_set_b%nset
722 nsgfb => basis_set_b%nsgf_set
723 rpgfb => basis_set_b%pgf_radius
724 set_radius_b => basis_set_b%set_radius
725 scon_b => basis_set_b%scon
726 zetb => basis_set_b%zet
727
728 IF (do_symmetric) THEN
729 IF (iatom <= jatom) THEN
730 irow = iatom
731 icol = jatom
732 ELSE
733 irow = jatom
734 icol = iatom
735 END IF
736 ELSE
737 irow = iatom
738 icol = jatom
739 END IF
740
741 NULLIFY (k_block)
742 CALL dbcsr_get_block_p(matrix, irow, icol, k_block, found)
743 cpassert(found)
744
745 IF (direction_or) THEN
746 IF (jatom /= lambda) k_block(:, :) = 0._dp
747 ELSE IF (.NOT. direction_or) THEN
748 IF (iatom /= lambda) k_block(:, :) = 0._dp
749 END IF
750 END DO
751!$OMP END PARALLEL
752 CALL neighbor_list_iterator_release(nl_iterator)
753
754 ! Release work storage
755 DEALLOCATE (basis_set_list)
756
757 CALL timestop(handle)
758
759 END SUBROUTINE hr_mult_by_delta_1d
760
761END MODULE qs_dcdr_ao
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_get_block_p(matrix, row, col, block, found, row_size, col_size)
...
subroutine, public dbcsr_set(matrix, alpha)
...
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_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
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
Provides Cartesian and spherical orbital pointers and indices.
integer, dimension(:), allocatable, public ncoset
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
functions related to the poisson solver on regular grids
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_drho_elec_dr(matrix_p, matrix_p_kp, drho, drho_gspace, qs_env, soft_valid, basis_type, beta, lambda)
Computes the gradient wrt. nuclear coordinates of a density on the grid The density is given in terms...
subroutine, public calculate_drho_core(drho_core, qs_env, beta, lambda)
Computes the derivative of the density of the core charges with respect to the nuclear coordinates on...
Calculation of the core Hamiltonian integral matrix <a|H|b> over Cartesian Gaussian-type functions.
subroutine, public core_matrices(qs_env, matrix_h, matrix_p, calculate_forces, nder, ec_env, dcdr_env, ec_env_matrices, basis_type, debug_forces, debug_stress, atcore)
...
Calculate the derivatives of the MO coefficients wrt nuclear coordinates.
Definition qs_dcdr_ao.F:13
subroutine, public core_dr(qs_env, dcdr_env)
Core Hamiltonian contributions to the operator (the pseudopotentials)
Definition qs_dcdr_ao.F:395
subroutine, public apply_op_constant_term(qs_env, dcdr_env, overlap1)
Build the perturbed density matrix correction depending on the overlap derivative.
Definition qs_dcdr_ao.F:110
subroutine, public hr_mult_by_delta_1d(matrix, qs_kind_set, basis_type, sab_nl, lambda, direction_or)
Enforce that one of the basis functions in < a | O | b > is centered on atom lambda.
Definition qs_dcdr_ao.F:637
subroutine, public vhxc_r_perturbed_basis_functions(qs_env, dcdr_env)
The derivatives of the basis functions over which the HXC potential is integrated,...
Definition qs_dcdr_ao.F:566
subroutine, public d_vhxc_dr(qs_env, dcdr_env)
The derivatives of the basis functions going into the HXC potential wrt nuclear positions.
Definition qs_dcdr_ao.F:435
subroutine, public d_core_charge_density_dr(qs_env, dcdr_env)
Calculate the derivative of the Hartree term due to the core charge density.
Definition qs_dcdr_ao.F:323
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, mimic, 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_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, 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, xcint_weights, 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, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_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, harris_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, eeq, rhs, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
Some utility functions for the calculation of integrals.
subroutine, public basis_set_list_setup(basis_set_list, basis_type, qs_kind_set)
Set up an easy accessible list of the basis sets for all kinds.
Integrate single or product functions over a potential on a RS grid.
Define the quickstep kind type and their sub types.
subroutine, public get_ks_env(ks_env, v_hartree_rspace, s_mstruct_changed, rho_changed, exc_accint, potential_changed, forces_up_to_date, complex_ks, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, kinetic, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_ks_im_kp, rho, rho_xc, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, vee, neighbor_list_id, sab_orb, sab_all, sac_ae, sac_ppl, sac_lri, sap_ppnl, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_vdw, sab_scp, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, task_list, task_list_soft, kpoints, do_kpoints, atomic_kind_set, qs_kind_set, cell, cell_ref, use_ref_cell, particle_set, energy, force, local_particles, local_molecules, molecule_kind_set, molecule_set, subsys, cp_subsys, virial, results, atprop, nkind, natom, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env, nelectron_total, nelectron_spin)
...
Type definitiona for linear response calculations.
Define the neighbor list data types and the corresponding functionality.
subroutine, public neighbor_list_iterator_create(iterator_set, nl, search, nthread)
Neighbor list iterator functions.
subroutine, public neighbor_list_iterator_release(iterator_set)
...
subroutine, public get_neighbor_list_set_p(neighbor_list_sets, nlist, symmetric)
Return the components of the first neighbor list set.
integer function, public neighbor_list_iterate(iterator_set, mepos)
...
subroutine, public get_iterator_info(iterator_set, mepos, ikind, jkind, nkind, ilist, nlist, inode, nnode, iatom, jatom, r, cell)
...
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_rebuild(rho, qs_env, rebuild_ao, rebuild_grids, admm, pw_env_external)
rebuilds rho (if necessary allocating and initializing it)
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 qs_rho_create(rho)
Allocates a new instance of rho.
subroutine, public qs_rho_release(rho_struct)
releases a rho_struct by decreasing the reference count by one and deallocating if it reaches 0 (to b...
subroutine, public qs_vxc_create(ks_env, rho_struct, xc_section, vxc_rho, vxc_tau, exc, just_energy, edisp, dispersion_env, adiabatic_rescale_factor, pw_env_external)
calculates and allocates the xc potential, already reducing it to the dependence on rho and the one o...
Definition qs_vxc.F:98
represent a group ofunctional derivatives
subroutine, public xc_dset_release(derivative_set)
releases a derivative set
contains the structure
subroutine, public xc_rho_set_release(rho_set, pw_pool)
releases the given rho_set
Exchange and Correlation functional calculations.
Definition xc.F:17
subroutine, public xc_prep_2nd_deriv(deriv_set, rho_set, rho_r, pw_pool, weights, xc_section, tau_r)
Prepare objects for the calculation of the 2nd derivatives of the density functional....
Definition xc.F:5533
subroutine, public xc_calc_2nd_deriv(v_xc, v_xc_tau, deriv_set, rho_set, rho1_r, rho1_g, tau1_r, pw_pool, weights, xc_section, gapw, vxg, do_excitations, do_sf, do_triplet, compute_virial, virial_xc)
Caller routine to calculate the second order potential in the direction of rho1_r.
Definition xc.F:924
represent a full matrix
type of a logger, at the moment it contains just a print level starting at which level it should be l...
contained for different pw related things
environment for the poisson solver
to create arrays of pools
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Provides all information about a quickstep kind.
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.
A derivative set contains the different derivatives of a xc-functional in form of a linked list.
represent a density, with all the representation and data needed to perform a functional evaluation