(git:e7e05ae)
qs_harmonics_atom.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 
11  gto_basis_set_type
12  USE kinds, ONLY: dp
13  USE lebedev, ONLY: lebedev_grid
14  USE memory_utilities, ONLY: reallocate
15  USE orbital_pointers, ONLY: indco,&
16  indso,&
17  nco,&
18  ncoset,&
19  nso,&
20  nsoset
22  USE spherical_harmonics, ONLY: dy_lm,&
23  y_lm
24 #include "./base/base_uses.f90"
25 
26  IMPLICIT NONE
27 
28  PRIVATE
29 
30  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_harmonics_atom'
31 
32  TYPE harmonics_atom_type
33  INTEGER :: max_s_harm, llmax, &
34  max_iso_not0, &
35  dmax_iso_not0, &
36  damax_iso_not0, &
37  ngrid
38  REAL(dp), DIMENSION(:, :), POINTER :: a, slm
39  REAL(dp), DIMENSION(:, :, :), POINTER :: dslm, dslm_dxyz
40  REAL(dp), DIMENSION(:, :, :), POINTER :: my_CG
41  REAL(dp), DIMENSION(:, :, :, :), POINTER :: my_CG_dxyz
42  REAL(dp), DIMENSION(:, :, :, :), POINTER :: my_CG_dxyz_asym
43  REAL(dp), DIMENSION(:), POINTER :: slm_int
44 
45  END TYPE harmonics_atom_type
46 
47  PUBLIC :: allocate_harmonics_atom, &
50  get_none0_cg_list
51 
52  PUBLIC :: harmonics_atom_type, get_maxl_cg
53 
54  INTERFACE get_none0_cg_list
55  MODULE PROCEDURE get_none0_cg_list3, get_none0_cg_list4
56  END INTERFACE
57 
58 CONTAINS
59 
60 ! **************************************************************************************************
61 !> \brief Allocate a spherical harmonics set for the atom grid.
62 !> \param harmonics ...
63 !> \version 1.0
64 ! **************************************************************************************************
65  SUBROUTINE allocate_harmonics_atom(harmonics)
66 
67  TYPE(harmonics_atom_type), POINTER :: harmonics
68 
69  IF (ASSOCIATED(harmonics)) CALL deallocate_harmonics_atom(harmonics)
70 
71  ALLOCATE (harmonics)
72 
73  harmonics%max_s_harm = 0
74  harmonics%llmax = 0
75  harmonics%max_iso_not0 = 0
76  harmonics%dmax_iso_not0 = 0
77  harmonics%damax_iso_not0 = 0
78  harmonics%ngrid = 0
79 
80  NULLIFY (harmonics%slm)
81  NULLIFY (harmonics%dslm)
82  NULLIFY (harmonics%dslm_dxyz)
83  NULLIFY (harmonics%slm_int)
84  NULLIFY (harmonics%my_CG)
85  NULLIFY (harmonics%my_CG_dxyz)
86  NULLIFY (harmonics%my_CG_dxyz_asym)
87  NULLIFY (harmonics%a)
88 
89  END SUBROUTINE allocate_harmonics_atom
90 
91 ! **************************************************************************************************
92 !> \brief Deallocate the spherical harmonics set for the atom grid.
93 !> \param harmonics ...
94 !> \version 1.0
95 ! **************************************************************************************************
96  SUBROUTINE deallocate_harmonics_atom(harmonics)
97 
98  TYPE(harmonics_atom_type), POINTER :: harmonics
99 
100  IF (ASSOCIATED(harmonics)) THEN
101 
102  IF (ASSOCIATED(harmonics%slm)) THEN
103  DEALLOCATE (harmonics%slm)
104  END IF
105 
106  IF (ASSOCIATED(harmonics%dslm)) THEN
107  DEALLOCATE (harmonics%dslm)
108  END IF
109 
110  IF (ASSOCIATED(harmonics%dslm_dxyz)) THEN
111  DEALLOCATE (harmonics%dslm_dxyz)
112  END IF
113 
114  IF (ASSOCIATED(harmonics%slm_int)) THEN
115  DEALLOCATE (harmonics%slm_int)
116  END IF
117 
118  IF (ASSOCIATED(harmonics%my_CG)) THEN
119  DEALLOCATE (harmonics%my_CG)
120  END IF
121 
122  IF (ASSOCIATED(harmonics%my_CG_dxyz)) THEN
123  DEALLOCATE (harmonics%my_CG_dxyz)
124  END IF
125 
126  IF (ASSOCIATED(harmonics%my_CG_dxyz_asym)) THEN
127  DEALLOCATE (harmonics%my_CG_dxyz_asym)
128  END IF
129 
130  IF (ASSOCIATED(harmonics%a)) THEN
131  DEALLOCATE (harmonics%a)
132  END IF
133 
134  DEALLOCATE (harmonics)
135  ELSE
136  CALL cp_abort(__location__, &
137  "The pointer harmonics is not associated and "// &
138  "cannot be deallocated")
139  END IF
140 
141  END SUBROUTINE deallocate_harmonics_atom
142 
143 ! **************************************************************************************************
144 !> \brief ...
145 !> \param harmonics ...
146 !> \param my_CG ...
147 !> \param na ...
148 !> \param llmax ...
149 !> \param maxs ...
150 !> \param max_s_harm ...
151 !> \param ll ...
152 !> \param wa ...
153 !> \param azi ...
154 !> \param pol ...
155 !> \note Slight refactoring + OMP parallelized (03.2020 A. Bussy)
156 ! **************************************************************************************************
157  SUBROUTINE create_harmonics_atom(harmonics, my_CG, na, llmax, maxs, max_s_harm, ll, wa, azi, pol)
158 
159  TYPE(harmonics_atom_type), POINTER :: harmonics
160  REAL(dp), DIMENSION(:, :, :), POINTER :: my_cg
161  INTEGER, INTENT(IN) :: na, llmax, maxs, max_s_harm, ll
162  REAL(dp), DIMENSION(:), INTENT(IN) :: wa, azi, pol
163 
164  CHARACTER(len=*), PARAMETER :: routinen = 'create_harmonics_atom'
165 
166  INTEGER :: handle, i, ia, ic, is, is1, is2, iso, &
167  iso1, iso2, l, l1, l2, lx, ly, lz, m, &
168  m1, m2, n
169  REAL(dp) :: drx, dry, drz, rx, ry, rz
170  REAL(dp), DIMENSION(2) :: cin, dylm
171  REAL(dp), DIMENSION(:), POINTER :: slm_int, y
172  REAL(dp), DIMENSION(:, :), POINTER :: dc, slm
173  REAL(dp), DIMENSION(:, :, :), POINTER :: dslm_dxyz
174 
175  CALL timeset(routinen, handle)
176 
177  NULLIFY (y, slm, dslm_dxyz, dc)
178 
179  cpassert(ASSOCIATED(harmonics))
180 
181  harmonics%max_s_harm = max_s_harm
182  harmonics%llmax = llmax
183  harmonics%ngrid = na
184 
185  NULLIFY (harmonics%my_CG, harmonics%my_CG_dxyz, harmonics%my_CG_dxyz_asym)
186  CALL reallocate(harmonics%my_CG, 1, maxs, 1, maxs, 1, max_s_harm)
187  CALL reallocate(harmonics%my_CG_dxyz, 1, 3, 1, maxs, 1, maxs, 1, max_s_harm)
188  CALL reallocate(harmonics%my_CG_dxyz_asym, 1, 3, 1, maxs, 1, maxs, 1, max_s_harm)
189 
190  DO i = 1, max_s_harm
191  DO is1 = 1, maxs
192  harmonics%my_CG(1:maxs, is1, i) = my_cg(1:maxs, is1, i)
193  END DO
194  END DO
195 
196  ! allocate and calculate the spherical harmonics LM for this grid
197  ! and their derivatives
198  NULLIFY (harmonics%slm, harmonics%dslm, harmonics%dslm_dxyz, harmonics%a, harmonics%slm_int)
199  CALL reallocate(harmonics%slm, 1, na, 1, max_s_harm)
200  CALL reallocate(harmonics%dslm, 1, 2, 1, na, 1, maxs)
201  CALL reallocate(harmonics%dslm_dxyz, 1, 3, 1, na, 1, max_s_harm)
202  CALL reallocate(harmonics%a, 1, 3, 1, na)
203  CALL reallocate(harmonics%slm_int, 1, max_s_harm)
204 
205  NULLIFY (slm, dslm_dxyz, slm_int)
206  slm => harmonics%slm
207  dslm_dxyz => harmonics%dslm_dxyz
208  dslm_dxyz = 0.0_dp
209  slm_int => harmonics%slm_int
210  slm_int = 0.0_dp
211 
212 !$OMP PARALLEL DEFAULT(NONE), &
213 !$OMP SHARED (slm,dslm_dxyz,slm_int,max_s_harm,ll,lebedev_grid,na,harmonics,wa,indco,orbtramat) &
214 !$OMP SHARED (nso,nsoset,nco,maxs,indso,ncoset,pol,azi,llmax) &
215 !$OMP PRIVATE(ia,iso,l,m,i,lx,ly,lz,rx,ry,rz,drx,dry,drz,ic,dc,iso1,iso2,cin,dylm) &
216 !$OMP PRIVATE(is1,l1,m1,is2,l2,m2,is,n,y)
217 
218  ALLOCATE (y(na))
219  ALLOCATE (dc(nco(llmax), 3))
220 
221 !$OMP DO
222  DO iso = 1, max_s_harm
223  l = indso(1, iso)
224  m = indso(2, iso)
225  CALL y_lm(lebedev_grid(ll)%r, y, l, m)
226 
227  DO ia = 1, na
228  slm(ia, iso) = y(ia)
229  slm_int(iso) = slm_int(iso) + slm(ia, iso)*wa(ia)
230  END DO ! ia
231  END DO ! iso
232 !$OMP END DO
233 
234 !$OMP DO
235  DO ia = 1, na
236  harmonics%a(:, ia) = lebedev_grid(ll)%r(:, ia)
237  END DO
238 !$OMP END DO
239 
240  !
241  ! The derivatives dslm_dxyz and its expansions my_CG_dxyz and my_CG_dxyz_asymm
242  ! are NOT the dSlm/dx but the scaled by r**(l-1) derivatives of the monomial
243  ! terms x^n1 y^n2 z^n3 transformed by spherical harmonics expansion coefficients
244  !
245 
246 !$OMP DO
247  DO ia = 1, na
248  DO l = 0, indso(1, max_s_harm)
249  DO ic = 1, nco(l)
250  lx = indco(1, ic + ncoset(l - 1))
251  ly = indco(2, ic + ncoset(l - 1))
252  lz = indco(3, ic + ncoset(l - 1))
253 
254  IF (lx == 0) THEN
255  rx = 1.0_dp
256  drx = 0.0_dp
257  ELSE IF (lx == 1) THEN
258  rx = lebedev_grid(ll)%r(1, ia)
259  drx = 1.0_dp
260  ELSE
261  rx = lebedev_grid(ll)%r(1, ia)**lx
262  drx = real(lx, dp)*lebedev_grid(ll)%r(1, ia)**(lx - 1)
263  END IF
264  IF (ly == 0) THEN
265  ry = 1.0_dp
266  dry = 0.0_dp
267  ELSE IF (ly == 1) THEN
268  ry = lebedev_grid(ll)%r(2, ia)
269  dry = 1.0_dp
270  ELSE
271  ry = lebedev_grid(ll)%r(2, ia)**ly
272  dry = real(ly, dp)*lebedev_grid(ll)%r(2, ia)**(ly - 1)
273  END IF
274  IF (lz == 0) THEN
275  rz = 1.0_dp
276  drz = 0.0_dp
277  ELSE IF (lz == 1) THEN
278  rz = lebedev_grid(ll)%r(3, ia)
279  drz = 1.0_dp
280  ELSE
281  rz = lebedev_grid(ll)%r(3, ia)**lz
282  drz = real(lz, dp)*lebedev_grid(ll)%r(3, ia)**(lz - 1)
283  END IF
284  dc(ic, 1) = drx*ry*rz
285  dc(ic, 2) = rx*dry*rz
286  dc(ic, 3) = rx*ry*drz
287  END DO
288  n = nsoset(l - 1)
289  DO is = 1, nso(l)
290  iso = is + n
291  DO ic = 1, nco(l)
292  dslm_dxyz(:, ia, iso) = dslm_dxyz(:, ia, iso) + &
293  orbtramat(l)%slm(is, ic)*dc(ic, :)
294  END DO
295  END DO
296  END DO ! l
297  END DO !ia
298 !$OMP END DO
299 
300  ! Expansion coefficients of the cartesian derivatives
301  ! of the product of two harmonics :
302  ! d(Y(l1m1) * Y(l2m2))/dx ; d(Y(l1m1) * Y(l2m2))/dy ; d(Y(l1m1) * Y(l2m2))/dz
303 
304 !$OMP DO COLLAPSE(3)
305  DO iso1 = 1, maxs
306  DO iso2 = 1, maxs
307 
308  DO iso = 1, max_s_harm
309  rx = 0.0_dp
310  ry = 0.0_dp
311  rz = 0.0_dp
312 
313  DO ia = 1, na
314  rx = rx + wa(ia)*slm(ia, iso)* &
315  (dslm_dxyz(1, ia, iso1)*slm(ia, iso2) + slm(ia, iso1)*dslm_dxyz(1, ia, iso2))
316  ry = ry + wa(ia)*slm(ia, iso)* &
317  (dslm_dxyz(2, ia, iso1)*slm(ia, iso2) + slm(ia, iso1)*dslm_dxyz(2, ia, iso2))
318  rz = rz + wa(ia)*slm(ia, iso)* &
319  (dslm_dxyz(3, ia, iso1)*slm(ia, iso2) + slm(ia, iso1)*dslm_dxyz(3, ia, iso2))
320  END DO
321 
322  harmonics%my_CG_dxyz(1, iso1, iso2, iso) = rx
323 
324  harmonics%my_CG_dxyz(2, iso1, iso2, iso) = ry
325 
326  harmonics%my_CG_dxyz(3, iso1, iso2, iso) = rz
327 
328  END DO
329  END DO
330  END DO
331 !$OMP END DO
332 
333  ! Expansion coefficients of the cartesian of the combinations
334  ! Y(l1m1) * d(Y(l2m2))/dx - d(Y(l1m1))/dx * Y(l2m2)
335  ! Y(l1m1) * d(Y(l2m2))/dy - d(Y(l1m1))/dy * Y(l2m2)
336  ! Y(l1m1) * d(Y(l2m2))/dz - d(Y(l1m1))/dz * Y(l2m2)
337 
338 !$OMP DO COLLAPSE(3)
339  DO iso1 = 1, maxs
340  DO iso2 = 1, maxs
341  DO iso = 1, max_s_harm
342  drx = 0.0_dp
343  dry = 0.0_dp
344  drz = 0.0_dp
345 
346  DO ia = 1, na
347  drx = drx + wa(ia)*slm(ia, iso)* &
348  (-dslm_dxyz(1, ia, iso1)*slm(ia, iso2) + &
349  slm(ia, iso1)*dslm_dxyz(1, ia, iso2))
350  dry = dry + wa(ia)*slm(ia, iso)* &
351  (-dslm_dxyz(2, ia, iso1)*slm(ia, iso2) + &
352  slm(ia, iso1)*dslm_dxyz(2, ia, iso2))
353  drz = drz + wa(ia)*slm(ia, iso)* &
354  (-dslm_dxyz(3, ia, iso1)*slm(ia, iso2) + &
355  slm(ia, iso1)*dslm_dxyz(3, ia, iso2))
356  END DO
357 
358  harmonics%my_CG_dxyz_asym(1, iso1, iso2, iso) = drx
359 
360  harmonics%my_CG_dxyz_asym(2, iso1, iso2, iso) = dry
361 
362  harmonics%my_CG_dxyz_asym(3, iso1, iso2, iso) = drz
363 
364  END DO ! iso
365  END DO ! iso2
366  END DO ! iso1
367 !$OMP END DO
368 
369  ! Calculate the derivatives of the harmonics with respect of the 2 angles
370  ! the first angle (polar) is acos(lebedev_grid(ll)%r(3))
371  ! the second angle (azimutal) is atan(lebedev_grid(ll)%r(2)/lebedev_grid(ll)%r(1))
372 !$OMP DO
373  DO iso = 1, maxs
374  l = indso(1, iso)
375  m = indso(2, iso)
376  DO ia = 1, na
377  cin(1) = pol(ia)
378  cin(2) = azi(ia)
379  CALL dy_lm(cin, dylm, l, m)
380  harmonics%dslm(:, ia, iso) = dylm(:)
381  END DO
382  END DO
383 !$OMP END DO
384 
385  ! expansion coefficients of product of polar angle derivatives (dslm(1...)) in
386  ! spherical harmonics (used for tau functionals)
387 
388  DEALLOCATE (y, dc)
389 !$OMP END PARALLEL
390 
391  CALL timestop(handle)
392 
393  END SUBROUTINE create_harmonics_atom
394 
395 ! **************************************************************************************************
396 !> \brief ...
397 !> \param harmonics ...
398 !> \param orb_basis ...
399 !> \param llmax ...
400 !> \param max_s_harm ...
401 ! **************************************************************************************************
402  SUBROUTINE get_maxl_cg(harmonics, orb_basis, llmax, max_s_harm)
403 
404  TYPE(harmonics_atom_type), POINTER :: harmonics
405  TYPE(gto_basis_set_type), POINTER :: orb_basis
406  INTEGER, INTENT(IN) :: llmax, max_s_harm
407 
408  CHARACTER(len=*), PARAMETER :: routinen = 'get_maxl_CG'
409 
410  INTEGER :: damax_iso_not0, dmax_iso_not0, handle, &
411  is1, is2, itmp, max_iso_not0, nset
412  INTEGER, DIMENSION(:), POINTER :: lmax, lmin
413 
414  CALL timeset(routinen, handle)
415 
416  cpassert(ASSOCIATED(harmonics))
417 
418  CALL get_gto_basis_set(gto_basis_set=orb_basis, lmax=lmax, lmin=lmin, nset=nset)
419 
420  ! *** Assign indexes for the non null CG coefficients ***
421  max_iso_not0 = 0
422  dmax_iso_not0 = 0
423  damax_iso_not0 = 0
424  DO is1 = 1, nset
425  DO is2 = 1, nset
426  CALL get_none0_cg_list(harmonics%my_CG, &
427  lmin(is1), lmax(is1), lmin(is2), lmax(is2), &
428  max_s_harm, llmax, max_iso_not0=itmp)
429  max_iso_not0 = max(max_iso_not0, itmp)
430  CALL get_none0_cg_list(harmonics%my_CG_dxyz, &
431  lmin(is1), lmax(is1), lmin(is2), lmax(is2), &
432  max_s_harm, llmax, max_iso_not0=itmp)
433  dmax_iso_not0 = max(dmax_iso_not0, itmp)
434  CALL get_none0_cg_list(harmonics%my_CG_dxyz_asym, &
435  lmin(is1), lmax(is1), lmin(is2), lmax(is2), &
436  max_s_harm, llmax, max_iso_not0=itmp)
437  damax_iso_not0 = max(damax_iso_not0, itmp)
438  END DO ! is2
439  END DO ! is1
440  harmonics%max_iso_not0 = max_iso_not0
441  harmonics%dmax_iso_not0 = dmax_iso_not0
442  harmonics%damax_iso_not0 = damax_iso_not0
443 
444  CALL timestop(handle)
445 
446  END SUBROUTINE get_maxl_cg
447 
448 ! **************************************************************************************************
449 !> \brief ...
450 !> \param cgc ...
451 !> \param lmin1 ...
452 !> \param lmax1 ...
453 !> \param lmin2 ...
454 !> \param lmax2 ...
455 !> \param max_s_harm ...
456 !> \param llmax ...
457 !> \param list ...
458 !> \param n_list ...
459 !> \param max_iso_not0 ...
460 ! **************************************************************************************************
461  SUBROUTINE get_none0_cg_list4(cgc, lmin1, lmax1, lmin2, lmax2, max_s_harm, llmax, &
462  list, n_list, max_iso_not0)
463 
464  REAL(dp), DIMENSION(:, :, :, :), INTENT(IN) :: cgc
465  INTEGER, INTENT(IN) :: lmin1, lmax1, lmin2, lmax2, max_s_harm, &
466  llmax
467  INTEGER, DIMENSION(:, :, :), INTENT(OUT), OPTIONAL :: list
468  INTEGER, DIMENSION(:), INTENT(OUT), OPTIONAL :: n_list
469  INTEGER, INTENT(OUT) :: max_iso_not0
470 
471  INTEGER :: iso, iso1, iso2, l1, l2, nlist
472 
473  cpassert(nsoset(lmax1) .LE. SIZE(cgc, 2))
474  cpassert(nsoset(lmax2) .LE. SIZE(cgc, 3))
475  cpassert(max_s_harm .LE. SIZE(cgc, 4))
476  IF (PRESENT(n_list) .AND. PRESENT(list)) THEN
477  cpassert(max_s_harm .LE. SIZE(list, 3))
478  END IF
479  max_iso_not0 = 0
480  IF (PRESENT(n_list) .AND. PRESENT(list)) n_list = 0
481  DO iso = 1, max_s_harm
482  nlist = 0
483  DO l1 = lmin1, lmax1
484  DO iso1 = nsoset(l1 - 1) + 1, nsoset(l1)
485  DO l2 = lmin2, lmax2
486  IF (l1 + l2 > llmax) cycle
487  DO iso2 = nsoset(l2 - 1) + 1, nsoset(l2)
488  IF (abs(cgc(1, iso1, iso2, iso)) + &
489  abs(cgc(2, iso1, iso2, iso)) + &
490  abs(cgc(3, iso1, iso2, iso)) > 1.e-8_dp) THEN
491  nlist = nlist + 1
492  IF (PRESENT(n_list) .AND. PRESENT(list)) THEN
493  list(1, nlist, iso) = iso1
494  list(2, nlist, iso) = iso2
495  END IF
496  max_iso_not0 = max(max_iso_not0, iso)
497  END IF
498  END DO
499  END DO
500  END DO
501  END DO
502  IF (PRESENT(n_list) .AND. PRESENT(list)) n_list(iso) = nlist
503  END DO
504  END SUBROUTINE get_none0_cg_list4
505 
506 ! **************************************************************************************************
507 !> \brief ...
508 !> \param cgc ...
509 !> \param lmin1 ...
510 !> \param lmax1 ...
511 !> \param lmin2 ...
512 !> \param lmax2 ...
513 !> \param max_s_harm ...
514 !> \param llmax ...
515 !> \param list ...
516 !> \param n_list ...
517 !> \param max_iso_not0 ...
518 ! **************************************************************************************************
519  SUBROUTINE get_none0_cg_list3(cgc, lmin1, lmax1, lmin2, lmax2, max_s_harm, llmax, &
520  list, n_list, max_iso_not0)
521 
522  REAL(dp), DIMENSION(:, :, :), INTENT(IN) :: cgc
523  INTEGER, INTENT(IN) :: lmin1, lmax1, lmin2, lmax2, max_s_harm, &
524  llmax
525  INTEGER, DIMENSION(:, :, :), INTENT(OUT), OPTIONAL :: list
526  INTEGER, DIMENSION(:), INTENT(OUT), OPTIONAL :: n_list
527  INTEGER, INTENT(OUT) :: max_iso_not0
528 
529  INTEGER :: iso, iso1, iso2, l1, l2, nlist
530 
531  cpassert(nsoset(lmax1) .LE. SIZE(cgc, 1))
532  cpassert(nsoset(lmax2) .LE. SIZE(cgc, 2))
533  cpassert(max_s_harm .LE. SIZE(cgc, 3))
534  IF (PRESENT(n_list) .AND. PRESENT(list)) THEN
535  cpassert(max_s_harm .LE. SIZE(list, 3))
536  END IF
537  max_iso_not0 = 0
538  IF (PRESENT(n_list) .AND. PRESENT(list)) n_list = 0
539  DO iso = 1, max_s_harm
540  nlist = 0
541  DO l1 = lmin1, lmax1
542  DO iso1 = nsoset(l1 - 1) + 1, nsoset(l1)
543  DO l2 = lmin2, lmax2
544  IF (l1 + l2 > llmax) cycle
545  DO iso2 = nsoset(l2 - 1) + 1, nsoset(l2)
546  IF (abs(cgc(iso1, iso2, iso)) > 1.e-8_dp) THEN
547  nlist = nlist + 1
548  IF (PRESENT(n_list) .AND. PRESENT(list)) THEN
549  list(1, nlist, iso) = iso1
550  list(2, nlist, iso) = iso2
551  END IF
552  max_iso_not0 = max(max_iso_not0, iso)
553  END IF
554  END DO
555  END DO
556  END DO
557  END DO
558  IF (PRESENT(n_list) .AND. PRESENT(list)) n_list(iso) = nlist
559  END DO
560  END SUBROUTINE get_none0_cg_list3
561 
562 END MODULE qs_harmonics_atom
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 the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
Generation of the spherical Lebedev grids. All Lebedev grids were generated with a precision of at le...
Definition: lebedev.F:57
type(oh_grid), dimension(nlg), target, public lebedev_grid
Definition: lebedev.F:85
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
Definition: list.F:24
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 allocate_harmonics_atom(harmonics)
Allocate a spherical harmonics set for the atom grid.
subroutine, public get_maxl_cg(harmonics, orb_basis, llmax, max_s_harm)
...
subroutine, public deallocate_harmonics_atom(harmonics)
Deallocate the spherical harmonics set for the atom grid.
subroutine, public create_harmonics_atom(harmonics, my_CG, na, llmax, maxs, max_s_harm, ll, wa, azi, pol)
...
Calculate spherical harmonics.