(git:374b731)
Loading...
Searching...
No Matches
qs_rho0_methods.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
10
11 USE ao_util, ONLY: exp_radius,&
27 USE kinds, ONLY: default_string_length,&
28 dp
29 USE mathconstants, ONLY: fourpi
31 USE orbital_pointers, ONLY: indco,&
32 indso,&
33 nco,&
34 ncoset,&
35 nso,&
36 nsoset
43 USE qs_kind_types, ONLY: get_qs_kind,&
49 rhoz_type,&
52 USE qs_rho0_types, ONLY: &
59#include "./base/base_uses.f90"
60
61 IMPLICIT NONE
62
63 PRIVATE
64
65 ! Global parameters (only in this module)
66
67 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_rho0_methods'
68
69 ! Public subroutines
70
72
73CONTAINS
74
75! **************************************************************************************************
76!> \brief ...
77!> \param mp_gau ...
78!> \param basis_1c ...
79!> \param harmonics ...
80!> \param nchannels ...
81!> \param nsotot ...
82! **************************************************************************************************
83 SUBROUTINE calculate_mpole_gau(mp_gau, basis_1c, harmonics, nchannels, nsotot)
84
85 TYPE(mpole_gau_overlap) :: mp_gau
86 TYPE(gto_basis_set_type), POINTER :: basis_1c
87 TYPE(harmonics_atom_type), POINTER :: harmonics
88 INTEGER, INTENT(IN) :: nchannels, nsotot
89
90 CHARACTER(len=*), PARAMETER :: routineN = 'calculate_mpole_gau'
91
92 INTEGER :: handle, icg, ig1, ig2, ipgf1, ipgf2, iset1, iset2, iso, iso1, iso2, l, l1, l2, &
93 llmax, m1, m2, max_iso_not0_local, max_s_harm, maxl, maxso, n1, n2, nset
94 INTEGER, ALLOCATABLE, DIMENSION(:) :: cg_n_list
95 INTEGER, ALLOCATABLE, DIMENSION(:, :, :) :: cg_list
96 INTEGER, DIMENSION(:), POINTER :: lmax, lmin, npgf
97 REAL(KIND=dp) :: zet1, zet2
98 REAL(KIND=dp), DIMENSION(:, :), POINTER :: zet
99 REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: my_cg
100
101 CALL timeset(routinen, handle)
102
103 NULLIFY (lmax, lmin, npgf, my_cg, zet)
104
105 CALL reallocate(mp_gau%Qlm_gg, 1, nsotot, 1, nsotot, 1, nchannels)
106
107 CALL get_gto_basis_set(gto_basis_set=basis_1c, &
108 lmax=lmax, lmin=lmin, maxso=maxso, &
109 npgf=npgf, nset=nset, zet=zet, maxl=maxl)
110
111 max_s_harm = harmonics%max_s_harm
112 llmax = harmonics%llmax
113
114 ALLOCATE (cg_list(2, nsoset(maxl)**2, max_s_harm), cg_n_list(max_s_harm))
115
116 my_cg => harmonics%my_CG
117
118 m1 = 0
119 DO iset1 = 1, nset
120 m2 = 0
121 DO iset2 = 1, nset
122
123 CALL get_none0_cg_list(my_cg, lmin(iset1), lmax(iset1), lmin(iset2), lmax(iset2), &
124 max_s_harm, llmax, cg_list, cg_n_list, max_iso_not0_local)
125
126 n1 = nsoset(lmax(iset1))
127 DO ipgf1 = 1, npgf(iset1)
128 zet1 = zet(ipgf1, iset1)
129
130 n2 = nsoset(lmax(iset2))
131 DO ipgf2 = 1, npgf(iset2)
132 zet2 = zet(ipgf2, iset2)
133
134 DO iso = 1, min(nchannels, max_iso_not0_local)
135 l = indso(1, iso)
136 DO icg = 1, cg_n_list(iso)
137 iso1 = cg_list(1, icg, iso)
138 iso2 = cg_list(2, icg, iso)
139
140 l1 = indso(1, iso1)
141 l2 = indso(1, iso2)
142 ig1 = iso1 + n1*(ipgf1 - 1) + m1
143 ig2 = iso2 + n2*(ipgf2 - 1) + m2
144
145 mp_gau%Qlm_gg(ig1, ig2, iso) = fourpi/(2._dp*l + 1._dp)* &
146 my_cg(iso1, iso2, iso)*gaussint_sph(zet1 + zet2, l + l1 + l2)
147 END DO ! icg
148 END DO ! iso
149
150 END DO ! ipgf2
151 END DO ! ipgf1
152 m2 = m2 + maxso
153 END DO ! iset2
154 m1 = m1 + maxso
155 END DO ! iset1
156
157 DEALLOCATE (cg_list, cg_n_list)
158
159 CALL timestop(handle)
160 END SUBROUTINE calculate_mpole_gau
161
162! **************************************************************************************************
163!> \brief ...
164!> \param gapw_control ...
165!> \param rho_atom_set ...
166!> \param rho0_atom_set ...
167!> \param rho0_mp ...
168!> \param a_list ...
169!> \param natom ...
170!> \param ikind ...
171!> \param qs_kind ...
172!> \param rho0_h_tot ...
173! **************************************************************************************************
174 SUBROUTINE calculate_rho0_atom(gapw_control, rho_atom_set, rho0_atom_set, &
175 rho0_mp, a_list, natom, ikind, qs_kind, rho0_h_tot)
176
177 TYPE(gapw_control_type), POINTER :: gapw_control
178 TYPE(rho_atom_type), DIMENSION(:), POINTER :: rho_atom_set
179 TYPE(rho0_atom_type), DIMENSION(:), POINTER :: rho0_atom_set
180 TYPE(rho0_mpole_type), POINTER :: rho0_mp
181 INTEGER, DIMENSION(:), INTENT(IN) :: a_list
182 INTEGER, INTENT(IN) :: natom, ikind
183 TYPE(qs_kind_type), INTENT(IN) :: qs_kind
184 REAL(kind=dp), INTENT(INOUT) :: rho0_h_tot
185
186 CHARACTER(len=*), PARAMETER :: routinen = 'calculate_rho0_atom'
187
188 INTEGER :: handle, iat, iatom, ic, ico, ir, is, &
189 iso, ispin, l, lmax0, lshell, lx, ly, &
190 lz, nr, nsotot, nspins
191 LOGICAL :: paw_atom
192 REAL(kind=dp) :: sum1
193 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: cpc_ah, cpc_as
194 REAL(kind=dp), DIMENSION(:), POINTER :: norm_g0l_h
195 REAL(kind=dp), DIMENSION(:, :), POINTER :: g0_h, vg0_h
196 TYPE(grid_atom_type), POINTER :: g_atom
197 TYPE(harmonics_atom_type), POINTER :: harmonics
198 TYPE(mpole_gau_overlap), POINTER :: mpole_gau
199 TYPE(mpole_rho_atom), POINTER :: mpole_rho
200 TYPE(rho_atom_coeff), DIMENSION(:), POINTER :: cpc_h, cpc_s
201 TYPE(rho_atom_type), POINTER :: rho_atom
202
203 CALL timeset(routinen, handle)
204
205 NULLIFY (mpole_gau)
206 NULLIFY (mpole_rho)
207 NULLIFY (g0_h, vg0_h, g_atom)
208 NULLIFY (norm_g0l_h, harmonics)
209
210 CALL get_rho0_mpole(rho0_mpole=rho0_mp, ikind=ikind, &
211 l0_ikind=lmax0, mp_gau_ikind=mpole_gau, &
212 g0_h=g0_h, &
213 vg0_h=vg0_h, &
214 norm_g0l_h=norm_g0l_h)
215
216 CALL get_qs_kind(qs_kind, harmonics=harmonics, paw_atom=paw_atom, grid_atom=g_atom)
217
218 nr = g_atom%nr
219
220 ! Set density coefficient to zero before the calculation
221 DO iat = 1, natom
222 iatom = a_list(iat)
223 rho0_atom_set(iatom)%rho0_rad_h%r_coef = 0.0_dp
224 rho0_mp%mp_rho(iatom)%Qlm_tot = 0.0_dp
225 rho0_mp%mp_rho(iatom)%Qlm_tot(1) = rho0_mp%mp_rho(iatom)%Qlm_z
226 rho0_mp%mp_rho(iatom)%Q0 = 0.0_dp
227 rho0_mp%mp_rho(iatom)%Qlm_car = 0.0_dp
228 END DO
229
230 IF (.NOT. (.NOT. paw_atom .AND. gapw_control%nopaw_as_gpw)) THEN
231 DO iat = 1, natom
232 iatom = a_list(iat)
233 mpole_rho => rho0_mp%mp_rho(iatom)
234 rho_atom => rho_atom_set(iatom)
235
236 IF (paw_atom) THEN
237 NULLIFY (cpc_h, cpc_s)
238 CALL get_rho_atom(rho_atom=rho_atom, cpc_h=cpc_h, cpc_s=cpc_s)
239 nspins = SIZE(cpc_h)
240 nsotot = SIZE(mpole_gau%Qlm_gg, 1)
241 ALLOCATE (cpc_ah(nsotot, nsotot, nspins))
242 cpc_ah = 0._dp
243 ALLOCATE (cpc_as(nsotot, nsotot, nspins))
244 cpc_as = 0._dp
245 DO ispin = 1, nspins
246 CALL prj_scatter(cpc_h(ispin)%r_coef, cpc_ah(:, :, ispin), qs_kind)
247 CALL prj_scatter(cpc_s(ispin)%r_coef, cpc_as(:, :, ispin), qs_kind)
248 END DO
249 END IF
250
251 ! Total charge (hard-soft) at atom
252 IF (paw_atom) THEN
253 DO ispin = 1, nspins
254 mpole_rho%Q0(ispin) = (trace_r_axb(mpole_gau%Qlm_gg(:, :, 1), nsotot, &
255 cpc_ah(:, :, ispin), nsotot, nsotot, nsotot) &
256 - trace_r_axb(mpole_gau%Qlm_gg(:, :, 1), nsotot, &
257 cpc_as(:, :, ispin), nsotot, nsotot, nsotot))/sqrt(fourpi)
258 END DO
259 END IF
260 ! Multipoles of local charge distribution
261 DO iso = 1, nsoset(lmax0)
262 l = indso(1, iso)
263 IF (paw_atom) THEN
264 mpole_rho%Qlm_h(iso) = 0.0_dp
265 mpole_rho%Qlm_s(iso) = 0.0_dp
266
267 DO ispin = 1, nspins
268 mpole_rho%Qlm_h(iso) = mpole_rho%Qlm_h(iso) + &
269 trace_r_axb(mpole_gau%Qlm_gg(:, :, iso), nsotot, &
270 cpc_ah(:, :, ispin), nsotot, nsotot, nsotot)
271 mpole_rho%Qlm_s(iso) = mpole_rho%Qlm_s(iso) + &
272 trace_r_axb(mpole_gau%Qlm_gg(:, :, iso), nsotot, &
273 cpc_as(:, :, ispin), nsotot, nsotot, nsotot)
274 END DO ! ispin
275
276 mpole_rho%Qlm_tot(iso) = mpole_rho%Qlm_tot(iso) + &
277 mpole_rho%Qlm_h(iso) - mpole_rho%Qlm_s(iso)
278 END IF
279
280 rho0_atom_set(iatom)%rho0_rad_h%r_coef(1:nr, iso) = &
281 g0_h(1:nr, l)*mpole_rho%Qlm_tot(iso)
282 rho0_atom_set(iatom)%vrho0_rad_h%r_coef(1:nr, iso) = &
283 vg0_h(1:nr, l)*mpole_rho%Qlm_tot(iso)
284
285 sum1 = 0.0_dp
286 DO ir = 1, nr
287 sum1 = sum1 + g_atom%wr(ir)* &
288 rho0_atom_set(iatom)%rho0_rad_h%r_coef(ir, iso)
289 END DO
290 rho0_h_tot = rho0_h_tot + sum1*harmonics%slm_int(iso)
291 END DO ! iso
292 IF (paw_atom) THEN
293 DEALLOCATE (cpc_ah, cpc_as)
294 END IF
295 END DO ! iat
296 END IF
297
298 ! Transform the coefficinets from spherical to Cartesian
299 IF (.NOT. paw_atom .AND. gapw_control%nopaw_as_gpw) THEN
300 DO iat = 1, natom
301 iatom = a_list(iat)
302 mpole_rho => rho0_mp%mp_rho(iatom)
303
304 DO lshell = 0, lmax0
305 DO ic = 1, nco(lshell)
306 ico = ic + ncoset(lshell - 1)
307 mpole_rho%Qlm_car(ico) = 0.0_dp
308 END DO
309 END DO
310 END DO
311 ELSE
312 DO iat = 1, natom
313 iatom = a_list(iat)
314 mpole_rho => rho0_mp%mp_rho(iatom)
315 DO lshell = 0, lmax0
316 DO ic = 1, nco(lshell)
317 ico = ic + ncoset(lshell - 1)
318 mpole_rho%Qlm_car(ico) = 0.0_dp
319 lx = indco(1, ico)
320 ly = indco(2, ico)
321 lz = indco(3, ico)
322 DO is = 1, nso(lshell)
323 iso = is + nsoset(lshell - 1)
324 mpole_rho%Qlm_car(ico) = mpole_rho%Qlm_car(ico) + &
325 norm_g0l_h(lshell)* &
326 orbtramat(lshell)%slm(is, ic)* &
327 mpole_rho%Qlm_tot(iso)
328
329 END DO
330 END DO
331 END DO ! lshell
332 END DO ! iat
333 END IF
334 !MI Get rid of full gapw
335
336 CALL timestop(handle)
337
338 END SUBROUTINE calculate_rho0_atom
339
340! **************************************************************************************************
341!> \brief ...
342!> \param local_rho_set ...
343!> \param qs_env ...
344!> \param gapw_control ...
345!> \param zcore ...
346! **************************************************************************************************
347 SUBROUTINE init_rho0(local_rho_set, qs_env, gapw_control, zcore)
348
349 TYPE(local_rho_type), POINTER :: local_rho_set
350 TYPE(qs_environment_type), POINTER :: qs_env
351 TYPE(gapw_control_type), POINTER :: gapw_control
352 REAL(kind=dp), INTENT(IN), OPTIONAL :: zcore
353
354 CHARACTER(len=*), PARAMETER :: routinen = 'init_rho0'
355
356 CHARACTER(LEN=default_string_length) :: unit_str
357 INTEGER :: handle, iat, iatom, ikind, l, l_rho1_max, laddg, lmaxg, maxl, maxnset, maxso, &
358 nat, natom, nchan_c, nchan_s, nkind, nr, nset, nsotot, output_unit
359 INTEGER, DIMENSION(:), POINTER :: atom_list
360 LOGICAL :: paw_atom
361 REAL(kind=dp) :: alpha_core, eps_vrho0, max_rpgf0_s, &
362 radius, rc_min, rc_orb, &
363 total_rho_core_rspace, zeff
364 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
365 TYPE(cp_logger_type), POINTER :: logger
366 TYPE(grid_atom_type), POINTER :: grid_atom
367 TYPE(gto_basis_set_type), POINTER :: basis_1c
368 TYPE(harmonics_atom_type), POINTER :: harmonics
369 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
370 TYPE(rho0_atom_type), DIMENSION(:), POINTER :: rho0_atom_set
371 TYPE(rho0_mpole_type), POINTER :: rho0_mpole
372 TYPE(rhoz_type), DIMENSION(:), POINTER :: rhoz_set
373 TYPE(section_vals_type), POINTER :: dft_section
374
375 CALL timeset(routinen, handle)
376
377 NULLIFY (logger)
378 logger => cp_get_default_logger()
379
380 NULLIFY (qs_kind_set)
381 NULLIFY (atomic_kind_set)
382 NULLIFY (harmonics)
383 NULLIFY (basis_1c)
384 NULLIFY (rho0_mpole)
385 NULLIFY (rho0_atom_set)
386 NULLIFY (rhoz_set)
387
388 CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set, &
389 atomic_kind_set=atomic_kind_set)
390
391 nkind = SIZE(atomic_kind_set)
392 eps_vrho0 = gapw_control%eps_Vrho0
393
394 ! Initialize rhoz total to zero
395 ! in gapw rhoz is calculated on local the lebedev grids
396 total_rho_core_rspace = 0.0_dp
397
398 CALL get_atomic_kind_set(atomic_kind_set, natom=natom)
399
400 ! Initialize the multipole and the compensation charge type
401 CALL allocate_rho0_mpole(rho0_mpole)
402 CALL allocate_rho0_atom(rho0_atom_set, natom)
403
404 ! Allocate the multipole set
405 CALL allocate_multipoles(rho0_mpole%mp_rho, natom, rho0_mpole%mp_gau, nkind)
406
407 ! Allocate the core density on the radial grid for each kind: rhoz_set
408 CALL allocate_rhoz(rhoz_set, nkind)
409
410 ! For each kind, determine the max l for the compensation charge density
411 lmaxg = gapw_control%lmax_rho0
412 laddg = gapw_control%ladd_rho0
413
414 CALL reallocate(rho0_mpole%lmax0_kind, 1, nkind)
415
416 rho0_mpole%lmax_0 = 0
417 rc_min = 100.0_dp
418 maxnset = 0
419 DO ikind = 1, nkind
420 CALL get_atomic_kind(atomic_kind_set(ikind), atom_list=atom_list, natom=nat)
421 CALL get_qs_kind(qs_kind_set(ikind), &
422 ngrid_rad=nr, &
423 grid_atom=grid_atom, &
424 harmonics=harmonics, &
425 paw_atom=paw_atom, &
426 hard0_radius=rc_orb, &
427 zeff=zeff, &
428 alpha_core_charge=alpha_core)
429 CALL get_qs_kind(qs_kind_set(ikind), &
430 basis_set=basis_1c, basis_type="GAPW_1C")
431
432 ! Set charge distribution of ionic cores to zero when computing the response-density
433 IF (PRESENT(zcore)) zeff = zcore
434
435 CALL get_gto_basis_set(gto_basis_set=basis_1c, &
436 maxl=maxl, &
437 maxso=maxso, nset=nset)
438
439 maxnset = max(maxnset, nset)
440
441 l_rho1_max = indso(1, harmonics%max_iso_not0)
442 IF (paw_atom) THEN
443 rho0_mpole%lmax0_kind(ikind) = min(2*maxl, l_rho1_max, maxl + laddg, lmaxg)
444 ELSE
445 rho0_mpole%lmax0_kind(ikind) = 0
446 END IF
447
448 CALL set_qs_kind(qs_kind_set(ikind), lmax_rho0=rho0_mpole%lmax0_kind(ikind))
449
450 IF (gapw_control%lrho1_eq_lrho0) harmonics%max_iso_not0 = &
451 nsoset(rho0_mpole%lmax0_kind(ikind))
452
453 rho0_mpole%lmax_0 = max(rho0_mpole%lmax_0, rho0_mpole%lmax0_kind(ikind))
454 rc_min = min(rc_min, rc_orb)
455
456 nchan_s = nsoset(rho0_mpole%lmax0_kind(ikind))
457 nchan_c = ncoset(rho0_mpole%lmax0_kind(ikind))
458 nsotot = maxso*nset
459
460 DO iat = 1, nat
461 iatom = atom_list(iat)
462 ! Allocate the multipole for rho1_h rho1_s and rho_z
463 CALL initialize_mpole_rho(rho0_mpole%mp_rho(iatom), nchan_s, nchan_c, zeff)
464 ! Allocate the radial part of rho0_h and rho0_s
465 ! This is calculated on the radial grid centered at the atomic position
466 CALL allocate_rho0_atom_rad(rho0_atom_set(iatom), nr, nchan_s)
467 END DO
468
469 IF (paw_atom) THEN
470 ! Calculate multipoles given by the product of 2 primitives Qlm_gg
471 CALL calculate_mpole_gau(rho0_mpole%mp_gau(ikind), &
472 basis_1c, harmonics, nchan_s, nsotot)
473 END IF
474
475 ! Calculate the core density rhoz
476 ! exp(-alpha_c**2 r**2)Z(alpha_c**2/pi)**(3/2)
477 ! on the logarithmic radial grid
478 ! WARNING: alpha_core_charge = alpha_c**2
479 CALL calculate_rhoz(rhoz_set(ikind), grid_atom, alpha_core, zeff, &
480 nat, total_rho_core_rspace, harmonics)
481 END DO ! ikind
482 total_rho_core_rspace = -total_rho_core_rspace
483
484 IF (gapw_control%alpha0_hard_from_input) THEN
485 ! The exponent for the compensation charge rho0_hard is read from input
486 rho0_mpole%zet0_h = gapw_control%alpha0_hard
487 ELSE
488 ! Calculate the exponent for the compensation charge rho0_hard
489 rho0_mpole%zet0_h = 0.1_dp
490 DO
491 radius = exp_radius(rho0_mpole%lmax_0, rho0_mpole%zet0_h, eps_vrho0, 1.0_dp)
492 IF (radius <= rc_min) EXIT
493 rho0_mpole%zet0_h = rho0_mpole%zet0_h + 0.1_dp
494 END DO
495
496 END IF
497
498 ! Allocate and calculate the normalization factors for g0_lm_h and g0_lm_s
499 CALL reallocate(rho0_mpole%norm_g0l_h, 0, rho0_mpole%lmax_0)
500 DO l = 0, rho0_mpole%lmax_0
501 rho0_mpole%norm_g0l_h(l) = (2._dp*l + 1._dp)/ &
502 (fourpi*gaussint_sph(rho0_mpole%zet0_h, 2*l))
503 END DO
504
505 ! Allocate and Initialize the g0 gaussians used to build the compensation density
506 ! and calculate the interaction radii
507 max_rpgf0_s = 0.0_dp
508 DO ikind = 1, nkind
509 CALL get_qs_kind(qs_kind_set(ikind), grid_atom=grid_atom)
510 CALL calculate_g0(rho0_mpole, grid_atom, ikind)
511 CALL interaction_radii_g0(rho0_mpole, ikind, eps_vrho0, max_rpgf0_s)
512 END DO
513 rho0_mpole%max_rpgf0_s = max_rpgf0_s
514
515 CALL set_local_rho(local_rho_set, rho0_atom_set=rho0_atom_set, rho0_mpole=rho0_mpole, rhoz_set=rhoz_set)
516 local_rho_set%rhoz_tot = total_rho_core_rspace
517
518 dft_section => section_vals_get_subs_vals(qs_env%input, "DFT")
519 output_unit = cp_print_key_unit_nr(logger, dft_section, "PRINT%GAPW%RHO0_INFORMATION", &
520 extension=".Log")
521 CALL section_vals_val_get(dft_section, "PRINT%GAPW%RHO0_INFORMATION%UNIT", c_val=unit_str)
522 IF (output_unit > 0) THEN
523 CALL write_rho0_info(rho0_mpole, unit_str, output_unit)
524 END IF
525 CALL cp_print_key_finished_output(output_unit, logger, dft_section, &
526 "PRINT%GAPW%RHO0_INFORMATION")
527
528 CALL timestop(handle)
529
530 END SUBROUTINE init_rho0
531
532! **************************************************************************************************
533!> \brief ...
534!> \param rho0_mpole ...
535!> \param ik ...
536!> \param eps_Vrho0 ...
537!> \param max_rpgf0_s ...
538! **************************************************************************************************
539 SUBROUTINE interaction_radii_g0(rho0_mpole, ik, eps_Vrho0, max_rpgf0_s)
540
541 TYPE(rho0_mpole_type), POINTER :: rho0_mpole
542 INTEGER, INTENT(IN) :: ik
543 REAL(kind=dp), INTENT(IN) :: eps_vrho0
544 REAL(kind=dp), INTENT(INOUT) :: max_rpgf0_s
545
546 INTEGER :: l, lmax
547 REAL(kind=dp) :: r_h, z0_h
548 REAL(kind=dp), DIMENSION(:), POINTER :: ng0_h
549
550 CALL get_rho0_mpole(rho0_mpole, ikind=ik, l0_ikind=lmax, &
551 zet0_h=z0_h, norm_g0l_h=ng0_h)
552 r_h = 0.0_dp
553 DO l = 0, lmax
554 r_h = max(r_h, exp_radius(l, z0_h, eps_vrho0, ng0_h(l), rlow=r_h))
555 END DO
556
557 rho0_mpole%mp_gau(ik)%rpgf0_h = r_h
558 rho0_mpole%mp_gau(ik)%rpgf0_s = r_h
559 max_rpgf0_s = max(max_rpgf0_s, r_h)
560
561 END SUBROUTINE interaction_radii_g0
562
563END MODULE qs_rho0_methods
All kind of helpful little routines.
Definition ao_util.F:14
real(dp) function, public gaussint_sph(alpha, l)
...
Definition ao_util.F:299
real(kind=dp) function, public exp_radius(l, alpha, threshold, prefactor, epsabs, epsrel, rlow)
The radius of a primitive Gaussian function for a given threshold is calculated. g(r) = prefactor*r**...
Definition ao_util.F:96
pure real(dp) function, public trace_r_axb(a, lda, b, ldb, m, n)
...
Definition ao_util.F:331
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind_set(atomic_kind_set, atom_of_kind, kind_of, natom_of_kind, maxatom, natom, nshell, fist_potential_present, shell_present, shell_adiabatic, shell_check_distance, damping_present)
Get attributes of an atomic kind set.
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.
subroutine, public get_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, lmin, lx, ly, lz, m, ncgf_set, npgf, nsgf_set, nshell, cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, last_cgf, last_sgf, n, gcc, maxco, maxl, maxpgf, maxsgf_set, maxshell, maxso, nco_sum, npgf_sum, nshell_sum, maxder, short_kind_radius)
...
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
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
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
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
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
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
Definition of mathematical constants and functions.
real(kind=dp), parameter, public fourpi
Utility routines for the memory handling.
Provides Cartesian and spherical orbital pointers and indices.
integer, dimension(:), allocatable, public nco
integer, dimension(:), allocatable, public nsoset
integer, dimension(:, :), allocatable, public indso
integer, dimension(:), allocatable, public ncoset
integer, dimension(:, :), allocatable, public indco
integer, dimension(:), allocatable, public nso
Calculation of the spherical harmonics and the corresponding orbital transformation matrices.
type(orbtramat_type), dimension(:), pointer, public orbtramat
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, rhs)
Get the QUICKSTEP environment.
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, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, 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_r3d_rs_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, floating, name, element_symbol, pao_basis_size, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
subroutine, public set_qs_kind(qs_kind, paw_atom, ghost, floating, hard_radius, hard0_radius, covalent_radius, vdw_radius, lmax_rho0, zeff, no_optimize, dispersion, u_minus_j, reltmat, dftb_parameter, xtb_parameter, elec_conf, pao_basis_size)
Set the components of an atomic kind data set.
subroutine, public allocate_rhoz(rhoz_set, nkind)
...
subroutine, public set_local_rho(local_rho_set, rho_atom_set, rho0_atom_set, rho0_mpole, rhoz_set)
...
subroutine, public calculate_rhoz(rhoz, grid_atom, alpha, zeff, natom, rhoz_tot, harmonics)
...
Routines for the construction of the coefficients for the expansion of the atomic densities rho1_hard...
subroutine, public prj_scatter(ain, aout, atom)
...
subroutine, public init_rho0(local_rho_set, qs_env, gapw_control, zcore)
...
subroutine, public calculate_rho0_atom(gapw_control, rho_atom_set, rho0_atom_set, rho0_mp, a_list, natom, ikind, qs_kind, rho0_h_tot)
...
subroutine, public allocate_multipoles(mp_rho, natom, mp_gau, nkind)
...
subroutine, public initialize_mpole_rho(mp_rho, nchan_s, nchan_c, zeff)
...
subroutine, public allocate_rho0_atom_rad(rho0_atom, nr, nchannels)
...
subroutine, public calculate_g0(rho0_mpole, grid_atom, ik)
...
subroutine, public write_rho0_info(rho0_mpole, unit_str, output_unit)
...
subroutine, public allocate_rho0_atom(rho0_set, natom)
...
subroutine, public allocate_rho0_mpole(rho0)
...
subroutine, public get_rho0_mpole(rho0_mpole, g0_h, vg0_h, iat, ikind, lmax_0, l0_ikind, mp_gau_ikind, mp_rho, norm_g0l_h, qlm_gg, qlm_car, qlm_tot, zet0_h, igrid_zet0_s, rpgf0_h, rpgf0_s, max_rpgf0_s, rho0_s_rs, rho0_s_gs)
...
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)
...
Provides all information about an atomic kind.
type of a logger, at the moment it contains just a print level starting at which level it should be l...
Provides all information about a quickstep kind.