24 #include "./base/base_uses.f90"
30 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_harmonics_atom'
32 TYPE harmonics_atom_type
33 INTEGER :: max_s_harm, llmax, &
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
45 END TYPE harmonics_atom_type
54 INTERFACE get_none0_cg_list
55 MODULE PROCEDURE get_none0_cg_list3, get_none0_cg_list4
67 TYPE(harmonics_atom_type),
POINTER :: harmonics
73 harmonics%max_s_harm = 0
75 harmonics%max_iso_not0 = 0
76 harmonics%dmax_iso_not0 = 0
77 harmonics%damax_iso_not0 = 0
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)
98 TYPE(harmonics_atom_type),
POINTER :: harmonics
100 IF (
ASSOCIATED(harmonics))
THEN
102 IF (
ASSOCIATED(harmonics%slm))
THEN
103 DEALLOCATE (harmonics%slm)
106 IF (
ASSOCIATED(harmonics%dslm))
THEN
107 DEALLOCATE (harmonics%dslm)
110 IF (
ASSOCIATED(harmonics%dslm_dxyz))
THEN
111 DEALLOCATE (harmonics%dslm_dxyz)
114 IF (
ASSOCIATED(harmonics%slm_int))
THEN
115 DEALLOCATE (harmonics%slm_int)
118 IF (
ASSOCIATED(harmonics%my_CG))
THEN
119 DEALLOCATE (harmonics%my_CG)
122 IF (
ASSOCIATED(harmonics%my_CG_dxyz))
THEN
123 DEALLOCATE (harmonics%my_CG_dxyz)
126 IF (
ASSOCIATED(harmonics%my_CG_dxyz_asym))
THEN
127 DEALLOCATE (harmonics%my_CG_dxyz_asym)
130 IF (
ASSOCIATED(harmonics%a))
THEN
131 DEALLOCATE (harmonics%a)
134 DEALLOCATE (harmonics)
136 CALL cp_abort(__location__, &
137 "The pointer harmonics is not associated and "// &
138 "cannot be deallocated")
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
164 CHARACTER(len=*),
PARAMETER :: routinen =
'create_harmonics_atom'
166 INTEGER :: handle, i, ia, ic, is, is1, is2, iso, &
167 iso1, iso2, l, l1, l2, lx, ly, lz, m, &
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
175 CALL timeset(routinen, handle)
177 NULLIFY (y, slm, dslm_dxyz, dc)
179 cpassert(
ASSOCIATED(harmonics))
181 harmonics%max_s_harm = max_s_harm
182 harmonics%llmax = llmax
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)
192 harmonics%my_CG(1:maxs, is1, i) = my_cg(1:maxs, is1, i)
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)
205 NULLIFY (slm, dslm_dxyz, slm_int)
207 dslm_dxyz => harmonics%dslm_dxyz
209 slm_int => harmonics%slm_int
219 ALLOCATE (dc(
nco(llmax), 3))
222 DO iso = 1, max_s_harm
229 slm_int(iso) = slm_int(iso) + slm(ia, iso)*wa(ia)
248 DO l = 0,
indso(1, max_s_harm)
257 ELSE IF (lx == 1)
THEN
267 ELSE IF (ly == 1)
THEN
277 ELSE IF (lz == 1)
THEN
284 dc(ic, 1) = drx*ry*rz
285 dc(ic, 2) = rx*dry*rz
286 dc(ic, 3) = rx*ry*drz
292 dslm_dxyz(:, ia, iso) = dslm_dxyz(:, ia, iso) + &
308 DO iso = 1, max_s_harm
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))
322 harmonics%my_CG_dxyz(1, iso1, iso2, iso) = rx
324 harmonics%my_CG_dxyz(2, iso1, iso2, iso) = ry
326 harmonics%my_CG_dxyz(3, iso1, iso2, iso) = rz
341 DO iso = 1, max_s_harm
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))
358 harmonics%my_CG_dxyz_asym(1, iso1, iso2, iso) = drx
360 harmonics%my_CG_dxyz_asym(2, iso1, iso2, iso) = dry
362 harmonics%my_CG_dxyz_asym(3, iso1, iso2, iso) = drz
379 CALL dy_lm(cin, dylm, l, m)
380 harmonics%dslm(:, ia, iso) = dylm(:)
391 CALL timestop(handle)
404 TYPE(harmonics_atom_type),
POINTER :: harmonics
405 TYPE(gto_basis_set_type),
POINTER :: orb_basis
406 INTEGER,
INTENT(IN) :: llmax, max_s_harm
408 CHARACTER(len=*),
PARAMETER :: routinen =
'get_maxl_CG'
410 INTEGER :: damax_iso_not0, dmax_iso_not0, handle, &
411 is1, is2, itmp, max_iso_not0, nset
412 INTEGER,
DIMENSION(:),
POINTER :: lmax, lmin
414 CALL timeset(routinen, handle)
416 cpassert(
ASSOCIATED(harmonics))
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)
440 harmonics%max_iso_not0 = max_iso_not0
441 harmonics%dmax_iso_not0 = dmax_iso_not0
442 harmonics%damax_iso_not0 = damax_iso_not0
444 CALL timestop(handle)
461 SUBROUTINE get_none0_cg_list4(cgc, lmin1, lmax1, lmin2, lmax2, max_s_harm, llmax, &
462 list, n_list, max_iso_not0)
464 REAL(
dp),
DIMENSION(:, :, :, :),
INTENT(IN) :: cgc
465 INTEGER,
INTENT(IN) :: lmin1, lmax1, lmin2, lmax2, max_s_harm, &
467 INTEGER,
DIMENSION(:, :, :),
INTENT(OUT),
OPTIONAL ::
list
468 INTEGER,
DIMENSION(:),
INTENT(OUT),
OPTIONAL :: n_list
469 INTEGER,
INTENT(OUT) :: max_iso_not0
471 INTEGER :: iso, iso1, iso2, l1, l2, nlist
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))
480 IF (
PRESENT(n_list) .AND.
PRESENT(
list)) n_list = 0
481 DO iso = 1, max_s_harm
486 IF (l1 + l2 > llmax) cycle
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
492 IF (
PRESENT(n_list) .AND.
PRESENT(
list))
THEN
493 list(1, nlist, iso) = iso1
494 list(2, nlist, iso) = iso2
496 max_iso_not0 = max(max_iso_not0, iso)
502 IF (
PRESENT(n_list) .AND.
PRESENT(
list)) n_list(iso) = nlist
504 END SUBROUTINE get_none0_cg_list4
519 SUBROUTINE get_none0_cg_list3(cgc, lmin1, lmax1, lmin2, lmax2, max_s_harm, llmax, &
520 list, n_list, max_iso_not0)
522 REAL(
dp),
DIMENSION(:, :, :),
INTENT(IN) :: cgc
523 INTEGER,
INTENT(IN) :: lmin1, lmax1, lmin2, lmax2, max_s_harm, &
525 INTEGER,
DIMENSION(:, :, :),
INTENT(OUT),
OPTIONAL ::
list
526 INTEGER,
DIMENSION(:),
INTENT(OUT),
OPTIONAL :: n_list
527 INTEGER,
INTENT(OUT) :: max_iso_not0
529 INTEGER :: iso, iso1, iso2, l1, l2, nlist
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))
538 IF (
PRESENT(n_list) .AND.
PRESENT(
list)) n_list = 0
539 DO iso = 1, max_s_harm
544 IF (l1 + l2 > llmax) cycle
546 IF (abs(cgc(iso1, iso2, iso)) > 1.e-8_dp)
THEN
548 IF (
PRESENT(n_list) .AND.
PRESENT(
list))
THEN
549 list(1, nlist, iso) = iso1
550 list(2, nlist, iso) = iso2
552 max_iso_not0 = max(max_iso_not0, iso)
558 IF (
PRESENT(n_list) .AND.
PRESENT(
list)) n_list(iso) = nlist
560 END SUBROUTINE get_none0_cg_list3
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.
integer, parameter, public dp
Generation of the spherical Lebedev grids. All Lebedev grids were generated with a precision of at le...
type(oh_grid), dimension(nlg), target, public lebedev_grid
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
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
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.