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