(git:1d3eaad)
Loading...
Searching...
No Matches
xc_pot_saop.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 saop potential
10! **************************************************************************************************
17 USE cp_dbcsr_api, ONLY: dbcsr_copy,&
24 USE cp_fm_types, ONLY: cp_fm_create,&
33 oe_gllb,&
34 oe_lb,&
35 oe_saop,&
37 USE input_section_types, ONLY: &
41 USE kinds, ONLY: dp
42 USE mathconstants, ONLY: pi
44 USE pw_env_types, ONLY: pw_env_get,&
46 USE pw_methods, ONLY: pw_axpy,&
47 pw_copy,&
48 pw_scale,&
51 USE pw_types, ONLY: pw_c1d_gs_type,&
59 USE qs_integrate_potential, ONLY: integrate_v_rspace
60 USE qs_kind_types, ONLY: get_qs_kind,&
62 USE qs_ks_atom, ONLY: update_ks_atom
67 USE qs_mo_types, ONLY: get_mo_set,&
76 USE qs_rho_types, ONLY: qs_rho_get,&
78 USE qs_vxc_atom, ONLY: calc_rho_angular,&
80 USE util, ONLY: get_limit
81 USE virial_types, ONLY: virial_type
82 USE xc, ONLY: xc_vxc_pw_create
83 USE xc_atom, ONLY: fill_rho_set,&
100 USE xc_xbecke88, ONLY: xb88_lda_info,&
102#include "./base/base_uses.f90"
103
104 IMPLICIT NONE
105
106 PRIVATE
107
108 PUBLIC :: add_saop_pot
109
110 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'xc_pot_saop'
111
112 ! should be eliminated
113 REAL(KIND=dp), PARAMETER :: alpha = 1.19_dp, beta = 0.01_dp, k_rho = 0.42_dp
114 REAL(KIND=dp), PARAMETER :: kappa = 0.804_dp, mu = 0.21951_dp, &
115 beta_ec = 0.066725_dp, gamma_saop = 0.031091_dp
116
117CONTAINS
118
119! **************************************************************************************************
120!> \brief ...
121!> \param ks_matrix ...
122!> \param qs_env ...
123!> \param oe_corr ...
124! **************************************************************************************************
125 SUBROUTINE add_saop_pot(ks_matrix, qs_env, oe_corr)
126
127 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: ks_matrix
128 TYPE(qs_environment_type), POINTER :: qs_env
129 INTEGER, INTENT(IN) :: oe_corr
130
131 INTEGER :: dft_method_id, homo, i, ispin, j, k, &
132 nspins, orb, xc_deriv_method_id, &
133 xc_rho_smooth_id
134 INTEGER, DIMENSION(2) :: ncol, nrow
135 INTEGER, DIMENSION(2, 3) :: bo
136 LOGICAL :: compute_virial, gapw, lsd
137 REAL(kind=dp) :: density_cut, efac, gradient_cut, &
138 tau_cut, we_gllb, we_lb, xc_energy
139 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: coeff_col
140 REAL(kind=dp), DIMENSION(3, 3) :: virial_xc_tmp
141 REAL(kind=dp), DIMENSION(:), POINTER :: mo_eigenvalues
142 REAL(kind=dp), DIMENSION(:, :, :), POINTER :: e_uniform
143 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:) :: single_mo_coeff
144 TYPE(cp_fm_type), POINTER :: mo_coeff
145 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: orbital_density_matrix, rho_struct_ao
146 TYPE(mo_set_type), DIMENSION(:), POINTER :: molecular_orbitals
147 TYPE(pw_c1d_gs_type) :: orbital_g
148 TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER :: rho_g
149 TYPE(pw_env_type), POINTER :: pw_env
150 TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
151 TYPE(pw_r3d_rs_type) :: orbital
152 TYPE(pw_r3d_rs_type), ALLOCATABLE, DIMENSION(:) :: vxc_gllb, vxc_saop
153 TYPE(pw_r3d_rs_type), DIMENSION(:), POINTER :: rho_r, rho_struct_r, tau, vxc_lb, &
154 vxc_tau, vxc_tmp
155 TYPE(pw_r3d_rs_type), POINTER :: weights
156 TYPE(qs_ks_env_type), POINTER :: ks_env
157 TYPE(qs_rho_type), POINTER :: rho_struct
158 TYPE(section_vals_type), POINTER :: input, xc_fun_section_orig, &
159 xc_fun_section_tmp, xc_section_orig, &
160 xc_section_tmp
161 TYPE(virial_type), POINTER :: virial
162 TYPE(xc_derivative_set_type) :: deriv_set
163 TYPE(xc_derivative_type), POINTER :: deriv
164 TYPE(xc_rho_cflags_type) :: needs
165 TYPE(xc_rho_set_type) :: rho_set
166
167 NULLIFY (ks_env, pw_env, auxbas_pw_pool, input)
168 NULLIFY (rho_g, rho_r, tau, rho_struct, e_uniform)
169 NULLIFY (vxc_lb, vxc_tmp, vxc_tau)
170 NULLIFY (mo_eigenvalues, deriv, rho_struct_r, rho_struct_ao)
171 NULLIFY (orbital_density_matrix, xc_section_tmp, xc_fun_section_tmp)
172
173 CALL get_qs_env(qs_env, &
174 ks_env=ks_env, &
175 rho=rho_struct, &
176 xcint_weights=weights, &
177 pw_env=pw_env, &
178 input=input, &
179 virial=virial, &
180 mos=molecular_orbitals)
181 compute_virial = virial%pv_calculate .AND. (.NOT. virial%pv_numer)
182 CALL section_vals_val_get(input, "DFT%QS%METHOD", i_val=dft_method_id)
183 gapw = (dft_method_id == do_method_gapw)
184
185 xc_section_orig => section_vals_get_subs_vals(input, "DFT%XC")
186 CALL section_vals_retain(xc_section_orig)
187 CALL section_vals_duplicate(xc_section_orig, xc_section_tmp)
188
189 CALL section_vals_val_get(xc_section_orig, "DENSITY_CUTOFF", &
190 r_val=density_cut)
191 CALL section_vals_val_get(xc_section_orig, "GRADIENT_CUTOFF", &
192 r_val=gradient_cut)
193 CALL section_vals_val_get(xc_section_orig, "TAU_CUTOFF", &
194 r_val=tau_cut)
195
196 CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
197
198 CALL section_vals_val_get(input, "DFT%LSD", l_val=lsd)
199 IF (lsd) THEN
200 nspins = 2
201 ELSE
202 nspins = 1
203 END IF
204
205 ALLOCATE (single_mo_coeff(nspins))
206 CALL dbcsr_allocate_matrix_set(orbital_density_matrix, nspins)
207 CALL qs_rho_get(rho_struct, rho_r=rho_struct_r, rho_ao=rho_struct_ao)
208 rho_r => rho_struct_r
209 DO ispin = 1, nspins
210 ALLOCATE (orbital_density_matrix(ispin)%matrix)
211 CALL dbcsr_copy(orbital_density_matrix(ispin)%matrix, &
212 rho_struct_ao(ispin)%matrix, "orbital density")
213 END DO
214 bo = rho_r(1)%pw_grid%bounds_local
215
216 !---------------------------!
217 ! create the density needed !
218 !---------------------------!
219 CALL xc_rho_set_create(rho_set, bo, &
220 density_cut, &
221 gradient_cut, &
222 tau_cut)
223 CALL xc_rho_cflags_setall(needs, .false.)
224 IF (lsd) THEN
225 CALL xb88_lsd_info(needs=needs)
226 needs%norm_drho = .true.
227 ELSE
228 CALL xb88_lda_info(needs=needs)
229 END IF
230 CALL section_vals_val_get(xc_section_orig, "XC_GRID%XC_DERIV", &
231 i_val=xc_deriv_method_id)
232 CALL section_vals_val_get(xc_section_orig, "XC_GRID%XC_SMOOTH_RHO", &
233 i_val=xc_rho_smooth_id)
234 CALL xc_rho_set_update(rho_set, rho_r, rho_g, tau, needs, &
235 xc_deriv_method_id, &
236 xc_rho_smooth_id, &
237 auxbas_pw_pool)
238
239 !----------------------------------------!
240 ! Construct the LB94 potential in vxc_LB !
241 !----------------------------------------!
242 xc_fun_section_orig => section_vals_get_subs_vals(xc_section_orig, &
243 "XC_FUNCTIONAL")
244 CALL section_vals_create(xc_fun_section_tmp, xc_fun_section_orig%section)
245 CALL section_vals_val_set(xc_fun_section_tmp, "_SECTION_PARAMETERS_", &
247 CALL section_vals_val_set(xc_fun_section_tmp, "XALPHA%_SECTION_PARAMETERS_", &
248 l_val=.true.)
249 CALL section_vals_set_subs_vals(xc_section_tmp, "XC_FUNCTIONAL", &
250 xc_fun_section_tmp)
251
252 cpassert(.NOT. compute_virial)
253 CALL xc_vxc_pw_create(vxc_tmp, vxc_tau, xc_energy, rho_r, rho_g, tau, &
254 xc_section_tmp, weights, auxbas_pw_pool, &
255 compute_virial=.false., virial_xc=virial_xc_tmp)
256
257 CALL section_vals_val_set(xc_fun_section_tmp, "XALPHA%_SECTION_PARAMETERS_", &
258 l_val=.false.)
259 CALL section_vals_val_set(xc_fun_section_tmp, "PZ81%_SECTION_PARAMETERS_", &
260 l_val=.true.)
261
262 cpassert(.NOT. compute_virial)
263 CALL xc_vxc_pw_create(vxc_lb, vxc_tau, xc_energy, rho_r, rho_g, tau, &
264 xc_section_tmp, weights, auxbas_pw_pool, &
265 compute_virial=.false., virial_xc=virial_xc_tmp)
266
267 DO ispin = 1, nspins
268 CALL pw_axpy(vxc_tmp(ispin), vxc_lb(ispin), alpha)
269 END DO
270
271 DO ispin = 1, nspins
272 CALL add_lb_pot(vxc_tmp(ispin)%array, rho_set, lsd, ispin)
273 CALL pw_axpy(vxc_tmp(ispin), vxc_lb(ispin), -1.0_dp)
274 END DO
275
276 !-----------------------------------------------------------------------------------!
277 ! Construct 2 times PBE one particle density from the PZ correlation energy density !
278 !-----------------------------------------------------------------------------------!
279 CALL xc_dset_create(deriv_set, local_bounds=bo)
280 CALL xc_functionals_eval(xc_fun_section_tmp, &
281 lsd=lsd, &
282 rho_set=rho_set, &
283 deriv_set=deriv_set, &
284 deriv_order=0)
285
286 deriv => xc_dset_get_derivative(deriv_set, [INTEGER::])
287 CALL xc_derivative_get(deriv, deriv_data=e_uniform)
288
289 ALLOCATE (vxc_gllb(nspins))
290 DO ispin = 1, nspins
291 CALL auxbas_pw_pool%create_pw(vxc_gllb(ispin))
292 END DO
293
294 DO ispin = 1, nspins
295 CALL calc_2excpbe(vxc_gllb(ispin)%array, rho_set, e_uniform, lsd)
296 END DO
297
298 CALL xc_dset_release(deriv_set)
299
300 CALL auxbas_pw_pool%create_pw(orbital)
301 CALL auxbas_pw_pool%create_pw(orbital_g)
302
303 DO ispin = 1, nspins
304
305 CALL get_mo_set(molecular_orbitals(ispin), &
306 mo_coeff=mo_coeff, &
307 eigenvalues=mo_eigenvalues, &
308 homo=homo)
309 CALL cp_fm_create(single_mo_coeff(ispin), &
310 mo_coeff%matrix_struct, &
311 "orbital density matrix")
312
313 CALL cp_fm_get_info(single_mo_coeff(ispin), &
314 nrow_global=nrow(ispin), ncol_global=ncol(ispin))
315 ALLOCATE (coeff_col(nrow(ispin), 1))
316
317 CALL pw_zero(vxc_tmp(ispin))
318
319 DO orb = 1, homo - 1
320
321 efac = k_rho*sqrt(mo_eigenvalues(homo) - mo_eigenvalues(orb))
322 IF (.NOT. lsd) efac = 2.0_dp*efac
323
324 CALL cp_fm_set_all(single_mo_coeff(ispin), 0.0_dp)
325 CALL cp_fm_get_submatrix(mo_coeff, coeff_col, &
326 1, orb, nrow(ispin), 1)
327 CALL cp_fm_set_submatrix(single_mo_coeff(ispin), coeff_col, &
328 1, orb)
329 CALL dbcsr_set(orbital_density_matrix(ispin)%matrix, 0.0_dp)
330 CALL cp_dbcsr_plus_fm_fm_t(orbital_density_matrix(ispin)%matrix, &
331 matrix_v=single_mo_coeff(ispin), &
332 ncol=ncol(ispin), &
333 alpha=1.0_dp)
334 CALL pw_zero(orbital)
335 CALL pw_zero(orbital_g)
336 CALL calculate_rho_elec(matrix_p=orbital_density_matrix(ispin)%matrix, &
337 rho=orbital, rho_gspace=orbital_g, &
338 ks_env=ks_env)
339
340 CALL pw_axpy(orbital, vxc_tmp(ispin), efac)
341
342 END DO
343 DEALLOCATE (coeff_col)
344
345 DO k = bo(1, 3), bo(2, 3)
346 DO j = bo(1, 2), bo(2, 2)
347 DO i = bo(1, 1), bo(2, 1)
348 IF (rho_r(ispin)%array(i, j, k) > density_cut) THEN
349 vxc_tmp(ispin)%array(i, j, k) = vxc_tmp(ispin)%array(i, j, k)/ &
350 rho_r(ispin)%array(i, j, k)
351 ELSE
352 vxc_tmp(ispin)%array(i, j, k) = 0.0_dp
353 END IF
354 END DO
355 END DO
356 END DO
357
358 CALL pw_axpy(vxc_tmp(ispin), vxc_gllb(ispin), 1.0_dp)
359
360 END DO
361
362 !---------------!
363 ! Assemble SAOP !
364 !---------------!
365 ALLOCATE (vxc_saop(nspins))
366
367 DO ispin = 1, nspins
368
369 CALL get_mo_set(molecular_orbitals(ispin), &
370 mo_coeff=mo_coeff, &
371 eigenvalues=mo_eigenvalues, &
372 homo=homo)
373 CALL auxbas_pw_pool%create_pw(vxc_saop(ispin))
374 CALL pw_zero(vxc_saop(ispin))
375
376 ALLOCATE (coeff_col(nrow(ispin), 1))
377
378 DO orb = 1, homo
379
380 we_lb = exp(-2.0_dp*(mo_eigenvalues(homo) - mo_eigenvalues(orb))**2)
381 we_gllb = 1.0_dp - we_lb
382 IF (.NOT. lsd) THEN
383 we_lb = 2.0_dp*we_lb
384 we_gllb = 2.0_dp*we_gllb
385 END IF
386
387 vxc_tmp(ispin)%array = we_lb*vxc_lb(ispin)%array + &
388 we_gllb*vxc_gllb(ispin)%array
389
390 CALL cp_fm_set_all(single_mo_coeff(ispin), 0.0_dp)
391 CALL cp_fm_get_submatrix(mo_coeff, coeff_col, &
392 1, orb, nrow(ispin), 1)
393 CALL cp_fm_set_submatrix(single_mo_coeff(ispin), coeff_col, &
394 1, orb)
395 CALL dbcsr_set(orbital_density_matrix(ispin)%matrix, 0.0_dp)
396 CALL cp_dbcsr_plus_fm_fm_t(orbital_density_matrix(ispin)%matrix, &
397 matrix_v=single_mo_coeff(ispin), &
398 ncol=ncol(ispin), &
399 alpha=1.0_dp)
400 CALL pw_zero(orbital)
401 CALL pw_zero(orbital_g)
402 CALL calculate_rho_elec(matrix_p=orbital_density_matrix(ispin)%matrix, &
403 rho=orbital, rho_gspace=orbital_g, &
404 ks_env=ks_env)
405
406 vxc_saop(ispin)%array = vxc_saop(ispin)%array + &
407 orbital%array*vxc_tmp(ispin)%array
408
409 END DO
410
411 CALL dbcsr_deallocate_matrix(orbital_density_matrix(ispin)%matrix)
412
413 DEALLOCATE (coeff_col)
414
415 DO k = bo(1, 3), bo(2, 3)
416 DO j = bo(1, 2), bo(2, 2)
417 DO i = bo(1, 1), bo(2, 1)
418 IF (rho_r(ispin)%array(i, j, k) > density_cut) THEN
419 vxc_saop(ispin)%array(i, j, k) = vxc_saop(ispin)%array(i, j, k)/ &
420 rho_r(ispin)%array(i, j, k)
421 ELSE
422 vxc_saop(ispin)%array(i, j, k) = 0.0_dp
423 END IF
424 END DO
425 END DO
426 END DO
427
428 END DO
429
430 CALL cp_fm_release(single_mo_coeff)
431
432 CALL xc_rho_set_release(rho_set, auxbas_pw_pool)
433 CALL auxbas_pw_pool%give_back_pw(orbital)
434 CALL auxbas_pw_pool%give_back_pw(orbital_g)
435
436 !--------------------!
437 ! Do the integration !
438 !--------------------!
439 DO ispin = 1, nspins
440
441 IF (oe_corr == oe_lb) THEN
442 CALL pw_copy(vxc_lb(ispin), vxc_saop(ispin))
443 ELSE IF (oe_corr == oe_gllb) THEN
444 CALL pw_copy(vxc_gllb(ispin), vxc_saop(ispin))
445 END IF
446 CALL pw_scale(vxc_saop(ispin), vxc_saop(ispin)%pw_grid%dvol)
447
448 CALL integrate_v_rspace(v_rspace=vxc_saop(ispin), pmat=rho_struct_ao(ispin), &
449 hmat=ks_matrix(ispin), qs_env=qs_env, &
450 calculate_forces=.false., &
451 gapw=gapw)
452
453 END DO
454
455 DO ispin = 1, nspins
456 CALL auxbas_pw_pool%give_back_pw(vxc_saop(ispin))
457 CALL auxbas_pw_pool%give_back_pw(vxc_gllb(ispin))
458 CALL vxc_lb(ispin)%release()
459 CALL vxc_tmp(ispin)%release()
460 END DO
461 DEALLOCATE (vxc_gllb, vxc_lb, vxc_tmp, orbital_density_matrix)
462
463 DEALLOCATE (vxc_saop)
464
465 CALL section_vals_release(xc_fun_section_tmp)
466 CALL section_vals_release(xc_section_tmp)
467 CALL section_vals_release(xc_section_orig)
468
469 !-----------------------!
470 ! Call the GAPW routine !
471 !-----------------------!
472 IF (gapw) THEN
473 CALL gapw_add_atomic_saop_pot(qs_env, oe_corr)
474 END IF
475
476 END SUBROUTINE add_saop_pot
477
478! **************************************************************************************************
479!> \brief ...
480!> \param qs_env ...
481!> \param oe_corr ...
482! **************************************************************************************************
483 SUBROUTINE gapw_add_atomic_saop_pot(qs_env, oe_corr)
484
485 TYPE(qs_environment_type), POINTER :: qs_env
486 INTEGER, INTENT(IN) :: oe_corr
487
488 INTEGER :: ia, iat, iatom, ikind, ir, ispin, na, &
489 natom, nr, ns, nspins, orb
490 INTEGER, DIMENSION(2) :: bo, homo, ncol, nrow
491 INTEGER, DIMENSION(2, 3) :: bounds
492 INTEGER, DIMENSION(:), POINTER :: atom_list
493 LOGICAL :: accint, lsd, paw_atom
494 REAL(dp), DIMENSION(:, :, :), POINTER :: tau
495 REAL(kind=dp) :: agr, alpha, density_cut, efac, exc, &
496 gradient_cut, tau_cut, we_gllb, we_lb
497 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: coeff_col, weight_h, weight_s
498 REAL(kind=dp), DIMENSION(:, :, :), POINTER :: dummy, e_uniform, rho_h, rho_s, vtau, &
499 vxc_gllb_h, vxc_gllb_s, vxc_lb_h, vxc_lb_s, vxc_saop_h, vxc_saop_s, vxc_tmp_h, vxc_tmp_s
500 REAL(kind=dp), DIMENSION(:, :, :, :), POINTER :: drho_h, drho_s, vxg
501 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
502 TYPE(cp_1d_r_p_type), DIMENSION(:), POINTER :: mo_eigenvalues
503 TYPE(cp_fm_p_type), ALLOCATABLE, DIMENSION(:) :: mo_coeff
504 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:) :: single_mo_coeff
505 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_ks, orbital_density_matrix, &
506 rho_struct_ao
507 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: ksmat, psmat
508 TYPE(dft_control_type), POINTER :: dft_control
509 TYPE(grid_atom_type), POINTER :: atomic_grid, grid_atom
510 TYPE(gto_basis_set_type), POINTER :: orb_basis
511 TYPE(harmonics_atom_type), POINTER :: harmonics
512 TYPE(local_rho_type), POINTER :: local_rho_set
513 TYPE(mo_set_type), DIMENSION(:), POINTER :: molecular_orbitals
514 TYPE(mp_para_env_type), POINTER :: para_env
515 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
516 POINTER :: sab
517 TYPE(oce_matrix_type), POINTER :: oce
518 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
519 TYPE(qs_rho_type), POINTER :: rho_structure
520 TYPE(rho_atom_coeff), DIMENSION(:), POINTER :: dr_h, dr_s, int_hh, int_ss, r_h, r_s
521 TYPE(rho_atom_coeff), DIMENSION(:, :), POINTER :: r_h_d, r_s_d
522 TYPE(rho_atom_type), DIMENSION(:), POINTER :: rho_atom_set
523 TYPE(rho_atom_type), POINTER :: rho_atom
524 TYPE(section_vals_type), POINTER :: input, xc_fun_section_orig, &
525 xc_fun_section_tmp, xc_section_orig, &
526 xc_section_tmp
527 TYPE(xc_derivative_set_type) :: deriv_set
528 TYPE(xc_derivative_type), POINTER :: deriv
529 TYPE(xc_rho_cflags_type) :: needs, needs_orbs
530 TYPE(xc_rho_set_type) :: orb_rho_set_h, orb_rho_set_s, rho_set_h, &
531 rho_set_s
532
533 NULLIFY (rho_h, rho_s, vxc_lb_h, vxc_lb_s, vxc_gllb_h, vxc_gllb_s, &
534 vxc_tmp_h, vxc_tmp_s, vtau, dummy, e_uniform, drho_h, drho_s, vxg, atom_list, &
535 atomic_kind_set, qs_kind_set, deriv, atomic_grid, rho_struct_ao, &
536 harmonics, molecular_orbitals, rho_structure, r_h, r_s, dr_h, dr_s, &
537 r_h_d, r_s_d, rho_atom_set, rho_atom, para_env, &
538 mo_eigenvalues, local_rho_set, matrix_ks, &
539 orbital_density_matrix, vxc_saop_h, vxc_saop_s)
540
541 ! tau is needed for fill_rho_set, but should never be used
542 NULLIFY (tau)
543 NULLIFY (dft_control, oce, sab)
544
545 CALL get_qs_env(qs_env, input=input, &
546 rho=rho_structure, &
547 mos=molecular_orbitals, &
548 atomic_kind_set=atomic_kind_set, &
549 qs_kind_set=qs_kind_set, &
550 rho_atom_set=rho_atom_set, &
551 matrix_ks=matrix_ks, &
552 dft_control=dft_control, &
553 para_env=para_env, &
554 oce=oce, sab_orb=sab)
555
556 CALL qs_rho_get(rho_structure, rho_ao=rho_struct_ao)
557
558 xc_section_orig => section_vals_get_subs_vals(input, "DFT%XC")
559 CALL section_vals_retain(xc_section_orig)
560 CALL section_vals_duplicate(xc_section_orig, xc_section_tmp)
561
562 accint = dft_control%qs_control%gapw_control%accurate_xcint
563
564 ! [SC] the following code can be traced back to SVN rev. 4296 (git:f97138b) that
565 ! has removed the component 'nspins' from the derived type 'dft_control_type'.
566 ! Is it worth to remove the code below in favour of 'dft_control%nspins'
567 ! since its reintroduction? Note that in case of ROKS calculations,
568 ! 'lsd == .FALSE.' but 'dft_control%nspins == 2'.
569 CALL section_vals_val_get(input, "DFT%LSD", l_val=lsd)
570 IF (lsd) THEN
571 nspins = 2
572 ELSE
573 nspins = 1
574 END IF
575
576 CALL section_vals_val_get(xc_section_orig, "DENSITY_CUTOFF", &
577 r_val=density_cut)
578 CALL section_vals_val_get(xc_section_orig, "GRADIENT_CUTOFF", &
579 r_val=gradient_cut)
580 CALL section_vals_val_get(xc_section_orig, "TAU_CUTOFF", &
581 r_val=tau_cut)
582
583 ! remap pointer
584 ns = SIZE(rho_struct_ao)
585 psmat(1:ns, 1:1) => rho_struct_ao(1:ns)
586 CALL calculate_rho_atom_coeff(qs_env, psmat, rho_atom_set, qs_kind_set, oce, sab, para_env)
587 CALL prepare_gapw_den(qs_env)
588
589 ALLOCATE (mo_coeff(nspins), single_mo_coeff(nspins), mo_eigenvalues(nspins))
590
591 CALL dbcsr_allocate_matrix_set(orbital_density_matrix, nspins)
592
593 DO ispin = 1, nspins
594 CALL get_mo_set(molecular_orbitals(ispin), &
595 mo_coeff=mo_coeff(ispin)%matrix, &
596 eigenvalues=mo_eigenvalues(ispin)%array, &
597 homo=homo(ispin))
598 CALL cp_fm_create(single_mo_coeff(ispin), &
599 mo_coeff(ispin)%matrix%matrix_struct, &
600 "orbital density matrix")
601 CALL cp_fm_get_info(single_mo_coeff(ispin), &
602 nrow_global=nrow(ispin), ncol_global=ncol(ispin))
603 ALLOCATE (orbital_density_matrix(ispin)%matrix)
604 CALL dbcsr_copy(orbital_density_matrix(ispin)%matrix, &
605 rho_struct_ao(ispin)%matrix, &
606 "orbital density")
607 END DO
608 CALL local_rho_set_create(local_rho_set)
609 CALL allocate_rho_atom_internals(local_rho_set%rho_atom_set, atomic_kind_set, &
610 qs_kind_set, dft_control, para_env)
611
612 DO ikind = 1, SIZE(atomic_kind_set)
613 CALL get_atomic_kind(atomic_kind_set(ikind), atom_list=atom_list, natom=natom)
614
615 CALL get_qs_kind(qs_kind_set(ikind), paw_atom=paw_atom, &
616 harmonics=harmonics, grid_atom=atomic_grid)
617 IF (.NOT. paw_atom) cycle
618
619 nr = atomic_grid%nr
620 na = atomic_grid%ng_sphere
621 bounds(1:2, 1:3) = 1
622 bounds(2, 1) = na
623 bounds(2, 2) = nr
624
625 CALL xc_dset_create(deriv_set, local_bounds=bounds)
626
627 CALL xc_rho_set_create(rho_set_h, bounds, density_cut, &
628 gradient_cut, tau_cut)
629 CALL xc_rho_set_create(rho_set_s, bounds, density_cut, &
630 gradient_cut, tau_cut)
631 CALL xc_rho_set_create(orb_rho_set_h, bounds, density_cut, &
632 gradient_cut, tau_cut)
633 CALL xc_rho_set_create(orb_rho_set_s, bounds, density_cut, &
634 gradient_cut, tau_cut)
635
636 CALL xc_rho_cflags_setall(needs, .false.)
637 IF (lsd) THEN
638 CALL xb88_lsd_info(needs=needs)
639 needs%norm_drho = .true.
640 ELSE
641 CALL xb88_lda_info(needs=needs)
642 END IF
643 CALL xc_rho_set_atom_update(rho_set_h, needs, nspins, bounds)
644 CALL xc_rho_set_atom_update(rho_set_s, needs, nspins, bounds)
645 CALL xc_rho_cflags_setall(needs_orbs, .false.)
646 needs_orbs%rho = .true.
647 IF (lsd) needs_orbs%rho_spin = .true.
648 CALL xc_rho_set_atom_update(orb_rho_set_h, needs, nspins, bounds)
649 CALL xc_rho_set_atom_update(orb_rho_set_s, needs, nspins, bounds)
650
651 ALLOCATE (rho_h(1:na, 1:nr, 1:nspins), rho_s(1:na, 1:nr, 1:nspins))
652 ALLOCATE (weight_h(1:na, 1:nr), weight_s(1:na, 1:nr))
653 ALLOCATE (vxc_lb_h(1:na, 1:nr, 1:nspins), vxc_lb_s(1:na, 1:nr, 1:nspins))
654 ALLOCATE (vxc_gllb_h(1:na, 1:nr, 1:nspins), vxc_gllb_s(1:na, 1:nr, 1:nspins))
655 ALLOCATE (vxc_tmp_h(1:na, 1:nr, 1:nspins), vxc_tmp_s(1:na, 1:nr, 1:nspins))
656 ALLOCATE (vxc_saop_h(1:na, 1:nr, 1:nspins), vxc_saop_s(1:na, 1:nr, 1:nspins))
657 ALLOCATE (drho_h(1:4, 1:na, 1:nr, 1:nspins), drho_s(1:4, 1:na, 1:nr, 1:nspins))
658
659 ! Distribute the atoms of this kind
660 bo = get_limit(natom, para_env%num_pe, para_env%mepos)
661
662 DO ir = 1, nr
663 DO ia = 1, na
664 weight_h(ia, ir) = atomic_grid%wr(ir)*atomic_grid%wa(ia)
665 END DO
666 IF (accint) THEN
667 alpha = dft_control%qs_control%gapw_control%aw(ikind)
668 agr = 1.0_dp - exp(-alpha*atomic_grid%rad2(ir))
669 DO ia = 1, na
670 weight_s(ia, ir) = atomic_grid%wr(ir)*atomic_grid%wa(ia)*agr
671 END DO
672 ELSE
673 DO ia = 1, na
674 weight_s(ia, ir) = atomic_grid%wr(ir)*atomic_grid%wa(ia)
675 END DO
676 END IF
677 END DO
678
679 DO iat = 1, natom !bo(1),bo(2)
680 iatom = atom_list(iat)
681
682 rho_atom => rho_atom_set(iatom)
683 NULLIFY (r_h, r_s, dr_h, dr_s, r_h_d, r_s_d)
684 CALL get_rho_atom(rho_atom=rho_atom, rho_rad_h=r_h, &
685 rho_rad_s=r_s, drho_rad_h=dr_h, &
686 drho_rad_s=dr_s, rho_rad_h_d=r_h_d, &
687 rho_rad_s_d=r_s_d)
688 rho_h = 0.0_dp
689 rho_s = 0.0_dp
690 drho_h = 0.0_dp
691 drho_s = 0.0_dp
692 DO ir = 1, nr
693 CALL calc_rho_angular(atomic_grid, harmonics, nspins, .true., &
694 ir, r_h, r_s, rho_h, rho_s, &
695 dr_h, dr_s, r_h_d, r_s_d, drho_h, drho_s)
696 END DO
697 DO ir = 1, nr
698 CALL fill_rho_set(rho_set_h, lsd, nspins, needs, rho_h, drho_h, tau, na, ir)
699 CALL fill_rho_set(rho_set_s, lsd, nspins, needs, rho_s, drho_s, tau, na, ir)
700 END DO
701
702 !-----------------------------!
703 ! 1. Slater exchange for LB94 !
704 !-----------------------------!
705 xc_fun_section_orig => section_vals_get_subs_vals(xc_section_orig, &
706 "XC_FUNCTIONAL")
707 CALL section_vals_create(xc_fun_section_tmp, xc_fun_section_orig%section)
708 CALL section_vals_val_set(xc_fun_section_tmp, "_SECTION_PARAMETERS_", &
710 CALL section_vals_val_set(xc_fun_section_tmp, "XALPHA%_SECTION_PARAMETERS_", &
711 l_val=.true.)
712 CALL section_vals_set_subs_vals(xc_section_tmp, "XC_FUNCTIONAL", &
713 xc_fun_section_tmp)
714
715 !---------------------!
716 ! Both: hard and soft !
717 !---------------------!
718 CALL xc_dset_zero_all(deriv_set)
719 CALL vxc_of_r_new(xc_fun_section_tmp, rho_set_h, deriv_set, 1, needs, &
720 weight_h, lsd, na, nr, exc, vxc_tmp_h, vxg, vtau)
721 CALL xc_dset_zero_all(deriv_set)
722 CALL vxc_of_r_new(xc_fun_section_tmp, rho_set_s, deriv_set, 1, needs, &
723 weight_s, lsd, na, nr, exc, vxc_tmp_s, vxg, vtau)
724
725 !-------------------------------------------!
726 ! 2. PZ correlation for LB94 and ec_uniform !
727 !-------------------------------------------!
728 CALL section_vals_val_set(xc_fun_section_tmp, "XALPHA%_SECTION_PARAMETERS_", &
729 l_val=.false.)
730 CALL section_vals_val_set(xc_fun_section_tmp, "PZ81%_SECTION_PARAMETERS_", &
731 l_val=.true.)
732
733 !------!
734 ! Hard !
735 !------!
736 CALL xc_dset_zero_all(deriv_set)
737 CALL vxc_of_r_new(xc_fun_section_tmp, rho_set_h, deriv_set, 1, needs, &
738 weight_h, lsd, na, nr, exc, vxc_lb_h, vxg, vtau)
739 vxc_lb_h = vxc_lb_h + alpha*vxc_tmp_h
740 DO ispin = 1, nspins
741 dummy => vxc_tmp_h(:, :, ispin:ispin)
742 CALL add_lb_pot(dummy, rho_set_h, lsd, ispin)
743 vxc_lb_h(:, :, ispin) = vxc_lb_h(:, :, ispin) - weight_h(:, :)*vxc_tmp_h(:, :, ispin)
744 END DO
745 NULLIFY (dummy)
746
747 vxc_gllb_h = 0.0_dp
748 deriv => xc_dset_get_derivative(deriv_set, [INTEGER::])
749 cpassert(ASSOCIATED(deriv))
750 CALL xc_derivative_get(deriv, deriv_data=e_uniform)
751 DO ispin = 1, nspins
752 dummy => vxc_gllb_h(:, :, ispin:ispin)
753 CALL calc_2excpbe(dummy, rho_set_h, e_uniform, lsd)
754 vxc_gllb_h(:, :, ispin) = vxc_gllb_h(:, :, ispin)*weight_h(:, :)
755 END DO
756 NULLIFY (deriv, dummy, e_uniform)
757
758 !------!
759 ! Soft !
760 !------!
761 CALL xc_dset_zero_all(deriv_set)
762 CALL vxc_of_r_new(xc_fun_section_tmp, rho_set_s, deriv_set, 1, needs, &
763 weight_s, lsd, na, nr, exc, vxc_lb_s, vxg, vtau)
764
765 vxc_lb_s = vxc_lb_s + alpha*vxc_tmp_s
766 DO ispin = 1, nspins
767 dummy => vxc_tmp_s(:, :, ispin:ispin)
768 CALL add_lb_pot(dummy, rho_set_s, lsd, ispin)
769 vxc_lb_s(:, :, ispin) = vxc_lb_s(:, :, ispin) - weight_s(:, :)*vxc_tmp_s(:, :, ispin)
770 END DO
771 NULLIFY (dummy)
772
773 vxc_gllb_s = 0.0_dp
774 deriv => xc_dset_get_derivative(deriv_set, [INTEGER::])
775 cpassert(ASSOCIATED(deriv))
776 CALL xc_derivative_get(deriv, deriv_data=e_uniform)
777 DO ispin = 1, nspins
778 dummy => vxc_gllb_s(:, :, ispin:ispin)
779 CALL calc_2excpbe(dummy, rho_set_s, e_uniform, lsd)
780 vxc_gllb_s(:, :, ispin) = vxc_gllb_s(:, :, ispin)*weight_s(:, :)
781 END DO
782 NULLIFY (deriv, dummy, e_uniform)
783
784 !------------------!
785 ! Now the orbitals !
786 !------------------!
787 vxc_tmp_h = 0.0_dp; vxc_tmp_s = 0.0_dp
788
789 DO ispin = 1, nspins
790
791 DO orb = 1, homo(ispin) - 1
792
793 ALLOCATE (coeff_col(nrow(ispin), 1))
794
795 efac = k_rho*sqrt(mo_eigenvalues(ispin)%array(homo(ispin)) - &
796 mo_eigenvalues(ispin)%array(orb))
797 IF (.NOT. lsd) efac = 2.0_dp*efac
798
799 CALL cp_fm_set_all(single_mo_coeff(ispin), 0.0_dp)
800 CALL cp_fm_get_submatrix(mo_coeff(ispin)%matrix, coeff_col, &
801 1, orb, nrow(ispin), 1)
802 CALL cp_fm_set_submatrix(single_mo_coeff(ispin), coeff_col, &
803 1, orb)
804 CALL dbcsr_set(orbital_density_matrix(ispin)%matrix, 0.0_dp)
805 CALL cp_dbcsr_plus_fm_fm_t(orbital_density_matrix(ispin)%matrix, &
806 matrix_v=single_mo_coeff(ispin), &
807 ncol=ncol(ispin), &
808 alpha=1.0_dp)
809
810 DEALLOCATE (coeff_col)
811
812 ! This calculates the CPC and density on the grids for every atom even though
813 ! we need it only for iatom at the moment. It seems that to circumvent this,
814 ! the routines must be adapted to calculate just iatom
815 ! remap pointer
816 ns = SIZE(orbital_density_matrix)
817 psmat(1:ns, 1:1) => orbital_density_matrix(1:ns)
818 CALL calculate_rho_atom_coeff(qs_env, psmat, local_rho_set%rho_atom_set, qs_kind_set, oce, sab, para_env)
819 CALL prepare_gapw_den(qs_env, local_rho_set, .false.)
820
821 rho_atom => local_rho_set%rho_atom_set(iatom)
822 NULLIFY (r_h, r_s, dr_h, dr_s, r_h_d, r_s_d)
823 CALL get_rho_atom(rho_atom=rho_atom, rho_rad_h=r_h, rho_rad_s=r_s)
824 rho_h = 0.0_dp
825 rho_s = 0.0_dp
826 drho_h = 0.0_dp
827 drho_s = 0.0_dp
828 DO ir = 1, nr
829 CALL calc_rho_angular(atomic_grid, harmonics, nspins, .false., &
830 ir, r_h, r_s, rho_h, rho_s, &
831 dr_h, dr_s, r_h_d, r_s_d, drho_h, drho_s)
832 END DO
833 DO ir = 1, nr
834 CALL fill_rho_set(orb_rho_set_h, lsd, nspins, needs_orbs, rho_h, drho_h, tau, na, ir)
835 CALL fill_rho_set(orb_rho_set_s, lsd, nspins, needs_orbs, rho_s, drho_s, tau, na, ir)
836 END DO
837
838 IF (lsd) THEN
839 IF (ispin == 1) THEN
840 vxc_tmp_h(:, :, 1) = vxc_tmp_h(:, :, 1) + efac*orb_rho_set_h%rhoa(:, :, 1)
841 vxc_tmp_s(:, :, 1) = vxc_tmp_s(:, :, 1) + efac*orb_rho_set_s%rhoa(:, :, 1)
842 ELSE
843 vxc_tmp_h(:, :, 2) = vxc_tmp_h(:, :, 2) + efac*orb_rho_set_h%rhob(:, :, 1)
844 vxc_tmp_s(:, :, 2) = vxc_tmp_s(:, :, 2) + efac*orb_rho_set_s%rhob(:, :, 1)
845 END IF
846 ELSE
847 vxc_tmp_h(:, :, 1) = vxc_tmp_h(:, :, 1) + efac*orb_rho_set_h%rho(:, :, 1)
848 vxc_tmp_s(:, :, 1) = vxc_tmp_s(:, :, 1) + efac*orb_rho_set_s%rho(:, :, 1)
849 END IF
850
851 END DO ! orb
852
853 END DO ! ispin
854
855 IF (lsd) THEN
856 DO ir = 1, nr
857 DO ia = 1, na
858 IF (rho_set_h%rhoa(ia, ir, 1) > rho_set_h%rho_cutoff) &
859 vxc_gllb_h(ia, ir, 1) = vxc_gllb_h(ia, ir, 1) + &
860 weight_h(ia, ir)*vxc_tmp_h(ia, ir, 1)/rho_set_h%rhoa(ia, ir, 1)
861 IF (rho_set_h%rhob(ia, ir, 1) > rho_set_h%rho_cutoff) &
862 vxc_gllb_h(ia, ir, 2) = vxc_gllb_h(ia, ir, 2) + &
863 weight_h(ia, ir)*vxc_tmp_h(ia, ir, 2)/rho_set_h%rhob(ia, ir, 1)
864 IF (rho_set_s%rhoa(ia, ir, 1) > rho_set_s%rho_cutoff) &
865 vxc_gllb_s(ia, ir, 1) = vxc_gllb_s(ia, ir, 1) + &
866 weight_s(ia, ir)*vxc_tmp_s(ia, ir, 1)/rho_set_s%rhoa(ia, ir, 1)
867 IF (rho_set_s%rhob(ia, ir, 1) > rho_set_s%rho_cutoff) &
868 vxc_gllb_s(ia, ir, 2) = vxc_gllb_s(ia, ir, 2) + &
869 weight_s(ia, ir)*vxc_tmp_s(ia, ir, 2)/rho_set_s%rhob(ia, ir, 1)
870 END DO
871 END DO
872 ELSE
873 DO ir = 1, nr
874 DO ia = 1, na
875 IF (rho_set_h%rho(ia, ir, 1) > rho_set_h%rho_cutoff) &
876 vxc_gllb_h(ia, ir, 1) = vxc_gllb_h(ia, ir, 1) + &
877 weight_h(ia, ir)*vxc_tmp_h(ia, ir, 1)/rho_set_h%rho(ia, ir, 1)
878 IF (rho_set_s%rho(ia, ir, 1) > rho_set_s%rho_cutoff) &
879 vxc_gllb_s(ia, ir, 1) = vxc_gllb_s(ia, ir, 1) + &
880 weight_s(ia, ir)*vxc_tmp_s(ia, ir, 1)/rho_set_s%rho(ia, ir, 1)
881 END DO
882 END DO
883 END IF
884
885 vxc_saop_h = 0.0_dp; vxc_saop_s = 0.0_dp
886
887 DO ispin = 1, nspins
888
889 DO orb = 1, homo(ispin)
890
891 ALLOCATE (coeff_col(nrow(ispin), 1))
892
893 we_lb = exp(-2.0_dp*(mo_eigenvalues(ispin)%array(homo(ispin)) - &
894 mo_eigenvalues(ispin)%array(orb))**2)
895 we_gllb = 1.0_dp - we_lb
896 IF (.NOT. lsd) THEN
897 we_lb = 2.0_dp*we_lb
898 we_gllb = 2.0_dp*we_gllb
899 END IF
900
901 vxc_tmp_h(:, :, ispin) = we_lb*vxc_lb_h(:, :, ispin) + &
902 we_gllb*vxc_gllb_h(:, :, ispin)
903 vxc_tmp_s(:, :, ispin) = we_lb*vxc_lb_s(:, :, ispin) + &
904 we_gllb*vxc_gllb_s(:, :, ispin)
905
906 CALL cp_fm_set_all(single_mo_coeff(ispin), 0.0_dp)
907 CALL cp_fm_get_submatrix(mo_coeff(ispin)%matrix, coeff_col, &
908 1, orb, nrow(ispin), 1)
909 CALL cp_fm_set_submatrix(single_mo_coeff(ispin), coeff_col, &
910 1, orb)
911 CALL dbcsr_set(orbital_density_matrix(ispin)%matrix, 0.0_dp)
912 CALL cp_dbcsr_plus_fm_fm_t(orbital_density_matrix(ispin)%matrix, &
913 matrix_v=single_mo_coeff(ispin), &
914 ncol=ncol(ispin), &
915 alpha=1.0_dp)
916
917 DEALLOCATE (coeff_col)
918
919 ! This calculates the CPC and density on the grids for every atom even though
920 ! we need it only for iatom at the moment. It seems that to circumvent this,
921 ! the routines must be adapted to calculate just iatom
922 ! remap pointer
923 ns = SIZE(orbital_density_matrix)
924 psmat(1:ns, 1:1) => orbital_density_matrix(1:ns)
925 CALL calculate_rho_atom_coeff(qs_env, psmat, local_rho_set%rho_atom_set, qs_kind_set, oce, sab, para_env)
926 CALL prepare_gapw_den(qs_env, local_rho_set, .false.)
927
928 rho_atom => local_rho_set%rho_atom_set(iatom)
929 NULLIFY (r_h, r_s, dr_h, dr_s, r_h_d, r_s_d)
930 CALL get_rho_atom(rho_atom=rho_atom, rho_rad_h=r_h, rho_rad_s=r_s)
931 rho_h = 0.0_dp
932 rho_s = 0.0_dp
933 drho_h = 0.0_dp
934 drho_s = 0.0_dp
935 DO ir = 1, nr
936 CALL calc_rho_angular(atomic_grid, harmonics, nspins, .false., &
937 ir, r_h, r_s, rho_h, rho_s, &
938 dr_h, dr_s, r_h_d, r_s_d, drho_h, drho_s)
939 END DO
940 DO ir = 1, nr
941 CALL fill_rho_set(orb_rho_set_h, lsd, nspins, needs_orbs, rho_h, drho_h, tau, na, ir)
942 CALL fill_rho_set(orb_rho_set_s, lsd, nspins, needs_orbs, rho_s, drho_s, tau, na, ir)
943 END DO
944
945 IF (lsd) THEN
946 IF (ispin == 1) THEN
947 vxc_saop_h(:, :, 1) = vxc_saop_h(:, :, 1) + vxc_tmp_h(:, :, 1)*orb_rho_set_h%rhoa(:, :, 1)
948 vxc_saop_s(:, :, 1) = vxc_saop_s(:, :, 1) + vxc_tmp_s(:, :, 1)*orb_rho_set_s%rhoa(:, :, 1)
949 ELSE
950 vxc_saop_h(:, :, 2) = vxc_saop_h(:, :, 2) + vxc_tmp_h(:, :, 2)*orb_rho_set_h%rhob(:, :, 1)
951 vxc_saop_s(:, :, 2) = vxc_saop_s(:, :, 2) + vxc_tmp_s(:, :, 2)*orb_rho_set_s%rhob(:, :, 1)
952 END IF
953 ELSE
954 vxc_saop_h(:, :, 1) = vxc_saop_h(:, :, 1) + vxc_tmp_h(:, :, 1)*orb_rho_set_h%rho(:, :, 1)
955 vxc_saop_s(:, :, 1) = vxc_saop_s(:, :, 1) + vxc_tmp_s(:, :, 1)*orb_rho_set_s%rho(:, :, 1)
956 END IF
957
958 END DO ! orb
959
960 END DO ! ispin
961
962 IF (lsd) THEN
963 DO ir = 1, nr
964 DO ia = 1, na
965 IF (rho_set_h%rhoa(ia, ir, 1) > rho_set_h%rho_cutoff) THEN
966 vxc_saop_h(ia, ir, 1) = vxc_saop_h(ia, ir, 1)/rho_set_h%rhoa(ia, ir, 1)
967 ELSE
968 vxc_saop_h(ia, ir, 1) = 0.0_dp
969 END IF
970 IF (rho_set_h%rhob(ia, ir, 1) > rho_set_h%rho_cutoff) THEN
971 vxc_saop_h(ia, ir, 2) = vxc_saop_h(ia, ir, 2)/rho_set_h%rhob(ia, ir, 1)
972 ELSE
973 vxc_saop_h(ia, ir, 2) = 0.0_dp
974 END IF
975 IF (rho_set_s%rhoa(ia, ir, 1) > rho_set_s%rho_cutoff) THEN
976 vxc_saop_s(ia, ir, 1) = vxc_saop_s(ia, ir, 1)/rho_set_s%rhoa(ia, ir, 1)
977 ELSE
978 vxc_saop_s(ia, ir, 1) = 0.0_dp
979 END IF
980 IF (rho_set_s%rhob(ia, ir, 1) > rho_set_s%rho_cutoff) THEN
981 vxc_saop_s(ia, ir, 2) = vxc_saop_s(ia, ir, 2)/rho_set_s%rhob(ia, ir, 1)
982 ELSE
983 vxc_saop_s(ia, ir, 2) = 0.0_dp
984 END IF
985 END DO
986 END DO
987 ELSE
988 DO ir = 1, nr
989 DO ia = 1, na
990 IF (rho_set_h%rho(ia, ir, 1) > rho_set_h%rho_cutoff) THEN
991 vxc_saop_h(ia, ir, 1) = vxc_saop_h(ia, ir, 1)/rho_set_h%rho(ia, ir, 1)
992 ELSE
993 vxc_saop_h(ia, ir, 1) = 0.0_dp
994 END IF
995 IF (rho_set_s%rho(ia, ir, 1) > rho_set_s%rho_cutoff) THEN
996 vxc_saop_s(ia, ir, 1) = vxc_saop_s(ia, ir, 1)/rho_set_s%rho(ia, ir, 1)
997 ELSE
998 vxc_saop_s(ia, ir, 1) = 0.0_dp
999 END IF
1000 END DO
1001 END DO
1002 END IF
1003
1004 rho_atom => rho_atom_set(iatom)
1005 CALL get_rho_atom(rho_atom=rho_atom, ga_vlocal_gb_h=int_hh, ga_vlocal_gb_s=int_ss)
1006 CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis, &
1007 harmonics=harmonics, grid_atom=grid_atom)
1008 SELECT CASE (oe_corr)
1009 CASE (oe_lb)
1010 CALL gavxcgb_nogc(vxc_lb_h, vxc_lb_s, int_hh, int_ss, grid_atom, orb_basis, harmonics, nspins)
1011 CASE (oe_gllb)
1012 CALL gavxcgb_nogc(vxc_gllb_h, vxc_gllb_s, int_hh, int_ss, grid_atom, orb_basis, harmonics, nspins)
1013 CASE (oe_saop)
1014 CALL gavxcgb_nogc(vxc_saop_h, vxc_saop_s, int_hh, int_ss, grid_atom, orb_basis, harmonics, nspins)
1015 CASE default
1016 cpabort("Unknown correction!")
1017 END SELECT
1018
1019 END DO
1020
1021 DEALLOCATE (rho_h, rho_s, weight_h, weight_s)
1022 DEALLOCATE (vxc_lb_h, vxc_lb_s)
1023 DEALLOCATE (vxc_gllb_h, vxc_gllb_s)
1024 DEALLOCATE (vxc_tmp_h, vxc_tmp_s)
1025 DEALLOCATE (vxc_saop_h, vxc_saop_s)
1026 DEALLOCATE (drho_h, drho_s)
1027
1028 CALL xc_dset_release(deriv_set)
1029 CALL xc_rho_set_release(rho_set_h)
1030 CALL xc_rho_set_release(rho_set_s)
1031 CALL xc_rho_set_release(orb_rho_set_h)
1032 CALL xc_rho_set_release(orb_rho_set_s)
1033
1034 END DO
1035
1036 ! remap pointer
1037 ns = SIZE(matrix_ks)
1038 ksmat(1:ns, 1:1) => matrix_ks(1:ns)
1039 ns = SIZE(rho_struct_ao)
1040 psmat(1:ns, 1:1) => rho_struct_ao(1:ns)
1041
1042 CALL update_ks_atom(qs_env, ksmat, psmat, forces=.false.)
1043
1044 !---------!
1045 ! Cleanup !
1046 !---------!
1047 CALL section_vals_release(xc_fun_section_tmp)
1048 CALL section_vals_release(xc_section_tmp)
1049 CALL section_vals_release(xc_section_orig)
1050
1051 CALL local_rho_set_release(local_rho_set)
1052 CALL cp_fm_release(single_mo_coeff)
1053 DEALLOCATE (mo_coeff, mo_eigenvalues)
1054 CALL dbcsr_deallocate_matrix_set(orbital_density_matrix)
1055
1056 END SUBROUTINE gapw_add_atomic_saop_pot
1057
1058! **************************************************************************************************
1059!> \brief ...
1060!> \param pot ...
1061!> \param rho_set ...
1062!> \param lsd ...
1063!> \param spin ...
1064! **************************************************************************************************
1065 SUBROUTINE add_lb_pot(pot, rho_set, lsd, spin)
1066
1067 REAL(kind=dp), DIMENSION(:, :, :), POINTER :: pot
1068 TYPE(xc_rho_set_type), INTENT(IN) :: rho_set
1069 LOGICAL, INTENT(IN) :: lsd
1070 INTEGER, INTENT(IN) :: spin
1071
1072 REAL(kind=dp), PARAMETER :: ob3 = 1.0_dp/3.0_dp
1073
1074 INTEGER :: i, j, k
1075 INTEGER, DIMENSION(2, 3) :: bo
1076 REAL(kind=dp) :: n, n_13, x, x2
1077
1078 bo = rho_set%local_bounds
1079
1080 DO k = bo(1, 3), bo(2, 3)
1081 DO j = bo(1, 2), bo(2, 2)
1082 DO i = bo(1, 1), bo(2, 1)
1083 IF (.NOT. lsd) THEN
1084 IF (rho_set%rho(i, j, k) > rho_set%rho_cutoff) THEN
1085 n = rho_set%rho(i, j, k)/2.0_dp
1086 n_13 = n**ob3
1087 x = (rho_set%norm_drho(i, j, k)/2.0_dp)/(n*n_13)
1088 x2 = x*x
1089 pot(i, j, k) = beta*x2*n_13/(1.0_dp + 3.0_dp*beta*x*log(x + sqrt(x2 + 1.0_dp)))
1090 END IF
1091 ELSE
1092 IF (spin == 1) THEN
1093 IF (rho_set%rhoa(i, j, k) > rho_set%rho_cutoff) THEN
1094 n_13 = rho_set%rhoa_1_3(i, j, k)
1095 x = rho_set%norm_drhoa(i, j, k)/(rho_set%rhoa(i, j, k)*n_13)
1096 x2 = x*x
1097 pot(i, j, k) = beta*x2*n_13/(1.0_dp + 3.0_dp*beta*x*log(sqrt(x2 + 1.0_dp) + x))
1098 END IF
1099 ELSE IF (spin == 2) THEN
1100 IF (rho_set%rhob(i, j, k) > rho_set%rho_cutoff) THEN
1101 n_13 = rho_set%rhob_1_3(i, j, k)
1102 x = rho_set%norm_drhob(i, j, k)/(rho_set%rhob(i, j, k)*n_13)
1103 x2 = x*x
1104 pot(i, j, k) = beta*x2*n_13/(1.0_dp + 3.0_dp*beta*x*log(sqrt(x2 + 1.0_dp) + x))
1105 END IF
1106 END IF
1107 END IF
1108 END DO
1109 END DO
1110 END DO
1111
1112 END SUBROUTINE add_lb_pot
1113
1114! **************************************************************************************************
1115!> \brief ...
1116!> \param pot ...
1117!> \param rho_set ...
1118!> \param e_uniform ...
1119!> \param lsd ...
1120! **************************************************************************************************
1121 SUBROUTINE calc_2excpbe(pot, rho_set, e_uniform, lsd)
1122
1123 REAL(kind=dp), DIMENSION(:, :, :), POINTER :: pot
1124 TYPE(xc_rho_set_type), INTENT(IN) :: rho_set
1125 REAL(kind=dp), DIMENSION(:, :, :), POINTER :: e_uniform
1126 LOGICAL, INTENT(IN) :: lsd
1127
1128 INTEGER :: i, j, k
1129 INTEGER, DIMENSION(2, 3) :: bo
1130 REAL(kind=dp) :: e_unif, rho
1131
1132 bo = rho_set%local_bounds
1133
1134 DO k = bo(1, 3), bo(2, 3)
1135 DO j = bo(1, 2), bo(2, 2)
1136 DO i = bo(1, 1), bo(2, 1)
1137 IF (.NOT. lsd) THEN
1138 IF (rho_set%rho(i, j, k) > rho_set%rho_cutoff) THEN
1139 e_unif = e_uniform(i, j, k)/rho_set%rho(i, j, k)
1140 ELSE
1141 e_unif = 0.0_dp
1142 END IF
1143 pot(i, j, k) = &
1144 2.0_dp* &
1145 calc_ecpbe_r(rho_set%rho(i, j, k), rho_set%norm_drho(i, j, k), &
1146 e_unif, rho_set%rho_cutoff, rho_set%drho_cutoff) + &
1147 2.0_dp* &
1148 calc_expbe_r(rho_set%rho(i, j, k), rho_set%norm_drho(i, j, k), &
1149 rho_set%rho_cutoff, rho_set%drho_cutoff)
1150 ELSE
1151 rho = rho_set%rhoa(i, j, k) + rho_set%rhob(i, j, k)
1152 IF (rho > rho_set%rho_cutoff) THEN
1153 e_unif = e_uniform(i, j, k)/rho
1154 ELSE
1155 e_unif = 0.0_dp
1156 END IF
1157 pot(i, j, k) = &
1158 2.0_dp* &
1159 calc_ecpbe_u(rho_set%rhoa(i, j, k), rho_set%rhob(i, j, k), rho_set%norm_drho(i, j, k), &
1160 e_unif, &
1161 rho_set%rho_cutoff, rho_set%drho_cutoff) + &
1162 2.0_dp* &
1163 calc_expbe_u(rho_set%rhoa(i, j, k), rho_set%rhob(i, j, k), rho_set%norm_drho(i, j, k), &
1164 rho_set%rho_cutoff, rho_set%drho_cutoff)
1165 END IF
1166 END DO
1167 END DO
1168 END DO
1169
1170 END SUBROUTINE calc_2excpbe
1171
1172! **************************************************************************************************
1173!> \brief ...
1174!> \param ra ...
1175!> \param rb ...
1176!> \param ngr ...
1177!> \param ec_unif ...
1178!> \param rc ...
1179!> \param ngrc ...
1180!> \return ...
1181! **************************************************************************************************
1182 FUNCTION calc_ecpbe_u(ra, rb, ngr, ec_unif, rc, ngrc) RESULT(res)
1183
1184 REAL(kind=dp), INTENT(in) :: ra, rb, ngr, ec_unif, rc, ngrc
1185 REAL(kind=dp) :: res
1186
1187 REAL(kind=dp), PARAMETER :: ob3 = 1.0_dp/3.0_dp, tb3 = 2.0_dp/3.0_dp
1188
1189 REAL(kind=dp) :: a, at2, h, kf, kl, ks, phi, phi3, r, t2, &
1190 zeta
1191
1192 r = ra + rb
1193 h = 0.0_dp
1194 IF (r > rc .AND. ngr > ngrc) THEN
1195 zeta = (ra - rb)/r
1196 IF (zeta > 1.0_dp) zeta = 1.0_dp ! machine precision problem
1197 IF (zeta < -1.0_dp) zeta = -1.0_dp ! machine precision problem
1198 phi = ((1.0_dp + zeta)**tb3 + (1.0_dp - zeta)**tb3)/2.0_dp
1199 phi3 = phi*phi*phi
1200 kf = (3.0_dp*r*pi*pi)**ob3
1201 ks = sqrt(4.0_dp*kf/pi)
1202 t2 = (ngr/(2.0_dp*phi*ks*r))**2
1203 a = beta_ec/gamma_saop/(exp(-ec_unif/(gamma_saop*phi3)) - 1.0_dp)
1204 at2 = a*t2
1205 kl = (1.0_dp + at2)/(1.0_dp + at2 + at2*at2)
1206 h = gamma_saop*log(1.0_dp + beta_ec/gamma_saop*t2*kl)
1207 END IF
1208 res = ec_unif + h
1209
1210 END FUNCTION calc_ecpbe_u
1211
1212! **************************************************************************************************
1213!> \brief ...
1214!> \param r ...
1215!> \param ngr ...
1216!> \param ec_unif ...
1217!> \param rc ...
1218!> \param ngrc ...
1219!> \return ...
1220! **************************************************************************************************
1221 FUNCTION calc_ecpbe_r(r, ngr, ec_unif, rc, ngrc) RESULT(res)
1222
1223 REAL(kind=dp), INTENT(in) :: r, ngr, ec_unif, rc, ngrc
1224 REAL(kind=dp) :: res
1225
1226 REAL(kind=dp) :: a, at2, h, kf, kl, ks, t2
1227
1228 h = 0.0_dp
1229 IF (r > rc .AND. ngr > ngrc) THEN
1230 kf = (3.0_dp*r*pi*pi)**(1.0_dp/3.0_dp)
1231 ks = sqrt(4.0_dp*kf/pi)
1232 t2 = (ngr/(2.0_dp*ks*r))**2
1233 a = beta_ec/gamma_saop/(exp(-ec_unif/gamma_saop) - 1.0_dp)
1234 at2 = a*t2
1235 kl = (1.0_dp + at2)/(1.0_dp + at2 + at2*at2)
1236 h = gamma_saop*log(1.0_dp + beta_ec/gamma_saop*t2*kl)
1237 END IF
1238 res = ec_unif + h
1239
1240 END FUNCTION calc_ecpbe_r
1241
1242! **************************************************************************************************
1243!> \brief ...
1244!> \param ra ...
1245!> \param rb ...
1246!> \param ngr ...
1247!> \param rc ...
1248!> \param ngrc ...
1249!> \return ...
1250! **************************************************************************************************
1251 FUNCTION calc_expbe_u(ra, rb, ngr, rc, ngrc) RESULT(res)
1252
1253 REAL(kind=dp), INTENT(in) :: ra, rb, ngr, rc, ngrc
1254 REAL(kind=dp) :: res
1255
1256 REAL(kind=dp) :: r
1257
1258 r = ra + rb
1259 res = calc_expbe_r(r, ngr, rc, ngrc)
1260
1261 END FUNCTION calc_expbe_u
1262
1263! **************************************************************************************************
1264!> \brief ...
1265!> \param r ...
1266!> \param ngr ...
1267!> \param rc ...
1268!> \param ngrc ...
1269!> \return ...
1270! **************************************************************************************************
1271 FUNCTION calc_expbe_r(r, ngr, rc, ngrc) RESULT(res)
1272
1273 REAL(kind=dp), INTENT(in) :: r, ngr, rc, ngrc
1274 REAL(kind=dp) :: res
1275
1276 REAL(kind=dp) :: ex_unif, fx, kf, s
1277
1278 IF (r > rc) THEN
1279 kf = (3.0_dp*r*pi*pi)**(1.0_dp/3.0_dp)
1280 ex_unif = -3.0_dp*kf/(4.0_dp*pi)
1281 fx = 1.0_dp
1282 IF (ngr > ngrc) THEN
1283 s = ngr/(2.0_dp*kf*r)
1284 fx = fx + kappa - kappa/(1.0_dp + mu*s*s/kappa)
1285 END IF
1286 res = ex_unif*fx
1287 ELSE
1288 res = 0.0_dp
1289 END IF
1290
1291 END FUNCTION calc_expbe_r
1292
1293END MODULE xc_pot_saop
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
various utilities that regard array of different kinds: output, allocation,... maybe it is not a good...
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_deallocate_matrix(matrix)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_set(matrix, alpha)
...
DBCSR operations in CP2K.
subroutine, public cp_dbcsr_plus_fm_fm_t(sparse_matrix, matrix_v, matrix_g, ncol, alpha, keep_sparsity, symmetry_mode)
performs the multiplication sparse_matrix+dense_mat*dens_mat^T if matrix_g is not explicitly given,...
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
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
subroutine, public cp_fm_set_submatrix(fm, new_values, start_row, start_col, n_rows, n_cols, alpha, beta, transpose)
sets a submatrix of a full matrix fm(start_row:start_row+n_rows,start_col:start_col+n_cols) = alpha*o...
subroutine, public cp_fm_set_all(matrix, alpha, beta)
set all elements of a matrix to the same value, and optionally the diagonal to a different one
subroutine, public cp_fm_get_submatrix(fm, target_m, start_row, start_col, n_rows, n_cols, transpose)
gets a submatrix of a full matrix op(target_m)(1:n_rows,1:n_cols) =fm(start_row:start_row+n_rows,...
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public oe_saop
integer, parameter, public xc_funct_no_shortcut
integer, parameter, public do_method_gapw
integer, parameter, public oe_lb
integer, parameter, public oe_gllb
objects that represent the structure of input sections and the data contained in an input section
recursive subroutine, public section_vals_create(section_vals, section)
creates a object where to store the values of a section
subroutine, public section_vals_val_set(section_vals, keyword_name, i_rep_section, i_rep_val, val, l_val, i_val, r_val, c_val, l_vals_ptr, i_vals_ptr, r_vals_ptr, c_vals_ptr)
sets the requested value
subroutine, public section_vals_retain(section_vals)
retains the given section values (see doc/ReferenceCounting.html)
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
subroutine, public section_vals_set_subs_vals(section_vals, subsection_name, new_section_vals, i_rep_section)
replaces of the requested subsection with the one given
subroutine, public section_vals_duplicate(section_vals_in, section_vals_out, i_rep_start, i_rep_end)
creates a deep copy from section_vals_in to section_vals_out
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
recursive subroutine, public section_vals_release(section_vals)
releases the given object
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
Interface to the message passing library MPI.
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, 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.
subroutine, public prepare_gapw_den(qs_env, local_rho_set, do_rho0, kind_set_external, pw_env_sub)
...
Integrate single or product functions over a potential on a RS grid.
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, cneo_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zatom, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, monovalent, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
routines that build the Kohn-Sham matrix contributions coming from local atomic densities
Definition qs_ks_atom.F:12
subroutine, public update_ks_atom(qs_env, ksmat, pmat, forces, tddft, rho_atom_external, kind_set_external, oce_external, sab_external, kscale, kintegral, kforce, fscale)
The correction to the KS matrix due to the GAPW local terms to the hartree and XC contributions is he...
Definition qs_ks_atom.F:110
subroutine, public local_rho_set_create(local_rho_set)
...
subroutine, public local_rho_set_release(local_rho_set)
...
Definition and initialisation of the mo data type.
Definition qs_mo_types.F:22
subroutine, public get_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, mo_coeff, mo_coeff_b, uniform_occupation, kts, mu, flexible_electron_count)
Get the components of a MO set data structure.
Define the neighbor list data types and the corresponding functionality.
subroutine, public allocate_rho_atom_internals(rho_atom_set, atomic_kind_set, qs_kind_set, dft_control, para_env)
...
subroutine, public calculate_rho_atom_coeff(qs_env, rho_ao, rho_atom_set, qs_kind_set, oce, sab, para_env)
...
subroutine, public get_rho_atom(rho_atom, cpc_h, cpc_s, rho_rad_h, rho_rad_s, drho_rad_h, drho_rad_s, vrho_rad_h, vrho_rad_s, rho_rad_h_d, rho_rad_s_d, ga_vlocal_gb_h, ga_vlocal_gb_s, int_scr_h, int_scr_s)
...
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...
routines that build the integrals of the Vxc potential calculated for the atomic density in the basis...
Definition qs_vxc_atom.F:12
subroutine, public calc_rho_angular(grid_atom, harmonics, nspins, grad_func, ir, r_h, r_s, rho_h, rho_s, dr_h, dr_s, r_h_d, r_s_d, drho_h, drho_s)
...
subroutine, public gavxcgb_nogc(vxc_h, vxc_s, int_hh, int_ss, grid_atom, basis_1c, harmonics, nspins)
...
All kind of helpful little routines.
Definition util.F:14
pure integer function, dimension(2), public get_limit(m, n, me)
divide m entries into n parts, return size of part me
Definition util.F:333
subroutine, public vxc_of_r_new(xc_fun_section, rho_set, deriv_set, deriv_order, needs, w, lsd, na, nr, exc, vxc, vxg, vtau, energy_only, adiabatic_rescale_factor)
...
Definition xc_atom.F:64
subroutine, public xc_rho_set_atom_update(rho_set, needs, nspins, bo)
...
Definition xc_atom.F:495
subroutine, public fill_rho_set(rho_set, lsd, nspins, needs, rho, drho, tau, na, ir)
...
Definition xc_atom.F:655
represent a group ofunctional derivatives
subroutine, public xc_dset_zero_all(deriv_set)
...
type(xc_derivative_type) function, pointer, public xc_dset_get_derivative(derivative_set, description, allocate_deriv)
returns the requested xc_derivative
subroutine, public xc_dset_release(derivative_set)
releases a derivative set
subroutine, public xc_dset_create(derivative_set, pw_pool, local_bounds)
creates a derivative set object
Provides types for the management of the xc-functionals and their derivatives.
subroutine, public xc_derivative_get(deriv, split_desc, order, deriv_data, accept_null_data)
returns various information on the given derivative
subroutine, public xc_functionals_eval(functionals, lsd, rho_set, deriv_set, deriv_order)
...
Calculate the saop potential.
Definition xc_pot_saop.F:11
subroutine, public add_saop_pot(ks_matrix, qs_env, oe_corr)
...
contains the structure
elemental subroutine, public xc_rho_cflags_setall(cflags, value)
sets all the flags to the given value
contains the structure
subroutine, public xc_rho_set_create(rho_set, local_bounds, rho_cutoff, drho_cutoff, tau_cutoff)
allocates and does (minimal) initialization of a rho_set
subroutine, public xc_rho_set_release(rho_set, pw_pool)
releases the given rho_set
subroutine, public xc_rho_set_update(rho_set, rho_r, rho_g, tau, needs, xc_deriv_method_id, xc_rho_smooth_id, pw_pool, spinflip)
updates the given rho set with the density given by rho_r (and rho_g). The rho set will contain the c...
calculates the Becke 88 exchange functional
Definition xc_xbecke88.F:14
subroutine, public xb88_lsd_info(reference, shortform, needs, max_deriv)
return various information on the functional
Definition xc_xbecke88.F:90
subroutine, public xb88_lda_info(reference, shortform, needs, max_deriv)
return various information on the functional
Definition xc_xbecke88.F:59
Exchange and Correlation functional calculations.
Definition xc.F:17
subroutine, public xc_vxc_pw_create(vxc_rho, vxc_tau, exc, rho_r, rho_g, tau, xc_section, weights, pw_pool, compute_virial, virial_xc, exc_r)
Exchange and Correlation functional calculations.
Definition xc.F:470
Provides all information about an atomic kind.
represent a pointer to a 1d array
just to build arrays of pointers to matrices
represent a full matrix
stores all the informations relevant to an mpi environment
Orbital angular momentum.
contained for different pw related things
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 derivative of a functional
contains a flag for each component of xc_rho_set, so that you can use it to tell which components you...
represent a density, with all the representation and data needed to perform a functional evaluation