(git:495eafe)
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-2026 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
48 USE qs_kind_types, ONLY: get_qs_kind,&
55 rhoz_type,&
58 USE qs_rho0_types, ONLY: &
65#include "./base/base_uses.f90"
66
67 IMPLICIT NONE
68
69 PRIVATE
70
71 ! Global parameters (only in this module)
72
73 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_rho0_methods'
74
75 ! Public subroutines
76
78
79CONTAINS
80
81! **************************************************************************************************
82!> \brief ...
83!> \param Qlm_gg ...
84!> \param basis_1c ...
85!> \param harmonics ...
86!> \param nchannels ...
87!> \param nsotot ...
88! **************************************************************************************************
89 SUBROUTINE calculate_mpole_gau(Qlm_gg, basis_1c, harmonics, nchannels, nsotot)
90
91 REAL(dp), DIMENSION(:, :, :), POINTER :: Qlm_gg
92 TYPE(gto_basis_set_type), POINTER :: basis_1c
93 TYPE(harmonics_atom_type), POINTER :: harmonics
94 INTEGER, INTENT(IN) :: nchannels, nsotot
95
96 CHARACTER(len=*), PARAMETER :: routineN = 'calculate_mpole_gau'
97
98 INTEGER :: handle, icg, ig1, ig2, ipgf1, ipgf2, iset1, iset2, iso, iso1, iso2, l, l1, l2, &
99 llmax, m1, m2, max_iso_not0_local, max_s_harm, maxl, maxso, n1, n2, nset
100 INTEGER, ALLOCATABLE, DIMENSION(:) :: cg_n_list
101 INTEGER, ALLOCATABLE, DIMENSION(:, :, :) :: cg_list
102 INTEGER, DIMENSION(:), POINTER :: lmax, lmin, npgf
103 REAL(KIND=dp) :: zet1, zet2
104 REAL(KIND=dp), DIMENSION(:, :), POINTER :: zet
105 REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: my_cg
106
107 CALL timeset(routinen, handle)
108
109 NULLIFY (lmax, lmin, npgf, my_cg, zet)
110
111 CALL reallocate(qlm_gg, 1, nsotot, 1, nsotot, 1, &
112 min(nchannels, harmonics%max_iso_not0))
113
114 CALL get_gto_basis_set(gto_basis_set=basis_1c, &
115 lmax=lmax, lmin=lmin, maxso=maxso, &
116 npgf=npgf, nset=nset, zet=zet, maxl=maxl)
117
118 max_s_harm = harmonics%max_s_harm
119 llmax = harmonics%llmax
120
121 ALLOCATE (cg_list(2, nsoset(maxl)**2, max_s_harm), cg_n_list(max_s_harm))
122
123 my_cg => harmonics%my_CG
124
125 m1 = 0
126 DO iset1 = 1, nset
127 m2 = 0
128 DO iset2 = 1, nset
129
130 CALL get_none0_cg_list(my_cg, lmin(iset1), lmax(iset1), lmin(iset2), lmax(iset2), &
131 max_s_harm, llmax, cg_list, cg_n_list, max_iso_not0_local)
132
133 n1 = nsoset(lmax(iset1))
134 DO ipgf1 = 1, npgf(iset1)
135 zet1 = zet(ipgf1, iset1)
136
137 n2 = nsoset(lmax(iset2))
138 DO ipgf2 = 1, npgf(iset2)
139 zet2 = zet(ipgf2, iset2)
140
141 DO iso = 1, min(nchannels, max_iso_not0_local)
142 l = indso(1, iso)
143 DO icg = 1, cg_n_list(iso)
144 iso1 = cg_list(1, icg, iso)
145 iso2 = cg_list(2, icg, iso)
146
147 l1 = indso(1, iso1)
148 l2 = indso(1, iso2)
149 ig1 = iso1 + n1*(ipgf1 - 1) + m1
150 ig2 = iso2 + n2*(ipgf2 - 1) + m2
151
152 qlm_gg(ig1, ig2, iso) = fourpi/(2._dp*l + 1._dp)* &
153 my_cg(iso1, iso2, iso)*gaussint_sph(zet1 + zet2, l + l1 + l2)
154 END DO ! icg
155 END DO ! iso
156
157 END DO ! ipgf2
158 END DO ! ipgf1
159 m2 = m2 + maxso
160 END DO ! iset2
161 m1 = m1 + maxso
162 END DO ! iset1
163
164 DEALLOCATE (cg_list, cg_n_list)
165
166 CALL timestop(handle)
167 END SUBROUTINE calculate_mpole_gau
168
169! **************************************************************************************************
170!> \brief ...
171!> \param gapw_control ...
172!> \param rho_atom_set ...
173!> \param rhoz_cneo_set ...
174!> \param rho0_atom_set ...
175!> \param rho0_mp ...
176!> \param a_list ...
177!> \param natom ...
178!> \param ikind ...
179!> \param qs_kind ...
180!> \param rho0_h_tot ...
181! **************************************************************************************************
182 SUBROUTINE calculate_rho0_atom(gapw_control, rho_atom_set, rhoz_cneo_set, rho0_atom_set, &
183 rho0_mp, a_list, natom, ikind, qs_kind, rho0_h_tot)
184
185 TYPE(gapw_control_type), POINTER :: gapw_control
186 TYPE(rho_atom_type), DIMENSION(:), POINTER :: rho_atom_set
187 TYPE(rhoz_cneo_type), DIMENSION(:), POINTER :: rhoz_cneo_set
188 TYPE(rho0_atom_type), DIMENSION(:), POINTER :: rho0_atom_set
189 TYPE(rho0_mpole_type), POINTER :: rho0_mp
190 INTEGER, DIMENSION(:), INTENT(IN) :: a_list
191 INTEGER, INTENT(IN) :: natom, ikind
192 TYPE(qs_kind_type), INTENT(IN) :: qs_kind
193 REAL(kind=dp), INTENT(INOUT) :: rho0_h_tot
194
195 CHARACTER(len=*), PARAMETER :: routinen = 'calculate_rho0_atom'
196
197 INTEGER :: handle, iat, iatom, ic, ico, ir, is, &
198 iso, ispin, l, lmax0, lshell, lx, ly, &
199 lz, nr, nsotot, nsotot_nuc, nspins
200 LOGICAL :: paw_atom
201 REAL(kind=dp) :: sum1
202 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: cpc_h_nuc, cpc_s_nuc
203 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: cpc_ah, cpc_as
204 REAL(kind=dp), DIMENSION(:), POINTER :: norm_g0l_h
205 REAL(kind=dp), DIMENSION(:, :), POINTER :: g0_h, vg0_h
206 TYPE(cneo_potential_type), POINTER :: cneo_potential
207 TYPE(grid_atom_type), POINTER :: g_atom
208 TYPE(harmonics_atom_type), POINTER :: harmonics
209 TYPE(mpole_gau_overlap), POINTER :: mpole_gau
210 TYPE(mpole_rho_atom), POINTER :: mpole_rho
211 TYPE(rho_atom_coeff), DIMENSION(:), POINTER :: cpc_h, cpc_s
212 TYPE(rho_atom_type), POINTER :: rho_atom
213
214 CALL timeset(routinen, handle)
215
216 NULLIFY (mpole_gau)
217 NULLIFY (mpole_rho)
218 NULLIFY (g0_h, vg0_h, g_atom)
219 NULLIFY (norm_g0l_h, harmonics)
220 NULLIFY (cneo_potential)
221
222 CALL get_rho0_mpole(rho0_mpole=rho0_mp, ikind=ikind, &
223 l0_ikind=lmax0, mp_gau_ikind=mpole_gau, &
224 g0_h=g0_h, &
225 vg0_h=vg0_h, &
226 norm_g0l_h=norm_g0l_h)
227
228 CALL get_qs_kind(qs_kind, harmonics=harmonics, paw_atom=paw_atom, grid_atom=g_atom, &
229 cneo_potential=cneo_potential)
230
231 nr = g_atom%nr
232
233 ! Set density coefficient to zero before the calculation
234 DO iat = 1, natom
235 iatom = a_list(iat)
236 rho0_atom_set(iatom)%rho0_rad_h%r_coef = 0.0_dp
237 rho0_mp%mp_rho(iatom)%Qlm_tot = 0.0_dp
238 ! When no nuclear density matrix is available, use spherical zeff
239 IF (.NOT. ASSOCIATED(cneo_potential)) THEN
240 rho0_mp%mp_rho(iatom)%Qlm_tot(1) = rho0_mp%mp_rho(iatom)%Qlm_z
241 ELSE
242 IF (.NOT. rhoz_cneo_set(iatom)%ready) THEN
243 rho0_mp%mp_rho(iatom)%Qlm_tot(1) = rho0_mp%mp_rho(iatom)%Qlm_z
244 END IF
245 END IF
246 rho0_mp%mp_rho(iatom)%Q0 = 0.0_dp
247 rho0_mp%mp_rho(iatom)%Qlm_car = 0.0_dp
248 END DO
249
250 IF (.NOT. (.NOT. paw_atom .AND. gapw_control%nopaw_as_gpw)) THEN
251 DO iat = 1, natom
252 iatom = a_list(iat)
253 mpole_rho => rho0_mp%mp_rho(iatom)
254 rho_atom => rho_atom_set(iatom)
255
256 IF (paw_atom) THEN
257 NULLIFY (cpc_h, cpc_s)
258 CALL get_rho_atom(rho_atom=rho_atom, cpc_h=cpc_h, cpc_s=cpc_s)
259 nspins = SIZE(cpc_h)
260 nsotot = SIZE(mpole_gau%Qlm_gg, 1)
261 IF (ALLOCATED(cpc_ah)) THEN
262 IF (SIZE(cpc_ah, 1) /= nsotot .OR. SIZE(cpc_ah, 2) /= nsotot .OR. &
263 SIZE(cpc_ah, 3) /= nspins) DEALLOCATE (cpc_ah)
264 END IF
265 IF (.NOT. ALLOCATED(cpc_ah)) ALLOCATE (cpc_ah(nsotot, nsotot, nspins))
266 cpc_ah = 0._dp
267 IF (ALLOCATED(cpc_as)) THEN
268 IF (SIZE(cpc_as, 1) /= nsotot .OR. SIZE(cpc_as, 2) /= nsotot .OR. &
269 SIZE(cpc_as, 3) /= nspins) DEALLOCATE (cpc_as)
270 END IF
271 IF (.NOT. ALLOCATED(cpc_as)) ALLOCATE (cpc_as(nsotot, nsotot, nspins))
272 cpc_as = 0._dp
273 DO ispin = 1, nspins
274 CALL prj_scatter(cpc_h(ispin)%r_coef, cpc_ah(:, :, ispin), qs_kind)
275 CALL prj_scatter(cpc_s(ispin)%r_coef, cpc_as(:, :, ispin), qs_kind)
276 END DO
277 nsotot_nuc = 0
278 IF (ASSOCIATED(cneo_potential)) THEN
279 IF (rhoz_cneo_set(iatom)%ready) THEN
280 nsotot_nuc = SIZE(cneo_potential%Qlm_gg, 1)
281 IF (ALLOCATED(cpc_h_nuc)) THEN
282 IF (SIZE(cpc_h_nuc, 1) /= nsotot_nuc .OR. &
283 SIZE(cpc_h_nuc, 2) /= nsotot_nuc) DEALLOCATE (cpc_h_nuc)
284 END IF
285 IF (.NOT. ALLOCATED(cpc_h_nuc)) ALLOCATE (cpc_h_nuc(nsotot_nuc, nsotot_nuc))
286 IF (ALLOCATED(cpc_s_nuc)) THEN
287 IF (SIZE(cpc_s_nuc, 1) /= nsotot_nuc .OR. &
288 SIZE(cpc_s_nuc, 2) /= nsotot_nuc) DEALLOCATE (cpc_s_nuc)
289 END IF
290 IF (.NOT. ALLOCATED(cpc_s_nuc)) ALLOCATE (cpc_s_nuc(nsotot_nuc, nsotot_nuc))
291 cpc_h_nuc = 0._dp
292 cpc_s_nuc = 0._dp
293 CALL cneo_scatter(rhoz_cneo_set(iatom)%cpc_h, cpc_h_nuc, cneo_potential%npsgf, &
294 cneo_potential%n2oindex)
295 CALL cneo_scatter(rhoz_cneo_set(iatom)%cpc_s, cpc_s_nuc, cneo_potential%npsgf, &
296 cneo_potential%n2oindex)
297 END IF
298 END IF
299
300 ! Total charge (hard-soft) at atom
301 IF (paw_atom) THEN
302 DO ispin = 1, nspins
303 mpole_rho%Q0(ispin) = (trace_r_axb(mpole_gau%Qlm_gg(:, :, 1), nsotot, &
304 cpc_ah(:, :, ispin), nsotot, nsotot, nsotot) &
305 - trace_r_axb(mpole_gau%Qlm_gg(:, :, 1), nsotot, &
306 cpc_as(:, :, ispin), nsotot, nsotot, nsotot))/sqrt(fourpi)
307 END DO
308 END IF
309 ! Multipoles of local charge distribution
310 DO iso = 1, min(nsoset(lmax0), harmonics%max_iso_not0)
311 IF (paw_atom) THEN
312 mpole_rho%Qlm_h(iso) = 0.0_dp
313 mpole_rho%Qlm_s(iso) = 0.0_dp
314
315 DO ispin = 1, nspins
316 mpole_rho%Qlm_h(iso) = mpole_rho%Qlm_h(iso) + &
317 trace_r_axb(mpole_gau%Qlm_gg(:, :, iso), nsotot, &
318 cpc_ah(:, :, ispin), nsotot, nsotot, nsotot)
319 mpole_rho%Qlm_s(iso) = mpole_rho%Qlm_s(iso) + &
320 trace_r_axb(mpole_gau%Qlm_gg(:, :, iso), nsotot, &
321 cpc_as(:, :, ispin), nsotot, nsotot, nsotot)
322 END DO ! ispin
323
324 mpole_rho%Qlm_tot(iso) = mpole_rho%Qlm_tot(iso) + &
325 mpole_rho%Qlm_h(iso) - mpole_rho%Qlm_s(iso)
326 END IF
327 END DO ! iso
328
329 ! Multipoles of CNEO quantum nuclear charge distribuition
330 IF (ASSOCIATED(cneo_potential)) THEN
331 IF (rhoz_cneo_set(iatom)%ready) THEN
332 DO iso = 1, min(nsoset(lmax0), cneo_potential%harmonics%max_iso_not0)
333 mpole_rho%Qlm_tot(iso) = mpole_rho%Qlm_tot(iso) - cneo_potential%zeff* &
334 trace_r_axb(cneo_potential%Qlm_gg(:, :, iso), nsotot_nuc, &
335 cpc_h_nuc - cpc_s_nuc, nsotot_nuc, &
336 nsotot_nuc, nsotot_nuc)
337 END DO ! iso
338 END IF
339 END IF
340 END IF
341
342 DO iso = 1, nsoset(lmax0)
343 l = indso(1, iso)
344 rho0_atom_set(iatom)%rho0_rad_h%r_coef(1:nr, iso) = &
345 g0_h(1:nr, l)*mpole_rho%Qlm_tot(iso)
346 rho0_atom_set(iatom)%vrho0_rad_h%r_coef(1:nr, iso) = &
347 vg0_h(1:nr, l)*mpole_rho%Qlm_tot(iso)
348
349 ! When CNEO is enabled, it is possible for rho0 to have a higher angualr momentum
350 ! than that of the electronic density. In that case, the nuclear density must have
351 ! a higher angular momentum, but cneo_potential%harmonics%slm_int is not initialized.
352 ! For higher angular momenta, simply use the fact that slm_int(iso>1)=0
353 IF (iso <= harmonics%max_iso_not0) THEN
354 sum1 = 0.0_dp
355 DO ir = 1, nr
356 sum1 = sum1 + g_atom%wr(ir)* &
357 rho0_atom_set(iatom)%rho0_rad_h%r_coef(ir, iso)
358 END DO
359 rho0_h_tot = rho0_h_tot + sum1*harmonics%slm_int(iso)
360 END IF
361 END DO ! iso
362 END DO ! iat
363 END IF
364
365 IF (ALLOCATED(cpc_ah)) DEALLOCATE (cpc_ah)
366 IF (ALLOCATED(cpc_as)) DEALLOCATE (cpc_as)
367 IF (ALLOCATED(cpc_h_nuc)) DEALLOCATE (cpc_h_nuc)
368 IF (ALLOCATED(cpc_s_nuc)) DEALLOCATE (cpc_s_nuc)
369
370 ! Transform the coefficinets from spherical to Cartesian
371 IF (.NOT. paw_atom .AND. gapw_control%nopaw_as_gpw) THEN
372 DO iat = 1, natom
373 iatom = a_list(iat)
374 mpole_rho => rho0_mp%mp_rho(iatom)
375
376 DO lshell = 0, lmax0
377 DO ic = 1, nco(lshell)
378 ico = ic + ncoset(lshell - 1)
379 mpole_rho%Qlm_car(ico) = 0.0_dp
380 END DO
381 END DO
382 END DO
383 ELSE
384 DO iat = 1, natom
385 iatom = a_list(iat)
386 mpole_rho => rho0_mp%mp_rho(iatom)
387 DO lshell = 0, lmax0
388 DO ic = 1, nco(lshell)
389 ico = ic + ncoset(lshell - 1)
390 mpole_rho%Qlm_car(ico) = 0.0_dp
391 lx = indco(1, ico)
392 ly = indco(2, ico)
393 lz = indco(3, ico)
394 DO is = 1, nso(lshell)
395 iso = is + nsoset(lshell - 1)
396 mpole_rho%Qlm_car(ico) = mpole_rho%Qlm_car(ico) + &
397 norm_g0l_h(lshell)* &
398 orbtramat(lshell)%slm(is, ic)* &
399 mpole_rho%Qlm_tot(iso)
400
401 END DO
402 END DO
403 END DO ! lshell
404 END DO ! iat
405 END IF
406 !MI Get rid of full gapw
407
408 CALL timestop(handle)
409
410 END SUBROUTINE calculate_rho0_atom
411
412! **************************************************************************************************
413!> \brief ...
414!> \param local_rho_set ...
415!> \param qs_env ...
416!> \param gapw_control ...
417!> \param zcore ...
418! **************************************************************************************************
419 SUBROUTINE init_rho0(local_rho_set, qs_env, gapw_control, zcore)
420
421 TYPE(local_rho_type), POINTER :: local_rho_set
422 TYPE(qs_environment_type), POINTER :: qs_env
423 TYPE(gapw_control_type), POINTER :: gapw_control
424 REAL(kind=dp), INTENT(IN), OPTIONAL :: zcore
425
426 CHARACTER(len=*), PARAMETER :: routinen = 'init_rho0'
427
428 CHARACTER(LEN=default_string_length) :: unit_str
429 INTEGER :: handle, iat, iatom, ikind, l, l_rho1_max, l_rho1_max_nuc, laddg, lmaxg, maxl, &
430 maxl_nuc, maxnset, maxso, maxso_nuc, nat, natom, nchan_c, nchan_s, nkind, nr, nset, &
431 nset_nuc, nsotot, nsotot_nuc, output_unit
432 INTEGER, DIMENSION(:), POINTER :: atom_list
433 LOGICAL :: cneo_potential_present, paw_atom
434 REAL(kind=dp) :: alpha_core, eps_vrho0, max_rpgf0_s, radius, rc_min, rc_orb, &
435 total_rho_core_rspace, total_rho_nuc_cneo_rspace, zeff
436 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
437 TYPE(cneo_potential_type), POINTER :: cneo_potential
438 TYPE(cp_logger_type), POINTER :: logger
439 TYPE(grid_atom_type), POINTER :: grid_atom
440 TYPE(gto_basis_set_type), POINTER :: basis_1c, nuc_basis, nuc_soft_basis
441 TYPE(harmonics_atom_type), POINTER :: harmonics
442 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
443 TYPE(rho0_atom_type), DIMENSION(:), POINTER :: rho0_atom_set
444 TYPE(rho0_mpole_type), POINTER :: rho0_mpole
445 TYPE(rhoz_cneo_type), DIMENSION(:), POINTER :: rhoz_cneo_set
446 TYPE(rhoz_type), DIMENSION(:), POINTER :: rhoz_set
447 TYPE(section_vals_type), POINTER :: dft_section
448
449 CALL timeset(routinen, handle)
450
451 NULLIFY (logger)
452 logger => cp_get_default_logger()
453
454 NULLIFY (qs_kind_set)
455 NULLIFY (atomic_kind_set)
456 NULLIFY (harmonics)
457 NULLIFY (basis_1c)
458 NULLIFY (rho0_mpole)
459 NULLIFY (rho0_atom_set)
460 NULLIFY (rhoz_set)
461 NULLIFY (cneo_potential)
462 NULLIFY (nuc_basis)
463 NULLIFY (nuc_soft_basis)
464 NULLIFY (rhoz_cneo_set)
465
466 CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set, &
467 atomic_kind_set=atomic_kind_set)
468
469 nkind = SIZE(atomic_kind_set)
470 eps_vrho0 = gapw_control%eps_Vrho0
471
472 ! Initialize rhoz total to zero
473 ! in gapw rhoz is calculated on local the lebedev grids
474 total_rho_core_rspace = 0.0_dp
475 total_rho_nuc_cneo_rspace = 0.0_dp
476
477 CALL get_atomic_kind_set(atomic_kind_set, natom=natom)
478
479 ! Initialize the multipole and the compensation charge type
480 CALL allocate_rho0_mpole(rho0_mpole)
481 CALL allocate_rho0_atom(rho0_atom_set, natom)
482
483 ! Allocate the multipole set
484 CALL allocate_multipoles(rho0_mpole%mp_rho, natom, rho0_mpole%mp_gau, nkind)
485
486 ! Allocate the core density on the radial grid for each kind: rhoz_set
487 CALL allocate_rhoz(rhoz_set, nkind)
488
489 ! For each kind, determine the max l for the compensation charge density
490 lmaxg = gapw_control%lmax_rho0
491 laddg = gapw_control%ladd_rho0
492
493 CALL reallocate(rho0_mpole%lmax0_kind, 1, nkind)
494
495 rho0_mpole%lmax_0 = 0
496 rc_min = 100.0_dp
497 maxnset = 0
498 DO ikind = 1, nkind
499 CALL get_atomic_kind(atomic_kind_set(ikind), atom_list=atom_list, natom=nat)
500 CALL get_qs_kind(qs_kind_set(ikind), &
501 ngrid_rad=nr, &
502 grid_atom=grid_atom, &
503 harmonics=harmonics, &
504 paw_atom=paw_atom, &
505 hard0_radius=rc_orb, &
506 zeff=zeff, &
507 cneo_potential=cneo_potential)
508 CALL get_qs_kind(qs_kind_set(ikind), &
509 basis_set=basis_1c, basis_type="GAPW_1C")
510
511 IF (ASSOCIATED(cneo_potential)) THEN
512 IF (PRESENT(zcore)) THEN
513 IF (zcore == 0.0_dp) THEN
514 cpabort("Electronic TDDFT with CNEO quantum nuclei is not implemented.")
515 END IF
516 END IF
517 cpassert(paw_atom)
518 NULLIFY (nuc_basis, nuc_soft_basis)
519 CALL get_qs_kind(qs_kind_set(ikind), &
520 basis_set=nuc_basis, basis_type="NUC")
521 CALL get_qs_kind(qs_kind_set(ikind), &
522 basis_set=nuc_soft_basis, basis_type="NUC_SOFT")
523 alpha_core = 1.0_dp
524 ELSE
525 CALL get_qs_kind(qs_kind_set(ikind), alpha_core_charge=alpha_core)
526 END IF
527
528 ! Set charge distribution of ionic cores to zero when computing the response-density
529 IF (PRESENT(zcore)) zeff = zcore
530
531 CALL get_gto_basis_set(gto_basis_set=basis_1c, &
532 maxl=maxl, &
533 maxso=maxso, nset=nset)
534
535 maxnset = max(maxnset, nset)
536
537 l_rho1_max = indso(1, harmonics%max_iso_not0)
538
539 maxl_nuc = -1
540 maxso_nuc = 0
541 nset_nuc = 0
542 l_rho1_max_nuc = -1
543 IF (ASSOCIATED(cneo_potential)) THEN
544 CALL get_gto_basis_set(gto_basis_set=nuc_basis, &
545 maxl=maxl_nuc, &
546 maxso=maxso_nuc, nset=nset_nuc)
547 ! Initialize CNEO potential internals
548 CALL init_cneo_potential_internals(cneo_potential, nuc_basis, nuc_soft_basis, &
549 gapw_control, grid_atom)
550 l_rho1_max_nuc = indso(1, cneo_potential%harmonics%max_iso_not0)
551 END IF
552
553 IF (paw_atom) THEN
554 rho0_mpole%lmax0_kind(ikind) = min(2*max(maxl, maxl_nuc), &
555 max(l_rho1_max, l_rho1_max_nuc), &
556 max(maxl, maxl_nuc) + laddg, lmaxg)
557 ELSE
558 rho0_mpole%lmax0_kind(ikind) = 0
559 END IF
560
561 CALL set_qs_kind(qs_kind_set(ikind), lmax_rho0=rho0_mpole%lmax0_kind(ikind))
562
563 rho0_mpole%lmax_0 = max(rho0_mpole%lmax_0, rho0_mpole%lmax0_kind(ikind))
564 rc_min = min(rc_min, rc_orb)
565
566 nchan_s = nsoset(rho0_mpole%lmax0_kind(ikind))
567 nchan_c = ncoset(rho0_mpole%lmax0_kind(ikind))
568 nsotot = maxso*nset
569 nsotot_nuc = maxso_nuc*nset_nuc
570
571 DO iat = 1, nat
572 iatom = atom_list(iat)
573 ! Allocate the multipole for rho1_h rho1_s and rho_z
574 CALL initialize_mpole_rho(rho0_mpole%mp_rho(iatom), nchan_s, nchan_c, zeff)
575 ! Allocate the radial part of rho0_h and rho0_s
576 ! This is calculated on the radial grid centered at the atomic position
577 CALL allocate_rho0_atom_rad(rho0_atom_set(iatom), nr, nchan_s)
578 END DO
579
580 IF (paw_atom) THEN
581 ! Calculate multipoles given by the product of 2 primitives Qlm_gg
582 CALL calculate_mpole_gau(rho0_mpole%mp_gau(ikind)%Qlm_gg, &
583 basis_1c, harmonics, nchan_s, nsotot)
584 END IF
585
586 IF (ASSOCIATED(cneo_potential)) THEN
587 rho0_mpole%do_cneo = .true.
588 ! Calculate multipoles given by the product of two nuclear primitives Qlm_gg
589 CALL calculate_mpole_gau(cneo_potential%Qlm_gg, nuc_basis, &
590 cneo_potential%harmonics, nchan_s, nsotot_nuc)
591 ! initial CNEO quantum nuclear charge density is a simple Zeff sum,
592 ! but it will be calculated from numerical integration during SCF
593 total_rho_nuc_cneo_rspace = total_rho_nuc_cneo_rspace - zeff*nat
594 ELSE
595 ! Calculate the core density rhoz
596 ! exp(-alpha_c**2 r**2)Z(alpha_c**2/pi)**(3/2)
597 ! on the logarithmic radial grid
598 ! WARNING: alpha_core_charge = alpha_c**2
599 CALL calculate_rhoz(rhoz_set(ikind), grid_atom, alpha_core, zeff, &
600 nat, total_rho_core_rspace, harmonics)
601 END IF
602 END DO ! ikind
603 total_rho_core_rspace = -total_rho_core_rspace
604 total_rho_nuc_cneo_rspace = -total_rho_nuc_cneo_rspace
605
606 ! Allocate internals for quantum nuclear densities, if requested
607 CALL get_qs_kind_set(qs_kind_set, cneo_potential_present=cneo_potential_present)
608 IF (cneo_potential_present) THEN
609 CALL allocate_rhoz_cneo_internals(rhoz_cneo_set, atomic_kind_set, &
610 qs_kind_set, qs_env)
611 END IF
612
613 IF (gapw_control%alpha0_hard_from_input) THEN
614 ! The exponent for the compensation charge rho0_hard is read from input
615 rho0_mpole%zet0_h = gapw_control%alpha0_hard
616 ELSE
617 ! Calculate the exponent for the compensation charge rho0_hard
618 rho0_mpole%zet0_h = 0.1_dp
619 DO
620 radius = exp_radius(rho0_mpole%lmax_0, rho0_mpole%zet0_h, eps_vrho0, 1.0_dp)
621 IF (radius <= rc_min) EXIT
622 rho0_mpole%zet0_h = rho0_mpole%zet0_h + 0.1_dp
623 END DO
624
625 END IF
626
627 ! Allocate and calculate the normalization factors for g0_lm_h and g0_lm_s
628 CALL reallocate(rho0_mpole%norm_g0l_h, 0, rho0_mpole%lmax_0)
629 DO l = 0, rho0_mpole%lmax_0
630 rho0_mpole%norm_g0l_h(l) = (2._dp*l + 1._dp)/ &
631 (fourpi*gaussint_sph(rho0_mpole%zet0_h, 2*l))
632 END DO
633
634 ! Allocate and Initialize the g0 gaussians used to build the compensation density
635 ! and calculate the interaction radii
636 max_rpgf0_s = 0.0_dp
637 DO ikind = 1, nkind
638 CALL get_qs_kind(qs_kind_set(ikind), grid_atom=grid_atom)
639 CALL calculate_g0(rho0_mpole, grid_atom, ikind)
640 CALL interaction_radii_g0(rho0_mpole, ikind, eps_vrho0, max_rpgf0_s)
641 END DO
642 rho0_mpole%max_rpgf0_s = max_rpgf0_s
643
644 CALL set_local_rho(local_rho_set, rho0_atom_set=rho0_atom_set, rho0_mpole=rho0_mpole, &
645 rhoz_set=rhoz_set, rhoz_cneo_set=rhoz_cneo_set)
646 local_rho_set%rhoz_tot = total_rho_core_rspace
647 local_rho_set%rhoz_cneo_tot = total_rho_nuc_cneo_rspace
648
649 dft_section => section_vals_get_subs_vals(qs_env%input, "DFT")
650 output_unit = cp_print_key_unit_nr(logger, dft_section, "PRINT%GAPW%RHO0_INFORMATION", &
651 extension=".Log")
652 CALL section_vals_val_get(dft_section, "PRINT%GAPW%RHO0_INFORMATION%UNIT", c_val=unit_str)
653 IF (output_unit > 0) THEN
654 CALL write_rho0_info(rho0_mpole, unit_str, output_unit)
655 END IF
656 CALL cp_print_key_finished_output(output_unit, logger, dft_section, &
657 "PRINT%GAPW%RHO0_INFORMATION")
658
659 CALL timestop(handle)
660
661 END SUBROUTINE init_rho0
662
663! **************************************************************************************************
664!> \brief ...
665!> \param rho0_mpole ...
666!> \param ik ...
667!> \param eps_Vrho0 ...
668!> \param max_rpgf0_s ...
669! **************************************************************************************************
670 SUBROUTINE interaction_radii_g0(rho0_mpole, ik, eps_Vrho0, max_rpgf0_s)
671
672 TYPE(rho0_mpole_type), POINTER :: rho0_mpole
673 INTEGER, INTENT(IN) :: ik
674 REAL(kind=dp), INTENT(IN) :: eps_vrho0
675 REAL(kind=dp), INTENT(INOUT) :: max_rpgf0_s
676
677 INTEGER :: l, lmax
678 REAL(kind=dp) :: r_h, z0_h
679 REAL(kind=dp), DIMENSION(:), POINTER :: ng0_h
680
681 CALL get_rho0_mpole(rho0_mpole, ikind=ik, l0_ikind=lmax, &
682 zet0_h=z0_h, norm_g0l_h=ng0_h)
683 r_h = 0.0_dp
684 DO l = 0, lmax
685 r_h = max(r_h, exp_radius(l, z0_h, eps_vrho0, ng0_h(l), rlow=r_h))
686 END DO
687
688 rho0_mpole%mp_gau(ik)%rpgf0_h = r_h
689 rho0_mpole%mp_gau(ik)%rpgf0_s = r_h
690 max_rpgf0_s = max(max_rpgf0_s, r_h)
691
692 END SUBROUTINE interaction_radii_g0
693
694END 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, npgf_seg_sum)
...
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
A collection of functions used by CNEO-DFT (see J. Chem. Theory Comput. 2025, 21, 16,...
subroutine, public init_cneo_potential_internals(potential, nuc_basis, nuc_soft_basis, gapw_control, grid_atom)
...
subroutine, public allocate_rhoz_cneo_internals(rhoz_cneo_set, atomic_kind_set, qs_kind_set, qs_env)
...
Types used by CNEO-DFT (see J. Chem. Theory Comput. 2025, 21, 16, 7865–7877)
Utility functions for CNEO-DFT (see J. Chem. Theory Comput. 2025, 21, 16, 7865–7877)
subroutine, public cneo_scatter(ain, aout, nbas, n2oindex)
Mostly copied from qs_oce_methods::prj_scatter.
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.
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.
subroutine, public get_qs_kind_set(qs_kind_set, all_potential_present, tnadd_potential_present, gth_potential_present, sgp_potential_present, paw_atom_present, dft_plus_u_atom_present, maxcgf, maxsgf, maxco, maxco_proj, maxgtops, maxlgto, maxlprj, maxnset, maxsgf_set, ncgf, npgf, nset, nsgf, nshell, maxpol, maxlppl, maxlppnl, maxppnl, nelectron, maxder, max_ngrid_rad, max_sph_harm, maxg_iso_not0, lmax_rho0, basis_rcut, basis_type, total_zeff_corr, npgf_seg, cneo_potential_present, nkind_q, natom_q)
Get attributes of an atomic kind set.
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 calculate_rhoz(rhoz, grid_atom, alpha, zeff, natom, rhoz_tot, harmonics)
...
subroutine, public set_local_rho(local_rho_set, rho_atom_set, rho0_atom_set, rho0_mpole, rhoz_set, rhoz_cneo_set)
...
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, rhoz_cneo_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, rhoz_cneo_s_rs, rhoz_cneo_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.