(git:e7e05ae)
qs_rho0_types.F
Go to the documentation of this file.
1 !--------------------------------------------------------------------------------------------------!
2 ! CP2K: A general program to perform molecular dynamics simulations !
3 ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 ! !
5 ! SPDX-License-Identifier: GPL-2.0-or-later !
6 !--------------------------------------------------------------------------------------------------!
7 ! **************************************************************************************************
9 
10  USE cp_units, ONLY: cp_unit_from_cp2k
11  USE kinds, ONLY: dp
12  USE mathconstants, ONLY: fourpi,&
13  pi,&
14  rootpi
15  USE memory_utilities, ONLY: reallocate
16  USE pw_types, ONLY: pw_c1d_gs_type,&
17  pw_r3d_rs_type
18  USE qs_grid_atom, ONLY: grid_atom_type
19  USE qs_rho_atom_types, ONLY: rho_atom_coeff
20  USE whittaker, ONLY: whittaker_c0a,&
22 #include "./base/base_uses.f90"
23 
24  IMPLICIT NONE
25 
26  PRIVATE
27 
28 ! *** Global parameters (only in this module)
29 
30  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_rho0_types'
31 
32 ! *** Define multipole type ***
33 
34 ! **************************************************************************************************
35  TYPE mpole_rho_atom
36  REAL(dp), DIMENSION(:), POINTER :: Qlm_h, &
37  Qlm_s, &
38  Qlm_tot, &
39  Qlm_car
40  REAL(dp) :: Qlm_z
41  REAL(dp), DIMENSION(2) :: Q0
42  END TYPE mpole_rho_atom
43 
44 ! **************************************************************************************************
45  TYPE mpole_gau_overlap
46  REAL(dp), DIMENSION(:, :, :), POINTER :: Qlm_gg
47  REAL(dp), DIMENSION(:, :), POINTER :: g0_h, vg0_h
48  REAL(dp) :: rpgf0_h, rpgf0_s
49  END TYPE mpole_gau_overlap
50 
51 ! **************************************************************************************************
52  TYPE rho0_mpole_type
53  TYPE(mpole_rho_atom), DIMENSION(:), POINTER :: mp_rho
54  TYPE(mpole_gau_overlap), DIMENSION(:), &
55  POINTER :: mp_gau
56  REAL(dp) :: zet0_h, &
57  total_rho0_h
58  REAL(dp) :: max_rpgf0_s
59  REAL(dp), DIMENSION(:), POINTER :: norm_g0l_h
60  INTEGER, DIMENSION(:), POINTER :: lmax0_kind
61  INTEGER :: lmax_0, igrid_zet0_s
62  TYPE(pw_r3d_rs_type), POINTER :: rho0_s_rs
63  TYPE(pw_c1d_gs_type), POINTER :: rho0_s_gs
64  END TYPE rho0_mpole_type
65 
66 ! **************************************************************************************************
67  TYPE rho0_atom_type
68  TYPE(rho_atom_coeff), POINTER :: rho0_rad_h, &
69  vrho0_rad_h
70  END TYPE rho0_atom_type
71 
72 ! Public Types
73 
74  PUBLIC :: mpole_rho_atom, mpole_gau_overlap, &
75  rho0_atom_type, rho0_mpole_type
76 
77 ! Public Subroutine
78 
84 
85 CONTAINS
86 
87 ! **************************************************************************************************
88 !> \brief ...
89 !> \param mp_rho ...
90 !> \param natom ...
91 !> \param mp_gau ...
92 !> \param nkind ...
93 ! **************************************************************************************************
94  SUBROUTINE allocate_multipoles(mp_rho, natom, mp_gau, nkind)
95 
96  TYPE(mpole_rho_atom), DIMENSION(:), POINTER :: mp_rho
97  INTEGER, INTENT(IN) :: natom
98  TYPE(mpole_gau_overlap), DIMENSION(:), POINTER :: mp_gau
99  INTEGER, INTENT(IN) :: nkind
100 
101  INTEGER :: iat, ikind
102 
103  IF (ASSOCIATED(mp_rho)) THEN
104  CALL deallocate_mpole_rho(mp_rho)
105  END IF
106 
107  ALLOCATE (mp_rho(natom))
108 
109  DO iat = 1, natom
110  NULLIFY (mp_rho(iat)%Qlm_h)
111  NULLIFY (mp_rho(iat)%Qlm_s)
112  NULLIFY (mp_rho(iat)%Qlm_tot)
113  NULLIFY (mp_rho(iat)%Qlm_car)
114  END DO
115 
116  IF (ASSOCIATED(mp_gau)) THEN
117  CALL deallocate_mpole_gau(mp_gau)
118  END IF
119 
120  ALLOCATE (mp_gau(nkind))
121 
122  DO ikind = 1, nkind
123  NULLIFY (mp_gau(ikind)%Qlm_gg)
124  NULLIFY (mp_gau(ikind)%g0_h)
125  NULLIFY (mp_gau(ikind)%vg0_h)
126  mp_gau(ikind)%rpgf0_h = 0.0_dp
127  mp_gau(ikind)%rpgf0_s = 0.0_dp
128  END DO
129 
130  END SUBROUTINE allocate_multipoles
131 
132 ! **************************************************************************************************
133 !> \brief ...
134 !> \param rho0_set ...
135 !> \param natom ...
136 ! **************************************************************************************************
137  SUBROUTINE allocate_rho0_atom(rho0_set, natom)
138 
139  TYPE(rho0_atom_type), DIMENSION(:), POINTER :: rho0_set
140  INTEGER, INTENT(IN) :: natom
141 
142  INTEGER :: iat
143 
144  IF (ASSOCIATED(rho0_set)) THEN
145  CALL deallocate_rho0_atom(rho0_set)
146  END IF
147 
148  ALLOCATE (rho0_set(natom))
149 
150  DO iat = 1, natom
151  NULLIFY (rho0_set(iat)%rho0_rad_h)
152  NULLIFY (rho0_set(iat)%vrho0_rad_h)
153  END DO
154 
155  END SUBROUTINE allocate_rho0_atom
156 
157 ! **************************************************************************************************
158 !> \brief ...
159 !> \param rho0_atom ...
160 !> \param nr ...
161 !> \param nchannels ...
162 ! **************************************************************************************************
163  SUBROUTINE allocate_rho0_atom_rad(rho0_atom, nr, nchannels)
164 
165  TYPE(rho0_atom_type), INTENT(OUT) :: rho0_atom
166  INTEGER, INTENT(IN) :: nr, nchannels
167 
168  ALLOCATE (rho0_atom%rho0_rad_h)
169 
170  NULLIFY (rho0_atom%rho0_rad_h%r_coef)
171  ALLOCATE (rho0_atom%rho0_rad_h%r_coef(1:nr, 1:nchannels))
172  rho0_atom%rho0_rad_h%r_coef = 0.0_dp
173 
174  ALLOCATE (rho0_atom%vrho0_rad_h)
175 
176  NULLIFY (rho0_atom%vrho0_rad_h%r_coef)
177  ALLOCATE (rho0_atom%vrho0_rad_h%r_coef(1:nr, 1:nchannels))
178  rho0_atom%vrho0_rad_h%r_coef = 0.0_dp
179 
180  END SUBROUTINE allocate_rho0_atom_rad
181 
182 ! **************************************************************************************************
183 !> \brief ...
184 !> \param rho0 ...
185 ! **************************************************************************************************
186  SUBROUTINE allocate_rho0_mpole(rho0)
187 
188  TYPE(rho0_mpole_type), POINTER :: rho0
189 
190  IF (ASSOCIATED(rho0)) THEN
191  CALL deallocate_rho0_mpole(rho0)
192  END IF
193 
194  ALLOCATE (rho0)
195 
196  NULLIFY (rho0%mp_rho)
197  NULLIFY (rho0%mp_gau)
198  NULLIFY (rho0%norm_g0l_h)
199  NULLIFY (rho0%lmax0_kind)
200  NULLIFY (rho0%rho0_s_rs)
201  NULLIFY (rho0%rho0_s_gs)
202 
203  END SUBROUTINE allocate_rho0_mpole
204 
205 ! **************************************************************************************************
206 !> \brief ...
207 !> \param rho0_mpole ...
208 !> \param grid_atom ...
209 !> \param ik ...
210 ! **************************************************************************************************
211  SUBROUTINE calculate_g0(rho0_mpole, grid_atom, ik)
212 
213  TYPE(rho0_mpole_type), POINTER :: rho0_mpole
214  TYPE(grid_atom_type), POINTER :: grid_atom
215  INTEGER, INTENT(IN) :: ik
216 
217  INTEGER :: ir, l, lmax, nr
218  REAL(dp) :: c1, prefactor, root_z_h, z_h
219  REAL(dp), ALLOCATABLE, DIMENSION(:) :: erf_z_h, gexp, gh_tmp, int1, int2
220 
221  nr = grid_atom%nr
222  lmax = rho0_mpole%lmax0_kind(ik)
223  z_h = rho0_mpole%zet0_h
224  root_z_h = sqrt(z_h)
225 
226 ! Allocate g0
227  CALL reallocate(rho0_mpole%mp_gau(ik)%g0_h, 1, nr, 0, lmax)
228  CALL reallocate(rho0_mpole%mp_gau(ik)%vg0_h, 1, nr, 0, lmax)
229 
230  ALLOCATE (gexp(nr), gh_tmp(nr), erf_z_h(nr), int1(nr), int2(nr))
231 
232  gh_tmp(1:nr) = exp(-z_h*grid_atom%rad2(1:nr))
233 
234  DO ir = 1, nr
235  erf_z_h(ir) = erf(grid_atom%rad(ir)*root_z_h)
236  END DO
237 
238  DO ir = 1, nr
239  IF (gh_tmp(ir) < 1.0e-30_dp) gh_tmp(ir) = 0.0_dp
240  END DO
241 
242  gexp(1:nr) = gh_tmp(1:nr)
243  rho0_mpole%mp_gau(ik)%g0_h(1:nr, 0) = gh_tmp(1:nr)* &
244  rho0_mpole%norm_g0l_h(0)
245  CALL whittaker_c0a(int1, grid_atom%rad, gh_tmp, erf_z_h, z_h, 0, 0, nr)
246  CALL whittaker_ci(int2, grid_atom%rad, gh_tmp, z_h, 0, nr)
247 
248  prefactor = fourpi*rho0_mpole%norm_g0l_h(0)
249 
250  c1 = sqrt(pi*pi*pi/(z_h*z_h*z_h))*rho0_mpole%norm_g0l_h(0)
251 
252  DO ir = 1, nr
253  rho0_mpole%mp_gau(ik)%vg0_h(ir, 0) = c1*erf_z_h(ir)*grid_atom%oorad2l(ir, 1)
254  END DO
255 
256  DO l = 1, lmax
257  gh_tmp(1:nr) = gh_tmp(1:nr)*grid_atom%rad(1:nr)
258  rho0_mpole%mp_gau(ik)%g0_h(1:nr, l) = gh_tmp(1:nr)* &
259  rho0_mpole%norm_g0l_h(l)
260 
261  prefactor = fourpi/(2.0_dp*l + 1.0_dp)*rho0_mpole%norm_g0l_h(l)
262  CALL whittaker_c0a(int1, grid_atom%rad, gexp, erf_z_h, z_h, l, l, nr)
263  DO ir = 1, nr
264  rho0_mpole%mp_gau(ik)%vg0_h(ir, l) = prefactor*(int1(ir) + &
265  int2(ir)*grid_atom%rad2l(ir, l))
266  END DO
267 
268  END DO ! l
269 
270  DEALLOCATE (gexp, erf_z_h, gh_tmp, int1, int2)
271  END SUBROUTINE calculate_g0
272 
273 ! **************************************************************************************************
274 !> \brief ...
275 !> \param mp_gau ...
276 ! **************************************************************************************************
277  SUBROUTINE deallocate_mpole_gau(mp_gau)
278 
279  TYPE(mpole_gau_overlap), DIMENSION(:), POINTER :: mp_gau
280 
281  INTEGER :: ikind, nkind
282 
283  nkind = SIZE(mp_gau)
284 
285  DO ikind = 1, nkind
286 
287  IF (ASSOCIATED(mp_gau(ikind)%Qlm_gg)) THEN
288  DEALLOCATE (mp_gau(ikind)%Qlm_gg)
289  END IF
290 
291  DEALLOCATE (mp_gau(ikind)%g0_h)
292 
293  DEALLOCATE (mp_gau(ikind)%vg0_h)
294  END DO
295 
296  DEALLOCATE (mp_gau)
297 
298  END SUBROUTINE deallocate_mpole_gau
299 
300 ! **************************************************************************************************
301 !> \brief ...
302 !> \param mp_rho ...
303 ! **************************************************************************************************
304  SUBROUTINE deallocate_mpole_rho(mp_rho)
305 
306  TYPE(mpole_rho_atom), DIMENSION(:), POINTER :: mp_rho
307 
308  INTEGER :: iat, natom
309 
310  natom = SIZE(mp_rho)
311 
312  DO iat = 1, natom
313  DEALLOCATE (mp_rho(iat)%Qlm_h)
314  DEALLOCATE (mp_rho(iat)%Qlm_s)
315  DEALLOCATE (mp_rho(iat)%Qlm_tot)
316  DEALLOCATE (mp_rho(iat)%Qlm_car)
317  END DO
318 
319  DEALLOCATE (mp_rho)
320 
321  END SUBROUTINE deallocate_mpole_rho
322 
323 ! **************************************************************************************************
324 !> \brief ...
325 !> \param rho0_atom_set ...
326 ! **************************************************************************************************
327  SUBROUTINE deallocate_rho0_atom(rho0_atom_set)
328 
329  TYPE(rho0_atom_type), DIMENSION(:), POINTER :: rho0_atom_set
330 
331  INTEGER :: iat, natom
332 
333  IF (ASSOCIATED(rho0_atom_set)) THEN
334 
335  natom = SIZE(rho0_atom_set)
336 
337  DO iat = 1, natom
338  IF (ASSOCIATED(rho0_atom_set(iat)%rho0_rad_h)) THEN
339  DEALLOCATE (rho0_atom_set(iat)%rho0_rad_h%r_coef)
340  DEALLOCATE (rho0_atom_set(iat)%rho0_rad_h)
341  END IF
342  IF (ASSOCIATED(rho0_atom_set(iat)%vrho0_rad_h)) THEN
343  DEALLOCATE (rho0_atom_set(iat)%vrho0_rad_h%r_coef)
344  DEALLOCATE (rho0_atom_set(iat)%vrho0_rad_h)
345  END IF
346  END DO
347 
348  DEALLOCATE (rho0_atom_set)
349  ELSE
350  CALL cp_abort(__location__, &
351  "The pointer rho0_atom_set is not associated and "// &
352  "cannot be deallocated")
353  END IF
354 
355  END SUBROUTINE deallocate_rho0_atom
356 ! **************************************************************************************************
357 !> \brief ...
358 !> \param rho0 ...
359 ! **************************************************************************************************
360  SUBROUTINE deallocate_rho0_mpole(rho0)
361 
362  TYPE(rho0_mpole_type), POINTER :: rho0
363 
364  IF (ASSOCIATED(rho0)) THEN
365 
366  IF (ASSOCIATED(rho0%mp_gau)) CALL deallocate_mpole_gau(rho0%mp_gau)
367 
368  IF (ASSOCIATED(rho0%mp_rho)) CALL deallocate_mpole_rho(rho0%mp_rho)
369 
370  IF (ASSOCIATED(rho0%lmax0_kind)) THEN
371  DEALLOCATE (rho0%lmax0_kind)
372  END IF
373 
374  IF (ASSOCIATED(rho0%norm_g0l_h)) THEN
375  DEALLOCATE (rho0%norm_g0l_h)
376  END IF
377 
378  IF (ASSOCIATED(rho0%rho0_s_rs)) THEN
379  CALL rho0%rho0_s_rs%release()
380  DEALLOCATE (rho0%rho0_s_rs)
381  END IF
382 
383  IF (ASSOCIATED(rho0%rho0_s_gs)) THEN
384  CALL rho0%rho0_s_gs%release()
385  DEALLOCATE (rho0%rho0_s_gs)
386 
387  END IF
388  DEALLOCATE (rho0)
389  ELSE
390  CALL cp_abort(__location__, &
391  "The pointer rho0 is not associated and "// &
392  "cannot be deallocated")
393  END IF
394 
395  END SUBROUTINE deallocate_rho0_mpole
396 
397 ! **************************************************************************************************
398 !> \brief ...
399 !> \param rho0_mpole ...
400 !> \param g0_h ...
401 !> \param vg0_h ...
402 !> \param iat ...
403 !> \param ikind ...
404 !> \param lmax_0 ...
405 !> \param l0_ikind ...
406 !> \param mp_gau_ikind ...
407 !> \param mp_rho ...
408 !> \param norm_g0l_h ...
409 !> \param Qlm_gg ...
410 !> \param Qlm_car ...
411 !> \param Qlm_tot ...
412 !> \param zet0_h ...
413 !> \param igrid_zet0_s ...
414 !> \param rpgf0_h ...
415 !> \param rpgf0_s ...
416 !> \param max_rpgf0_s ...
417 !> \param rho0_s_rs ...
418 !> \param rho0_s_gs ...
419 ! **************************************************************************************************
420  SUBROUTINE get_rho0_mpole(rho0_mpole, g0_h, vg0_h, iat, ikind, lmax_0, l0_ikind, &
421  mp_gau_ikind, mp_rho, norm_g0l_h, &
422  Qlm_gg, Qlm_car, Qlm_tot, &
423  zet0_h, igrid_zet0_s, rpgf0_h, rpgf0_s, &
424  max_rpgf0_s, rho0_s_rs, rho0_s_gs)
425 
426  TYPE(rho0_mpole_type), POINTER :: rho0_mpole
427  REAL(dp), DIMENSION(:, :), OPTIONAL, POINTER :: g0_h, vg0_h
428  INTEGER, INTENT(IN), OPTIONAL :: iat, ikind
429  INTEGER, INTENT(OUT), OPTIONAL :: lmax_0, l0_ikind
430  TYPE(mpole_gau_overlap), OPTIONAL, POINTER :: mp_gau_ikind
431  TYPE(mpole_rho_atom), DIMENSION(:), OPTIONAL, &
432  POINTER :: mp_rho
433  REAL(dp), DIMENSION(:), OPTIONAL, POINTER :: norm_g0l_h
434  REAL(dp), DIMENSION(:, :, :), OPTIONAL, POINTER :: qlm_gg
435  REAL(dp), DIMENSION(:), OPTIONAL, POINTER :: qlm_car, qlm_tot
436  REAL(dp), INTENT(OUT), OPTIONAL :: zet0_h
437  INTEGER, INTENT(OUT), OPTIONAL :: igrid_zet0_s
438  REAL(dp), INTENT(OUT), OPTIONAL :: rpgf0_h, rpgf0_s, max_rpgf0_s
439  TYPE(pw_r3d_rs_type), OPTIONAL, POINTER :: rho0_s_rs
440  TYPE(pw_c1d_gs_type), OPTIONAL, POINTER :: rho0_s_gs
441 
442  IF (ASSOCIATED(rho0_mpole)) THEN
443 
444  IF (PRESENT(lmax_0)) lmax_0 = rho0_mpole%lmax_0
445  IF (PRESENT(mp_rho)) mp_rho => rho0_mpole%mp_rho
446  IF (PRESENT(norm_g0l_h)) norm_g0l_h => rho0_mpole%norm_g0l_h
447  IF (PRESENT(zet0_h)) zet0_h = rho0_mpole%zet0_h
448  IF (PRESENT(igrid_zet0_s)) igrid_zet0_s = rho0_mpole%igrid_zet0_s
449  IF (PRESENT(max_rpgf0_s)) max_rpgf0_s = rho0_mpole%max_rpgf0_s
450  IF (PRESENT(rho0_s_rs)) rho0_s_rs => rho0_mpole%rho0_s_rs
451  IF (PRESENT(rho0_s_gs)) rho0_s_gs => rho0_mpole%rho0_s_gs
452 
453  IF (PRESENT(ikind)) THEN
454  IF (PRESENT(l0_ikind)) l0_ikind = rho0_mpole%lmax0_kind(ikind)
455  IF (PRESENT(mp_gau_ikind)) mp_gau_ikind => rho0_mpole%mp_gau(ikind)
456  IF (PRESENT(g0_h)) g0_h => rho0_mpole%mp_gau(ikind)%g0_h
457  IF (PRESENT(vg0_h)) vg0_h => rho0_mpole%mp_gau(ikind)%vg0_h
458  IF (PRESENT(qlm_gg)) qlm_gg => rho0_mpole%mp_gau(ikind)%Qlm_gg
459  IF (PRESENT(rpgf0_h)) rpgf0_h = rho0_mpole%mp_gau(ikind)%rpgf0_h
460  IF (PRESENT(rpgf0_s)) rpgf0_s = rho0_mpole%mp_gau(ikind)%rpgf0_s
461  END IF
462  IF (PRESENT(iat)) THEN
463  IF (PRESENT(qlm_car)) qlm_car => rho0_mpole%mp_rho(iat)%Qlm_car
464  IF (PRESENT(qlm_tot)) qlm_tot => rho0_mpole%mp_rho(iat)%Qlm_tot
465  END IF
466 
467  ELSE
468  cpabort("The pointer rho0_mpole is not associated")
469  END IF
470 
471  END SUBROUTINE get_rho0_mpole
472 
473 ! **************************************************************************************************
474 !> \brief ...
475 !> \param mp_rho ...
476 !> \param nchan_s ...
477 !> \param nchan_c ...
478 !> \param zeff ...
479 ! **************************************************************************************************
480  SUBROUTINE initialize_mpole_rho(mp_rho, nchan_s, nchan_c, zeff)
481 
482  TYPE(mpole_rho_atom) :: mp_rho
483  INTEGER, INTENT(IN) :: nchan_s, nchan_c
484  REAL(kind=dp), INTENT(IN) :: zeff
485 
486  CALL reallocate(mp_rho%Qlm_h, 1, nchan_s)
487  CALL reallocate(mp_rho%Qlm_s, 1, nchan_s)
488  CALL reallocate(mp_rho%Qlm_tot, 1, nchan_s)
489  CALL reallocate(mp_rho%Qlm_car, 1, nchan_c)
490 
491  mp_rho%Qlm_h = 0.0_dp
492  mp_rho%Qlm_s = 0.0_dp
493  mp_rho%Qlm_tot = 0.0_dp
494  mp_rho%Qlm_car = 0.0_dp
495  mp_rho%Qlm_z = -2.0_dp*rootpi*zeff
496  mp_rho%Q0 = 0.0_dp
497 
498  END SUBROUTINE initialize_mpole_rho
499 
500 ! **************************************************************************************************
501 !> \brief ...
502 !> \param rho0_mpole ...
503 !> \param unit_str ...
504 !> \param output_unit ...
505 ! **************************************************************************************************
506  SUBROUTINE write_rho0_info(rho0_mpole, unit_str, output_unit)
507 
508  TYPE(rho0_mpole_type), POINTER :: rho0_mpole
509  CHARACTER(LEN=*), INTENT(IN) :: unit_str
510  INTEGER, INTENT(in) :: output_unit
511 
512  INTEGER :: ikind, l, nkind
513  REAL(dp) :: conv
514 
515  IF (ASSOCIATED(rho0_mpole)) THEN
516  conv = cp_unit_from_cp2k(1.0_dp, trim(unit_str))
517 
518  WRITE (unit=output_unit, fmt="(/,T5,A,/)") &
519  "*** Compensation density charges data set ***"
520  WRITE (unit=output_unit, fmt="(T2,A,T35,f16.10)") &
521  "- Rho0 exponent :", rho0_mpole%zet0_h
522  WRITE (unit=output_unit, fmt="(T2,A,T35,I5)") &
523  "- Global max l :", rho0_mpole%lmax_0
524 
525  WRITE (unit=output_unit, fmt="(T2,A)") &
526  "- Normalization constants for g0"
527  DO l = 0, rho0_mpole%lmax_0
528  WRITE (unit=output_unit, fmt="(T20,A,T31,I2,T38,A,f15.5)") &
529  "ang. mom.= ", l, " hard= ", rho0_mpole%norm_g0l_h(l)
530  END DO
531 
532  nkind = SIZE(rho0_mpole%lmax0_kind, 1)
533  DO ikind = 1, nkind
534  WRITE (unit=output_unit, fmt="(/,T2,A,T55,I2)") &
535  "- rho0 max L and radii in "//trim(unit_str)// &
536  " for the atom kind :", ikind
537 
538  WRITE (unit=output_unit, fmt="(T2,T20,A,T55,I5)") &
539  "=> l max :", rho0_mpole%lmax0_kind(ikind)
540 
541  WRITE (unit=output_unit, fmt="(T2,T20,A,T55,f20.10)") &
542  "=> max radius of g0: ", &
543  rho0_mpole%mp_gau(ikind)%rpgf0_h*conv
544  END DO ! ikind
545 
546  ELSE
547  WRITE (unit=output_unit, fmt="(/,T5,A,/)") &
548  ' WARNING: I cannot print rho0, it is not associated'
549  END IF
550 
551  END SUBROUTINE write_rho0_info
552 END MODULE qs_rho0_types
unit conversion facility
Definition: cp_units.F:30
real(kind=dp) function, public cp_unit_from_cp2k(value, unit_str, defaults, power)
converts from the internal cp2k units to the given unit
Definition: cp_units.F:1179
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
Definition of mathematical constants and functions.
Definition: mathconstants.F:16
real(kind=dp), parameter, public pi
real(kind=dp), parameter, public rootpi
real(kind=dp), parameter, public fourpi
Utility routines for the memory handling.
subroutine, public allocate_multipoles(mp_rho, natom, mp_gau, nkind)
...
Definition: qs_rho0_types.F:95
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 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 deallocate_rho0_mpole(rho0)
...
subroutine, public deallocate_rho0_atom(rho0_atom_set)
...
Calculates special integrals.
Definition: whittaker.F:12
subroutine, public whittaker_c0a(wc, r, expa, erfa, alpha, l1, l2, n)
int(y^(2+l1+l2) * exp(-alpha*y*y),y=0..x) / x^(l2+1); wc(:) :: output r(:) :: coordinate expa(:) :: e...
Definition: whittaker.F:52
subroutine, public whittaker_ci(wc, r, expa, alpha, l, n)
int(y^(l+1) * exp(-alpha*y*y),y=x..infinity);
Definition: whittaker.F:340