42#include "./base/base_uses.f90"
48 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'kpoint_coulomb_2c'
55 REAL(KIND=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: block
56 END TYPE two_d_util_type
75 atomic_kind_set, size_lattice_sum, operator_type, ikp_start, ikp_end)
77 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_v_kp
79 CHARACTER(LEN=*),
INTENT(IN) :: basis_type
84 INTEGER :: size_lattice_sum, operator_type, &
87 CHARACTER(LEN=*),
PARAMETER :: routinen =
'build_2c_coulomb_matrix_kp'
92 CALL timeset(routinen, handle)
94 CALL allocate_tmp(matrix_v_l_tmp, matrix_v_kp, ikp_start)
96 CALL lattice_sum(matrix_v_kp, kpoints, basis_type, cell, particle_set, &
97 qs_kind_set, atomic_kind_set, size_lattice_sum, matrix_v_l_tmp, &
98 operator_type, ikp_start, ikp_end)
100 CALL deallocate_tmp(matrix_v_l_tmp)
102 CALL timestop(handle)
121 SUBROUTINE lattice_sum(matrix_v_kp, kpoints, basis_type, cell, particle_set, &
122 qs_kind_set, atomic_kind_set, size_lattice_sum, matrix_v_L_tmp, &
123 operator_type, ikp_start, ikp_end)
125 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_v_kp
127 CHARACTER(LEN=*),
INTENT(IN) :: basis_type
130 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
132 INTEGER :: size_lattice_sum
134 INTEGER :: operator_type, ikp_start, ikp_end
136 CHARACTER(LEN=*),
PARAMETER :: routinen =
'lattice_sum'
138 INTEGER :: factor, handle, handle2, i_block, i_x, i_x_inner, i_x_outer, ik, j_y, j_y_inner, &
139 j_y_outer, k_z, k_z_inner, k_z_outer, x_max, x_min, y_max, y_min, z_max, z_min
140 INTEGER,
DIMENSION(3) :: nkp_grid
141 REAL(kind=
dp) :: coskl, sinkl
142 REAL(kind=
dp),
DIMENSION(3) :: vec_l, vec_s
143 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat
144 TYPE(two_d_util_type),
ALLOCATABLE,
DIMENSION(:) :: blocks_v_l, blocks_v_l_store
145 TYPE(two_d_util_type),
ALLOCATABLE, &
146 DIMENSION(:, :, :) :: blocks_v_kp
148 CALL timeset(routinen, handle)
150 CALL get_factor_and_xyz_min_max(cell, kpoints, size_lattice_sum, factor, hmat, &
151 x_min, x_max, y_min, y_max, z_min, z_max, nkp_grid)
153 CALL allocate_blocks_v_kp(blocks_v_kp, matrix_v_kp, ikp_start, ikp_end)
154 CALL allocate_blocks_v_l(blocks_v_l, matrix_v_l_tmp)
155 CALL allocate_blocks_v_l(blocks_v_l_store, matrix_v_l_tmp)
157 DO i_x_inner = 0, 2*nkp_grid(1) - 1
158 DO j_y_inner = 0, 2*nkp_grid(2) - 1
159 DO k_z_inner = 0, 2*nkp_grid(3) - 1
161 DO i_x_outer = x_min, x_max + nkp_grid(1), 2*nkp_grid(1)
162 DO j_y_outer = y_min, y_max + nkp_grid(2), 2*nkp_grid(2)
163 DO k_z_outer = z_min, z_max + nkp_grid(3), 2*nkp_grid(3)
165 i_x = i_x_inner + i_x_outer
166 j_y = j_y_inner + j_y_outer
167 k_z = k_z_inner + k_z_outer
169 IF (i_x > x_max .OR. i_x < x_min .OR. &
170 j_y > y_max .OR. j_y < y_min .OR. &
171 k_z > z_max .OR. k_z < z_min) cycle
173 vec_s = [real(i_x,
dp), real(j_y,
dp), real(k_z,
dp)]
175 vec_l = matmul(hmat, vec_s)
178 CALL compute_v_transl(matrix_v_l_tmp, blocks_v_l, vec_l, particle_set, &
179 qs_kind_set, atomic_kind_set, basis_type, cell, &
185 DO i_block = 1,
SIZE(blocks_v_l)
186 blocks_v_l_store(i_block)%block(:, :) = blocks_v_l_store(i_block)%block(:, :) &
187 + blocks_v_l(i_block)%block(:, :)
195 CALL timeset(routinen//
"_R_to_k", handle2)
198 DO ik = ikp_start, ikp_end
201 coskl = cos(
twopi*dot_product(vec_s(1:3), kpoints%xkp(1:3, ik)))
202 sinkl = sin(
twopi*dot_product(vec_s(1:3), kpoints%xkp(1:3, ik)))
204 DO i_block = 1,
SIZE(blocks_v_l)
206 blocks_v_kp(ik, 1, i_block)%block(:, :) = blocks_v_kp(ik, 1, i_block)%block(:, :) &
207 + coskl*blocks_v_l_store(i_block)%block(:, :)
208 blocks_v_kp(ik, 2, i_block)%block(:, :) = blocks_v_kp(ik, 2, i_block)%block(:, :) &
209 + sinkl*blocks_v_l_store(i_block)%block(:, :)
215 DO i_block = 1,
SIZE(blocks_v_l)
217 blocks_v_l_store(i_block)%block(:, :) = 0.0_dp
221 CALL timestop(handle2)
227 CALL set_blocks_to_matrix_v_kp(matrix_v_kp, blocks_v_kp, ikp_start, ikp_end)
229 CALL deallocate_blocks_v_kp(blocks_v_kp)
230 CALL deallocate_blocks_v_l(blocks_v_l)
231 CALL deallocate_blocks_v_l(blocks_v_l_store)
233 CALL timestop(handle)
235 END SUBROUTINE lattice_sum
244 SUBROUTINE set_blocks_to_matrix_v_kp(matrix_v_kp, blocks_v_kp, ikp_start, ikp_end)
246 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_v_kp
247 TYPE(two_d_util_type),
ALLOCATABLE, &
248 DIMENSION(:, :, :) :: blocks_v_kp
249 INTEGER :: ikp_start, ikp_end
251 CHARACTER(LEN=*),
PARAMETER :: routinen =
'set_blocks_to_matrix_v_kp'
253 INTEGER :: col, handle, i_block, i_real_im, ik, row
254 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: data_block
257 CALL timeset(routinen, handle)
259 DO ik = ikp_start, ikp_end
271 data_block(:, :) = blocks_v_kp(ik, i_real_im, i_block)%block(:, :)
273 i_block = i_block + 1
283 CALL timestop(handle)
285 END SUBROUTINE set_blocks_to_matrix_v_kp
299 SUBROUTINE compute_v_transl(matrix_v_L_tmp, blocks_v_L, vec_L, particle_set, &
300 qs_kind_set, atomic_kind_set, basis_type, cell, operator_type)
303 TYPE(two_d_util_type),
ALLOCATABLE,
DIMENSION(:) :: blocks_v_l
304 REAL(kind=
dp),
DIMENSION(3) :: vec_l
306 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
308 CHARACTER(LEN=*),
INTENT(IN) :: basis_type
310 INTEGER :: operator_type
312 CHARACTER(LEN=*),
PARAMETER :: routinen =
'compute_v_transl'
314 INTEGER :: col, handle, i_block, kind_a, kind_b, row
315 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: kind_of
316 REAL(
dp),
DIMENSION(3) :: ra, rab_l, rb
317 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: data_block
318 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: contr_a, contr_b
322 CALL timeset(routinen, handle)
324 NULLIFY (basis_set_a, basis_set_b, data_block)
338 kind_a = kind_of(row)
339 kind_b = kind_of(col)
341 CALL get_qs_kind(qs_kind=qs_kind_set(kind_a), basis_set=basis_set_a, basis_type=basis_type)
342 CALL get_qs_kind(qs_kind=qs_kind_set(kind_b), basis_set=basis_set_b, basis_type=basis_type)
344 ra(1:3) =
pbc(particle_set(row)%r(1:3), cell)
345 rb(1:3) =
pbc(particle_set(col)%r(1:3), cell)
347 rab_l(1:3) = rb(1:3) - ra(1:3) + vec_l(1:3)
352 blocks_v_l(i_block)%block = 0.0_dp
355 fba=basis_set_a, fbb=basis_set_b, scona_shg=contr_a, sconb_shg=contr_b, &
356 calculate_forces=.false.)
358 i_block = i_block + 1
360 DEALLOCATE (contr_a, contr_b)
368 CALL timestop(handle)
370 END SUBROUTINE compute_v_transl
376 SUBROUTINE deallocate_blocks_v_kp(blocks_v_kp)
377 TYPE(two_d_util_type),
ALLOCATABLE, &
378 DIMENSION(:, :, :) :: blocks_v_kp
380 CHARACTER(LEN=*),
PARAMETER :: routinen =
'deallocate_blocks_v_kp'
382 INTEGER :: handle, i_block, i_real_img, ik
384 CALL timeset(routinen, handle)
386 DO ik = lbound(blocks_v_kp, 1), ubound(blocks_v_kp, 1)
387 DO i_real_img = 1,
SIZE(blocks_v_kp, 2)
388 DO i_block = 1,
SIZE(blocks_v_kp, 3)
389 DEALLOCATE (blocks_v_kp(ik, i_real_img, i_block)%block)
394 DEALLOCATE (blocks_v_kp)
396 CALL timestop(handle)
398 END SUBROUTINE deallocate_blocks_v_kp
404 SUBROUTINE deallocate_blocks_v_l(blocks_v_L)
405 TYPE(two_d_util_type),
ALLOCATABLE,
DIMENSION(:) :: blocks_v_l
407 CHARACTER(LEN=*),
PARAMETER :: routinen =
'deallocate_blocks_v_L'
409 INTEGER :: handle, i_block
411 CALL timeset(routinen, handle)
413 DO i_block = 1,
SIZE(blocks_v_l, 1)
414 DEALLOCATE (blocks_v_l(i_block)%block)
417 DEALLOCATE (blocks_v_l)
419 CALL timestop(handle)
421 END SUBROUTINE deallocate_blocks_v_l
428 SUBROUTINE allocate_blocks_v_l(blocks_v_L, matrix_v_L_tmp)
429 TYPE(two_d_util_type),
ALLOCATABLE,
DIMENSION(:) :: blocks_v_l
432 CHARACTER(LEN=*),
PARAMETER :: routinen =
'allocate_blocks_v_L'
434 INTEGER :: col, handle, i_block, nblocks, row
435 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: data_block
438 CALL timeset(routinen, handle)
448 nblocks = nblocks + 1
454 ALLOCATE (blocks_v_l(nblocks))
464 ALLOCATE (blocks_v_l(i_block)%block(
SIZE(data_block, 1),
SIZE(data_block, 2)))
465 blocks_v_l(i_block)%block = 0.0_dp
467 i_block = i_block + 1
473 CALL timestop(handle)
475 END SUBROUTINE allocate_blocks_v_l
484 SUBROUTINE allocate_blocks_v_kp(blocks_v_kp, matrix_v_kp, ikp_start, ikp_end)
485 TYPE(two_d_util_type),
ALLOCATABLE, &
486 DIMENSION(:, :, :) :: blocks_v_kp
487 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_v_kp
488 INTEGER :: ikp_start, ikp_end
490 CHARACTER(LEN=*),
PARAMETER :: routinen =
'allocate_blocks_v_kp'
492 INTEGER :: col, handle, i_block, i_real_img, ik, &
494 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: data_block
497 CALL timeset(routinen, handle)
507 nblocks = nblocks + 1
513 ALLOCATE (blocks_v_kp(ikp_start:ikp_end,
SIZE(matrix_v_kp, 2), nblocks))
515 DO ik = ikp_start, ikp_end
517 DO i_real_img = 1,
SIZE(matrix_v_kp, 2)
527 ALLOCATE (blocks_v_kp(ik, i_real_img, i_block)%block(
SIZE(data_block, 1), &
528 SIZE(data_block, 2)))
529 blocks_v_kp(ik, i_real_img, i_block)%block = 0.0_dp
531 i_block = i_block + 1
541 CALL timestop(handle)
543 END SUBROUTINE allocate_blocks_v_kp
560 SUBROUTINE get_factor_and_xyz_min_max(cell, kpoints, size_lattice_sum, factor, hmat, &
561 x_min, x_max, y_min, y_max, z_min, z_max, nkp_grid)
565 INTEGER :: size_lattice_sum, factor
566 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat
567 INTEGER :: x_min, x_max, y_min, y_max, z_min, z_max
568 INTEGER,
DIMENSION(3) :: nkp_grid
570 CHARACTER(LEN=*),
PARAMETER :: routinen =
'get_factor_and_xyz_min_max'
572 INTEGER :: handle, nkp
573 INTEGER,
DIMENSION(3) :: periodic
575 CALL timeset(routinen, handle)
578 CALL get_cell(cell=cell, h=hmat, periodic=periodic)
580 IF (periodic(1) == 0)
THEN
581 cpassert(nkp_grid(1) == 1)
583 IF (periodic(2) == 0)
THEN
584 cpassert(nkp_grid(2) == 1)
586 IF (periodic(3) == 0)
THEN
587 cpassert(nkp_grid(3) == 1)
590 IF (
modulo(nkp_grid(1), 2) == 1)
THEN
591 factor = 3**(size_lattice_sum - 1)
592 ELSE IF (
modulo(nkp_grid(1), 2) == 0)
THEN
593 factor = 2**(size_lattice_sum - 1)
596 IF (
modulo(nkp_grid(1), 2) == 1)
THEN
597 x_min = -(factor*nkp_grid(1) - 1)/2
598 x_max = (factor*nkp_grid(1) - 1)/2
599 ELSE IF (
modulo(nkp_grid(1), 2) == 0)
THEN
600 x_min = -factor*nkp_grid(1)/2
601 x_max = factor*nkp_grid(1)/2 - 1
603 IF (periodic(1) == 0)
THEN
608 IF (
modulo(nkp_grid(2), 2) == 1)
THEN
609 y_min = -(factor*nkp_grid(2) - 1)/2
610 y_max = (factor*nkp_grid(2) - 1)/2
611 ELSE IF (
modulo(nkp_grid(2), 2) == 0)
THEN
612 y_min = -factor*nkp_grid(2)/2
613 y_max = factor*nkp_grid(2)/2 - 1
615 IF (periodic(2) == 0)
THEN
620 IF (
modulo(nkp_grid(3), 2) == 1)
THEN
621 z_min = -(factor*nkp_grid(3) - 1)/2
622 z_max = (factor*nkp_grid(3) - 1)/2
623 ELSE IF (
modulo(nkp_grid(3), 2) == 0)
THEN
624 z_min = -factor*nkp_grid(3)/2
625 z_max = factor*nkp_grid(3)/2 - 1
627 IF (periodic(3) == 0)
THEN
632 CALL timestop(handle)
634 END SUBROUTINE get_factor_and_xyz_min_max
642 SUBROUTINE allocate_tmp(matrix_v_L_tmp, matrix_v_kp, ikp_start)
645 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_v_kp
648 CHARACTER(LEN=*),
PARAMETER :: routinen =
'allocate_tmp'
652 CALL timeset(routinen, handle)
654 NULLIFY (matrix_v_l_tmp)
657 template=matrix_v_kp(ikp_start, 1)%matrix, &
658 matrix_type=dbcsr_type_no_symmetry)
662 CALL timestop(handle)
664 END SUBROUTINE allocate_tmp
670 SUBROUTINE deallocate_tmp(matrix_v_L_tmp)
674 CHARACTER(LEN=*),
PARAMETER :: routinen =
'deallocate_tmp'
678 CALL timeset(routinen, handle)
682 CALL timestop(handle)
684 END SUBROUTINE deallocate_tmp
697 basis_type, ikp_start, ikp_end)
698 COMPLEX(KIND=dp),
DIMENSION(:, :, :) :: v_k
701 INTEGER :: size_lattice_sum
702 CHARACTER(LEN=*),
INTENT(IN) :: basis_type
703 INTEGER :: ikp_start, ikp_end
705 CHARACTER(LEN=*),
PARAMETER :: routinen =
'build_2c_coulomb_matrix_kp_small_cell'
707 INTEGER :: factor, handle, handle2, i_cell, i_x, i_x_inner, i_x_outer, ik, ikp_local, j_y, &
708 j_y_inner, j_y_outer, k_z, k_z_inner, k_z_outer, n_atom, n_bf, x_max, x_min, y_max, &
710 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: bf_end_from_atom, bf_start_from_atom
711 INTEGER,
DIMENSION(3) :: nkp_grid
712 REAL(kind=
dp) :: coskl, sinkl
713 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: v_l
714 REAL(kind=
dp),
DIMENSION(3) :: vec_l, vec_s
715 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat
720 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
722 CALL timeset(routinen, handle)
726 particle_set=particle_set, &
728 qs_kind_set=qs_kind_set, &
729 atomic_kind_set=atomic_kind_set)
731 CALL get_factor_and_xyz_min_max(cell, kpoints, size_lattice_sum, factor, hmat, &
732 x_min, x_max, y_min, y_max, z_min, z_max, nkp_grid)
734 CALL get_basis_sizes(qs_env, n_atom, basis_type, bf_start_from_atom, bf_end_from_atom, n_bf)
736 ALLOCATE (v_l(n_bf, n_bf))
738 DO i_x_inner = 0, 2*nkp_grid(1) - 1
739 DO j_y_inner = 0, 2*nkp_grid(2) - 1
740 DO k_z_inner = 0, 2*nkp_grid(3) - 1
745 DO i_x_outer = x_min, x_max + nkp_grid(1), 2*nkp_grid(1)
746 DO j_y_outer = y_min, y_max + nkp_grid(2), 2*nkp_grid(2)
747 DO k_z_outer = z_min, z_max + nkp_grid(3), 2*nkp_grid(3)
749 i_x = i_x_inner + i_x_outer
750 j_y = j_y_inner + j_y_outer
751 k_z = k_z_inner + k_z_outer
753 IF (i_x > x_max .OR. i_x < x_min .OR. &
754 j_y > y_max .OR. j_y < y_min .OR. &
755 k_z > z_max .OR. k_z < z_min) cycle
759 vec_s = [real(i_x,
dp), real(j_y,
dp), real(k_z,
dp)]
761 IF (
modulo(i_cell, para_env%num_pe) /= para_env%mepos) cycle
763 vec_l = matmul(hmat, vec_s)
766 CALL add_v_l(v_l, vec_l, n_atom, bf_start_from_atom, bf_end_from_atom, &
767 particle_set, qs_kind_set, atomic_kind_set, basis_type, cell)
774 CALL para_env%sum(v_l)
776 CALL timeset(routinen//
"_R_to_k", handle2)
783 IF (
modulo(ik, para_env%num_pe) /= para_env%mepos) cycle
785 ikp_local = ikp_local + 1
787 IF (ik < ikp_start) cycle
790 coskl = cos(
twopi*dot_product(vec_s(1:3), kpoints%xkp(1:3, ik)))
791 sinkl = sin(
twopi*dot_product(vec_s(1:3), kpoints%xkp(1:3, ik)))
793 v_k(:, :, ikp_local) = v_k(:, :, ikp_local) +
z_one*coskl*v_l(:, :) + &
798 CALL timestop(handle2)
804 CALL timestop(handle)
817 SUBROUTINE get_basis_sizes(qs_env, n_atom, basis_type, bf_start_from_atom, bf_end_from_atom, n_bf)
821 CHARACTER(LEN=*),
INTENT(IN) :: basis_type
822 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: bf_start_from_atom, bf_end_from_atom
825 CHARACTER(LEN=*),
PARAMETER :: routinen =
'get_basis_sizes'
827 INTEGER :: handle, iatom, ikind, n_kind, nsgf
828 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: kind_of
832 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
834 CALL timeset(routinen, handle)
836 CALL get_qs_env(qs_env=qs_env, particle_set=particle_set, &
837 qs_kind_set=qs_kind_set, atomic_kind_set=atomic_kind_set)
840 n_atom =
SIZE(particle_set)
841 n_kind =
SIZE(qs_kind_set)
844 CALL get_qs_kind(qs_kind=qs_kind_set(ikind), basis_set=basis_set_a, &
845 basis_type=basis_type)
846 cpassert(
ASSOCIATED(basis_set_a))
849 ALLOCATE (bf_start_from_atom(n_atom), bf_end_from_atom(n_atom))
853 bf_start_from_atom(iatom) = n_bf + 1
854 ikind = kind_of(iatom)
855 CALL get_qs_kind(qs_kind=qs_kind_set(ikind), nsgf=nsgf, basis_type=basis_type)
857 bf_end_from_atom(iatom) = n_bf
860 CALL timestop(handle)
862 END SUBROUTINE get_basis_sizes
877 SUBROUTINE add_v_l(V_L, vec_L, n_atom, bf_start_from_atom, bf_end_from_atom, &
878 particle_set, qs_kind_set, atomic_kind_set, basis_type, cell)
880 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: v_l
881 REAL(kind=
dp),
DIMENSION(3) :: vec_l
883 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: bf_start_from_atom, bf_end_from_atom
885 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
887 CHARACTER(LEN=*),
INTENT(IN) :: basis_type
890 CHARACTER(LEN=*),
PARAMETER :: routinen =
'add_V_L'
892 INTEGER :: a_1, a_2, atom_a, atom_b, b_1, b_2, &
893 handle, kind_a, kind_b
894 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: kind_of
895 REAL(
dp),
DIMENSION(3) :: ra, rab_l, rb
896 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: v_l_ab
897 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: contr_a, contr_b
900 CALL timeset(routinen, handle)
902 NULLIFY (basis_set_a, basis_set_b)
906 DO atom_a = 1, n_atom
908 DO atom_b = 1, n_atom
910 kind_a = kind_of(atom_a)
911 kind_b = kind_of(atom_b)
913 CALL get_qs_kind(qs_kind=qs_kind_set(kind_a), basis_set=basis_set_a, &
914 basis_type=basis_type)
915 CALL get_qs_kind(qs_kind=qs_kind_set(kind_b), basis_set=basis_set_b, &
916 basis_type=basis_type)
918 ra(1:3) =
pbc(particle_set(atom_a)%r(1:3), cell)
919 rb(1:3) =
pbc(particle_set(atom_b)%r(1:3), cell)
921 rab_l(1:3) = rb(1:3) - ra(1:3) + vec_l(1:3)
926 a_1 = bf_start_from_atom(atom_a)
927 a_2 = bf_end_from_atom(atom_a)
928 b_1 = bf_start_from_atom(atom_b)
929 b_2 = bf_end_from_atom(atom_b)
931 ALLOCATE (v_l_ab(a_2 - a_1 + 1, b_2 - b_1 + 1))
934 fba=basis_set_a, fbb=basis_set_b, &
935 scona_shg=contr_a, sconb_shg=contr_b, &
936 calculate_forces=.false.)
938 v_l(a_1:a_2, b_1:b_2) = v_l(a_1:a_2, b_1:b_2) + v_l_ab(:, :)
940 DEALLOCATE (contr_a, contr_b, v_l_ab)
948 CALL timestop(handle)
950 END SUBROUTINE add_v_l
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
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.
Handles all functions related to the CELL.
subroutine, public get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, h, h_inv, symmetry_id, tag)
Get informations about a simulation cell.
constants for the different operators of the 2c-integrals
integer, parameter, public operator_coulomb
subroutine, public dbcsr_release_p(matrix)
...
subroutine, public dbcsr_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset)
...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_init_p(matrix)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_reserve_all_blocks(matrix)
Reserves all blocks.
Initialization for solid harmonic Gaussian (SHG) integral scheme. Scheme for calculation of contracte...
subroutine, public contraction_matrix_shg(basis, scon_shg)
contraction matrix for SHG integrals
Calculation of contracted, spherical Gaussian integrals using the solid harmonic Gaussian (SHG) integ...
subroutine, public int_operators_r12_ab_shg(r12_operator, vab, dvab, rab, fba, fbb, scona_shg, sconb_shg, omega, calculate_forces)
Calcululates the two-center integrals of the type (a|O(r12)|b) using the SHG scheme.
Defines the basic variable types.
integer, parameter, public dp
Routines to compute the Coulomb integral V_(alpha beta)(k) for a k-point k using lattice summation in...
subroutine, public build_2c_coulomb_matrix_kp(matrix_v_kp, kpoints, basis_type, cell, particle_set, qs_kind_set, atomic_kind_set, size_lattice_sum, operator_type, ikp_start, ikp_end)
...
subroutine, public build_2c_coulomb_matrix_kp_small_cell(v_k, qs_env, kpoints, size_lattice_sum, basis_type, ikp_start, ikp_end)
...
Types and basic routines needed for a kpoint calculation.
subroutine, public get_kpoint_info(kpoint, kp_scheme, nkp_grid, kp_shift, symmetry, verbose, full_grid, use_real_wfn, eps_geo, parallel_group_size, kp_range, nkp, xkp, wkp, para_env, blacs_env_all, para_env_kp, para_env_inter_kp, blacs_env, kp_env, kp_aux_env, mpools, iogrp, nkp_groups, kp_dist, cell_to_index, index_to_cell, sab_nl, sab_nl_nosym)
Retrieve information from a kpoint environment.
Definition of mathematical constants and functions.
complex(kind=dp), parameter, public z_one
complex(kind=dp), parameter, public gaussi
real(kind=dp), parameter, public twopi
Interface to the message passing library MPI.
Define the data structure for the particle information.
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, sab_cneo, 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, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_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, do_rixs, tb_tblite)
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, cneo_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, monovalent, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Contains information about kpoints.
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.