33#include "./base/base_uses.f90"
38 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'dkh_main'
101 SUBROUTINE dkh_full_transformation(qs_env, matrix_s, matrix_v, matrix_t, matrix_pVp, n, dkh_order)
102 TYPE(qs_environment_type),
POINTER :: qs_env
103 TYPE(cp_fm_type),
INTENT(IN) :: matrix_s, matrix_v, matrix_t, matrix_pVp
104 INTEGER,
INTENT(IN) :: n, dkh_order
106 CHARACTER(LEN=*),
PARAMETER :: routineN =
'DKH_full_transformation'
109 REAL(KIND=
dp),
ALLOCATABLE,
DIMENSION(:) :: aa, e, ev0t, rr, tt
110 TYPE(cp_blacs_env_type),
POINTER :: blacs_env
111 TYPE(cp_fm_struct_type),
POINTER :: matrix_full
112 TYPE(cp_fm_type) :: matrix_aux, matrix_aux2, matrix_eig, matrix_ev1, matrix_ev2, matrix_ev3, &
113 matrix_ev4, matrix_pe1p, matrix_rev, matrix_se, matrix_sinv
115 CALL timeset(routinen, handle)
168 CALL kintegral(n, ev0t, tt, e)
183 CALL parallel_gemm(
"N",
"T", n, n, n, 1.0_dp, matrix_rev, matrix_aux, 0.0_dp, matrix_t)
191 aa(i) = sqrt((
c_light_au**2 + e(i))/(2.0_dp*e(i)))
201 CALL parallel_gemm(
"T",
"N", n, n, n, 1.0_dp, matrix_eig, matrix_v, 0.0_dp, matrix_aux)
202 CALL parallel_gemm(
"N",
"N", n, n, n, 1.0_dp, matrix_aux, matrix_eig, 0.0_dp, matrix_v)
210 CALL parallel_gemm(
"T",
"N", n, n, n, 1.0_dp, matrix_eig, matrix_pvp, 0.0_dp, matrix_aux)
211 CALL parallel_gemm(
"N",
"N", n, n, n, 1.0_dp, matrix_aux, matrix_eig, 0.0_dp, matrix_pvp)
217 CALL even1(matrix_ev1, matrix_v, matrix_pvp, aa, rr, matrix_aux, matrix_aux2)
223 CALL even2c(n, matrix_ev2, matrix_v, matrix_pvp, aa, rr, tt, e, matrix_aux)
229 IF (dkh_order .GE. 3)
THEN
230 CALL peven1p(n, matrix_pe1p, matrix_v, matrix_pvp, matrix_aux, matrix_aux2, aa, rr, tt)
231 CALL even3b(n, matrix_ev3, matrix_ev1, matrix_pe1p, matrix_v, matrix_pvp, aa, rr, tt, e, matrix_aux)
237 CALL parallel_gemm(
"N",
"N", n, n, n, 1.0_dp, matrix_rev, matrix_ev3, 0.0_dp, matrix_aux)
238 CALL parallel_gemm(
"N",
"T", n, n, n, 1.0_dp, matrix_aux, matrix_rev, 0.0_dp, matrix_ev3)
244 IF (dkh_order .GE. 4)
THEN
245 cpabort(
"DKH order greater than 3 not yet available")
263 CALL parallel_gemm(
"N",
"N", n, n, n, 1.0_dp, matrix_rev, matrix_ev1, 0.0_dp, matrix_aux)
264 CALL parallel_gemm(
"N",
"T", n, n, n, 1.0_dp, matrix_aux, matrix_rev, 0.0_dp, matrix_ev1)
270 CALL parallel_gemm(
"N",
"N", n, n, n, 1.0_dp, matrix_rev, matrix_ev2, 0.0_dp, matrix_aux)
271 CALL parallel_gemm(
"N",
"T", n, n, n, 1.0_dp, matrix_aux, matrix_rev, 0.0_dp, matrix_ev2)
280 IF (dkh_order .GE. 3)
THEN
282 IF (dkh_order .GE. 4)
THEN
301 DEALLOCATE (ev0t, e, aa, rr, tt)
303 CALL timestop(handle)
305 END SUBROUTINE dkh_full_transformation
314 SUBROUTINE kintegral(n, ev0t, tt, e)
315 INTEGER,
INTENT(IN) :: n
316 REAL(KIND=
dp),
DIMENSION(:),
INTENT(OUT) :: ev0t
317 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN) :: tt
318 REAL(KIND=
dp),
DIMENSION(:),
INTENT(OUT) :: e
321 REAL(KIND=
dp) :: con, con2, prea, ratio, tv1, tv2, tv3, &
329 IF (tt(i) .LT. 0.0_dp)
THEN
330 WRITE (*, *)
' dkh_main.F | tt(', i,
') = ', tt(i)
338 IF (ratio .LE. 0.02_dp)
THEN
340 tv2 = -tv1*tt(i)*prea*0.5_dp
341 tv3 = -tv2*tt(i)*prea
342 tv4 = -tv3*tt(i)*prea*1.25_dp
343 ev0t(i) = tv1 + tv2 + tv3 + tv4
345 ev0t(i) = con*(sqrt(1.0_dp + con2*tt(i)) - 1.0_dp)
350 END SUBROUTINE kintegral
362 SUBROUTINE even1(matrix_ev1, matrix_v, matrix_pvp, aa, rr, matrix_aux, matrix_aux2)
363 TYPE(cp_fm_type),
INTENT(IN) :: matrix_ev1, matrix_v, matrix_pVp
364 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN) :: aa, rr
365 TYPE(cp_fm_type),
INTENT(IN) :: matrix_aux, matrix_aux2
395 SUBROUTINE peven1p(n, matrix_pe1p, matrix_v, matrix_pvp, matrix_aux, matrix_aux2, aa, rr, tt)
397 INTEGER,
INTENT(IN) :: n
398 TYPE(cp_fm_type),
INTENT(IN) :: matrix_pe1p, matrix_v, matrix_pvp, &
399 matrix_aux, matrix_aux2
400 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN) :: aa, rr, tt
402 INTEGER :: i, nrow_local
403 INTEGER,
DIMENSION(:),
POINTER :: row_indices
404 REAL(KIND=
dp),
DIMENSION(n) :: vec_ar, vec_arrt
405 TYPE(cp_blacs_env_type),
POINTER :: context
406 TYPE(cp_fm_struct_type),
POINTER :: vec_full
407 TYPE(cp_fm_type) :: vec_a
410 vec_ar(i) = aa(i)*rr(i)
411 vec_arrt(i) = vec_ar(i)*rr(i)*tt(i)
423 row_indices=row_indices)
426 vec_a%local_data(i, 1) = vec_arrt(row_indices(i))
429 CALL cp_fm_syrk(
'U',
'N', 1, 1.0_dp, vec_a, 1, 1, 0.0_dp, matrix_aux)
434 vec_a%local_data(i, 1) = vec_ar(row_indices(i))
437 CALL cp_fm_syrk(
'U',
'N', 1, 1.0_dp, vec_a, 1, 1, 0.0_dp, matrix_aux)
446 END SUBROUTINE peven1p
460 SUBROUTINE even2c(n, matrix_ev2, matrix_v, matrix_pVp, aa, rr, tt, e, matrix_aux)
477 INTEGER,
INTENT(IN) :: n
478 TYPE(cp_fm_type),
INTENT(IN) :: matrix_ev2, matrix_v, matrix_pVp
479 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN) :: aa, rr, tt, e
480 TYPE(cp_fm_type),
INTENT(IN) :: matrix_aux
482 TYPE(cp_blacs_env_type),
POINTER :: context
483 TYPE(cp_fm_struct_type),
POINTER :: matrix_full
484 TYPE(cp_fm_type) :: matrix_apVpa, matrix_apVVpa, &
485 matrix_aux2, matrix_ava, matrix_avva
511 CALL mat_axa(matrix_v, matrix_ava, n, aa, matrix_aux)
515 CALL mat_arxra(matrix_pvp, matrix_apvpa, n, aa, rr, matrix_aux)
519 CALL mat_1_over_h(matrix_v, matrix_avva, e, matrix_aux)
521 CALL mat_axa(matrix_aux2, matrix_avva, n, aa, matrix_aux)
525 CALL mat_1_over_h(matrix_pvp, matrix_apvvpa, e, matrix_aux)
527 CALL mat_arxra(matrix_aux2, matrix_apvvpa, n, aa, rr, matrix_aux)
531 CALL parallel_gemm(
"N",
"N", n, n, n, -1.0_dp, matrix_apvvpa, matrix_ava, 0.0_dp, matrix_aux2)
532 CALL mat_muld(matrix_aux2, matrix_apvvpa, matrix_apvpa, n, 1.0_dp, 1.0_dp, tt, rr, matrix_aux)
533 CALL mat_mulm(matrix_aux2, matrix_avva, matrix_ava, n, 1.0_dp, 1.0_dp, tt, rr, matrix_aux)
534 CALL parallel_gemm(
"N",
"N", n, n, n, -1.0_dp, matrix_avva, matrix_apvpa, 1.0_dp, matrix_aux2)
538 CALL parallel_gemm(
"N",
"N", n, n, n, 1.0_dp, matrix_apvpa, matrix_avva, 0.0_dp, matrix_ev2)
539 CALL mat_muld(matrix_ev2, matrix_apvpa, matrix_apvvpa, n, -1.0_dp, 1.0_dp, tt, rr, matrix_aux)
540 CALL mat_mulm(matrix_ev2, matrix_ava, matrix_avva, n, -1.0_dp, 1.0_dp, tt, rr, matrix_aux)
541 CALL parallel_gemm(
"N",
"N", n, n, n, 1.0_dp, matrix_ava, matrix_apvvpa, 1.0_dp, matrix_ev2)
564 END SUBROUTINE even2c
580 SUBROUTINE even3b(n, matrix_ev3, matrix_ev1, matrix_pe1p, matrix_v, matrix_pVp, aa, rr, tt, e, matrix_aux)
619 INTEGER,
INTENT(IN) :: n
620 TYPE(cp_fm_type),
INTENT(IN) :: matrix_ev3, matrix_ev1, matrix_pe1p, &
622 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN) :: aa, rr, tt, e
623 TYPE(cp_fm_type),
INTENT(IN) :: matrix_aux
625 TYPE(cp_blacs_env_type),
POINTER :: context
626 TYPE(cp_fm_struct_type),
POINTER :: matrix_full
627 TYPE(cp_fm_type) :: matrix_apVVpa, matrix_aux2, matrix_avva, &
628 matrix_w1e1w1, matrix_w1w1
651 CALL mat_1_over_h(matrix_v, matrix_avva, e, matrix_aux)
653 CALL mat_axa(matrix_aux2, matrix_avva, n, aa, matrix_aux)
657 CALL mat_1_over_h(matrix_pvp, matrix_apvvpa, e, matrix_aux)
659 CALL mat_arxra(matrix_aux2, matrix_apvvpa, n, aa, rr, matrix_aux)
663 CALL parallel_gemm(
"N",
"N", n, n, n, 1.0_dp, matrix_apvvpa, matrix_avva, 0.0_dp, matrix_w1w1)
664 CALL mat_muld(matrix_w1w1, matrix_apvvpa, matrix_apvvpa, n, -1.0_dp, 1.0_dp, tt, rr, matrix_aux2)
665 CALL mat_mulm(matrix_w1w1, matrix_avva, matrix_avva, n, -1.0_dp, 1.0_dp, tt, rr, matrix_aux2)
666 CALL parallel_gemm(
"N",
"N", n, n, n, 1.0_dp, matrix_avva, matrix_apvvpa, 1.0_dp, matrix_w1w1)
670 CALL mat_muld(matrix_aux, matrix_apvvpa, matrix_pe1p, n, 1.0_dp, 0.0_dp, tt, rr, matrix_aux2)
671 CALL parallel_gemm(
"N",
"N", n, n, n, 1.0_dp, matrix_avva, matrix_pe1p, 0.0_dp, matrix_aux2)
672 CALL parallel_gemm(
"N",
"N", n, n, n, 1.0_dp, matrix_aux, matrix_avva, 0.0_dp, matrix_w1e1w1)
673 CALL mat_muld(matrix_w1e1w1, matrix_aux, matrix_apvvpa, n, -1.0_dp, 1.0_dp, tt, rr, matrix_ev3)
674 CALL parallel_gemm(
"N",
"N", n, n, n, -1.0_dp, matrix_aux2, matrix_avva, 1.0_dp, matrix_w1e1w1)
675 CALL mat_muld(matrix_w1e1w1, matrix_aux2, matrix_apvvpa, n, 1.0_dp, 1.0_dp, tt, rr, matrix_ev3)
681 CALL parallel_gemm(
"N",
"N", n, n, n, 0.5_dp, matrix_w1w1, matrix_ev1, 0.0_dp, matrix_ev3)
682 CALL parallel_gemm(
"N",
"N", n, n, n, 0.5_dp, matrix_ev1, matrix_w1w1, 1.0_dp, matrix_ev3)
700 END SUBROUTINE even3b
713 SUBROUTINE even4a(n, ev4, e1, pe1p, vv, gg)
762 INTEGER,
INTENT(IN) :: n
763 REAL(KIND=
dp),
DIMENSION(:, :),
INTENT(OUT) :: ev4
764 REAL(KIND=
dp),
DIMENSION(:, :),
INTENT(IN) :: e1, pe1p, vv, gg
766 REAL(KIND=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: o1w1, pvp, pvph, scr_1, scr_2, scr_3, &
767 scr_4, scrh_1, scrh_2, scrh_3, scrh_4, &
768 sum_1, sum_2, v, vh, w1o1, w1w1
777 ALLOCATE (pvph(n, n))
782 v(1:n, 1:n) = vv(1:n, 1:n)
783 vh(1:n, 1:n) = vv(1:n, 1:n)
784 pvp(1:n, 1:n) = gg(1:n, 1:n)
785 pvph(1:n, 1:n) = gg(1:n, 1:n)
807 ALLOCATE (w1w1(n, n))
809 ALLOCATE (w1o1(n, n))
811 ALLOCATE (o1w1(n, n))
813 ALLOCATE (sum_1(n, n))
815 ALLOCATE (sum_2(n, n))
817 ALLOCATE (scr_1(n, n))
819 ALLOCATE (scr_2(n, n))
821 ALLOCATE (scr_3(n, n))
823 ALLOCATE (scr_4(n, n))
825 ALLOCATE (scrh_1(n, n))
827 ALLOCATE (scrh_2(n, n))
829 ALLOCATE (scrh_3(n, n))
831 ALLOCATE (scrh_4(n, n))
835 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pvph, n, vh, n, 0.0_dp, w1w1, n)
838 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, vh, n, pvph, n, 1.0_dp, w1w1, n)
841 CALL dgemm(
"N",
"N", n, n, n, -1.0_dp, pvph, n, v, n, 0.0_dp, w1o1, n)
844 CALL dgemm(
"N",
"N", n, n, n, -1.0_dp, vh, n, pvp, n, 1.0_dp, w1o1, n)
846 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pvp, n, vh, n, 0.0_dp, o1w1, n)
849 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, v, n, pvph, n, 1.0_dp, o1w1, n)
857 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, vh, n, e1, n, 0.0_dp, scr_1, n)
858 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pvph, n, e1, n, 0.0_dp, scr_2, n)
859 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pe1p, n, vh, n, 0.0_dp, scr_3, n)
864 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, vh, n, pe1p, n, 0.0_dp, scrh_2, n)
866 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, e1, n, pvph, n, 0.0_dp, scrh_3, n)
868 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, e1, n, vh, n, 0.0_dp, scrh_4, n)
873 CALL dgemm(
"N",
"N", n, n, n, 0.5_dp, scrh_1, n, scr_1, n, 0.0_dp, sum_1, n)
875 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scrh_2, n, scr_1, n, 1.0_dp, sum_1, n)
877 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scrh_3, n, scr_1, n, 1.0_dp, sum_1, n)
880 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scrh_4, n, scr_2, n, 1.0_dp, sum_1, n)
890 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scrh_4, n, scr_3, n, 1.0_dp, sum_1, n)
891 CALL dgemm(
"N",
"N", n, n, n, 0.5_dp, scrh_4, n, scr_4, n, 1.0_dp, sum_1, n)
896 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, vh, n, pe1p, n, 0.0_dp, scr_2, n)
897 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, e1, n, pvph, n, 0.0_dp, scr_3, n)
898 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, e1, n, vh, n, 0.0_dp, scr_4, n)
900 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, vh, n, e1, n, 0.0_dp, scrh_1, n)
902 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pvph, n, e1, n, 0.0_dp, scrh_2, n)
904 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pe1p, n, vh, n, 0.0_dp, scr_3, n)
911 CALL dgemm(
"N",
"N", n, n, n, 0.5_dp, scr_1, n, scrh_1, n, 0.0_dp, sum_1, n)
913 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scr_2, n, scrh_1, n, 1.0_dp, sum_1, n)
922 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scr_3, n, scrh_1, n, 0.0_dp, sum_1, n)
925 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scr_4, n, scrh_2, n, 1.0_dp, sum_1, n)
928 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scr_4, n, scrh_3, n, 1.0_dp, sum_1, n)
929 CALL dgemm(
"N",
"N", n, n, n, 0.5_dp, scr_4, n, scrh_4, n, 1.0_dp, sum_1, n)
937 CALL dgemm(
"N",
"N", n, n, n, 0.125_dp, w1w1, n, w1o1, n, 0.0_dp, sum_2, n)
938 CALL dgemm(
"N",
"N", n, n, n, -0.375_dp, w1w1, n, o1w1, n, 1.0_dp, sum_2, n)
939 CALL dgemm(
"N",
"N", n, n, n, 0.375_dp, w1o1, n, w1w1, n, 1.0_dp, sum_2, n)
940 CALL dgemm(
"N",
"N", n, n, n, -0.125_dp, o1w1, n, w1w1, n, 1.0_dp, sum_2, n)
946 CALL mat_add(ev4, 1.0_dp, sum_1, 1.0_dp, sum_2, n)
952 DEALLOCATE (v, pvp, vh, pvph, w1w1, w1o1, o1w1, sum_1, sum_2)
953 DEALLOCATE (scr_1, scr_2, scr_3, scr_4, scrh_1, scrh_2, scrh_3, scrh_4)
958 END SUBROUTINE even4a
983 SUBROUTINE mat_1_over_h(matrix_p, matrix_pp, e, matrix_aux)
994 TYPE(cp_fm_type),
INTENT(IN) :: matrix_p, matrix_pp
995 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN) :: e
996 TYPE(cp_fm_type),
INTENT(IN) :: matrix_aux
998 INTEGER :: i, j, ncol_local, nrow_local
999 INTEGER,
DIMENSION(:),
POINTER :: col_indices, row_indices
1001 CALL cp_fm_get_info(matrix_aux, nrow_local=nrow_local, ncol_local=ncol_local, &
1002 row_indices=row_indices, col_indices=col_indices)
1004 DO j = 1, ncol_local
1005 DO i = 1, nrow_local
1006 matrix_aux%local_data(i, j) = 1/(e(row_indices(i)) + e(col_indices(j)))
1012 END SUBROUTINE mat_1_over_h
1022 SUBROUTINE mat_axa(matrix_x, matrix_axa, n, a, matrix_aux)
1034 TYPE(cp_fm_type),
INTENT(IN) :: matrix_x, matrix_axa
1035 INTEGER,
INTENT(IN) :: n
1036 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN) :: a
1037 TYPE(cp_fm_type),
INTENT(IN) :: matrix_aux
1039 INTEGER :: i, nrow_local
1040 INTEGER,
DIMENSION(:),
POINTER :: row_indices
1041 TYPE(cp_blacs_env_type),
POINTER :: context
1042 TYPE(cp_fm_struct_type),
POINTER :: vec_full
1043 TYPE(cp_fm_type) :: vec_a
1054 row_indices=row_indices)
1056 DO i = 1, nrow_local
1057 vec_a%local_data(i, 1) = a(row_indices(i))
1060 CALL cp_fm_syrk(
'U',
'N', 1, 1.0_dp, vec_a, 1, 1, 0.0_dp, matrix_aux)
1073 END SUBROUTINE mat_axa
1084 SUBROUTINE mat_arxra(matrix_x, matrix_axa, n, a, r, matrix_aux)
1097 TYPE(cp_fm_type),
INTENT(IN) :: matrix_x, matrix_axa
1098 INTEGER,
INTENT(IN) :: n
1099 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN) :: a, r
1100 TYPE(cp_fm_type),
INTENT(IN) :: matrix_aux
1102 INTEGER :: i, nrow_local
1103 INTEGER,
DIMENSION(:),
POINTER :: row_indices
1104 TYPE(cp_blacs_env_type),
POINTER :: context
1105 TYPE(cp_fm_struct_type),
POINTER :: vec_full
1106 TYPE(cp_fm_type) :: vec_a
1115 row_indices=row_indices)
1119 DO i = 1, nrow_local
1120 vec_a%local_data(i, 1) = a(row_indices(i))*r(row_indices(i))
1123 CALL cp_fm_syrk(
'U',
'N', 1, 1.0_dp, vec_a, 1, 1, 0.0_dp, matrix_aux)
1130 END SUBROUTINE mat_arxra
1144 SUBROUTINE mat_mulm(matrix_p, matrix_q, matrix_r, n, alpha, beta, t, rr, matrix_aux)
1163 TYPE(cp_fm_type),
INTENT(IN) :: matrix_p, matrix_q, matrix_r
1164 INTEGER,
INTENT(IN) :: n
1165 REAL(KIND=
dp),
INTENT(IN) :: alpha, beta
1166 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN) :: t, rr
1167 TYPE(cp_fm_type),
INTENT(IN) :: matrix_aux
1170 REAL(KIND=
dp),
DIMENSION(n) :: vec
1175 vec(i) = 2.0_dp*t(i)*rr(i)*rr(i)
1179 CALL parallel_gemm(
"N",
"N", n, n, n, alpha, matrix_aux, matrix_r, beta, matrix_p)
1181 END SUBROUTINE mat_mulm
1195 SUBROUTINE mat_muld(matrix_p, matrix_q, matrix_r, n, alpha, beta, t, rr, matrix_aux)
1214 TYPE(cp_fm_type),
INTENT(IN) :: matrix_p, matrix_q, matrix_r
1215 INTEGER,
INTENT(IN) :: n
1216 REAL(KIND=
dp),
INTENT(IN) :: alpha, beta
1217 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN) :: t, rr
1218 TYPE(cp_fm_type),
INTENT(IN) :: matrix_aux
1221 REAL(KIND=
dp),
DIMENSION(n) :: vec
1226 vec(i) = 0.5_dp/(t(i)*rr(i)*rr(i))
1231 CALL parallel_gemm(
"N",
"N", n, n, n, alpha, matrix_aux, matrix_r, beta, matrix_p)
1233 END SUBROUTINE mat_muld
1282 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(INOUT) :: s, v, h, pvp
1283 INTEGER,
INTENT(IN) :: n, dkh_order
1286 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: aa, e, ev0t, rr, tt
1287 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: aux, eig, ev1, ev1t, ev2, ev2t, ev3, &
1288 ev3t, ev4, ev4t, ove, pev1tp, pvpt, &
1291 IF (dkh_order < 0)
RETURN
1299 ALLOCATE (eig(n, n))
1300 ALLOCATE (sinv(n, n))
1301 ALLOCATE (revt(n, n))
1302 ALLOCATE (aux(n, n))
1303 ALLOCATE (ove(n, n))
1309 ALLOCATE (ev1t(n, n))
1310 ALLOCATE (ev2t(n, n))
1311 ALLOCATE (ev3t(n, n))
1312 ALLOCATE (ev4t(n, n))
1314 ALLOCATE (pvpt(n, n))
1315 ALLOCATE (pev1tp(n, n))
1316 ALLOCATE (ev1(n, n))
1317 ALLOCATE (ev2(n, n))
1318 ALLOCATE (ev3(n, n))
1319 ALLOCATE (ev4(n, n))
1325 CALL sog(n, s, sinv)
1331 CALL dkh_diag(h, n, eig, tt, sinv, aux, 0)
1337 CALL kintegral_a(n, ev0t, tt, e)
1343 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, sinv, n, eig, n, 0.0_dp, aux, n)
1344 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, s, n, aux, n, 0.0_dp, revt, n)
1354 h(i, j) = h(i, j) + revt(i, k)*revt(j, k)*ev0t(k)
1365 aa(i) = sqrt((
c_light_au**2 + e(i))/(2.0_dp*e(i)))
1373 CALL trsm(v, sinv, ove, n, aux)
1374 CALL trsm(ove, eig, vt, n, aux)
1380 CALL trsm(pvp, sinv, ove, n, aux)
1381 CALL trsm(ove, eig, pvpt, n, aux)
1387 IF (dkh_order .GE. 1)
THEN
1388 CALL even1_a(n, ev1t, vt, pvpt, aa, rr)
1394 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, revt, n, ev1t, n, 0.0_dp, aux, n)
1395 CALL dgemm(
"N",
"T", n, n, n, 1.0_dp, aux, n, revt, n, 0.0_dp, ev1, n)
1402 IF (dkh_order .GE. 2)
THEN
1403 CALL even2c_a(n, ev2t, vt, pvpt, aa, rr, tt, e)
1410 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, revt, n, ev2t, n, 0.0_dp, aux, n)
1411 CALL dgemm(
"N",
"T", n, n, n, 1.0_dp, aux, n, revt, n, 0.0_dp, ev2, n)
1418 IF (dkh_order .GE. 3)
THEN
1419 CALL peven1p_a(n, pev1tp, vt, pvpt, aa, rr, tt)
1420 CALL even3b_a(n, ev3t, ev1t, pev1tp, vt, pvpt, aa, rr, tt, e)
1426 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, revt, n, ev3t, n, 0.0_dp, aux, n)
1427 CALL dgemm(
"N",
"T", n, n, n, 1.0_dp, aux, n, revt, n, 0.0_dp, ev3, n)
1433 IF (dkh_order .GE. 4)
THEN
1434 CALL even4a_a(n, ev4t, ev1t, pev1tp, vt, pvpt, aa, rr, tt, e)
1440 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, revt, n, ev4t, n, 0.0_dp, aux, n)
1441 CALL dgemm(
"N",
"T", n, n, n, 1.0_dp, aux, n, revt, n, 0.0_dp, ev4, n)
1445 IF (dkh_order .GE. 4)
THEN
1452 IF (dkh_order .GE. 1)
THEN
1453 CALL mat_add2(v, 0.0_dp, 1.0_dp, ev1, n)
1455 IF (dkh_order .GE. 2)
THEN
1456 CALL mat_add2(v, 1.0_dp, 1.0_dp, ev2, n)
1458 IF (dkh_order .GE. 3)
THEN
1459 CALL mat_add2(v, 1.0_dp, 1.0_dp, ev3, n)
1461 IF (dkh_order .GE. 4)
THEN
1462 CALL mat_add2(v, 1.0_dp, 1.0_dp, ev4, n)
1467 DEALLOCATE (eig, sinv, revt, ove, aux, vt, pvpt, ev1, ev2, ev3, ev4, ev1t, ev2t, ev3t, ev4t, pev1tp)
1468 DEALLOCATE (ev0t, e, aa, rr, tt)
1479 SUBROUTINE kintegral_a(n, ev0t, tt, e)
1481 INTEGER,
INTENT(IN) :: n
1482 REAL(kind=
dp),
DIMENSION(:),
INTENT(OUT) :: ev0t
1483 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: tt
1484 REAL(kind=
dp),
DIMENSION(:),
INTENT(OUT) :: e
1487 REAL(kind=
dp) :: con, con2, prea, ratio, tv1, tv2, tv3, &
1491 IF (tt(i) .LT. 0.0_dp)
THEN
1492 WRITE (*, *)
' dkh_main.F | tt(', i,
') = ', tt(i)
1506 IF (ratio .LE. 0.02_dp)
THEN
1508 tv2 = -tv1*tt(i)*prea/2.0_dp
1509 tv3 = -tv2*tt(i)*prea
1510 tv4 = -tv3*tt(i)*prea*1.25_dp
1511 ev0t(i) = tv1 + tv2 + tv3 + tv4
1513 ev0t(i) = con*(sqrt(1.0_dp + con2*tt(i)) - 1.0_dp)
1515 e(i) = ev0t(i) + con
1518 END SUBROUTINE kintegral_a
1529 SUBROUTINE even1_a(n, ev1t, vt, pvpt, aa, rr)
1544 INTEGER,
INTENT(IN) :: n
1545 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(OUT) :: ev1t
1546 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: vt, pvpt
1547 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: aa, rr
1555 ev1t(i, j) = vt(i, j)*aa(i)*aa(j) + pvpt(i, j)*aa(i)*rr(i)*aa(j)*rr(j)
1556 ev1t(j, i) = ev1t(i, j)
1560 END SUBROUTINE even1_a
1572 SUBROUTINE peven1p_a(n, pev1tp, vt, pvpt, aa, rr, tt)
1587 INTEGER,
INTENT(IN) :: n
1588 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(OUT) :: pev1tp
1589 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: vt, pvpt
1590 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: aa, rr, tt
1598 pev1tp(i, j) = 4.0_dp*vt(i, j)*aa(i)*aa(j)*rr(i)*rr(i)*rr(j)*rr(j)*tt(i)*tt(j) + &
1599 pvpt(i, j)*aa(i)*rr(i)*aa(j)*rr(j)
1600 pev1tp(j, i) = pev1tp(i, j)
1604 END SUBROUTINE peven1p_a
1617 SUBROUTINE even2c_a(n, ev2, vv, gg, aa, rr, tt, e)
1653 INTEGER,
INTENT(IN) :: n
1654 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(OUT) :: ev2
1655 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: vv, gg
1656 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: aa, rr, tt, e
1658 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: o1w1, pvp, pvph, v, vh, w1o1
1665 ALLOCATE (pvp(n, n))
1667 ALLOCATE (pvph(n, n))
1672 v(1:n, 1:n) = vv(1:n, 1:n)
1673 vh(1:n, 1:n) = vv(1:n, 1:n)
1674 pvp(1:n, 1:n) = gg(1:n, 1:n)
1675 pvph(1:n, 1:n) = gg(1:n, 1:n)
1680 CALL mat_axa_a(v, n, aa)
1684 CALL mat_arxra_a(pvp, n, aa, rr)
1688 CALL mat_1_over_h_a(vh, n, e)
1689 CALL mat_axa_a(vh, n, aa)
1693 CALL mat_1_over_h_a(pvph, n, e)
1694 CALL mat_arxra_a(pvph, n, aa, rr)
1697 ALLOCATE (w1o1(n, n))
1698 ALLOCATE (o1w1(n, n))
1703 CALL dgemm(
"N",
"N", n, n, n, -1.0_dp, pvph, n, v, n, 0.0_dp, w1o1, n)
1704 CALL mat_muld_a(w1o1, pvph, pvp, n, 1.0_dp, 1.0_dp, tt, rr)
1705 CALL mat_mulm_a(w1o1, vh, v, n, 1.0_dp, 1.0_dp, tt, rr)
1706 CALL dgemm(
"N",
"N", n, n, n, -1.0_dp, vh, n, pvp, n, 1.0_dp, w1o1, n)
1708 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pvp, n, vh, n, 0.0_dp, o1w1, n)
1709 CALL mat_muld_a(o1w1, pvp, pvph, n, -1.0_dp, 1.0_dp, tt, rr)
1710 CALL mat_mulm_a(o1w1, v, vh, n, -1.0_dp, 1.0_dp, tt, rr)
1711 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, v, n, pvph, n, 1.0_dp, o1w1, n)
1718 CALL mat_add(ev2, 0.5_dp, w1o1, -0.5_dp, o1w1, n)
1724 DEALLOCATE (v, vh, pvp, pvph, w1o1, o1w1)
1729 END SUBROUTINE even2c_a
1744 SUBROUTINE even3b_a(n, ev3, e1, pe1p, vv, gg, aa, rr, tt, e)
1783 INTEGER,
INTENT(IN) :: n
1784 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(OUT) :: ev3
1785 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: e1, pe1p, vv, gg
1786 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: aa, rr, tt, e
1788 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: pvph, scr_1, scr_2, vh, w1e1w1, w1w1
1795 ALLOCATE (pvph(n, n))
1798 vh(1:n, 1:n) = vv(1:n, 1:n)
1799 pvph(1:n, 1:n) = gg(1:n, 1:n)
1805 CALL mat_1_over_h_a(vh, n, e)
1806 CALL mat_axa_a(vh, n, aa)
1810 CALL mat_1_over_h_a(pvph, n, e)
1811 CALL mat_arxra_a(pvph, n, aa, rr)
1814 ALLOCATE (w1w1(n, n))
1815 ALLOCATE (w1e1w1(n, n))
1816 ALLOCATE (scr_1(n, n))
1817 ALLOCATE (scr_2(n, n))
1824 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pvph, n, vh, n, 0.0_dp, w1w1, n)
1825 CALL mat_muld_a(w1w1, pvph, pvph, n, -1.0_dp, 1.0_dp, tt, rr)
1826 CALL mat_mulm_a(w1w1, vh, vh, n, -1.0_dp, 1.0_dp, tt, rr)
1827 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, vh, n, pvph, n, 1.0_dp, w1w1, n)
1830 CALL mat_muld_a(scr_1, pvph, pe1p, n, 1.0_dp, 0.0_dp, tt, rr)
1831 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, vh, n, pe1p, n, 0.0_dp, scr_2, n)
1832 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, scr_1, n, vh, n, 0.0_dp, w1e1w1, n)
1833 CALL mat_muld_a(w1e1w1, scr_1, pvph, n, -1.0_dp, 1.0_dp, tt, rr)
1834 CALL dgemm(
"N",
"N", n, n, n, -1.0_dp, scr_2, n, vh, n, 1.0_dp, w1e1w1, n)
1835 CALL mat_muld_a(w1e1w1, scr_2, pvph, n, 1.0_dp, 1.0_dp, tt, rr)
1841 CALL dgemm(
"N",
"N", n, n, n, 0.5_dp, w1w1, n, e1, n, 0.0_dp, ev3, n)
1842 CALL dgemm(
"N",
"N", n, n, n, 0.5_dp, e1, n, w1w1, n, 1.0_dp, ev3, n)
1843 CALL mat_add2(ev3, 1.0_dp, -1.0_dp, w1e1w1, n)
1849 DEALLOCATE (vh, pvph, w1w1, w1e1w1, scr_1, scr_2)
1854 END SUBROUTINE even3b_a
1869 SUBROUTINE even4a_a(n, ev4, e1, pe1p, vv, gg, aa, rr, tt, e)
1917 INTEGER,
INTENT(IN) :: n
1918 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(OUT) :: ev4
1919 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: e1, pe1p, vv, gg
1920 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: aa, rr, tt, e
1922 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: o1w1, pvp, pvph, scr_1, scr_2, scr_3, &
1923 scr_4, scrh_1, scrh_2, scrh_3, scrh_4, &
1924 sum_1, sum_2, v, vh, w1o1, w1w1
1931 ALLOCATE (pvp(n, n))
1933 ALLOCATE (pvph(n, n))
1938 v(1:n, 1:n) = vv(1:n, 1:n)
1939 vh(1:n, 1:n) = vv(1:n, 1:n)
1940 pvp(1:n, 1:n) = gg(1:n, 1:n)
1941 pvph(1:n, 1:n) = gg(1:n, 1:n)
1946 CALL mat_axa_a(v, n, aa)
1950 CALL mat_arxra_a(pvp, n, aa, rr)
1954 CALL mat_1_over_h_a(vh, n, e)
1955 CALL mat_axa_a(vh, n, aa)
1959 CALL mat_1_over_h_a(pvph, n, e)
1960 CALL mat_arxra_a(pvph, n, aa, rr)
1963 ALLOCATE (w1w1(n, n))
1965 ALLOCATE (w1o1(n, n))
1967 ALLOCATE (o1w1(n, n))
1969 ALLOCATE (sum_1(n, n))
1971 ALLOCATE (sum_2(n, n))
1973 ALLOCATE (scr_1(n, n))
1975 ALLOCATE (scr_2(n, n))
1977 ALLOCATE (scr_3(n, n))
1979 ALLOCATE (scr_4(n, n))
1981 ALLOCATE (scrh_1(n, n))
1983 ALLOCATE (scrh_2(n, n))
1985 ALLOCATE (scrh_3(n, n))
1987 ALLOCATE (scrh_4(n, n))
1991 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pvph, n, vh, n, 0.0_dp, w1w1, n)
1992 CALL mat_muld_a(w1w1, pvph, pvph, n, -1.0_dp, 1.0_dp, tt, rr)
1993 CALL mat_mulm_a(w1w1, vh, vh, n, -1.0_dp, 1.0_dp, tt, rr)
1994 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, vh, n, pvph, n, 1.0_dp, w1w1, n)
1997 CALL dgemm(
"N",
"N", n, n, n, -1.0_dp, pvph, n, v, n, 0.0_dp, w1o1, n)
1998 CALL mat_muld_a(w1o1, pvph, pvp, n, 1.0_dp, 1.0_dp, tt, rr)
1999 CALL mat_mulm_a(w1o1, vh, v, n, 1.0_dp, 1.0_dp, tt, rr)
2000 CALL dgemm(
"N",
"N", n, n, n, -1.0_dp, vh, n, pvp, n, 1.0_dp, w1o1, n)
2002 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pvp, n, vh, n, 0.0_dp, o1w1, n)
2003 CALL mat_muld_a(o1w1, pvp, pvph, n, -1.0_dp, 1.0_dp, tt, rr)
2004 CALL mat_mulm_a(o1w1, v, vh, n, -1.0_dp, 1.0_dp, tt, rr)
2005 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, v, n, pvph, n, 1.0_dp, o1w1, n)
2013 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, vh, n, e1, n, 0.0_dp, scr_1, n)
2014 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pvph, n, e1, n, 0.0_dp, scr_2, n)
2015 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pe1p, n, vh, n, 0.0_dp, scr_3, n)
2016 CALL mat_muld_a(scr_4, pe1p, pvph, n, 1.0_dp, 0.0_dp, tt, rr)
2018 CALL mat_muld_a(scrh_1, pvph, pe1p, n, 1.0_dp, 0.0_dp, tt, rr)
2019 CALL mat_1_over_h_a(scrh_1, n, e)
2020 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, vh, n, pe1p, n, 0.0_dp, scrh_2, n)
2021 CALL mat_1_over_h_a(scrh_2, n, e)
2022 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, e1, n, pvph, n, 0.0_dp, scrh_3, n)
2023 CALL mat_1_over_h_a(scrh_3, n, e)
2024 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, e1, n, vh, n, 0.0_dp, scrh_4, n)
2025 CALL mat_1_over_h_a(scrh_4, n, e)
2029 CALL dgemm(
"N",
"N", n, n, n, 0.5_dp, scrh_1, n, scr_1, n, 0.0_dp, sum_1, n)
2030 CALL mat_muld_a(sum_1, scrh_1, scr_2, n, -0.5_dp, 1.0_dp, tt, rr)
2031 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scrh_2, n, scr_1, n, 1.0_dp, sum_1, n)
2032 CALL mat_muld_a(sum_1, scrh_2, scr_2, n, 0.5_dp, 1.0_dp, tt, rr)
2033 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scrh_3, n, scr_1, n, 1.0_dp, sum_1, n)
2034 CALL mat_muld_a(sum_1, scrh_3, scr_2, n, 0.5_dp, 1.0_dp, tt, rr)
2035 CALL mat_mulm_a(sum_1, scrh_4, scr_1, n, 0.5_dp, 1.0_dp, tt, rr)
2036 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scrh_4, n, scr_2, n, 1.0_dp, sum_1, n)
2040 CALL mat_muld_a(sum_1, scrh_1, scr_3, n, -0.5_dp, 1.0_dp, tt, rr)
2041 CALL mat_muld_a(sum_1, scrh_1, scr_4, n, 0.5_dp, 1.0_dp, tt, rr)
2042 CALL mat_muld_a(sum_1, scrh_2, scr_3, n, 0.5_dp, 1.0_dp, tt, rr)
2043 CALL mat_muld_a(sum_1, scrh_2, scr_4, n, -0.5_dp, 1.0_dp, tt, rr)
2044 CALL mat_muld_a(sum_1, scrh_3, scr_3, n, 0.5_dp, 1.0_dp, tt, rr)
2045 CALL mat_muld_a(sum_1, scrh_3, scr_4, n, -0.5_dp, 1.0_dp, tt, rr)
2046 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scrh_4, n, scr_3, n, 1.0_dp, sum_1, n)
2047 CALL dgemm(
"N",
"N", n, n, n, 0.5_dp, scrh_4, n, scr_4, n, 1.0_dp, sum_1, n)
2051 CALL mat_muld_a(scr_1, pvph, pe1p, n, 1.0_dp, 0.0_dp, tt, rr)
2052 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, vh, n, pe1p, n, 0.0_dp, scr_2, n)
2053 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, e1, n, pvph, n, 0.0_dp, scr_3, n)
2054 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, e1, n, vh, n, 0.0_dp, scr_4, n)
2056 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, vh, n, e1, n, 0.0_dp, scrh_1, n)
2057 CALL mat_1_over_h_a(scrh_1, n, e)
2058 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pvph, n, e1, n, 0.0_dp, scrh_2, n)
2059 CALL mat_1_over_h_a(scrh_2, n, e)
2060 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, pe1p, n, vh, n, 0.0_dp, scr_3, n)
2061 CALL mat_1_over_h_a(scrh_3, n, e)
2062 CALL mat_muld_a(scrh_4, pe1p, pvph, n, 1.0_dp, 0.0_dp, tt, rr)
2063 CALL mat_1_over_h_a(scrh_4, n, e)
2067 CALL dgemm(
"N",
"N", n, n, n, 0.5_dp, scr_1, n, scrh_1, n, 0.0_dp, sum_1, n)
2068 CALL mat_muld_a(sum_1, scr_1, scrh_2, n, -0.5_dp, 1.0_dp, tt, rr)
2069 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scr_2, n, scrh_1, n, 1.0_dp, sum_1, n)
2070 CALL mat_muld_a(sum_1, scr_2, scrh_2, n, 0.5_dp, 1.0_dp, tt, rr)
2071 CALL mat_muld_a(sum_1, scr_1, scrh_3, n, -0.5_dp, 1.0_dp, tt, rr)
2072 CALL mat_muld_a(sum_1, scr_1, scrh_4, n, 0.5_dp, 1.0_dp, tt, rr)
2073 CALL mat_muld_a(sum_1, scr_2, scrh_3, n, 0.5_dp, 1.0_dp, tt, rr)
2074 CALL mat_muld_a(sum_1, scr_2, scrh_4, n, -0.5_dp, 1.0_dp, tt, rr)
2078 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scr_3, n, scrh_1, n, 0.0_dp, sum_1, n)
2079 CALL mat_muld_a(sum_1, scr_3, scrh_2, n, 0.5_dp, 1.0_dp, tt, rr)
2080 CALL mat_mulm_a(sum_1, scr_4, scrh_1, n, 0.5_dp, 1.0_dp, tt, rr)
2081 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scr_4, n, scrh_2, n, 1.0_dp, sum_1, n)
2082 CALL mat_muld_a(sum_1, scr_3, scrh_3, n, 0.5_dp, 1.0_dp, tt, rr)
2083 CALL mat_muld_a(sum_1, scr_3, scrh_4, n, -0.5_dp, 1.0_dp, tt, rr)
2084 CALL dgemm(
"N",
"N", n, n, n, -0.5_dp, scr_4, n, scrh_3, n, 1.0_dp, sum_1, n)
2085 CALL dgemm(
"N",
"N", n, n, n, 0.5_dp, scr_4, n, scrh_4, n, 1.0_dp, sum_1, n)
2093 CALL dgemm(
"N",
"N", n, n, n, 0.125_dp, w1w1, n, w1o1, n, 0.0_dp, sum_2, n)
2094 CALL dgemm(
"N",
"N", n, n, n, -0.375_dp, w1w1, n, o1w1, n, 1.0_dp, sum_2, n)
2095 CALL dgemm(
"N",
"N", n, n, n, 0.375_dp, w1o1, n, w1w1, n, 1.0_dp, sum_2, n)
2096 CALL dgemm(
"N",
"N", n, n, n, -0.125_dp, o1w1, n, w1w1, n, 1.0_dp, sum_2, n)
2102 CALL mat_add(ev4, 1.0_dp, sum_1, 1.0_dp, sum_2, n)
2108 DEALLOCATE (v, pvp, vh, pvph, w1w1, w1o1, o1w1, sum_1, sum_2)
2109 DEALLOCATE (scr_1, scr_2, scr_3, scr_4, scrh_1, scrh_2, scrh_3, scrh_4)
2114 END SUBROUTINE even4a_a
2138 SUBROUTINE mat_1_over_h_a(p, n, e)
2150 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(INOUT) :: p
2151 INTEGER,
INTENT(IN) :: n
2152 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: e
2158 p(i, j) = p(i, j)/(e(i) + e(j))
2162 END SUBROUTINE mat_1_over_h_a
2170 SUBROUTINE mat_axa_a(p, n, a)
2182 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(INOUT) :: p
2183 INTEGER,
INTENT(IN) :: n
2184 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: a
2190 p(i, j) = p(i, j)*a(i)*a(j)
2194 END SUBROUTINE mat_axa_a
2203 SUBROUTINE mat_arxra_a(p, n, a, r)
2216 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(INOUT) :: p
2217 INTEGER,
INTENT(IN) :: n
2218 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: a, r
2224 p(i, j) = p(i, j)*a(i)*a(j)*r(i)*r(j)
2228 END SUBROUTINE mat_arxra_a
2241 SUBROUTINE mat_mulm_a(p, q, r, n, alpha, beta, t, rr)
2260 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: p, q, r
2261 INTEGER,
INTENT(IN) :: n
2262 REAL(kind=
dp),
INTENT(IN) :: alpha, beta
2263 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: t, rr
2266 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: qtemp
2268 ALLOCATE (qtemp(n, n))
2272 qtemp(i, j) = q(i, j)*2.0_dp*t(j)*rr(j)*rr(j)
2276 CALL dgemm(
"N",
"N", n, n, n, alpha, qtemp, n, r, n, beta, p, n)
2280 END SUBROUTINE mat_mulm_a
2293 SUBROUTINE mat_muld_a(p, q, r, n, alpha, beta, t, rr)
2312 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: p, q, r
2313 INTEGER,
INTENT(IN) :: n
2314 REAL(kind=
dp),
INTENT(IN) :: alpha, beta
2315 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: t, rr
2318 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: qtemp
2320 ALLOCATE (qtemp(n, n))
2324 qtemp(i, j) = q(i, j)*0.5_dp/(t(j)*rr(j)*rr(j))
2328 CALL dgemm(
"N",
"N", n, n, n, alpha, qtemp, n, r, n, beta, p, n)
2332 END SUBROUTINE mat_muld_a
2342 SUBROUTINE mat_add2(p, alpha, beta, r, n)
2362 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(INOUT) :: p
2363 REAL(kind=
dp),
INTENT(IN) :: alpha, beta
2364 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: r
2365 INTEGER,
INTENT(IN) :: n
2373 p(i, j) = alpha*p(i, j) + beta*r(i, j)
2377 END SUBROUTINE mat_add2
2388 SUBROUTINE mat_add(p, alpha, q, beta, r, n)
2408 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(OUT) :: p
2409 REAL(kind=
dp),
INTENT(IN) :: alpha
2410 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: q
2411 REAL(kind=
dp),
INTENT(IN) :: beta
2412 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: r
2413 INTEGER,
INTENT(IN) :: n
2421 p(i, j) = alpha*q(i, j) + beta*r(i, j)
2425 END SUBROUTINE mat_add
2435 SUBROUTINE trsm(W, B, C, N, H)
2437 REAL(kind=
dp),
DIMENSION(:, :) :: w, b, c
2439 REAL(kind=
dp),
DIMENSION(:, :) :: h
2441 INTEGER :: i, ij, j, k, l
2462 h(i, l) = b(k, i)*w(k, l) + h(i, l)
2472 c(i, j) = h(i, l)*b(l, j) + c(i, j)
2490 SUBROUTINE dkh_diag(matrix_t_pgf, n, eig, ew, matrix_sinv_pgf, aux, ic)
2492 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(INOUT) :: matrix_t_pgf
2493 INTEGER,
INTENT(IN) :: n
2494 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(INOUT) :: eig
2495 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: ew
2496 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: matrix_sinv_pgf
2497 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(INOUT) :: aux
2505 CALL dgemm(
"N",
"N", n, n, n, 1.0_dp, matrix_t_pgf, n, matrix_sinv_pgf, n, 0.0_dp, eig, n)
2509 CALL dgemm(
"T",
"N", n, n, n, 1.0_dp, matrix_sinv_pgf, n, eig, n, 0.0_dp, aux, n)
2513 CALL jacob2(aux, eig, ew, n, ic)
2515 END SUBROUTINE dkh_diag
2525 SUBROUTINE jacob2(sogt, eigv, eigw, n, ic)
2527 INTEGER,
INTENT(IN) :: n
2528 REAL(kind=
dp),
DIMENSION(n),
INTENT(OUT) :: eigw
2529 REAL(kind=
dp),
DIMENSION(n, n),
INTENT(OUT) :: eigv
2530 REAL(kind=
dp),
DIMENSION(n, n),
INTENT(INOUT) :: sogt
2531 INTEGER,
INTENT(IN) :: ic
2533 INTEGER :: i, il, im, ind, j, k, l, ll, m, mm
2534 REAL(kind=
dp) :: cost, cost2, ext_norm, sincs, sint, &
2535 sint2, thr, thr_min, tol, u1, x, xy, y
2539 u1 = real(n, kind=
dp)
2542 eigw(i) = sogt(i, i)
2547 ext_norm = ext_norm + sogt(i, j)*sogt(i, j)
2552 IF (ext_norm .GT. 0.0_dp)
THEN
2553 ext_norm = sqrt(2.0_dp*ext_norm)
2554 thr_min = ext_norm*tol/u1
2565 IF ((abs(sogt(m, l)) - thr) .GE. 0.0_dp)
THEN
2567 x = 0.5_dp*(eigw(l) - eigw(m))
2568 y = -sogt(m, l)/sqrt(sogt(m, l)*sogt(m, l) + x*x)
2569 IF (x .LT. 0.0_dp) y = -y
2571 IF (y .GT. 1.0_dp) y = 1.0_dp
2572 IF (y .LT. -1.0_dp) y = -1.0_dp
2574 sint = y/sqrt(2.0_dp*(1.0_dp + sqrt(xy)))
2576 cost2 = 1.0_dp - sint2
2581 IF ((i - m) .NE. 0)
THEN
2582 IF ((i - m) .LT. 0)
THEN
2589 IF ((i - l) .NE. 0)
THEN
2590 IF ((i - l) .LT. 0)
THEN
2597 x = sogt(il, ll)*cost - sogt(im, mm)*sint
2598 sogt(im, mm) = sogt(il, ll)*sint + sogt(im, mm)*cost
2603 x = eigv(i, l)*cost - eigv(i, m)*sint
2604 eigv(i, m) = eigv(i, l)*sint + eigv(i, m)*cost
2608 x = 2.0_dp*sogt(m, l)*sincs
2609 y = eigw(l)*cost2 + eigw(m)*sint2 - x
2610 x = eigw(l)*sint2 + eigw(m)*cost2 + x
2611 sogt(m, l) = (eigw(l) - eigw(m))*sincs + sogt(m, l)*(cost2 - sint2)
2615 IF ((m - n) .EQ. 0)
EXIT
2618 IF ((l - m + 1) .EQ. 0)
EXIT
2621 IF ((ind - 1) .NE. 0.0_dp)
EXIT
2624 IF ((thr - thr_min) .LE. 0.0_dp)
EXIT
2631 IF ((eigw(i) - eigw(j)) .GT. 0.0_dp)
THEN
2637 eigv(k, i) = eigv(k, j)
2646 END SUBROUTINE jacob2
2654 SUBROUTINE sog(n, matrix_s_pgf, matrix_sinv_pgf)
2657 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: matrix_s_pgf
2658 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(INOUT) :: matrix_sinv_pgf
2660 INTEGER :: i, j, jn, k
2661 REAL(kind=
dp) :: diag_s, row_sum, scalar
2662 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: a
2663 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: g
2678 diag_s = matrix_s_pgf(jn, jn)
2685 scalar = scalar + matrix_s_pgf(i, jn)*g(i, j)
2687 diag_s = diag_s - scalar*scalar
2694 row_sum = row_sum + a(k)*g(j, k)
2700 diag_s = 1.0_dp/sqrt(diag_s)
2702 g(i, jn) = g(i, jn)*diag_s
2708 matrix_sinv_pgf(j, i) = 0.0_dp
2709 matrix_sinv_pgf(i, j) = g(i, j)
static void dgemm(const char transa, const char transb, const int m, const int n, const int k, const double alpha, const double *a, const int lda, const double *b, const int ldb, const double beta, double *c, const int ldc)
Convenient wrapper to hide Fortran nature of dgemm_, swapping a and b.
methods related to the blacs parallel environment
basic linear algebra operations for full matrices
subroutine, public cp_fm_column_scale(matrixa, scaling)
scales column i of matrix a with scaling(i)
subroutine, public cp_fm_transpose(matrix, matrixt)
transposes a matrix matrixt = matrix ^ T
subroutine, public cp_fm_syrk(uplo, trans, k, alpha, matrix_a, ia, ja, beta, matrix_c)
performs a rank-k update of a symmetric matrix_c matrix_c = beta * matrix_c + alpha * matrix_a * tran...
subroutine, public cp_fm_schur_product(matrix_a, matrix_b, matrix_c)
computes the schur product of two matrices c_ij = a_ij * b_ij
subroutine, public cp_fm_scale_and_add(alpha, matrix_a, beta, matrix_b)
calc A <- alpha*A + beta*B optimized for alpha == 1.0 (just add beta*B) and beta == 0....
subroutine, public cp_fm_uplo_to_full(matrix, work, uplo)
given a triangular matrix according to uplo, computes the corresponding full matrix
subroutine, public cp_fm_triangular_multiply(triangular_matrix, matrix_b, side, transpose_tr, invert_tr, uplo_tr, unit_diag_tr, n_rows, n_cols, alpha)
multiplies in place by a triangular matrix: matrix_b = alpha op(triangular_matrix) matrix_b or (if si...
various cholesky decomposition related routines
subroutine, public cp_fm_cholesky_decompose(matrix, n, info_out)
used to replace a symmetric positive def. matrix M with its cholesky decomposition U: M = U^T * U,...
subroutine, public cp_fm_cholesky_reduce(matrix, matrixb, itype)
reduce a matrix pencil A,B to normal form B has to be cholesky decomposed with cp_fm_cholesky_decompo...
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
subroutine, public cp_fm_syevd(matrix, eigenvectors, eigenvalues, info)
Computes all eigenvalues and vectors of a real symmetric matrix significantly faster than syevx,...
represent the structure of a full matrix
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
subroutine, public cp_fm_struct_get(fmstruct, para_env, context, descriptor, ncol_block, nrow_block, nrow_global, ncol_global, first_p_pos, row_indices, col_indices, nrow_local, ncol_local, nrow_locals, ncol_locals, local_leading_dimension)
returns the values of various attributes of the matrix structure
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
represent a full matrix distributed on many processors
subroutine, public cp_fm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp)
creates a new full matrix with the given structure
subroutine, public dkh_atom_transformation(s, v, h, pvp, n, dkh_order)
...
Defines the basic variable types.
integer, parameter, public dp
basic linear algebra operations for full matrixes
Definition of physical constants:
real(kind=dp), parameter, public c_light_au
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.
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
keeps the information about the structure of a full matrix