102 calculate_forces, just_energy)
105 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: ks_matrix
107 REAL(
dp),
DIMENSION(:, :),
INTENT(in) :: charges
108 REAL(
dp),
DIMENSION(:),
INTENT(in) :: mcharge
110 LOGICAL,
INTENT(in) :: calculate_forces, just_energy
112 CHARACTER(len=*),
PARAMETER :: routinen =
'build_xtb_coulomb'
114 INTEGER :: atom_i, atom_j, ewald_type, handle, i, ia, iac, iatom, ic, icol, ikind, img, &
115 irow, is, j, jatom, jkind, la, lb, lmaxa, lmaxb, natom, natorb_a, natorb_b, ni, nimg, nj, &
117 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: atom_of_kind, kind_of
118 INTEGER,
DIMENSION(25) :: laoa, laob
119 INTEGER,
DIMENSION(3) :: cellind, periodic
120 INTEGER,
DIMENSION(5) :: occ
121 INTEGER,
DIMENSION(:, :, :),
POINTER :: cell_to_index
122 LOGICAL :: defined, do_ewald, do_gamma_stress, &
124 REAL(kind=
dp) :: alpha, deth, dr, ecsr, etaa, etab, f1, &
125 f2, fi, gmij, kg, rcut, rcuta, rcutb, &
127 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: xgamma, zeffk
128 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: gammab, gcij, gmcharge
129 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: gchrg
130 REAL(kind=
dp),
DIMENSION(25) :: gcint
131 REAL(kind=
dp),
DIMENSION(3) :: fij, rij
132 REAL(kind=
dp),
DIMENSION(5) :: kappaa, kappab
133 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: dsblock, ksblock, pblock, sblock
134 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: dsint
139 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_p, matrix_s
147 DIMENSION(:),
POINTER :: nl_iterator
152 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
155 TYPE(
xtb_atom_type),
POINTER :: xtb_atom_a, xtb_atom_b, xtb_kind
158 CALL timeset(routinen, handle)
160 NULLIFY (matrix_p, matrix_s, virial, atprop, dft_control)
163 qs_kind_set=qs_kind_set, &
164 particle_set=particle_set, &
168 dft_control=dft_control)
170 xtb_control => dft_control%qs_control%xtb_control
173 IF (calculate_forces)
THEN
174 use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
177 do_gamma_stress = .false.
178 IF (.NOT. just_energy .AND. use_virial)
THEN
179 IF (dft_control%nimages == 1) do_gamma_stress = .true.
182 IF (atprop%energy)
THEN
183 CALL get_qs_env(qs_env=qs_env, particle_set=particle_set)
184 natom =
SIZE(particle_set)
188 IF (calculate_forces)
THEN
194 CALL get_qs_env(qs_env, nkind=nkind, natom=natom)
195 ALLOCATE (gchrg(natom, 5, nmat))
197 ALLOCATE (gmcharge(natom, nmat))
204 CALL get_qs_env(qs_env=qs_env, sab_xtbe=n_list)
205 IF (.NOT.
ASSOCIATED(n_list))
THEN
206 cpabort(
"sab_xtbe neighbor list is not associated in build_xtb_coulomb")
211 iatom=iatom, jatom=jatom, r=rij, cell=cellind)
212 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
214 IF (.NOT. defined .OR. natorb_a < 1) cycle
215 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
217 IF (.NOT. defined .OR. natorb_b < 1) cycle
224 ALLOCATE (gammab(ni, nj))
226 dr = sqrt(sum(rij(:)**2))
227 CALL gamma_rab_sr(gammab, dr, ni, kappaa, etaa, nj, kappab, etab, kg, rcut)
228 gchrg(iatom, 1:ni, 1) = gchrg(iatom, 1:ni, 1) + matmul(gammab, charges(jatom, 1:nj))
229 IF (iatom /= jatom)
THEN
230 gchrg(jatom, 1:nj, 1) = gchrg(jatom, 1:nj, 1) + matmul(charges(iatom, 1:ni), gammab)
232 IF (calculate_forces)
THEN
233 IF (dr > 1.e-6_dp)
THEN
234 CALL dgamma_rab_sr(gammab, dr, ni, kappaa, etaa, nj, kappab, etab, kg, rcut)
236 gchrg(iatom, 1:ni, i + 1) = gchrg(iatom, 1:ni, i + 1) &
237 + matmul(gammab, charges(jatom, 1:nj))*rij(i)/dr
238 IF (iatom /= jatom)
THEN
239 gchrg(jatom, 1:nj, i + 1) = gchrg(jatom, 1:nj, i + 1) &
240 - matmul(charges(iatom, 1:ni), gammab)*rij(i)/dr
244 gcint(1:ni) = matmul(gammab, charges(jatom, 1:nj))
246 fij(i) = -sum(charges(iatom, 1:ni)*gcint(1:ni))*rij(i)/dr
249 IF (iatom == jatom) fi = 0.5_dp
260 IF (xtb_control%coulomb_lr)
THEN
261 do_ewald = xtb_control%do_ewald
264 NULLIFY (ewald_env, ewald_pw)
266 ewald_env=ewald_env, ewald_pw=ewald_pw)
267 CALL get_cell(cell=cell, periodic=periodic, deth=deth)
268 CALL ewald_env_get(ewald_env, alpha=alpha, ewald_type=ewald_type)
269 CALL get_qs_env(qs_env=qs_env, sab_tbe=n_list)
271 SELECT CASE (ewald_type)
273 cpabort(
"Invalid Ewald type")
275 cpabort(
"Not allowed with xTB/DFTB")
277 cpabort(
"Standard Ewald not implemented in xTB/DFTB")
279 cpabort(
"PME not implemented in xTB/DFTB")
282 gmcharge, mcharge, calculate_forces, virial, use_virial)
287 local_particles=local_particles)
288 DO ikind = 1,
SIZE(local_particles%n_el)
289 DO ia = 1, local_particles%n_el(ikind)
290 iatom = local_particles%list(ikind)%array(ia)
291 DO jatom = 1, iatom - 1
292 rij = particle_set(iatom)%r - particle_set(jatom)%r
294 dr = sqrt(sum(rij(:)**2))
295 IF (dr > 1.e-6_dp)
THEN
296 gmcharge(iatom, 1) = gmcharge(iatom, 1) + mcharge(jatom)/dr
297 gmcharge(jatom, 1) = gmcharge(jatom, 1) + mcharge(iatom)/dr
299 gmcharge(iatom, i) = gmcharge(iatom, i) + rij(i - 1)*mcharge(jatom)/dr**3
300 gmcharge(jatom, i) = gmcharge(jatom, i) - rij(i - 1)*mcharge(iatom)/dr**3
304 fij(i) = mcharge(iatom)*mcharge(jatom)*rij(i)/dr**3
317 atomic_kind_set=atomic_kind_set, &
318 force=force, para_env=para_env)
319 CALL para_env%sum(gmcharge(:, 1))
320 CALL para_env%sum(gchrg(:, :, 1))
322 IF (xtb_control%coulomb_lr)
THEN
325 gmcharge(:, 1) = gmcharge(:, 1) - 2._dp*alpha*
oorootpi*mcharge(:)
326 IF (any(periodic(:) == 1))
THEN
327 gmcharge(:, 1) = gmcharge(:, 1) -
pi/alpha**2/deth
335 atom_of_kind=atom_of_kind)
338 ikind = kind_of(iatom)
339 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_kind)
342 ecsr = ecsr + sum(charges(iatom, 1:ni)*gchrg(iatom, 1:ni, 1))
345 energy%hartree = energy%hartree + 0.5_dp*ecsr
346 energy%hartree = energy%hartree + 0.5_dp*sum(mcharge(:)*gmcharge(:, 1))
348 IF (atprop%energy)
THEN
349 CALL get_qs_env(qs_env=qs_env, local_particles=local_particles)
350 DO ikind = 1,
SIZE(local_particles%n_el)
351 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_kind)
354 zeff = sum(real(occ, kind=
dp))
355 DO ia = 1, local_particles%n_el(ikind)
356 iatom = local_particles%list(ikind)%array(ia)
357 atprop%atecoul(iatom) = atprop%atecoul(iatom) + &
358 0.5_dp*sum(real(occ(1:ni), kind=
dp)*gchrg(iatom, 1:ni, 1))
359 atprop%atecoul(iatom) = atprop%atecoul(iatom) + &
360 0.5_dp*zeff*gmcharge(iatom, 1)
365 IF (calculate_forces)
THEN
367 ikind = kind_of(iatom)
368 atom_i = atom_of_kind(iatom)
369 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_kind)
374 fij(i) = sum(charges(iatom, 1:ni)*gchrg(iatom, 1:ni, i + 1))
376 force(ikind)%rho_elec(1, atom_i) = force(ikind)%rho_elec(1, atom_i) - fij(1)
377 force(ikind)%rho_elec(2, atom_i) = force(ikind)%rho_elec(2, atom_i) - fij(2)
378 force(ikind)%rho_elec(3, atom_i) = force(ikind)%rho_elec(3, atom_i) - fij(3)
381 fij(i) = gmcharge(iatom, i + 1)*mcharge(iatom)
383 force(ikind)%rho_elec(1, atom_i) = force(ikind)%rho_elec(1, atom_i) - fij(1)
384 force(ikind)%rho_elec(2, atom_i) = force(ikind)%rho_elec(2, atom_i) - fij(2)
385 force(ikind)%rho_elec(3, atom_i) = force(ikind)%rho_elec(3, atom_i) - fij(3)
389 IF (.NOT. just_energy)
THEN
390 CALL get_qs_env(qs_env=qs_env, matrix_s_kp=matrix_s)
393 nimg = dft_control%nimages
394 NULLIFY (cell_to_index)
397 CALL get_qs_env(qs_env=qs_env, kpoints=kpoints)
401 IF (calculate_forces .AND.
SIZE(matrix_p, 1) == 2)
THEN
403 CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
404 alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
409 IF (do_gamma_stress)
THEN
419 ikind = kind_of(irow)
420 jkind = kind_of(icol)
423 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
424 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
430 ALLOCATE (gcij(ni, nj))
435 gcij(i, j) = 0.5_dp*(gchrg(irow, la, 1) + gchrg(icol, lb, 1))
438 gmij = 0.5_dp*(gmcharge(irow, 1) + gmcharge(icol, 1))
439 DO is = 1,
SIZE(ks_matrix, 1)
442 row=irow, col=icol, block=ksblock, found=found)
444 ksblock = ksblock - gcij*sblock
445 ksblock = ksblock - gmij*sblock
447 IF (calculate_forces)
THEN
448 atom_i = atom_of_kind(irow)
449 atom_j = atom_of_kind(icol)
452 row=irow, col=icol, block=pblock, found=found)
457 row=irow, col=icol, block=dsblock, found=found)
461 fi = -2.0_dp*sum(pblock*dsblock*gcij)
462 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
463 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
466 fi = -2.0_dp*gmij*sum(pblock*dsblock)
467 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
468 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
476 IF (do_gamma_stress)
THEN
479 iac = ikind + nkind*(jkind - 1)
480 IF (.NOT.
ASSOCIATED(sap_int(iac)%alist)) cycle
482 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
483 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
486 DO ia = 1, sap_int(iac)%nalist
487 IF (.NOT.
ASSOCIATED(sap_int(iac)%alist(ia)%clist)) cycle
488 iatom = sap_int(iac)%alist(ia)%aatom
489 DO ic = 1, sap_int(iac)%alist(ia)%nclist
490 jatom = sap_int(iac)%alist(ia)%clist(ic)%catom
491 rij = sap_int(iac)%alist(ia)%clist(ic)%rac
492 dr = sqrt(sum(rij(:)**2))
493 IF (dr > 1.e-6_dp)
THEN
494 dsint => sap_int(iac)%alist(ia)%clist(ic)%acint
495 ALLOCATE (gcij(ni, nj))
500 gcij(i, j) = 0.5_dp*(gchrg(iatom, la, 1) + gchrg(jatom, lb, 1))
503 gmij = 0.5_dp*(gmcharge(iatom, 1) + gmcharge(jatom, 1))
504 icol = max(iatom, jatom)
505 irow = min(iatom, jatom)
508 row=irow, col=icol, block=pblock, found=found)
513 IF (irow == iatom)
THEN
514 f1 = -2.0_dp*sum(pblock*dsint(:, :, i)*gcij)
515 f2 = -2.0_dp*gmij*sum(pblock*dsint(:, :, i))
517 f1 = -2.0_dp*sum(transpose(pblock)*dsint(:, :, i)*gcij)
518 f2 = -2.0_dp*gmij*sum(transpose(pblock)*dsint(:, :, i))
524 IF (iatom == jatom) fi = 0.5_dp
534 CALL get_qs_env(qs_env=qs_env, sab_orb=n_list)
538 iatom=iatom, jatom=jatom, r=rij, cell=cellind)
540 icol = max(iatom, jatom)
541 irow = min(iatom, jatom)
543 ic = cell_to_index(cellind(1), cellind(2), cellind(3))
548 row=irow, col=icol, block=sblock, found=found)
552 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
553 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
559 ALLOCATE (gcij(ni, nj))
562 IF (irow == iatom)
THEN
565 gcij(i, j) = 0.5_dp*(gchrg(iatom, la, 1) + gchrg(jatom, lb, 1))
569 gcij(i, j) = 0.5_dp*(gchrg(iatom, la, 1) + gchrg(jatom, lb, 1))
573 gmij = 0.5_dp*(gmcharge(iatom, 1) + gmcharge(jatom, 1))
574 DO is = 1,
SIZE(ks_matrix, 1)
577 row=irow, col=icol, block=ksblock, found=found)
579 ksblock = ksblock - gcij*sblock
580 ksblock = ksblock - gmij*sblock
583 IF (calculate_forces)
THEN
584 atom_i = atom_of_kind(iatom)
585 atom_j = atom_of_kind(jatom)
586 IF (irow /= iatom)
THEN
592 row=irow, col=icol, block=pblock, found=found)
597 row=irow, col=icol, block=dsblock, found=found)
601 fi = -2.0_dp*sum(pblock*dsblock*gcij)
602 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
603 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
606 fi = -2.0_dp*gmij*sum(pblock*dsblock)
607 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
608 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
612 dr = sqrt(sum(rij(:)**2))
613 IF (dr > 1.e-6_dp)
THEN
615 IF (iatom == jatom) fi = 0.5_dp
626 IF (calculate_forces .AND.
SIZE(matrix_p, 1) == 2)
THEN
628 CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
629 alpha_scalar=1.0_dp, beta_scalar=-1.0_dp)
634 IF (xtb_control%tb3_interaction)
THEN
636 ALLOCATE (zeffk(nkind), xgamma(nkind))
638 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_kind)
643 sap_int, calculate_forces, just_energy)
644 DEALLOCATE (zeffk, xgamma)
648 IF (qs_env%qmmm .AND. qs_env%qmmm_periodic)
THEN
650 calculate_forces, just_energy)
653 IF (do_gamma_stress)
THEN
657 CALL timestop(handle)
814 CHARACTER(LEN=*),
PARAMETER :: routinen =
'xtb_dsint_list'
816 INTEGER :: handle, i, iac, iatom, ikind, ilist, iset, jatom, jkind, jneighbor, jset, ldsab, &
817 n1, n2, natorb_a, natorb_b, ncoa, ncob, nkind, nlist, nneighbor, nseta, nsetb, sgfa, sgfb
818 INTEGER,
DIMENSION(3) :: cell
819 INTEGER,
DIMENSION(:),
POINTER :: la_max, la_min, lb_max, lb_min, npgfa, &
821 INTEGER,
DIMENSION(:, :),
POINTER :: first_sgfa, first_sgfb
824 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: owork
825 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: oint, sint
826 REAL(kind=
dp),
DIMENSION(3) :: rij
827 REAL(kind=
dp),
DIMENSION(:),
POINTER :: set_radius_a, set_radius_b
828 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: rpgfa, rpgfb, scon_a, scon_b, zeta, zetb
834 DIMENSION(:),
POINTER :: nl_iterator
837 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
840 CALL timeset(routinen, handle)
843 cpassert(.NOT.
ASSOCIATED(sap_int))
844 ALLOCATE (sap_int(nkind*nkind))
845 DO i = 1, nkind*nkind
846 NULLIFY (sap_int(i)%alist, sap_int(i)%asort, sap_int(i)%aindex)
847 sap_int(i)%nalist = 0
851 qs_kind_set=qs_kind_set, &
852 dft_control=dft_control, &
856 ALLOCATE (basis_set_list(nkind))
863 jatom=jatom, nlist=nlist, ilist=ilist, nnode=nneighbor, &
864 inode=jneighbor, cell=cell, r=rij)
865 iac = ikind + nkind*(jkind - 1)
867 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
869 IF (.NOT. defined .OR. natorb_a < 1) cycle
870 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
872 IF (.NOT. defined .OR. natorb_b < 1) cycle
874 dr = sqrt(sum(rij(:)**2))
877 IF (.NOT.
ASSOCIATED(sap_int(iac)%alist))
THEN
878 sap_int(iac)%a_kind = ikind
879 sap_int(iac)%p_kind = jkind
880 sap_int(iac)%nalist = nlist
881 ALLOCATE (sap_int(iac)%alist(nlist))
883 NULLIFY (sap_int(iac)%alist(i)%clist)
884 sap_int(iac)%alist(i)%aatom = 0
885 sap_int(iac)%alist(i)%nclist = 0
888 IF (.NOT.
ASSOCIATED(sap_int(iac)%alist(ilist)%clist))
THEN
889 sap_int(iac)%alist(ilist)%aatom = iatom
890 sap_int(iac)%alist(ilist)%nclist = nneighbor
891 ALLOCATE (sap_int(iac)%alist(ilist)%clist(nneighbor))
893 sap_int(iac)%alist(ilist)%clist(i)%catom = 0
896 clist => sap_int(iac)%alist(ilist)%clist(jneighbor)
900 ALLOCATE (clist%acint(natorb_a, natorb_b, 3))
901 NULLIFY (clist%achint)
904 NULLIFY (clist%sgf_list)
907 basis_set_a => basis_set_list(ikind)%gto_basis_set
908 IF (.NOT.
ASSOCIATED(basis_set_a)) cycle
909 basis_set_b => basis_set_list(jkind)%gto_basis_set
910 IF (.NOT.
ASSOCIATED(basis_set_b)) cycle
912 first_sgfa => basis_set_a%first_sgf
913 la_max => basis_set_a%lmax
914 la_min => basis_set_a%lmin
915 npgfa => basis_set_a%npgf
916 nseta = basis_set_a%nset
917 nsgfa => basis_set_a%nsgf_set
918 rpgfa => basis_set_a%pgf_radius
919 set_radius_a => basis_set_a%set_radius
920 scon_a => basis_set_a%scon
921 zeta => basis_set_a%zet
923 first_sgfb => basis_set_b%first_sgf
924 lb_max => basis_set_b%lmax
925 lb_min => basis_set_b%lmin
926 npgfb => basis_set_b%npgf
927 nsetb = basis_set_b%nset
928 nsgfb => basis_set_b%nsgf_set
929 rpgfb => basis_set_b%pgf_radius
930 set_radius_b => basis_set_b%set_radius
931 scon_b => basis_set_b%scon
932 zetb => basis_set_b%zet
935 ALLOCATE (oint(ldsab, ldsab, 4), owork(ldsab, ldsab))
936 ALLOCATE (sint(natorb_a, natorb_b, 4))
940 ncoa = npgfa(iset)*
ncoset(la_max(iset))
941 n1 = npgfa(iset)*(
ncoset(la_max(iset)) -
ncoset(la_min(iset) - 1))
942 sgfa = first_sgfa(1, iset)
944 IF (set_radius_a(iset) + set_radius_b(jset) < dr) cycle
945 ncob = npgfb(jset)*
ncoset(lb_max(jset))
946 n2 = npgfb(jset)*(
ncoset(lb_max(jset)) -
ncoset(lb_min(jset) - 1))
947 sgfb = first_sgfb(1, jset)
948 CALL overlap_ab(la_max(iset), la_min(iset), npgfa(iset), rpgfa(:, iset), zeta(:, iset), &
949 lb_max(jset), lb_min(jset), npgfb(jset), rpgfb(:, jset), zetb(:, jset), &
950 rij, sab=oint(:, :, 1), dab=oint(:, :, 2:4))
953 CALL contraction(oint(:, :, i), owork, ca=scon_a(:, sgfa:), na=n1, ma=nsgfa(iset), &
954 cb=scon_b(:, sgfb:), nb=n2, mb=nsgfb(jset), fscale=1.0_dp, trans=.false.)
955 CALL block_add(
"IN", owork, nsgfa(iset), nsgfb(jset), sint(:, :, i), &
956 sgfa, sgfb, trans=.false.)
961 clist%acint(1:natorb_a, 1:natorb_b, 1:3) = sint(1:natorb_a, 1:natorb_b, 2:4)
963 DEALLOCATE (oint, owork, sint)
968 DEALLOCATE (basis_set_list)
970 CALL timestop(handle)
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, mimic, 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, xcint_weights, 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.
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.