61#include "./base/base_uses.f90"
67 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'se_core_matrix'
83 LOGICAL,
INTENT(IN) :: calculate_forces
85 CHARACTER(LEN=*),
PARAMETER :: routinen =
'build_se_core_matrix'
87 INTEGER :: after, atom_a, atom_b, handle, i, iatom, icol, icor, ikind, inode, irow, itype, &
88 iw, j, jatom, jkind, natom, natorb_a, nkind, nr_a, nra, nrb
89 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: atom_of_kind, nrt
90 LOGICAL :: defined, found, omit_headers, use_virial
91 LOGICAL,
ALLOCATABLE,
DIMENSION(:) :: se_defined
92 REAL(kind=
dp) :: delta, dr, econst, eheat, eisol, kh, &
93 udd, uff, upp, uss, zpa, zpb, zsa, zsb
94 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: zpt, zst
95 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: hmt, umt
96 REAL(kind=
dp),
DIMENSION(16) :: ha, hb, ua
97 REAL(kind=
dp),
DIMENSION(3) :: force_ab, rij
98 REAL(kind=
dp),
DIMENSION(:),
POINTER :: beta_a, sto_exponents_a
99 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: dsmat, h_block, h_blocka, pabmat, pamat, &
103 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_h, matrix_p, matrix_s
104 TYPE(
dbcsr_type),
POINTER :: diagmat_h, diagmat_p
107 DIMENSION(:),
POINTER :: nl_iterator
113 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
121 CALL timeset(routinen, handle)
123 NULLIFY (logger, energy)
126 NULLIFY (rho, force, atomic_kind_set, qs_kind_set, sab_orb, &
127 diagmat_h, diagmat_p, particle_set, matrix_p, ks_env)
133 particle_set=particle_set, &
134 atomic_kind_set=atomic_kind_set, &
135 qs_kind_set=qs_kind_set, &
136 dft_control=dft_control, &
144 IF (calculate_forces)
THEN
146 matrix_name=
"OVERLAP", &
147 basis_type_a=
"ORB", &
148 basis_type_b=
"ORB", &
151 use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
154 matrix_name=
"OVERLAP", &
155 basis_type_a=
"ORB", &
156 basis_type_b=
"ORB", &
162 IF (calculate_forces)
THEN
165 IF (
SIZE(matrix_p) == 2)
THEN
166 CALL dbcsr_add(matrix_p(1)%matrix, matrix_p(2)%matrix, alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
168 delta = dft_control%qs_control%se_control%delta
170 atom_of_kind=atom_of_kind)
178 ALLOCATE (matrix_h(1)%matrix)
179 CALL dbcsr_copy(matrix_h(1)%matrix, matrix_s(1)%matrix,
"CORE HAMILTONIAN MATRIX")
180 CALL dbcsr_set(matrix_h(1)%matrix, 0.0_dp)
189 itype =
get_se_type(dft_control%qs_control%method_id)
192 nkind =
SIZE(atomic_kind_set)
194 ALLOCATE (se_defined(nkind))
195 ALLOCATE (hmt(16, nkind))
196 ALLOCATE (umt(16, nkind))
198 ALLOCATE (zst(nkind))
199 ALLOCATE (zpt(nkind))
200 ALLOCATE (nrt(nkind))
205 CALL get_qs_kind(qs_kind_set(ikind), se_parameter=se_kind_a)
206 CALL get_se_param(se_kind_a, defined=defined, natorb=natorb_a, &
207 beta=beta_a, uss=uss, upp=upp, udd=udd, uff=uff, eisol=eisol, eheat=eheat, &
208 nr=nr_a, sto_exponents=sto_exponents_a)
209 econst = econst - (eisol - eheat)*real(natom,
dp)
210 se_defined(ikind) = (defined .AND. natorb_a >= 1)
211 hmt(1, ikind) = beta_a(0)
212 hmt(2:4, ikind) = beta_a(1)
213 hmt(5:9, ikind) = beta_a(2)
214 hmt(10:16, ikind) = beta_a(3)
216 umt(2:4, ikind) = upp
217 umt(5:9, ikind) = udd
218 umt(10:16, ikind) = uff
220 zst(ikind) = sto_exponents_a(0)
221 zpt(ikind) = sto_exponents_a(1)
225 energy%core_self = econst
229 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, iatom=iatom, jatom=jatom, inode=inode, r=rij)
230 IF (.NOT. se_defined(ikind)) cycle
231 IF (.NOT. se_defined(jkind)) cycle
232 ha(1:16) = hmt(1:16, ikind)
233 ua(1:16) = umt(1:16, ikind)
234 hb(1:16) = hmt(1:16, jkind)
244 SELECT CASE (dft_control%qs_control%method_id)
249 cpassert(
ASSOCIATED(h_blocka))
250 IF (calculate_forces)
THEN
252 cpassert(
ASSOCIATED(pamat))
259 SELECT CASE (dft_control%qs_control%method_id)
264 DO i = 1,
SIZE(h_blocka, 1)
265 h_blocka(i, i) = h_blocka(i, i) + ua(i)
270 IF (iatom <= jatom)
THEN
279 irow, icol, h_block, found)
280 cpassert(
ASSOCIATED(h_block))
303 irow, icol, s_block, found)
304 cpassert(
ASSOCIATED(s_block))
305 IF (irow == iatom)
THEN
306 DO i = 1,
SIZE(h_block, 1)
307 DO j = 1,
SIZE(h_block, 2)
308 h_block(i, j) = h_block(i, j) + kh*(ha(i) + hb(j))*s_block(i, j)
312 DO i = 1,
SIZE(h_block, 1)
313 DO j = 1,
SIZE(h_block, 2)
314 h_block(i, j) = h_block(i, j) + kh*(ha(j) + hb(i))*s_block(i, j)
318 IF (calculate_forces)
THEN
319 atom_a = atom_of_kind(iatom)
320 atom_b = atom_of_kind(jatom)
331 cpassert(
ASSOCIATED(pabmat))
333 force_ab(icor) = 0._dp
345 cpassert(
ASSOCIATED(dsmat))
346 dsmat = 2._dp*kh*dsmat*pabmat
347 IF (irow == iatom)
THEN
348 DO i = 1,
SIZE(h_block, 1)
349 DO j = 1,
SIZE(h_block, 2)
350 force_ab(icor) = force_ab(icor) + (ha(i) + hb(j))*dsmat(i, j)
354 DO i = 1,
SIZE(h_block, 1)
355 DO j = 1,
SIZE(h_block, 2)
356 force_ab(icor) = force_ab(icor) + (ha(j) + hb(i))*dsmat(i, j)
365 IF (calculate_forces .AND. (iatom /= jatom .OR. dr >
rij_threshold))
THEN
366 IF (irow == iatom) force_ab = -force_ab
367 force(ikind)%all_potential(:, atom_a) = &
368 force(ikind)%all_potential(:, atom_a) - force_ab(:)
369 force(jkind)%all_potential(:, atom_b) = &
370 force(jkind)%all_potential(:, atom_b) + force_ab(:)
379 DEALLOCATE (se_defined, hmt, umt, zst, zpt, nrt)
383 CALL dbcsr_add(matrix_h(1)%matrix, diagmat_h, 1.0_dp, 1.0_dp)
387 qs_env%input,
"DFT%PRINT%AO_MATRICES/CORE_HAMILTONIAN"),
cp_p_file))
THEN
391 CALL section_vals_val_get(qs_env%input,
"DFT%PRINT%AO_MATRICES%OMIT_HEADERS", l_val=omit_headers)
392 after = min(max(after, 1), 16)
394 scale=
evolt, output_unit=iw, omit_headers=omit_headers)
396 "DFT%PRINT%AO_MATRICES/CORE_HAMILTONIAN")
399 IF (calculate_forces)
THEN
400 IF (
SIZE(matrix_p) == 2)
THEN
401 CALL dbcsr_add(matrix_p(1)%matrix, matrix_p(2)%matrix, alpha_scalar=1.0_dp, beta_scalar=-1.0_dp)
403 DEALLOCATE (atom_of_kind)
409 CALL timestop(handle)
424 SUBROUTINE makes(R, nra, nrb, ZSA, ZSB, ZPA, ZPB, S)
426 REAL(kind=
dp),
DIMENSION(3) :: r
428 REAL(kind=
dp) :: zsa, zsb, zpa, zpb
429 REAL(kind=
dp),
DIMENSION(4, 4) :: s
431 INTEGER,
DIMENSION(4, 4),
PARAMETER :: &
432 nc1 = reshape((/2, 4, 4, 6, 4, 3, 6, 7, 4, 6, 4, 8, 6, 7, 8, 5/), (/4, 4/)), &
433 nc2 = reshape((/4, 4, 8, 8, 6, 8, 6, 12, 8, 8, 12, 8, 10, 12, 14, 16/), (/4, 4/)), &
434 nc3 = reshape((/4, 6, 8, 10, 4, 8, 8, 12, 8, 6, 12, 14, 8, 12, 8, 16/), (/4, 4/)), &
435 nc4 = reshape((/4, 8, 11, 14, 8, 6, 12, 14, 11, 12, 10, 20, 14, 14, 20, 12/), (/4, 4/)), &
436 nc5 = reshape((/2, 4, 6, 8, 4, 4, 8, 8, 6, 8, 6, 12, 8, 8, 12, 8/), (/4, 4/))
437 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: c1 = reshape((/1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1 &
438 , 1, 1, 1, 1, -1, 1, 2, 3, -1, -2, 1, 2, -2, -1, -3, 1, -3, -2, -1, -4, 0, -1, -2, 2, -1, &
439 1, -2, -1, 2, -2, 3, -3, 2, -1, -3, 6, 0, -1, -1, -2, 1, 0, -2, -4, -1, 2, -1, -3, 2, 4, 3&
440 , -4, 0, 0, 0, -3, 0, 0, 1, -1, 0, 1, 0, 3, -3, -1, 3, 1, 0, 0, 0, -1, 0, 0, 1, 2, 0, -1, &
441 0, 3, 1, -2, -3, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, -1, 0, 1, -1, 0, 0, 0, 0, 0, 0, 0, 0 &
442 , 0, 0, 0, 0, -1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
443 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
444 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
445 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
446 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
447 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
448 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
449 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: c2 = reshape((/1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1 &
450 , 1, 1, 1, 1, -1, 1, 1, 2, -2, -1, 1, 1, -3, -2, -1, 1, -4, -3, -2, -1, 1, -1, 1, 1, 1, 1 &
451 , -2, 1, 1, 1, 1, -3, 1, 1, 1, 1, -1, -1, -1, 2, 1, -1, -2, -2, 3, -2, -2, -3, 6, 2, -1, &
452 -3, 0, 0, 1, -2, -2, -1, 1, 1, -3, 2, -1, 3, -4, -3, -2, -1, 0, 0, -1, -1, 1, 1, 1, -2, -1&
453 , -1, 2, 3, -4, 2, 4, 3, 0, 0, -1, -2, 0, -1, 0, -2, 3, 2, -2, -1, 6, 2, -1, -3, 0, 0, -1,&
454 -1, 0, 1, 0, 1, -1, -1, 1, -1, 1, -3, -1, 3, 0, 0, 0, 0, 0, 0, 0, -2, 0, 0, 2, 0, -4, 2, 4&
455 , 3, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, -1, 0, 1, 1, -2, -3, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0&
456 , 0, -3, -1, 3, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, -1, 0, 0, 1, 1, -1, 0, 0, 0, 0, 0, 0, 0, 0, &
457 0, 0, 0, 0, 0, 0, -2, -3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0&
458 , 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, &
459 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
460 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
462 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: c3 = reshape((/-1, -1, -1, -1, -1, -1, -1, -1, -1 &
463 , -1, -1, -1, -1, -1, -1, -1, -1, -2, -3, -4, 1, -1, -2, -3, 1, 1, -1, -2, 2, 1, 1, -1, 1 &
464 , 1, 1, 1, 1, 1, 1, 1, 1, 2, 1, 1, 1, 1, 3, 1, 1, -1, -3, -6, -1, 1, 2, -2, 1, -2, 2, 1, &
465 -2, 2, -3, 3, 0, 2, 3, 4, 0, 1, 2, 3, -1, -1, 1, 2, -2, -1, -3, 1, 0, 1, -1, -4, 0, 1, 1, &
466 2, -1, 1, 2, 4, 1, -2, 3, 3, 0, 0, 3, 6, 0, -1, -2, 2, -1, 0, -2, -1, 2, -2, 1, -3, 0, 0, &
467 1, -1, 0, -1, -1, 3, 1, 0, -1, 1, -1, -1, -1, -3, 0, 0, 0, 4, 0, 0, 0, -2, 0, 0, -2, -4, 0&
468 , 2, 0, -3, 0, 0, 0, 1, 0, 0, 0, 1, 0, 0, -1, -2, 0, 1, 0, -3, 0, 0, 0, 0, 0, 0, 0, -3, 0 &
469 , 0, 1, -1, 0, 1, 0, 3, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 1, -1, 0, -1, 0, 1, 0, 0, 0, 0, 0, &
470 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, &
471 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0&
472 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
473 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
474 , 0, 0, 0/), (/4, 4, 20/))
475 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: c4 = reshape((/-1, -1, -1, -1, -1, -1, -1, -1, -1 &
476 , -1, -1, -1, -1, -1, -1, -1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, -1, -2, &
477 -3, 1, 1, -1, -2, 2, 1, 2, -1, 3, 2, 1, 3, -1, 1, 2, 3, -1, -1, 1, 2, -2, -1, -1, 1, -3, &
478 -2, -1, -2, 0, 1, -1, -3, 1, -1, 1, 1, -1, 1, -1, 2, -3, 1, 2, -1, 0, -1, 2, 4, -1, 1, -1 &
479 , -1, 2, -1, -1, -1, 4, -1, -1, -3, 0, 1, -1, -1, -1, 0, 1, 2, -1, -1, -1, -1, -1, -2, -1 &
480 , 3, 0, -1, 2, -1, 1, 0, -1, -2, -2, 1, 2, 2, 1, 2, -2, 1, 0, 0, -2, 4, 0, 0, -1, 1, 2, -1&
481 , 1, -1, -4, 1, 1, 2, 0, 0, 1, -3, 0, 0, 1, -1, 1, 1, -1, -1, 3, -1, 1, -3, 0, 0, -1, 3, 0&
482 , 0, -1, -2, -1, 1, 0, -1, 3, 2, -1, -1, 0, 0, 0, -3, 0, 0, 1, 2, 0, -1, 0, -1, -3, -2, -1&
483 , 1, 0, 0, 0, 1, 0, 0, 0, -1, 0, 0, 0, 2, -1, -1, 2, 0, 0, 0, 0, -1, 0, 0, 0, 1, 0, 0, 0, &
484 -1, 1, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
485 , 0, 0, 2, 0, 0, -2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, &
486 0, 0, 0, 0, 0, -1, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, -1, 0, 0, 0, 0, &
487 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 1, 0/), (/4, 4, 20/))
488 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: c5 = reshape((/-1, -1, -1, -1, -1, -1, -1, -1, -1 &
489 , -1, -1, -1, -1, -1, -1, -1, 1, -1, -2, -3, 1, 1, -1, -2, 2, 1, 2, -1, 3, 2, 1, 3, 0, 1, &
490 -1, -3, 1, 1, 1, 1, -1, 1, 1, 2, -3, 1, 2, 1, 0, 1, 1, 1, -1, -1, 1, 2, 1, 1, -1, 1, 1, -2&
491 , 1, -3, 0, 0, 2, -1, 0, 0, 1, 2, -2, -1, -2, 2, 1, -2, -2, -3, 0, 0, 1, 3, 0, 0, 1, 1, 1,&
492 -1, 1, 1, -3, 1, -1, 1, 0, 0, 0, 3, 0, 0, -1, -2, 0, -1, 0, -1, 3, 2, -1, 3, 0, 0, 0, 1, 0&
493 , 0, -1, -1, 0, 1, 0, -2, -1, -1, -2, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 1, 0,&
494 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -2, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0,&
495 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 &
496 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
497 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
498 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
499 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, &
501 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: ma1 = reshape((/2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7&
502 , 5, 6, 7, 8, 0, 2, 3, 4, 2, 2, 4, 5, 3, 4, 4, 6, 4, 5, 6, 6, 0, 1, 1, 3, 1, 0, 3, 4, 1, 3&
503 , 2, 5, 3, 4, 5, 4, 0, 0, 0, 2, 0, 0, 2, 3, 0, 2, 0, 4, 2, 3, 4, 2, 0, 0, 0, 1, 0, 0, 1, 2&
504 , 0, 1, 0, 3, 1, 2, 3, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 2, 0, 1, 2, 0, 0, 0, 0, 0, 0, 0&
505 , 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
506 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
507 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
508 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
509 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
510 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
511 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
512 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
513 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: ma2 = reshape((/2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7&
514 , 5, 6, 7, 8, 1, 2, 3, 4, 2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7, 1, 1, 3, 4, 2, 3, 3, 5, 3, 4&
515 , 5, 5, 4, 5, 6, 7, 0, 0, 2, 3, 1, 2, 2, 4, 2, 3, 4, 4, 3, 4, 5, 6, 0, 0, 2, 2, 1, 2, 1, 4&
516 , 2, 2, 4, 3, 3, 4, 5, 6, 0, 0, 1, 1, 0, 1, 0, 3, 1, 1, 3, 2, 2, 3, 4, 5, 0, 0, 1, 1, 0, 1&
517 , 0, 3, 1, 1, 3, 1, 2, 3, 4, 5, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 2, 0, 1, 2, 3, 4, 0, 0, 0, 0&
518 , 0, 0, 0, 2, 0, 0, 2, 0, 1, 2, 3, 4, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 2, 3, 0, 0&
519 , 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 2, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 2&
520 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
521 , 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
522 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
523 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
524 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
525 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: ma3 = reshape((/2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7&
526 , 5, 6, 7, 8, 1, 2, 3, 4, 2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7, 1, 2, 3, 4, 1, 3, 4, 5, 3, 3&
527 , 5, 6, 4, 5, 5, 7, 0, 1, 2, 3, 0, 2, 3, 4, 2, 2, 4, 5, 3, 4, 4, 6, 0, 1, 2, 3, 0, 2, 2, 4&
528 , 2, 1, 4, 5, 2, 4, 3, 6, 0, 0, 1, 2, 0, 1, 1, 3, 1, 0, 3, 4, 1, 3, 2, 5, 0, 0, 1, 2, 0, 1&
529 , 1, 3, 1, 0, 3, 4, 1, 3, 1, 5, 0, 0, 0, 1, 0, 0, 0, 2, 0, 0, 2, 3, 0, 2, 0, 4, 0, 0, 0, 1&
530 , 0, 0, 0, 2, 0, 0, 2, 3, 0, 2, 0, 4, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 2, 0, 1, 0, 3, 0, 0&
531 , 0, 0, 0, 0, 0, 1, 0, 0, 1, 2, 0, 1, 0, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 2&
532 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
533 , 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
534 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
535 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
536 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
537 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: ma4 = reshape((/2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7&
538 , 5, 6, 7, 8, 2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7, 5, 6, 7, 8, 0, 2, 3, 4, 2, 2, 4, 5, 3, 4&
539 , 4, 6, 4, 5, 6, 6, 0, 2, 3, 4, 2, 2, 4, 5, 3, 4, 4, 6, 4, 5, 6, 6, 0, 1, 2, 3, 1, 0, 3, 4&
540 , 2, 3, 4, 5, 3, 4, 5, 6, 0, 1, 2, 3, 1, 0, 3, 4, 2, 3, 2, 5, 3, 4, 5, 4, 0, 0, 2, 3, 0, 0&
541 , 2, 3, 2, 2, 2, 5, 3, 3, 5, 4, 0, 0, 1, 2, 0, 0, 2, 3, 1, 2, 2, 4, 2, 3, 4, 2, 0, 0, 1, 2&
542 , 0, 0, 1, 2, 1, 1, 0, 4, 2, 2, 4, 2, 0, 0, 0, 2, 0, 0, 1, 2, 0, 1, 0, 4, 2, 2, 4, 2, 0, 0&
543 , 0, 1, 0, 0, 0, 1, 0, 0, 0, 3, 1, 1, 3, 0, 0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 3, 1, 1, 3, 0&
544 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, 0, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0&
545 , 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2&
546 , 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
547 , 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
548 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
549 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: ma5 = reshape((/2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7&
550 , 5, 6, 7, 8, 0, 2, 3, 4, 2, 2, 4, 5, 3, 4, 4, 6, 4, 5, 6, 6, 0, 1, 2, 3, 1, 2, 3, 4, 2, 3&
551 , 4, 5, 3, 4, 5, 6, 0, 0, 2, 3, 0, 0, 3, 3, 2, 3, 2, 5, 3, 3, 5, 4, 0, 0, 1, 2, 0, 0, 2, 3&
552 , 1, 2, 2, 4, 2, 3, 4, 4, 0, 0, 0, 2, 0, 0, 2, 2, 0, 2, 0, 4, 2, 2, 4, 2, 0, 0, 0, 1, 0, 0&
553 , 1, 1, 0, 1, 0, 3, 1, 1, 3, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, 0, 3, 0, 0, 0, 0, 0&
554 , 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 2, 0, 0, 0&
555 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
556 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
557 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
558 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
559 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
560 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
561 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: mb1 = reshape((/0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
562 , 0, 0, 0, 0, 2, 1, 1, 1, 1, 2, 1, 1, 1, 1, 2, 1, 1, 1, 1, 2, 0, 2, 3, 2, 2, 4, 2, 2, 3, 2&
563 , 4, 2, 2, 2, 2, 4, 0, 3, 4, 3, 3, 0, 3, 3, 4, 3, 6, 3, 3, 3, 3, 6, 0, 0, 0, 4, 0, 0, 4, 4&
564 , 0, 4, 0, 4, 4, 4, 4, 8, 0, 0, 0, 5, 0, 0, 5, 5, 0, 5, 0, 5, 5, 5, 5, 0, 0, 0, 0, 0, 0, 0&
565 , 0, 6, 0, 0, 0, 6, 0, 6, 6, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 7, 0, 0, 7, 0, 0, 0, 0, 0&
566 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
567 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
568 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
569 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
570 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
571 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
572 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
573 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: mb2 = reshape((/1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1&
574 , 1, 1, 1, 1, 2, 0, 2, 2, 2, 2, 0, 2, 2, 2, 2, 0, 2, 2, 2, 2, 0, 3, 0, 0, 0, 0, 3, 0, 0, 0&
575 , 0, 3, 0, 0, 0, 0, 1, 2, 3, 1, 3, 3, 2, 3, 3, 1, 3, 2, 3, 3, 3, 3, 0, 0, 1, 4, 1, 1, 5, 1&
576 , 1, 4, 1, 5, 1, 1, 1, 1, 0, 0, 4, 5, 2, 4, 4, 4, 4, 5, 4, 4, 4, 4, 4, 4, 0, 0, 2, 3, 0, 2&
577 , 0, 2, 2, 3, 2, 7, 2, 2, 2, 2, 0, 0, 3, 4, 0, 3, 0, 5, 3, 4, 5, 6, 5, 5, 5, 5, 0, 0, 0, 0&
578 , 0, 0, 0, 3, 0, 0, 3, 0, 3, 3, 3, 3, 0, 0, 0, 0, 0, 0, 0, 6, 0, 0, 6, 0, 4, 6, 6, 6, 0, 0&
579 , 0, 0, 0, 0, 0, 4, 0, 0, 4, 0, 0, 4, 4, 4, 0, 0, 0, 0, 0, 0, 0, 5, 0, 0, 5, 0, 0, 5, 7, 7&
580 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 5, 5, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
581 , 6, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 6, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
582 , 0, 0, 0, 7, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
583 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
584 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
585 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: mb3 = reshape((/1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1&
586 , 1, 1, 1, 1, 2, 2, 2, 2, 0, 2, 2, 2, 2, 0, 2, 2, 2, 2, 0, 2, 0, 0, 0, 0, 3, 0, 0, 0, 0, 3&
587 , 0, 0, 0, 0, 3, 0, 1, 3, 3, 3, 2, 3, 1, 3, 3, 2, 3, 3, 1, 3, 2, 3, 0, 1, 1, 1, 0, 1, 4, 1&
588 , 1, 5, 1, 1, 4, 1, 5, 1, 0, 2, 4, 4, 0, 4, 5, 4, 4, 4, 4, 4, 5, 4, 4, 4, 0, 0, 2, 2, 0, 2&
589 , 3, 2, 2, 0, 2, 2, 3, 2, 7, 2, 0, 0, 3, 5, 0, 3, 4, 5, 3, 0, 5, 5, 4, 5, 6, 5, 0, 0, 0, 3&
590 , 0, 0, 0, 3, 0, 0, 3, 3, 0, 3, 0, 3, 0, 0, 0, 4, 0, 0, 0, 6, 0, 0, 6, 6, 0, 6, 0, 6, 0, 0&
591 , 0, 0, 0, 0, 0, 4, 0, 0, 4, 4, 0, 4, 0, 4, 0, 0, 0, 0, 0, 0, 0, 5, 0, 0, 5, 7, 0, 5, 0, 7&
592 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 5, 0, 0, 0, 5, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 6, 0, 0&
593 , 0, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 6, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
594 , 0, 0, 0, 7, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
595 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
596 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
597 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: mb4 = reshape((/2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2&
598 , 2, 2, 2, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 3, 3, 3, 3, 4, 3, 3, 3, 3&
599 , 4, 3, 3, 3, 3, 4, 0, 1, 1, 1, 1, 0, 1, 1, 1, 1, 2, 1, 1, 1, 1, 2, 0, 2, 4, 4, 2, 4, 4, 2&
600 , 4, 4, 0, 4, 4, 2, 4, 0, 0, 0, 2, 2, 0, 2, 0, 0, 2, 0, 6, 2, 2, 0, 2, 6, 0, 3, 0, 0, 3, 0&
601 , 5, 5, 0, 5, 4, 0, 0, 5, 0, 2, 0, 1, 3, 5, 1, 0, 1, 1, 3, 1, 2, 5, 5, 1, 5, 8, 0, 0, 1, 3&
602 , 0, 0, 4, 6, 1, 4, 6, 3, 3, 6, 3, 6, 0, 0, 4, 1, 0, 0, 2, 4, 4, 2, 4, 1, 1, 4, 1, 4, 0, 0&
603 , 2, 4, 0, 0, 5, 5, 2, 5, 0, 6, 4, 5, 6, 8, 0, 0, 0, 2, 0, 0, 3, 3, 0, 3, 0, 4, 2, 3, 4, 6&
604 , 0, 0, 0, 5, 0, 0, 0, 6, 0, 0, 0, 2, 5, 6, 2, 0, 0, 0, 0, 3, 0, 0, 0, 4, 0, 0, 0, 7, 3, 4&
605 , 7, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 5, 0, 0, 5, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 3&
606 , 0, 0, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 6, 0, 0, 6, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
607 , 0, 4, 0, 0, 4, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 7, 0, 0, 7, 0, 0, 0, 0, 0, 0, 0, 0, 0&
608 , 0, 0, 0, 5, 0, 0, 5, 0/), (/4, 4, 20/))
609 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: mb5 = reshape((/2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2&
610 , 2, 2, 2, 2, 0, 3, 3, 3, 3, 4, 3, 3, 3, 3, 4, 3, 3, 3, 3, 4, 0, 0, 4, 4, 0, 0, 4, 0, 4, 4&
611 , 0, 4, 4, 0, 4, 0, 0, 1, 0, 0, 1, 2, 0, 5, 0, 0, 6, 0, 0, 5, 0, 6, 0, 0, 1, 5, 0, 0, 5, 1&
612 , 1, 5, 2, 5, 5, 1, 5, 2, 0, 0, 2, 1, 0, 0, 1, 6, 2, 1, 4, 1, 1, 6, 1, 8, 0, 0, 0, 2, 0, 0&
613 , 2, 3, 0, 2, 0, 6, 2, 3, 6, 4, 0, 0, 0, 3, 0, 0, 3, 4, 0, 3, 0, 2, 3, 4, 2, 6, 0, 0, 0, 0&
614 , 0, 0, 0, 0, 0, 0, 0, 7, 0, 0, 7, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, 0, 3, 0, 0, 0&
615 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 4, 0, 0, 4, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 5, 0, 0, 5, 0&
616 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
617 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
618 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
619 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
620 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
622 INTEGER :: k, k1, k2, mu
623 REAL(kind=
dp) :: cp, ct, fac1, fac2, j, jc, jcc, jss, rr, &
624 sp, st, xx, yy, za, zb
625 REAL(kind=
dp),
DIMENSION(3) :: v
626 REAL(kind=
dp),
DIMENSION(3, 3) :: arot
631 rr = sqrt(dot_product(v, v))
633 IF (rr < 1.0e-20_dp)
THEN
646 fac1 = fac1*sqrt(4.0_dp/3.0_dp)
649 fac1 = fac1*sqrt(8.0_dp/45.0_dp)
652 fac1 = fac1*sqrt(4.0_dp/315.0_dp)
654 WRITE (*, *)
'nra= ', nra
664 fac1 = fac1*sqrt(4.0_dp/3.0_dp)
667 fac1 = fac1*sqrt(8.0_dp/45.0_dp)
670 fac1 = fac1*sqrt(4.0_dp/315.0_dp)
672 WRITE (*, *)
'nrb= ', nrb
680 IF (abs(ct) < 1.0_dp)
THEN
681 st = sqrt(1.0_dp - ct**2)
707 fac2 = sqrt(za**(2*nra + 1)*zb**(2*nrb + 1))
708 xx = 0.5_dp*rr*(za + zb)
709 yy = 0.5_dp*rr*(za - zb)
712 DO k = 1, nc1(nra, nrb)
713 j = j + real(c1(nra, nrb, k),
dp)*aa(ma1(nra, nrb, k), xx)*bb(mb1(nra, nrb, k), yy)
715 j = j*rr**(nra + nrb + 1)
716 j = j/2.0_dp**(nra + nrb + 2)
718 s(1, 1) = s(1, 1) + fac1*fac2*j
722 fac2 = sqrt(za**(2*nra + 1)*zb**(2*nrb + 1))
723 xx = 0.5_dp*rr*(za + zb)
724 yy = 0.5_dp*rr*(za - zb)
727 DO k = 1, nc2(nra, nrb)
728 jc = jc + real(c2(nra, nrb, k),
dp)*aa(ma2(nra, nrb, k), xx)*bb(mb2(nra, nrb, k), yy)
730 jc = jc*rr**(nra + nrb + 1)
731 jc = jc/2.0_dp**(nra + nrb + 2)
734 s(k1 + 1, 1) = s(k1 + 1, 1) &
735 & + sqrt(3.0_dp)*arot(k1, 3)*fac1*fac2*jc
740 fac2 = sqrt(za**(2*nra + 1)*zb**(2*nrb + 1))
741 xx = 0.5_dp*rr*(za + zb)
742 yy = 0.5_dp*rr*(za - zb)
745 DO k = 1, nc3(nra, nrb)
746 jc = jc + real(c3(nra, nrb, k),
dp)*aa(ma3(nra, nrb, k), xx)*bb(mb3(nra, nrb, k), yy)
748 jc = jc*rr**(nra + nrb + 1)
749 jc = jc/2.0_dp**(nra + nrb + 2)
752 s(1, k1 + 1) = s(1, k1 + 1) &
753 & - sqrt(3.0_dp)*arot(k1, 3)*fac1*fac2*jc
758 fac2 = sqrt(za**(2*nra + 1)*zb**(2*nrb + 1))
759 xx = 0.5_dp*rr*(za + zb)
760 yy = 0.5_dp*rr*(za - zb)
763 DO k = 1, nc4(nra, nrb)
764 jss = jss + real(c4(nra, nrb, k),
dp)*aa(ma4(nra, nrb, k), xx)*bb(mb4(nra, nrb, k), yy)
766 jss = jss*rr**(nra + nrb + 1)
767 jss = jss/2.0_dp**(nra + nrb + 2)
770 DO k = 1, nc5(nra, nrb)
771 jcc = jcc + real(c5(nra, nrb, k),
dp)*aa(ma5(nra, nrb, k), xx)*bb(mb5(nra, nrb, k), yy)
773 jcc = jcc*rr**(nra + nrb + 1)
774 jcc = jcc/2.0_dp**(nra + nrb + 2)
778 s(k1 + 1, k2 + 1) = s(k1 + 1, k2 + 1) &
779 & + 1.5_dp*arot(k1, 1)*arot(k2, 1)*fac1*fac2*jss &
780 & + 1.5_dp*arot(k1, 2)*arot(k2, 2)*fac1*fac2*jss &
781 & - 3.0_dp*arot(k1, 3)*arot(k2, 3)*fac1*fac2*jcc
800 SUBROUTINE makeds(R, nra, nrb, ZSA, ZSB, ZPA, ZPB, dS)
802 REAL(kind=
dp),
DIMENSION(3) :: r
804 REAL(kind=
dp) :: zsa, zsb, zpa, zpb
805 REAL(kind=
dp),
DIMENSION(4, 4, 3) :: ds
807 INTEGER,
DIMENSION(4, 4),
PARAMETER :: &
808 nc1 = reshape((/2, 4, 4, 6, 4, 3, 6, 7, 4, 6, 4, 8, 6, 7, 8, 5/), (/4, 4/)), &
809 nc2 = reshape((/4, 4, 8, 8, 6, 8, 6, 12, 8, 8, 12, 8, 10, 12, 14, 16/), (/4, 4/)), &
810 nc3 = reshape((/4, 6, 8, 10, 4, 8, 8, 12, 8, 6, 12, 14, 8, 12, 8, 16/), (/4, 4/)), &
811 nc4 = reshape((/4, 8, 11, 14, 8, 6, 12, 14, 11, 12, 10, 20, 14, 14, 20, 12/), (/4, 4/)), &
812 nc5 = reshape((/2, 4, 6, 8, 4, 4, 8, 8, 6, 8, 6, 12, 8, 8, 12, 8/), (/4, 4/))
813 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: c1 = reshape((/1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1 &
814 , 1, 1, 1, 1, -1, 1, 2, 3, -1, -2, 1, 2, -2, -1, -3, 1, -3, -2, -1, -4, 0, -1, -2, 2, -1, &
815 1, -2, -1, 2, -2, 3, -3, 2, -1, -3, 6, 0, -1, -1, -2, 1, 0, -2, -4, -1, 2, -1, -3, 2, 4, 3&
816 , -4, 0, 0, 0, -3, 0, 0, 1, -1, 0, 1, 0, 3, -3, -1, 3, 1, 0, 0, 0, -1, 0, 0, 1, 2, 0, -1, &
817 0, 3, 1, -2, -3, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, -1, 0, 1, -1, 0, 0, 0, 0, 0, 0, 0, 0 &
818 , 0, 0, 0, 0, -1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
819 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
820 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
821 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
822 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
823 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
824 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
825 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: c2 = reshape((/1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1 &
826 , 1, 1, 1, 1, -1, 1, 1, 2, -2, -1, 1, 1, -3, -2, -1, 1, -4, -3, -2, -1, 1, -1, 1, 1, 1, 1 &
827 , -2, 1, 1, 1, 1, -3, 1, 1, 1, 1, -1, -1, -1, 2, 1, -1, -2, -2, 3, -2, -2, -3, 6, 2, -1, &
828 -3, 0, 0, 1, -2, -2, -1, 1, 1, -3, 2, -1, 3, -4, -3, -2, -1, 0, 0, -1, -1, 1, 1, 1, -2, -1&
829 , -1, 2, 3, -4, 2, 4, 3, 0, 0, -1, -2, 0, -1, 0, -2, 3, 2, -2, -1, 6, 2, -1, -3, 0, 0, -1,&
830 -1, 0, 1, 0, 1, -1, -1, 1, -1, 1, -3, -1, 3, 0, 0, 0, 0, 0, 0, 0, -2, 0, 0, 2, 0, -4, 2, 4&
831 , 3, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, -1, 0, 1, 1, -2, -3, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0&
832 , 0, -3, -1, 3, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, -1, 0, 0, 1, 1, -1, 0, 0, 0, 0, 0, 0, 0, 0, &
833 0, 0, 0, 0, 0, 0, -2, -3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0&
834 , 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, &
835 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
836 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
838 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: c3 = reshape((/-1, -1, -1, -1, -1, -1, -1, -1, -1 &
839 , -1, -1, -1, -1, -1, -1, -1, -1, -2, -3, -4, 1, -1, -2, -3, 1, 1, -1, -2, 2, 1, 1, -1, 1 &
840 , 1, 1, 1, 1, 1, 1, 1, 1, 2, 1, 1, 1, 1, 3, 1, 1, -1, -3, -6, -1, 1, 2, -2, 1, -2, 2, 1, &
841 -2, 2, -3, 3, 0, 2, 3, 4, 0, 1, 2, 3, -1, -1, 1, 2, -2, -1, -3, 1, 0, 1, -1, -4, 0, 1, 1, &
842 2, -1, 1, 2, 4, 1, -2, 3, 3, 0, 0, 3, 6, 0, -1, -2, 2, -1, 0, -2, -1, 2, -2, 1, -3, 0, 0, &
843 1, -1, 0, -1, -1, 3, 1, 0, -1, 1, -1, -1, -1, -3, 0, 0, 0, 4, 0, 0, 0, -2, 0, 0, -2, -4, 0&
844 , 2, 0, -3, 0, 0, 0, 1, 0, 0, 0, 1, 0, 0, -1, -2, 0, 1, 0, -3, 0, 0, 0, 0, 0, 0, 0, -3, 0 &
845 , 0, 1, -1, 0, 1, 0, 3, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 1, -1, 0, -1, 0, 1, 0, 0, 0, 0, 0, &
846 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, &
847 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0&
848 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
849 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
850 , 0, 0, 0/), (/4, 4, 20/))
851 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: c4 = reshape((/-1, -1, -1, -1, -1, -1, -1, -1, -1 &
852 , -1, -1, -1, -1, -1, -1, -1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, -1, -2, &
853 -3, 1, 1, -1, -2, 2, 1, 2, -1, 3, 2, 1, 3, -1, 1, 2, 3, -1, -1, 1, 2, -2, -1, -1, 1, -3, &
854 -2, -1, -2, 0, 1, -1, -3, 1, -1, 1, 1, -1, 1, -1, 2, -3, 1, 2, -1, 0, -1, 2, 4, -1, 1, -1 &
855 , -1, 2, -1, -1, -1, 4, -1, -1, -3, 0, 1, -1, -1, -1, 0, 1, 2, -1, -1, -1, -1, -1, -2, -1 &
856 , 3, 0, -1, 2, -1, 1, 0, -1, -2, -2, 1, 2, 2, 1, 2, -2, 1, 0, 0, -2, 4, 0, 0, -1, 1, 2, -1&
857 , 1, -1, -4, 1, 1, 2, 0, 0, 1, -3, 0, 0, 1, -1, 1, 1, -1, -1, 3, -1, 1, -3, 0, 0, -1, 3, 0&
858 , 0, -1, -2, -1, 1, 0, -1, 3, 2, -1, -1, 0, 0, 0, -3, 0, 0, 1, 2, 0, -1, 0, -1, -3, -2, -1&
859 , 1, 0, 0, 0, 1, 0, 0, 0, -1, 0, 0, 0, 2, -1, -1, 2, 0, 0, 0, 0, -1, 0, 0, 0, 1, 0, 0, 0, &
860 -1, 1, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
861 , 0, 0, 2, 0, 0, -2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, &
862 0, 0, 0, 0, 0, -1, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, -1, 0, 0, 0, 0, &
863 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 1, 0/), (/4, 4, 20/))
864 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: c5 = reshape((/-1, -1, -1, -1, -1, -1, -1, -1, -1 &
865 , -1, -1, -1, -1, -1, -1, -1, 1, -1, -2, -3, 1, 1, -1, -2, 2, 1, 2, -1, 3, 2, 1, 3, 0, 1, &
866 -1, -3, 1, 1, 1, 1, -1, 1, 1, 2, -3, 1, 2, 1, 0, 1, 1, 1, -1, -1, 1, 2, 1, 1, -1, 1, 1, -2&
867 , 1, -3, 0, 0, 2, -1, 0, 0, 1, 2, -2, -1, -2, 2, 1, -2, -2, -3, 0, 0, 1, 3, 0, 0, 1, 1, 1,&
868 -1, 1, 1, -3, 1, -1, 1, 0, 0, 0, 3, 0, 0, -1, -2, 0, -1, 0, -1, 3, 2, -1, 3, 0, 0, 0, 1, 0&
869 , 0, -1, -1, 0, 1, 0, -2, -1, -1, -2, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 1, 0,&
870 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -2, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0,&
871 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 &
872 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
873 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
874 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
875 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, &
877 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: ma1 = reshape((/2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7&
878 , 5, 6, 7, 8, 0, 2, 3, 4, 2, 2, 4, 5, 3, 4, 4, 6, 4, 5, 6, 6, 0, 1, 1, 3, 1, 0, 3, 4, 1, 3&
879 , 2, 5, 3, 4, 5, 4, 0, 0, 0, 2, 0, 0, 2, 3, 0, 2, 0, 4, 2, 3, 4, 2, 0, 0, 0, 1, 0, 0, 1, 2&
880 , 0, 1, 0, 3, 1, 2, 3, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 2, 0, 1, 2, 0, 0, 0, 0, 0, 0, 0&
881 , 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
882 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
883 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
884 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
885 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
886 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
887 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
888 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
889 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: ma2 = reshape((/2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7&
890 , 5, 6, 7, 8, 1, 2, 3, 4, 2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7, 1, 1, 3, 4, 2, 3, 3, 5, 3, 4&
891 , 5, 5, 4, 5, 6, 7, 0, 0, 2, 3, 1, 2, 2, 4, 2, 3, 4, 4, 3, 4, 5, 6, 0, 0, 2, 2, 1, 2, 1, 4&
892 , 2, 2, 4, 3, 3, 4, 5, 6, 0, 0, 1, 1, 0, 1, 0, 3, 1, 1, 3, 2, 2, 3, 4, 5, 0, 0, 1, 1, 0, 1&
893 , 0, 3, 1, 1, 3, 1, 2, 3, 4, 5, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 2, 0, 1, 2, 3, 4, 0, 0, 0, 0&
894 , 0, 0, 0, 2, 0, 0, 2, 0, 1, 2, 3, 4, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 2, 3, 0, 0&
895 , 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 2, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 2&
896 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
897 , 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
898 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
899 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
900 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
901 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: ma3 = reshape((/2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7&
902 , 5, 6, 7, 8, 1, 2, 3, 4, 2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7, 1, 2, 3, 4, 1, 3, 4, 5, 3, 3&
903 , 5, 6, 4, 5, 5, 7, 0, 1, 2, 3, 0, 2, 3, 4, 2, 2, 4, 5, 3, 4, 4, 6, 0, 1, 2, 3, 0, 2, 2, 4&
904 , 2, 1, 4, 5, 2, 4, 3, 6, 0, 0, 1, 2, 0, 1, 1, 3, 1, 0, 3, 4, 1, 3, 2, 5, 0, 0, 1, 2, 0, 1&
905 , 1, 3, 1, 0, 3, 4, 1, 3, 1, 5, 0, 0, 0, 1, 0, 0, 0, 2, 0, 0, 2, 3, 0, 2, 0, 4, 0, 0, 0, 1&
906 , 0, 0, 0, 2, 0, 0, 2, 3, 0, 2, 0, 4, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 2, 0, 1, 0, 3, 0, 0&
907 , 0, 0, 0, 0, 0, 1, 0, 0, 1, 2, 0, 1, 0, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 2&
908 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
909 , 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
910 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
911 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
912 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
913 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: ma4 = reshape((/2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7&
914 , 5, 6, 7, 8, 2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7, 5, 6, 7, 8, 0, 2, 3, 4, 2, 2, 4, 5, 3, 4&
915 , 4, 6, 4, 5, 6, 6, 0, 2, 3, 4, 2, 2, 4, 5, 3, 4, 4, 6, 4, 5, 6, 6, 0, 1, 2, 3, 1, 0, 3, 4&
916 , 2, 3, 4, 5, 3, 4, 5, 6, 0, 1, 2, 3, 1, 0, 3, 4, 2, 3, 2, 5, 3, 4, 5, 4, 0, 0, 2, 3, 0, 0&
917 , 2, 3, 2, 2, 2, 5, 3, 3, 5, 4, 0, 0, 1, 2, 0, 0, 2, 3, 1, 2, 2, 4, 2, 3, 4, 2, 0, 0, 1, 2&
918 , 0, 0, 1, 2, 1, 1, 0, 4, 2, 2, 4, 2, 0, 0, 0, 2, 0, 0, 1, 2, 0, 1, 0, 4, 2, 2, 4, 2, 0, 0&
919 , 0, 1, 0, 0, 0, 1, 0, 0, 0, 3, 1, 1, 3, 0, 0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 3, 1, 1, 3, 0&
920 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, 0, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0&
921 , 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2&
922 , 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
923 , 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
924 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
925 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: ma5 = reshape((/2, 3, 4, 5, 3, 4, 5, 6, 4, 5, 6, 7&
926 , 5, 6, 7, 8, 0, 2, 3, 4, 2, 2, 4, 5, 3, 4, 4, 6, 4, 5, 6, 6, 0, 1, 2, 3, 1, 2, 3, 4, 2, 3&
927 , 4, 5, 3, 4, 5, 6, 0, 0, 2, 3, 0, 0, 3, 3, 2, 3, 2, 5, 3, 3, 5, 4, 0, 0, 1, 2, 0, 0, 2, 3&
928 , 1, 2, 2, 4, 2, 3, 4, 4, 0, 0, 0, 2, 0, 0, 2, 2, 0, 2, 0, 4, 2, 2, 4, 2, 0, 0, 0, 1, 0, 0&
929 , 1, 1, 0, 1, 0, 3, 1, 1, 3, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, 0, 3, 0, 0, 0, 0, 0&
930 , 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 2, 0, 0, 0&
931 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
932 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
933 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
934 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
935 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
936 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
937 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: mb1 = reshape((/0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
938 , 0, 0, 0, 0, 2, 1, 1, 1, 1, 2, 1, 1, 1, 1, 2, 1, 1, 1, 1, 2, 0, 2, 3, 2, 2, 4, 2, 2, 3, 2&
939 , 4, 2, 2, 2, 2, 4, 0, 3, 4, 3, 3, 0, 3, 3, 4, 3, 6, 3, 3, 3, 3, 6, 0, 0, 0, 4, 0, 0, 4, 4&
940 , 0, 4, 0, 4, 4, 4, 4, 8, 0, 0, 0, 5, 0, 0, 5, 5, 0, 5, 0, 5, 5, 5, 5, 0, 0, 0, 0, 0, 0, 0&
941 , 0, 6, 0, 0, 0, 6, 0, 6, 6, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 7, 0, 0, 7, 0, 0, 0, 0, 0&
942 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
943 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
944 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
945 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
946 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
947 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
948 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
949 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: mb2 = reshape((/1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1&
950 , 1, 1, 1, 1, 2, 0, 2, 2, 2, 2, 0, 2, 2, 2, 2, 0, 2, 2, 2, 2, 0, 3, 0, 0, 0, 0, 3, 0, 0, 0&
951 , 0, 3, 0, 0, 0, 0, 1, 2, 3, 1, 3, 3, 2, 3, 3, 1, 3, 2, 3, 3, 3, 3, 0, 0, 1, 4, 1, 1, 5, 1&
952 , 1, 4, 1, 5, 1, 1, 1, 1, 0, 0, 4, 5, 2, 4, 4, 4, 4, 5, 4, 4, 4, 4, 4, 4, 0, 0, 2, 3, 0, 2&
953 , 0, 2, 2, 3, 2, 7, 2, 2, 2, 2, 0, 0, 3, 4, 0, 3, 0, 5, 3, 4, 5, 6, 5, 5, 5, 5, 0, 0, 0, 0&
954 , 0, 0, 0, 3, 0, 0, 3, 0, 3, 3, 3, 3, 0, 0, 0, 0, 0, 0, 0, 6, 0, 0, 6, 0, 4, 6, 6, 6, 0, 0&
955 , 0, 0, 0, 0, 0, 4, 0, 0, 4, 0, 0, 4, 4, 4, 0, 0, 0, 0, 0, 0, 0, 5, 0, 0, 5, 0, 0, 5, 7, 7&
956 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 5, 5, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
957 , 6, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 6, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
958 , 0, 0, 0, 7, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
959 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
960 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
961 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: mb3 = reshape((/1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1&
962 , 1, 1, 1, 1, 2, 2, 2, 2, 0, 2, 2, 2, 2, 0, 2, 2, 2, 2, 0, 2, 0, 0, 0, 0, 3, 0, 0, 0, 0, 3&
963 , 0, 0, 0, 0, 3, 0, 1, 3, 3, 3, 2, 3, 1, 3, 3, 2, 3, 3, 1, 3, 2, 3, 0, 1, 1, 1, 0, 1, 4, 1&
964 , 1, 5, 1, 1, 4, 1, 5, 1, 0, 2, 4, 4, 0, 4, 5, 4, 4, 4, 4, 4, 5, 4, 4, 4, 0, 0, 2, 2, 0, 2&
965 , 3, 2, 2, 0, 2, 2, 3, 2, 7, 2, 0, 0, 3, 5, 0, 3, 4, 5, 3, 0, 5, 5, 4, 5, 6, 5, 0, 0, 0, 3&
966 , 0, 0, 0, 3, 0, 0, 3, 3, 0, 3, 0, 3, 0, 0, 0, 4, 0, 0, 0, 6, 0, 0, 6, 6, 0, 6, 0, 6, 0, 0&
967 , 0, 0, 0, 0, 0, 4, 0, 0, 4, 4, 0, 4, 0, 4, 0, 0, 0, 0, 0, 0, 0, 5, 0, 0, 5, 7, 0, 5, 0, 7&
968 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 5, 0, 0, 0, 5, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 6, 0, 0&
969 , 0, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 6, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
970 , 0, 0, 0, 7, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
971 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
972 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
973 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: mb4 = reshape((/2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2&
974 , 2, 2, 2, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 3, 3, 3, 3, 4, 3, 3, 3, 3&
975 , 4, 3, 3, 3, 3, 4, 0, 1, 1, 1, 1, 0, 1, 1, 1, 1, 2, 1, 1, 1, 1, 2, 0, 2, 4, 4, 2, 4, 4, 2&
976 , 4, 4, 0, 4, 4, 2, 4, 0, 0, 0, 2, 2, 0, 2, 0, 0, 2, 0, 6, 2, 2, 0, 2, 6, 0, 3, 0, 0, 3, 0&
977 , 5, 5, 0, 5, 4, 0, 0, 5, 0, 2, 0, 1, 3, 5, 1, 0, 1, 1, 3, 1, 2, 5, 5, 1, 5, 8, 0, 0, 1, 3&
978 , 0, 0, 4, 6, 1, 4, 6, 3, 3, 6, 3, 6, 0, 0, 4, 1, 0, 0, 2, 4, 4, 2, 4, 1, 1, 4, 1, 4, 0, 0&
979 , 2, 4, 0, 0, 5, 5, 2, 5, 0, 6, 4, 5, 6, 8, 0, 0, 0, 2, 0, 0, 3, 3, 0, 3, 0, 4, 2, 3, 4, 6&
980 , 0, 0, 0, 5, 0, 0, 0, 6, 0, 0, 0, 2, 5, 6, 2, 0, 0, 0, 0, 3, 0, 0, 0, 4, 0, 0, 0, 7, 3, 4&
981 , 7, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 5, 0, 0, 5, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 3&
982 , 0, 0, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 6, 0, 0, 6, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
983 , 0, 4, 0, 0, 4, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 7, 0, 0, 7, 0, 0, 0, 0, 0, 0, 0, 0, 0&
984 , 0, 0, 0, 5, 0, 0, 5, 0/), (/4, 4, 20/))
985 INTEGER,
DIMENSION(4, 4, 20),
PARAMETER :: mb5 = reshape((/2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2&
986 , 2, 2, 2, 2, 0, 3, 3, 3, 3, 4, 3, 3, 3, 3, 4, 3, 3, 3, 3, 4, 0, 0, 4, 4, 0, 0, 4, 0, 4, 4&
987 , 0, 4, 4, 0, 4, 0, 0, 1, 0, 0, 1, 2, 0, 5, 0, 0, 6, 0, 0, 5, 0, 6, 0, 0, 1, 5, 0, 0, 5, 1&
988 , 1, 5, 2, 5, 5, 1, 5, 2, 0, 0, 2, 1, 0, 0, 1, 6, 2, 1, 4, 1, 1, 6, 1, 8, 0, 0, 0, 2, 0, 0&
989 , 2, 3, 0, 2, 0, 6, 2, 3, 6, 4, 0, 0, 0, 3, 0, 0, 3, 4, 0, 3, 0, 2, 3, 4, 2, 6, 0, 0, 0, 0&
990 , 0, 0, 0, 0, 0, 0, 0, 7, 0, 0, 7, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, 0, 3, 0, 0, 0&
991 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 4, 0, 0, 4, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 5, 0, 0, 5, 0&
992 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
993 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
994 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
995 , 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0&
996 , 0, 0, 0, 0, 0, 0, 0, 0/), (/4, 4, 20/))
998 INTEGER :: k, k1, k2, mu
999 REAL(kind=
dp) :: cp, ct, dj, djc, djcc, djss, dxx, dyy, &
1000 f, fac1, fac2, j, jc, jcc, jss, rr, &
1001 sp, st, w, w1, w2, xx, yy, za, zb
1002 REAL(kind=
dp),
DIMENSION(3) :: dcp,
dct, dsp, dst, v
1003 REAL(kind=
dp),
DIMENSION(3, 3) :: arot
1004 REAL(kind=
dp),
DIMENSION(3, 3, 3) :: darot
1006 ds(:, :, :) = 0.0_dp
1009 rr = sqrt(dot_product(v, v))
1011 IF (rr < 1.0e-20_dp)
THEN
1014 ds(mu, mu, :) = 0.0_dp
1024 fac1 = fac1*sqrt(4.0_dp/3.0_dp)
1027 fac1 = fac1*sqrt(8.0_dp/45.0_dp)
1030 fac1 = fac1*sqrt(4.0_dp/315.0_dp)
1032 WRITE (*, *)
'nra= ', nra
1042 fac1 = fac1*sqrt(4.0_dp/3.0_dp)
1045 fac1 = fac1*sqrt(8.0_dp/45.0_dp)
1048 fac1 = fac1*sqrt(4.0_dp/315.0_dp)
1050 WRITE (*, *)
'nrb= ', nrb
1058 IF (abs(ct) >= 1.0_dp)
THEN
1060 dct(:) = v(:)*v(3)/rr**3
1061 dct(3) =
dct(3) - 1.0_dp/rr
1073 darot(1, 1, :) =
dct(:)
1074 darot(1, 2, :) = 0.0_dp
1075 darot(1, 3, :) = 0.0_dp
1076 darot(2, 1, :) = 0.0_dp
1077 darot(2, 2, :) = 0.0_dp
1078 darot(2, 3, :) = 0.0_dp
1079 darot(3, 1, :) = 0.0_dp
1080 darot(3, 2, :) = 0.0_dp
1081 darot(3, 3, :) =
dct(:)
1085 xx = sqrt(v(1)**2 + v(2)**2)
1090 dct(:) = v(:)*v(3)/rr**3
1091 dct(3) =
dct(3) - 1.0_dp/rr
1092 dst(:) = -ct*
dct(:)/st
1093 dcp(:) = v(:)*v(1)/(rr**3*st)
1094 dcp(:) = dcp(:) + v(1)*dst(:)/(rr*st**2)
1095 dcp(1) = dcp(1) - 1.0_dp/(rr*st)
1096 dsp(:) = v(:)*v(2)/(rr**3*st)
1097 dsp(:) = dsp(:) + v(2)*dst(:)/(rr*st**2)
1098 dsp(2) = dsp(2) - 1.0_dp/(rr*st)
1110 darot(1, 1, :) =
dct(:)*cp + ct*dcp(:)
1111 darot(1, 2, :) = -dsp(:)
1112 darot(1, 3, :) = dst(:)*cp + st*dcp(:)
1113 darot(2, 1, :) =
dct(:)*
sp + ct*dsp(:)
1114 darot(2, 2, :) = dcp(:)
1115 darot(2, 3, :) = dst(:)*
sp + st*dsp(:)
1116 darot(3, 1, :) = -dst(:)
1117 darot(3, 2, :) = 0.0_dp
1118 darot(3, 3, :) =
dct(:)
1124 fac2 = sqrt(za**(2*nra + 1)*zb**(2*nrb + 1))
1125 xx = 0.5_dp*rr*(za + zb)
1126 yy = 0.5_dp*rr*(za - zb)
1127 dxx = 0.5_dp*(za + zb)
1128 dyy = 0.5_dp*(za - zb)
1133 f = rr**(nra + nrb + 1)/2.0_dp**(nra + nrb + 2)
1134 DO k = 1, nc1(nra, nrb)
1135 w = w + real(c1(nra, nrb, k),
dp)*aa(ma1(nra, nrb, k), xx)*bb(mb1(nra, nrb, k), yy)
1136 w1 = w1 + real(c1(nra, nrb, k),
dp)*aa(ma1(nra, nrb, k) + 1, xx)*bb(mb1(nra, nrb, k), yy)
1137 w2 = w2 + real(c1(nra, nrb, k),
dp)*aa(ma1(nra, nrb, k), xx)*bb(mb1(nra, nrb, k) + 1, yy)
1140 dj = f*real(nra + nrb + 1,
dp)*w/rr
1144 ds(1, 1, :) = ds(1, 1, :) + fac1*fac2*dj*v(:)/rr
1148 fac2 = sqrt(za**(2*nra + 1)*zb**(2*nrb + 1))
1149 xx = 0.5_dp*rr*(za + zb)
1150 yy = 0.5_dp*rr*(za - zb)
1151 dxx = 0.5_dp*(za + zb)
1152 dyy = 0.5_dp*(za - zb)
1157 f = rr**(nra + nrb + 1)/2.0_dp**(nra + nrb + 2)
1158 DO k = 1, nc2(nra, nrb)
1159 w = w + real(c2(nra, nrb, k),
dp)*aa(ma2(nra, nrb, k), xx)*bb(mb2(nra, nrb, k), yy)
1160 w1 = w1 + real(c2(nra, nrb, k),
dp)*aa(ma2(nra, nrb, k) + 1, xx)*bb(mb2(nra, nrb, k), yy)
1161 w2 = w2 + real(c2(nra, nrb, k),
dp)*aa(ma2(nra, nrb, k), xx)*bb(mb2(nra, nrb, k) + 1, yy)
1164 djc = f*real(nra + nrb + 1,
dp)*w/rr
1165 djc = djc - dxx*f*w1
1166 djc = djc - dyy*f*w2
1169 ds(k1 + 1, 1, :) = ds(k1 + 1, 1, :) &
1170 & + sqrt(3.0_dp)*arot(k1, 3)*fac1*fac2*djc*v(:)/rr &
1171 & + sqrt(3.0_dp)*darot(k1, 3, :)*fac1*fac2*jc
1176 fac2 = sqrt(za**(2*nra + 1)*zb**(2*nrb + 1))
1177 xx = 0.5_dp*rr*(za + zb)
1178 yy = 0.5_dp*rr*(za - zb)
1179 dxx = 0.5_dp*(za + zb)
1180 dyy = 0.5_dp*(za - zb)
1185 f = rr**(nra + nrb + 1)/2.0_dp**(nra + nrb + 2)
1186 DO k = 1, nc3(nra, nrb)
1187 w = w + real(c3(nra, nrb, k),
dp)*aa(ma3(nra, nrb, k), xx)*bb(mb3(nra, nrb, k), yy)
1188 w1 = w1 + real(c3(nra, nrb, k),
dp)*aa(ma3(nra, nrb, k) + 1, xx)*bb(mb3(nra, nrb, k), yy)
1189 w2 = w2 + real(c3(nra, nrb, k),
dp)*aa(ma3(nra, nrb, k), xx)*bb(mb3(nra, nrb, k) + 1, yy)
1192 djc = f*real(nra + nrb + 1,
dp)*w/rr
1193 djc = djc - dxx*f*w1
1194 djc = djc - dyy*f*w2
1197 ds(1, k1 + 1, :) = ds(1, k1 + 1, :) &
1198 & - sqrt(3.0_dp)*arot(k1, 3)*fac1*fac2*djc*v(:)/rr &
1199 & - sqrt(3.0_dp)*darot(k1, 3, :)*fac1*fac2*jc
1204 fac2 = sqrt(za**(2*nra + 1)*zb**(2*nrb + 1))
1205 xx = 0.5_dp*rr*(za + zb)
1206 yy = 0.5_dp*rr*(za - zb)
1207 dxx = 0.5_dp*(za + zb)
1208 dyy = 0.5_dp*(za - zb)
1213 f = rr**(nra + nrb + 1)/2.0_dp**(nra + nrb + 2)
1214 DO k = 1, nc4(nra, nrb)
1215 w = w + real(c4(nra, nrb, k),
dp)*aa(ma4(nra, nrb, k), xx)*bb(mb4(nra, nrb, k), yy)
1216 w1 = w1 + real(c4(nra, nrb, k),
dp)*aa(ma4(nra, nrb, k) + 1, xx)*bb(mb4(nra, nrb, k), yy)
1217 w2 = w2 + real(c4(nra, nrb, k),
dp)*aa(ma4(nra, nrb, k), xx)*bb(mb4(nra, nrb, k) + 1, yy)
1220 djss = f*real(nra + nrb + 1,
dp)*w/rr
1221 djss = djss - dxx*f*w1
1222 djss = djss - dyy*f*w2
1227 f = rr**(nra + nrb + 1)/2.0_dp**(nra + nrb + 2)
1228 DO k = 1, nc5(nra, nrb)
1229 w = w + real(c5(nra, nrb, k),
dp)*aa(ma5(nra, nrb, k), xx)*bb(mb5(nra, nrb, k), yy)
1230 w1 = w1 + real(c5(nra, nrb, k),
dp)*aa(ma5(nra, nrb, k) + 1, xx)*bb(mb5(nra, nrb, k), yy)
1231 w2 = w2 + real(c5(nra, nrb, k),
dp)*aa(ma5(nra, nrb, k), xx)*bb(mb5(nra, nrb, k) + 1, yy)
1234 djcc = f*real(nra + nrb + 1,
dp)*w/rr
1235 djcc = djcc - dxx*f*w1
1236 djcc = djcc - dyy*f*w2
1240 ds(k1 + 1, k2 + 1, :) = ds(k1 + 1, k2 + 1, :) &
1241 & + 1.5_dp*arot(k1, 1)*arot(k2, 1)*fac1*fac2*djss*v(:)/rr &
1242 & + 1.5_dp*darot(k1, 1, :)*arot(k2, 1)*fac1*fac2*jss &
1243 & + 1.5_dp*arot(k1, 1)*darot(k2, 1, :)*fac1*fac2*jss &
1244 & + 1.5_dp*arot(k1, 2)*arot(k2, 2)*fac1*fac2*djss*v(:)/rr &
1245 & + 1.5_dp*darot(k1, 2, :)*arot(k2, 2)*fac1*fac2*jss &
1246 & + 1.5_dp*arot(k1, 2)*darot(k2, 2, :)*fac1*fac2*jss &
1247 & - 3.0_dp*arot(k1, 3)*arot(k2, 3)*fac1*fac2*djcc*v(:)/rr &
1248 & - 3.0_dp*darot(k1, 3, :)*arot(k2, 3)*fac1*fac2*jcc &
1249 & - 3.0_dp*arot(k1, 3)*darot(k2, 3, :)*fac1*fac2*jcc
1255 END SUBROUTINE makeds
1266 REAL(kind=
dp) :: x, aa
1292 p = 120.0_dp + x*( &
1299 p = 720.0_dp + x*( &
1307 p = 5040.0_dp + x*( &
1316 p = 40320.0_dp + x*( &
1326 p = 362880.0_dp + x*( &
1337 p = 3628800.0_dp + x*( &
1338 3628800.0_dp + x*( &
1339 1814400.0_dp + x*( &
1346 10.0_dp + x)))))))))
1349 WRITE (*, *)
' n= ', n,
' in AA(n,x) '
1362 aa = exp(-x)*p/x**(n + 1)
1375 REAL(kind=
dp) :: y, bb
1377 IF (abs(y) > 1.0e-20_dp)
THEN
1378 bb = real((-1)**(n + 1),
dp)*aa(n, -y) - aa(n, y)
1380 IF (mod(n, 2) == 0)
THEN
1381 bb = 2.0_dp/real(n + 1,
dp)
1383 bb = -y*2.0_dp/real(n + 2,
dp)
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.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_deallocate_matrix(matrix)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_get_block_p(matrix, row, col, block, found, row_size, col_size)
...
subroutine, public dbcsr_replicate_all(matrix)
...
subroutine, public dbcsr_distribute(matrix)
...
subroutine, public dbcsr_sum_replicated(matrix)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_add(matrix_a, matrix_b, alpha_scalar, beta_scalar)
...
subroutine, public dbcsr_get_block_diag(matrix, diag)
Copies the diagonal blocks of matrix into diag.
DBCSR operations in CP2K.
subroutine, public cp_dbcsr_write_sparse_matrix(sparse_matrix, before, after, qs_env, para_env, first_row, last_row, first_col, last_col, scale, output_unit, omit_headers)
...
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,...
integer, parameter, public cp_p_file
integer function, public cp_print_key_should_output(iteration_info, basis_section, print_key_path, used_print_key, first_time)
returns what should be done with the given property if btest(res,cp_p_store) then the property should...
the type I Discrete Cosine Transform (DCT-I)
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public sp
Interface to the message passing library MPI.
Define the data structure for the particle information.
Definition of physical constants:
real(kind=dp), parameter, public evolt
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, 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, 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, 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, ecoul_1c, rho0_s_rs, rho0_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)
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, 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, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
subroutine, public set_ks_env(ks_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, complex_ks, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, kinetic, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_ks_im_kp, vppl, rho_core, rho_nlcc, rho_nlcc_g, vee, neighbor_list_id, kpoints, sab_orb, sab_all, sac_ae, sac_ppl, sac_lri, sap_ppnl, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_vdw, sab_scp, sab_almo, sab_kp, sab_kp_nosym, task_list, task_list_soft, subsys, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env)
...
Define the neighbor list data types and the corresponding functionality.
subroutine, public neighbor_list_iterator_create(iterator_set, nl, search, nthread)
Neighbor list iterator functions.
subroutine, public neighbor_list_iterator_release(iterator_set)
...
integer function, public neighbor_list_iterate(iterator_set, mepos)
...
subroutine, public get_iterator_info(iterator_set, mepos, ikind, jkind, nkind, ilist, nlist, inode, nnode, iatom, jatom, r, cell)
...
Calculation of overlap matrix, its derivatives and forces.
subroutine, public build_overlap_matrix(ks_env, matrix_s, matrixkp_s, matrix_name, nderivative, basis_type_a, basis_type_b, sab_nl, calculate_forces, matrix_p, matrixkp_p)
Calculation of the overlap matrix over Cartesian Gaussian functions.
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public qs_rho_get(rho_struct, rho_ao, rho_ao_im, rho_ao_kp, rho_ao_im_kp, rho_r, drho_r, rho_g, drho_g, tau_r, tau_g, rho_r_valid, drho_r_valid, rho_g_valid, drho_g_valid, tau_r_valid, tau_g_valid, tot_rho_r, tot_rho_g, rho_r_sccs, soft_valid, complex_rho_ao)
returns info about the density described by this object. If some representation is not available an e...
Calculation of the Hamiltonian integral matrix <a|H|b> for semi-empirical methods.
subroutine, public build_se_core_matrix(qs_env, para_env, calculate_forces)
...
Arrays of parameters used in the semi-empirical calculations \References Everywhere in this module TC...
real(kind=dp), parameter, public rij_threshold
Definition of the semi empirical parameter types.
subroutine, public get_se_param(sep, name, typ, defined, z, zeff, natorb, eheat, beta, sto_exponents, uss, upp, udd, uff, alp, eisol, gss, gsp, gpp, gp2, acoul, nr, de, ass, asp, app, hsp, gsd, gpd, gdd, ppddg, dpddg, ngauss)
Get info from the semi-empirical type.
Working with the semi empirical parameter types.
integer function, public get_se_type(se_method)
Gives back the unique semi_empirical METHOD type.
pure subroutine, public virial_pair_force(pv_virial, f0, force, rab)
Computes the contribution to the stress tensor from two-body pair-wise forces.
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...
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...
keeps the density in various representations, keeping track of which ones are valid.