56#include "./base/base_uses.f90"
64 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'lri_forces'
88 CHARACTER(LEN=*),
PARAMETER :: routinen =
'calculate_lri_forces'
90 INTEGER :: handle, iatom, ikind, ispin, natom, &
92 INTEGER,
DIMENSION(:, :, :),
POINTER :: cell_to_index
94 REAL(kind=
dp),
DIMENSION(:),
POINTER :: v_dadr, v_dfdr
102 CALL timeset(routinen, handle)
103 NULLIFY (atomic_kind, force, lri_coef, v_dadr, v_dfdr, virial)
105 IF (
ASSOCIATED(lri_env%soo_list))
THEN
107 nkind = lri_env%lri_ints%nkind
108 nspin =
SIZE(pmatrix, 1)
110 CALL get_qs_env(qs_env=qs_env, force=force, virial=virial)
111 use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
113 CALL get_ks_env(ks_env=ks_env, kpoints=kpoints)
117 CALL calculate_v_dadr_sr(lri_env, lri_density, pmatrix, cell_to_index, atomic_kind_set, &
119 IF (lri_env%distant_pair_approximation)
THEN
120 CALL calculate_v_dadr_ff(lri_env, lri_density, pmatrix, cell_to_index, atomic_kind_set, &
126 lri_coef => lri_density%lri_coefs(ispin)%lri_kinds
129 atomic_kind => atomic_kind_set(ikind)
132 v_dadr => lri_coef(ikind)%v_dadr(iatom, :)
133 v_dfdr => lri_coef(ikind)%v_dfdr(iatom, :)
135 force(ikind)%rho_lri_elec(:, iatom) = force(ikind)%rho_lri_elec(:, iatom) &
136 + v_dfdr(:) + v_dadr(:)
144 CALL timestop(handle)
159 SUBROUTINE calculate_v_dadr_sr(lri_env, lri_density, pmatrix, cell_to_index, atomic_kind_set, &
165 INTEGER,
DIMENSION(:, :, :),
POINTER :: cell_to_index
167 LOGICAL,
INTENT(IN) :: use_virial
170 CHARACTER(LEN=*),
PARAMETER :: routinen =
'calculate_v_dadr_sr'
172 INTEGER :: atom_a, atom_b, handle, i, iac, iatom, ic, ikind, ilist, ispin, jatom, jkind, &
173 jneighbor, k, mepos, nba, nbb, nfa, nfb, nkind, nlist, nn, nspin, nthread
174 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: atom_of_kind
175 INTEGER,
DIMENSION(3) :: cell
176 LOGICAL :: found, use_cell_mapping
177 REAL(kind=
dp) :: ai, dab, dfw, fw, isn, threshold
178 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: vint
179 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: pab
180 REAL(kind=
dp),
DIMENSION(3) :: dcharge, dlambda, force_a, force_b, &
181 idav, nsdssn, nsdsst, nsdt, rab
182 REAL(kind=
dp),
DIMENSION(:),
POINTER :: st, v_dadra, v_dadrb
183 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: dssn, dsst, dtvec, pbij, sdssn, sdsst, &
194 DIMENSION(:),
POINTER :: nl_iterator
198 CALL timeset(routinen, handle)
200 NULLIFY (lri_coef, lri_force, lrii, lri_rho, lrho, &
201 nl_iterator, pbij, pmat, soo_list, v_dadra, v_dadrb)
202 NULLIFY (dssn, dsst, dtvec, sdssn, sdsst, sdt, st)
204 IF (
ASSOCIATED(lri_env%soo_list))
THEN
205 soo_list => lri_env%soo_list
207 nkind = lri_env%lri_ints%nkind
208 nspin =
SIZE(pmatrix, 1)
211 use_cell_mapping = (
SIZE(pmatrix, 2) > 1)
217 lri_coef => lri_density%lri_coefs(ispin)%lri_kinds
218 lri_rho => lri_density%lri_rhos(ispin)%lri_list
235 iatom=iatom, jatom=jatom, nlist=nlist, ilist=ilist, &
236 inode=jneighbor, r=rab, cell=cell)
238 iac = ikind + nkind*(jkind - 1)
240 IF (.NOT.
ASSOCIATED(lri_env%lri_ints%lri_atom(iac)%lri_node)) cycle
242 lrii => lri_env%lri_ints%lri_atom(iac)%lri_node(ilist)%lri_int(jneighbor)
243 lrho => lri_rho%lri_atom(iac)%lri_node(ilist)%lri_rhoab(jneighbor)
246 IF (.NOT. lrii%lrisr) cycle
252 IF (.NOT. lrii%calc_force_pair) cycle
260 IF (use_cell_mapping)
THEN
261 ic = cell_to_index(cell(1), cell(2), cell(3))
266 pmat => pmatrix(ispin, ic)%matrix
269 ALLOCATE (pab(nba, nbb))
271 IF (iatom <= jatom)
THEN
272 CALL dbcsr_get_block_p(matrix=pmat, row=iatom, col=jatom, block=pbij, found=found)
274 pab(1:nba, 1:nbb) = pbij(1:nba, 1:nbb)
276 CALL dbcsr_get_block_p(matrix=pmat, row=jatom, col=iatom, block=pbij, found=found)
278 pab(1:nba, 1:nbb) = transpose(pbij(1:nbb, 1:nba))
281 obasa => lri_env%orb_basis(ikind)%gto_basis_set
282 obasb => lri_env%orb_basis(jkind)%gto_basis_set
283 fbasa => lri_env%ri_basis(ikind)%gto_basis_set
284 fbasb => lri_env%ri_basis(jkind)%gto_basis_set
288 CALL lri_dint2(lri_env, lrii, lridint, rab, obasa, obasb, fbasa, fbasb, &
289 iatom, jatom, ikind, jkind)
294 dsst => lri_force%dsst
295 dssn => lri_force%dssn
296 dtvec => lri_force%dtvec
300 threshold = 0.01_dp*lri_env%eps_o3_int/max(sum(abs(pab(1:nba, 1:nbb))), 1.0e-14_dp)
301 dcharge(k) = sum(pab(1:nba, 1:nbb)*lridint%dsooint(1:nba, 1:nbb, k))
303 IF (lrii%abascr(i) > threshold)
THEN
304 dtvec(i, k) = sum(pab(1:nba, 1:nbb)*lridint%dabaint(1:nba, 1:nbb, i, k))
308 IF (lrii%abbscr(i) > threshold)
THEN
309 dtvec(nfa + i, k) = sum(pab(1:nba, 1:nbb)*lridint%dabbint(1:nba, 1:nbb, i, k))
316 atom_a = atom_of_kind(iatom)
317 atom_b = atom_of_kind(jatom)
319 vint(1:nfa) = lri_coef(ikind)%v_int(atom_a, 1:nfa)
320 vint(nfa + 1:nn) = lri_coef(jkind)%v_int(atom_b, 1:nfb)
322 isn = sum(vint(1:nn)*lrii%sn(1:nn))
324 ai = isn/lrii%nsn*dcharge(k)
325 force_a(k) = 2.0_dp*fw*ai
326 force_b(k) = -2.0_dp*fw*ai
331 st(1:nn) = matmul(lrii%sinv(1:nn, 1:nn), lrho%tvec(1:nn))
333 dsst(1:nfa, k) = matmul(lridint%dsabint(1:nfa, 1:nfb, k), st(nfa + 1:nn))
334 dsst(nfa + 1:nn, k) = matmul(st(1:nfa), lridint%dsabint(1:nfa, 1:nfb, k))
335 nsdsst(k) = sum(lrii%sn(1:nn)*dsst(1:nn, k))
336 dssn(1:nfa, k) = matmul(lridint%dsabint(1:nfa, 1:nfb, k), lrii%sn(nfa + 1:nn))
337 dssn(nfa + 1:nn, k) = matmul(lrii%sn(1:nfa), lridint%dsabint(1:nfa, 1:nfb, k))
338 nsdssn(k) = sum(lrii%sn(1:nn)*dssn(1:nn, k))
339 nsdt(k) = sum(dtvec(1:nn, k)*lrii%sn(1:nn))
343 dlambda(k) = (nsdsst(k) - nsdt(k))/lrii%nsn &
344 + (lrho%charge - lrho%nst)*nsdssn(k)/(lrii%nsn*lrii%nsn)
347 force_a(k) = force_a(k) + 2.0_dp*fw*isn*dlambda(k)
348 force_b(k) = force_b(k) - 2.0_dp*fw*isn*dlambda(k)
351 st(1:nn) = dtvec(1:nn, k) - dsst(1:nn, k) - lrho%lambda*dssn(1:nn, k)
352 idav(k) = sum(vint(1:nn)*matmul(lrii%sinv(1:nn, 1:nn), st(1:nn)))
360 ai = 2.0_dp*fw*idav(k)
361 force_a(k) = force_a(k) + ai
362 force_b(k) = force_b(k) - ai
364 IF (abs(dfw) > 0.0_dp)
THEN
365 dab = sqrt(sum(rab(1:3)*rab(1:3)))
366 ai = 2.0_dp*dfw/dab*sum(lrho%avec(1:nn)*vint(1:nn))
368 force_a(k) = force_a(k) - ai*rab(k)
369 force_b(k) = force_b(k) + ai*rab(k)
375 v_dadra => lri_coef(ikind)%v_dadr(atom_a, :)
376 v_dadrb => lri_coef(jkind)%v_dadr(atom_b, :)
379 v_dadra(k) = v_dadra(k) + force_a(k)
380 v_dadrb(k) = v_dadrb(k) + force_b(k)
388 IF (iatom == jatom)
THEN
409 CALL timestop(handle)
411 END SUBROUTINE calculate_v_dadr_sr
424 SUBROUTINE calculate_v_dadr_ff(lri_env, lri_density, pmatrix, cell_to_index, atomic_kind_set, &
430 INTEGER,
DIMENSION(:, :, :),
POINTER :: cell_to_index
432 LOGICAL,
INTENT(IN) :: use_virial
435 CHARACTER(LEN=*),
PARAMETER :: routinen =
'calculate_v_dadr_ff'
437 INTEGER :: atom_a, atom_b, handle, i, iac, iatom, ic, ikind, ilist, ispin, jatom, jkind, &
438 jneighbor, k, mepos, nba, nbb, nfa, nfb, nkind, nlist, nn, nspin, nthread
439 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: atom_of_kind
440 INTEGER,
DIMENSION(3) :: cell
441 LOGICAL :: found, use_cell_mapping
442 REAL(kind=
dp) :: ai, dab, dfw, fw, isna, isnb, threshold
443 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: pab, wab, wbb
444 REAL(kind=
dp),
DIMENSION(3) :: dchargea, dchargeb, force_a, force_b, &
446 REAL(kind=
dp),
DIMENSION(:),
POINTER :: sta, stb, v_dadra, v_dadrb, vinta, vintb
447 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: dtveca, dtvecb, pbij
457 DIMENSION(:),
POINTER :: nl_iterator
461 CALL timeset(routinen, handle)
463 NULLIFY (lri_coef, lrii, lri_rho, lrho, &
464 nl_iterator, pbij, pmat, soo_list, v_dadra, v_dadrb)
466 IF (
ASSOCIATED(lri_env%soo_list))
THEN
467 soo_list => lri_env%soo_list
469 nkind = lri_env%lri_ints%nkind
470 nspin =
SIZE(pmatrix, 1)
473 use_cell_mapping = (
SIZE(pmatrix, 2) > 1)
479 lri_coef => lri_density%lri_coefs(ispin)%lri_kinds
480 lri_rho => lri_density%lri_rhos(ispin)%lri_list
497 iatom=iatom, jatom=jatom, nlist=nlist, ilist=ilist, &
498 inode=jneighbor, r=rab, cell=cell)
500 iac = ikind + nkind*(jkind - 1)
502 IF (.NOT.
ASSOCIATED(lri_env%lri_ints%lri_atom(iac)%lri_node)) cycle
504 lrii => lri_env%lri_ints%lri_atom(iac)%lri_node(ilist)%lri_int(jneighbor)
505 lrho => lri_rho%lri_atom(iac)%lri_node(ilist)%lri_rhoab(jneighbor)
508 IF (.NOT. lrii%lriff) cycle
514 IF (.NOT. lrii%calc_force_pair) cycle
522 IF (use_cell_mapping)
THEN
523 ic = cell_to_index(cell(1), cell(2), cell(3))
528 pmat => pmatrix(ispin, ic)%matrix
531 ALLOCATE (pab(nba, nbb))
533 IF (iatom <= jatom)
THEN
534 CALL dbcsr_get_block_p(matrix=pmat, row=iatom, col=jatom, block=pbij, found=found)
536 pab(1:nba, 1:nbb) = pbij(1:nba, 1:nbb)
538 CALL dbcsr_get_block_p(matrix=pmat, row=jatom, col=iatom, block=pbij, found=found)
540 pab(1:nba, 1:nbb) = transpose(pbij(1:nbb, 1:nba))
543 ALLOCATE (wab(nba, nbb), wbb(nba, nbb))
544 wab(1:nba, 1:nbb) = lri_env%wmat(ikind, jkind)%mat(1:nba, 1:nbb)
545 wbb(1:nba, 1:nbb) = 1.0_dp - lri_env%wmat(ikind, jkind)%mat(1:nba, 1:nbb)
547 obasa => lri_env%orb_basis(ikind)%gto_basis_set
548 obasb => lri_env%orb_basis(jkind)%gto_basis_set
549 fbasa => lri_env%ri_basis(ikind)%gto_basis_set
550 fbasb => lri_env%ri_basis(jkind)%gto_basis_set
555 CALL lri_dint2(lri_env, lrii, lridint, rab, obasa, obasb, fbasa, fbasb, &
556 iatom, jatom, ikind, jkind)
558 NULLIFY (lri_force_a, lri_force_b)
561 dtveca => lri_force_a%dtvec
562 dtvecb => lri_force_b%dtvec
563 sta => lri_force_a%st
564 stb => lri_force_b%st
567 threshold = 0.01_dp*lri_env%eps_o3_int/max(sum(abs(pab(1:nba, 1:nbb))), 1.0e-14_dp)
569 dchargea(k) = sum(wab(1:nba, 1:nbb)*pab(1:nba, 1:nbb)*lridint%dsooint(1:nba, 1:nbb, k))
571 IF (lrii%abascr(i) > threshold)
THEN
572 dtveca(i, k) = sum(wab(1:nba, 1:nbb)*pab(1:nba, 1:nbb)*lridint%dabaint(1:nba, 1:nbb, i, k))
575 dchargeb(k) = sum(wbb(1:nba, 1:nbb)*pab(1:nba, 1:nbb)*lridint%dsooint(1:nba, 1:nbb, k))
577 IF (lrii%abbscr(i) > threshold)
THEN
578 dtvecb(i, k) = sum(wbb(1:nba, 1:nbb)*pab(1:nba, 1:nbb)*lridint%dabbint(1:nba, 1:nbb, i, k))
583 DEALLOCATE (pab, wab, wbb)
585 atom_a = atom_of_kind(iatom)
586 atom_b = atom_of_kind(jatom)
587 vinta => lri_coef(ikind)%v_int(atom_a, :)
588 vintb => lri_coef(jkind)%v_int(atom_b, :)
590 isna = sum(vinta(1:nfa)*lrii%sna(1:nfa))
591 isnb = sum(vintb(1:nfb)*lrii%snb(1:nfb))
593 ai = isna/lrii%nsna*dchargea(k) + isnb/lrii%nsnb*dchargeb(k)
594 force_a(k) = 2.0_dp*fw*ai
595 force_b(k) = -2.0_dp*fw*ai
599 sta(1:nfa) = matmul(lrii%asinv(1:nfa, 1:nfa), dtveca(1:nfa, k))
600 idava(k) = sum(vinta(1:nfa)*sta(1:nfa)) - isna/lrii%nsna*sum(lrii%na(1:nfa)*sta(1:nfa))
601 stb(1:nfb) = matmul(lrii%bsinv(1:nfb, 1:nfb), dtvecb(1:nfb, k))
602 idavb(k) = sum(vintb(1:nfb)*stb(1:nfb)) - isnb/lrii%nsnb*sum(lrii%nb(1:nfb)*stb(1:nfb))
610 ai = 2.0_dp*fw*(idava(k) + idavb(k))
611 force_a(k) = force_a(k) + ai
612 force_b(k) = force_b(k) - ai
614 IF (abs(dfw) > 0.0_dp)
THEN
615 dab = sqrt(sum(rab(1:3)*rab(1:3)))
616 ai = 2.0_dp*dfw/dab* &
617 (sum(lrho%aveca(1:nfa)*vinta(1:nfa)) + &
618 sum(lrho%avecb(1:nfb)*vintb(1:nfb)))
620 force_a(k) = force_a(k) - ai*rab(k)
621 force_b(k) = force_b(k) + ai*rab(k)
624 v_dadra => lri_coef(ikind)%v_dadr(atom_a, :)
625 v_dadrb => lri_coef(jkind)%v_dadr(atom_b, :)
629 v_dadra(k) = v_dadra(k) + force_a(k)
630 v_dadrb(k) = v_dadrb(k) + force_b(k)
638 IF (iatom == jatom)
THEN
660 CALL timestop(handle)
662 END SUBROUTINE calculate_v_dadr_ff
680 CHARACTER(LEN=*),
PARAMETER :: routinen =
'calculate_ri_forces'
682 INTEGER :: handle, iatom, ikind, ispin, natom, &
684 LOGICAL :: use_virial
685 REAL(kind=
dp),
DIMENSION(:),
POINTER :: v_dadr, v_dfdr
691 CALL timeset(routinen, handle)
692 NULLIFY (atomic_kind, force, lri_coef, v_dadr, v_dfdr, virial)
694 nkind =
SIZE(atomic_kind_set)
695 nspin =
SIZE(pmatrix)
697 CALL get_qs_env(qs_env=qs_env, force=force, virial=virial)
698 use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
701 CALL calculate_v_dadr_ri(lri_env, lri_density, pmatrix, atomic_kind_set, &
706 lri_coef => lri_density%lri_coefs(ispin)%lri_kinds
709 atomic_kind => atomic_kind_set(ikind)
712 v_dadr => lri_coef(ikind)%v_dadr(iatom, :)
713 v_dfdr => lri_coef(ikind)%v_dfdr(iatom, :)
715 force(ikind)%rho_lri_elec(:, iatom) = force(ikind)%rho_lri_elec(:, iatom) &
716 + v_dfdr(:) + v_dadr(:)
722 CALL timestop(handle)
738 SUBROUTINE calculate_v_dadr_ri(lri_env, lri_density, pmatrix, atomic_kind_set, &
745 LOGICAL,
INTENT(IN) :: use_virial
748 CHARACTER(LEN=*),
PARAMETER :: routinen =
'calculate_v_dadr_ri'
750 INTEGER :: atom_a, atom_b, atom_c, handle, i, i1, &
751 i2, iatom, ikind, ispin, jatom, jkind, &
752 katom, kkind, m, mepos, nkind, nspin, &
754 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: atom_of_kind
755 INTEGER,
DIMENSION(:, :),
POINTER :: bas_ptr
756 REAL(kind=
dp) :: fscal
757 REAL(kind=
dp),
DIMENSION(3) :: force_a, force_b, force_c, rij, rik, rjk
758 REAL(kind=
dp),
DIMENSION(:),
POINTER :: v_dadra, v_dadrb, v_dadrc
759 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: fi, fj, fk, fo
763 CALL timeset(routinen, handle)
765 nspin =
SIZE(pmatrix)
766 nkind =
SIZE(atomic_kind_set)
768 lri_coef => lri_density%lri_coefs(ispin)%lri_kinds
770 lri_coef(ikind)%v_dadr(:, :) = 0.0_dp
774 bas_ptr => lri_env%ri_fit%bas_ptr
779 fo => lri_env%ri_fit%fout
794 DO WHILE (
o3c_iterate(o3c_iterator, mepos=mepos) == 0)
796 ikind=ikind, jkind=jkind, kkind=kkind, &
797 iatom=iatom, jatom=jatom, katom=katom, &
798 rij=rij, rik=rik, force_i=fi, force_j=fj, force_k=fk)
799 i1 = bas_ptr(1, katom)
800 i2 = bas_ptr(2, katom)
807 force_a(i) = force_a(i) + sum(fi(1:m, i)*fo(i1:i2, ispin))
808 force_b(i) = force_b(i) + sum(fj(1:m, i)*fo(i1:i2, ispin))
809 force_c(i) = force_c(i) + sum(fk(1:m, i)*fo(i1:i2, ispin))
812 atom_a = atom_of_kind(iatom)
813 atom_b = atom_of_kind(jatom)
814 atom_c = atom_of_kind(katom)
816 v_dadra => lri_coef(ikind)%v_dadr(atom_a, :)
817 v_dadrb => lri_coef(jkind)%v_dadr(atom_b, :)
818 v_dadrc => lri_coef(kkind)%v_dadr(atom_c, :)
821 v_dadra(1:3) = v_dadra(1:3) + force_a(1:3)
822 v_dadrb(1:3) = v_dadrb(1:3) + force_b(1:3)
823 v_dadrc(1:3) = v_dadrc(1:3) + force_c(1:3)
826 rjk(1:3) = rik(1:3) - rij(1:3)
829 IF (iatom == jatom) fscal = 1.0_dp
840 CALL timestop(handle)
842 END SUBROUTINE calculate_v_dadr_ri
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.
subroutine, public dbcsr_get_block_p(matrix, row, col, block, found, row_size, col_size)
...
Defines the basic variable types.
integer, parameter, public dp
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.
contains the types and subroutines for dealing with the lri_env lri : local resolution of the identit...
subroutine, public deallocate_lri_force_components(lri_force)
releases the given lri_force_type
subroutine, public allocate_lri_force_components(lri_force, nfa, nfb)
creates and initializes lri_force
Calculates forces for LRIGPW method lri : local resolution of the identity.
subroutine, public calculate_lri_forces(lri_env, lri_density, qs_env, pmatrix, atomic_kind_set)
calculates the lri forces
subroutine, public calculate_ri_forces(lri_env, lri_density, qs_env, pmatrix, atomic_kind_set)
calculates the ri forces
Calculates integrals for LRIGPW method lri : local resolution of the identity.
subroutine, public deallocate_int_type(lriint, lridint)
...
subroutine, public allocate_int_type(lriint, lridint, nba, nbb, nfa, nfb, skip_sab, skip_soo, skip_aba, skip_abb, skip_dsab, skip_dsoo, skip_daba, skip_dabb)
...
subroutine, public lri_dint2(lri_env, lrii, lridint, rab, obasa, obasb, fbasa, fbasb, iatom, jatom, ikind, jkind)
...
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.
subroutine, public get_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, rho, rho_xc, vppl, rho_core, rho_nlcc, rho_nlcc_g, vee, neighbor_list_id, 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, kpoints, do_kpoints, atomic_kind_set, qs_kind_set, cell, cell_ref, use_ref_cell, particle_set, energy, force, local_particles, local_molecules, molecule_kind_set, molecule_set, subsys, cp_subsys, virial, results, atprop, nkind, natom, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env, nelectron_total, nelectron_spin)
...
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)
...
3-center overlap type integrals containers
subroutine, public get_o3c_iterator_info(o3c_iterator, mepos, iatom, jatom, katom, ikind, jkind, kkind, rij, rik, cellj, cellk, integral, tvec, force_i, force_j, force_k)
...
subroutine, public o3c_iterator_create(o3c, o3c_iterator, nthread)
...
subroutine, public o3c_iterator_release(o3c_iterator)
...
integer function, public o3c_iterate(o3c_iterator, mepos)
...
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.
Contains information about kpoints.
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...