73#include "./base/base_uses.f90"
79 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_dftb_coulomb'
82 REAL(dp),
PARAMETER :: tol_gamma = 1.e-4_dp
84 REAL(dp),
PARAMETER :: rtiny = 1.e-10_dp
101 calculate_forces, just_energy)
104 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: ks_matrix
106 REAL(dp),
DIMENSION(:) :: mcharge
108 LOGICAL,
INTENT(in) :: calculate_forces, just_energy
110 CHARACTER(len=*),
PARAMETER :: routinen =
'build_dftb_coulomb'
112 INTEGER :: atom_i, atom_j, ewald_type, handle, i, ia, iac, iatom, ic, icol, ikind, img, &
113 irow, is, jatom, jkind, natom, natorb_a, natorb_b, nimg, nkind, nmat
114 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: atom_of_kind, kind_of
115 INTEGER,
DIMENSION(3) :: cellind, periodic
116 INTEGER,
DIMENSION(:, :, :),
POINTER :: cell_to_index
117 LOGICAL :: defined, do_ewald, do_gamma_stress, &
118 found, hb_sr_damp, use_virial
119 REAL(kind=dp) :: alpha, ddr, deth, dgam, dr, drm, drp, &
120 fi, ga, gb, gmat, gmij, hb_para, zeff
121 REAL(kind=dp),
ALLOCATABLE,
DIMENSION(:) :: xgamma, zeffk
122 REAL(kind=dp),
DIMENSION(0:3) :: eta_a, eta_b
123 REAL(kind=dp),
DIMENSION(3) :: fij, rij
124 REAL(kind=dp),
DIMENSION(:, :),
POINTER :: dsblock, gmcharge, ksblock, pblock, &
126 REAL(kind=dp),
DIMENSION(:, :, :),
POINTER :: dsint
131 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_p, matrix_s
139 DIMENSION(:),
POINTER :: nl_iterator
145 POINTER :: dftb_potential
147 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
151 CALL timeset(routinen, handle)
153 NULLIFY (matrix_p, matrix_s, virial, atprop, dft_control)
157 IF (calculate_forces)
THEN
163 CALL get_qs_env(qs_env, nkind=nkind, natom=natom)
164 ALLOCATE (gmcharge(natom, nmat))
168 particle_set=particle_set, &
172 dft_control=dft_control, &
173 atomic_kind_set=atomic_kind_set, &
174 qs_kind_set=qs_kind_set, &
175 force=force, para_env=para_env)
177 IF (calculate_forces)
THEN
178 use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
181 NULLIFY (dftb_potential)
182 CALL get_qs_env(qs_env=qs_env, dftb_potential=dftb_potential)
184 CALL get_qs_env(qs_env=qs_env, sab_orb=n_list)
188 iatom=iatom, jatom=jatom, r=rij, cell=cellind)
190 CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind_a)
192 defined=defined, eta=eta_a, natorb=natorb_a)
193 IF (.NOT. defined .OR. natorb_a < 1) cycle
194 CALL get_qs_kind(qs_kind_set(jkind), dftb_parameter=dftb_kind_b)
196 defined=defined, eta=eta_b, natorb=natorb_b)
197 IF (.NOT. defined .OR. natorb_b < 1) cycle
200 hb_sr_damp = dft_control%qs_control%dftb_control%hb_sr_damp
203 hb_sr_damp =
is_hydrogen(particle_set(iatom)%atomic_kind) .NEQV. &
207 hb_para = dft_control%qs_control%dftb_control%hb_sr_para
213 dr = sqrt(sum(rij(:)**2))
215 gmcharge(jatom, 1) = gmcharge(jatom, 1) + gmat*mcharge(iatom)
216 IF (iatom /= jatom)
THEN
217 gmcharge(iatom, 1) = gmcharge(iatom, 1) + gmat*mcharge(jatom)
219 IF (calculate_forces .AND. (iatom /= jatom .OR. dr > 0.001_dp))
THEN
220 ddr = 0.1_dp*dftb_potential(ikind, jkind)%dgrd
225 gmcharge(jatom, i + 1) = gmcharge(jatom, i + 1) - dgam*mcharge(iatom)*rij(i)/dr
226 IF (dr > 0.001_dp)
THEN
227 gmcharge(iatom, i + 1) = gmcharge(iatom, i + 1) + dgam*mcharge(jatom)*rij(i)/dr
230 fij(i) = -mcharge(iatom)*mcharge(jatom)*dgam*rij(i)/dr
235 IF (iatom == jatom) fi = 0.5_dp
242 IF (atprop%energy)
THEN
243 CALL get_qs_env(qs_env=qs_env, particle_set=particle_set)
244 natom =
SIZE(particle_set)
249 do_ewald = dft_control%qs_control%dftb_control%do_ewald
252 NULLIFY (ewald_env, ewald_pw)
254 ewald_env=ewald_env, ewald_pw=ewald_pw)
255 CALL get_cell(cell=cell, periodic=periodic, deth=deth)
256 CALL ewald_env_get(ewald_env, alpha=alpha, ewald_type=ewald_type)
257 CALL get_qs_env(qs_env=qs_env, sab_tbe=n_list)
259 SELECT CASE (ewald_type)
261 cpabort(
"Invalid Ewald type")
263 cpabort(
"Not allowed with DFTB")
265 cpabort(
"Standard Ewald not implemented in DFTB")
267 cpabort(
"PME not implemented in DFTB")
270 gmcharge, mcharge, calculate_forces, virial, use_virial)
275 local_particles=local_particles)
276 DO ikind = 1,
SIZE(local_particles%n_el)
277 DO ia = 1, local_particles%n_el(ikind)
278 iatom = local_particles%list(ikind)%array(ia)
279 DO jatom = 1, iatom - 1
280 rij = particle_set(iatom)%r - particle_set(jatom)%r
282 dr = sqrt(sum(rij(:)**2))
283 gmcharge(iatom, 1) = gmcharge(iatom, 1) + mcharge(jatom)/dr
284 gmcharge(jatom, 1) = gmcharge(jatom, 1) + mcharge(iatom)/dr
286 gmcharge(iatom, i) = gmcharge(iatom, i) + rij(i - 1)*mcharge(jatom)/dr**3
287 gmcharge(jatom, i) = gmcharge(jatom, i) - rij(i - 1)*mcharge(iatom)/dr**3
292 cpassert(.NOT. use_virial)
295 CALL para_env%sum(gmcharge(:, 1))
299 gmcharge(:, 1) = gmcharge(:, 1) - 2._dp*alpha*
oorootpi*mcharge(:)
300 IF (any(periodic(:) == 1))
THEN
301 gmcharge(:, 1) = gmcharge(:, 1) -
pi/alpha**2/deth
305 energy%hartree = energy%hartree + 0.5_dp*sum(mcharge(:)*gmcharge(:, 1))
306 IF (atprop%energy)
THEN
308 local_particles=local_particles)
309 DO ikind = 1,
SIZE(local_particles%n_el)
310 CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind)
312 DO ia = 1, local_particles%n_el(ikind)
313 iatom = local_particles%list(ikind)%array(ia)
314 atprop%atecoul(iatom) = atprop%atecoul(iatom) + &
315 0.5_dp*zeff*gmcharge(iatom, 1)
320 IF (calculate_forces)
THEN
323 atom_of_kind=atom_of_kind)
325 gmcharge(:, 2) = gmcharge(:, 2)*mcharge(:)
326 gmcharge(:, 3) = gmcharge(:, 3)*mcharge(:)
327 gmcharge(:, 4) = gmcharge(:, 4)*mcharge(:)
329 ikind = kind_of(iatom)
330 atom_i = atom_of_kind(iatom)
331 force(ikind)%rho_elec(1, atom_i) = force(ikind)%rho_elec(1, atom_i) - gmcharge(iatom, 2)
332 force(ikind)%rho_elec(2, atom_i) = force(ikind)%rho_elec(2, atom_i) - gmcharge(iatom, 3)
333 force(ikind)%rho_elec(3, atom_i) = force(ikind)%rho_elec(3, atom_i) - gmcharge(iatom, 4)
337 do_gamma_stress = .false.
338 IF (.NOT. just_energy .AND. use_virial)
THEN
339 IF (dft_control%nimages == 1) do_gamma_stress = .true.
342 IF (.NOT. just_energy)
THEN
343 CALL get_qs_env(qs_env=qs_env, matrix_s_kp=matrix_s)
346 nimg = dft_control%nimages
347 NULLIFY (cell_to_index)
350 CALL get_qs_env(qs_env=qs_env, kpoints=kpoints)
354 IF (calculate_forces .AND.
SIZE(matrix_p, 1) == 2)
THEN
356 CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
357 alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
362 IF (do_gamma_stress)
THEN
364 CALL dftb_dsint_list(qs_env, sap_int)
372 gmij = 0.5_dp*(gmcharge(irow, 1) + gmcharge(icol, 1))
373 DO is = 1,
SIZE(ks_matrix, 1)
376 row=irow, col=icol, block=ksblock, found=found)
378 ksblock = ksblock - gmij*sblock
380 IF (calculate_forces)
THEN
381 ikind = kind_of(irow)
382 atom_i = atom_of_kind(irow)
383 jkind = kind_of(icol)
384 atom_j = atom_of_kind(icol)
387 row=irow, col=icol, block=pblock, found=found)
392 row=irow, col=icol, block=dsblock, found=found)
394 fi = -2.0_dp*gmij*sum(pblock*dsblock)
395 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
396 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
403 IF (do_gamma_stress)
THEN
406 iac = ikind + nkind*(jkind - 1)
407 IF (.NOT.
ASSOCIATED(sap_int(iac)%alist)) cycle
408 DO ia = 1, sap_int(iac)%nalist
409 IF (.NOT.
ASSOCIATED(sap_int(iac)%alist(ia)%clist)) cycle
410 iatom = sap_int(iac)%alist(ia)%aatom
411 DO ic = 1, sap_int(iac)%alist(ia)%nclist
412 jatom = sap_int(iac)%alist(ia)%clist(ic)%catom
413 rij = sap_int(iac)%alist(ia)%clist(ic)%rac
414 dr = sqrt(sum(rij(:)**2))
415 IF (dr > 1.e-6_dp)
THEN
416 dsint => sap_int(iac)%alist(ia)%clist(ic)%acint
417 gmij = 0.5_dp*(gmcharge(iatom, 1) + gmcharge(jatom, 1))
418 icol = max(iatom, jatom)
419 irow = min(iatom, jatom)
422 row=irow, col=icol, block=pblock, found=found)
425 IF (irow == iatom)
THEN
426 fij(i) = -2.0_dp*gmij*sum(pblock*dsint(:, :, i))
428 fij(i) = -2.0_dp*gmij*sum(transpose(pblock)*dsint(:, :, i))
432 IF (iatom == jatom) fi = 0.5_dp
442 CALL get_qs_env(qs_env=qs_env, sab_orb=n_list)
446 iatom=iatom, jatom=jatom, r=rij, cell=cellind)
448 icol = max(iatom, jatom)
449 irow = min(iatom, jatom)
450 ic = cell_to_index(cellind(1), cellind(2), cellind(3))
453 gmij = 0.5_dp*(gmcharge(iatom, 1) + gmcharge(jatom, 1))
457 row=irow, col=icol, block=sblock, found=found)
459 DO is = 1,
SIZE(ks_matrix, 1)
462 row=irow, col=icol, block=ksblock, found=found)
464 ksblock = ksblock - gmij*sblock
467 IF (calculate_forces)
THEN
468 ikind = kind_of(iatom)
469 atom_i = atom_of_kind(iatom)
470 jkind = kind_of(jatom)
471 atom_j = atom_of_kind(jatom)
472 IF (irow == jatom) gmij = -gmij
475 row=irow, col=icol, block=pblock, found=found)
480 row=irow, col=icol, block=dsblock, found=found)
482 fi = -2.0_dp*gmij*sum(pblock*dsblock)
483 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
484 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
489 IF (iatom == jatom) fi = 0.5_dp
497 IF (calculate_forces .AND.
SIZE(matrix_p, 1) == 2)
THEN
499 CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
500 alpha_scalar=1.0_dp, beta_scalar=-1.0_dp)
505 IF (dft_control%qs_control%dftb_control%dftb3_diagonal)
THEN
507 ALLOCATE (zeffk(nkind), xgamma(nkind))
509 CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind)
514 sap_int, calculate_forces, just_energy)
515 DEALLOCATE (zeffk, xgamma)
519 IF (qs_env%qmmm .AND. qs_env%qmmm_periodic)
THEN
521 calculate_forces, just_energy)
524 IF (do_gamma_stress)
THEN
528 DEALLOCATE (gmcharge)
530 CALL timestop(handle)
550 REAL(dp),
INTENT(in) :: r, ga, gb, hb_para
553 REAL(dp) :: a, b,
fac, g_sum
559 IF (g_sum < tol_gamma)
RETURN
562 gamma = 0.5_dp*(a*b/g_sum + (a*b)**2/g_sum**3)
569 IF (abs(a - b) < rtiny)
THEN
570 fac = 1.6_dp*r*a*b/g_sum*(1.0_dp + a*b/g_sum**2)
573 gamma = -exp(-a*r)*(0.5_dp*a*b**4/(a**2 - b**2)**2 - &
574 (b**6 - 3._dp*a**2*b**4)/(r*(a**2 - b**2)**3)) - &
575 exp(-b*r)*(0.5_dp*b*a**4/(b**2 - a**2)**2 - &
576 (a**6 - 3._dp*b**2*a**4)/(r*(b**2 - a**2)**3))
583 IF (hb_para > 0.0_dp)
THEN
584 gamma =
gamma*exp(-(0.5_dp*(ga + gb))**hb_para*r*r)
593 SUBROUTINE dftb_dsint_list(qs_env, sap_int)
598 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dftb_dsint_list'
600 INTEGER :: handle, i, iac, iatom, ikind, ilist, jatom, jkind, jneighbor, llm, lmaxi, lmaxj, &
601 n1, n2, natorb_a, natorb_b, ngrd, ngrdcut, nkind, nlist, nneighbor
602 INTEGER,
DIMENSION(3) :: cell
604 REAL(kind=dp) :: ddr, dgrd, dr
605 REAL(kind=dp),
DIMENSION(3) :: drij, rij
606 REAL(kind=dp),
DIMENSION(:, :),
POINTER :: dsblock, dsblockm, smatij, smatji
611 DIMENSION(:),
POINTER :: nl_iterator
616 POINTER :: dftb_potential
618 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
620 CALL timeset(routinen, handle)
623 cpassert(.NOT.
ASSOCIATED(sap_int))
624 ALLOCATE (sap_int(nkind*nkind))
625 DO i = 1, nkind*nkind
626 NULLIFY (sap_int(i)%alist, sap_int(i)%asort, sap_int(i)%aindex)
627 sap_int(i)%nalist = 0
631 qs_kind_set=qs_kind_set, &
632 dft_control=dft_control, &
635 dftb_control => dft_control%qs_control%dftb_control
637 NULLIFY (dftb_potential)
639 dftb_potential=dftb_potential)
645 jatom=jatom, nlist=nlist, ilist=ilist, nnode=nneighbor, &
646 inode=jneighbor, cell=cell, r=rij)
647 iac = ikind + nkind*(jkind - 1)
649 CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind_a)
651 defined=defined, lmax=lmaxi, natorb=natorb_a)
652 IF (.NOT. defined .OR. natorb_a < 1) cycle
653 CALL get_qs_kind(qs_kind_set(jkind), dftb_parameter=dftb_kind_b)
655 defined=defined, lmax=lmaxj, natorb=natorb_b)
656 IF (.NOT. defined .OR. natorb_b < 1) cycle
658 dr = sqrt(sum(rij(:)**2))
661 IF (.NOT.
ASSOCIATED(sap_int(iac)%alist))
THEN
662 sap_int(iac)%a_kind = ikind
663 sap_int(iac)%p_kind = jkind
664 sap_int(iac)%nalist = nlist
665 ALLOCATE (sap_int(iac)%alist(nlist))
667 NULLIFY (sap_int(iac)%alist(i)%clist)
668 sap_int(iac)%alist(i)%aatom = 0
669 sap_int(iac)%alist(i)%nclist = 0
672 IF (.NOT.
ASSOCIATED(sap_int(iac)%alist(ilist)%clist))
THEN
673 sap_int(iac)%alist(ilist)%aatom = iatom
674 sap_int(iac)%alist(ilist)%nclist = nneighbor
675 ALLOCATE (sap_int(iac)%alist(ilist)%clist(nneighbor))
677 sap_int(iac)%alist(ilist)%clist(i)%catom = 0
680 clist => sap_int(iac)%alist(ilist)%clist(jneighbor)
684 ALLOCATE (clist%acint(natorb_a, natorb_b, 3))
685 NULLIFY (clist%achint)
688 NULLIFY (clist%sgf_list)
691 dftb_param_ij => dftb_potential(ikind, jkind)
692 dftb_param_ji => dftb_potential(jkind, ikind)
694 ngrd = dftb_param_ij%ngrd
695 ngrdcut = dftb_param_ij%ngrdcut
696 dgrd = dftb_param_ij%dgrd
698 cpassert(dftb_param_ij%llm == dftb_param_ji%llm)
699 llm = dftb_param_ij%llm
700 smatij => dftb_param_ij%smat
701 smatji => dftb_param_ji%smat
703 IF (nint(dr/dgrd) <= ngrdcut)
THEN
704 IF (iatom == jatom .AND. dr < 0.001_dp)
THEN
710 ALLOCATE (dsblock(n1, n2), dsblockm(n1, n2))
716 drij(i) = rij(i) - ddr
717 CALL compute_block_sk(dsblockm, smatij, smatji, drij, ngrd, ngrdcut, dgrd, &
718 llm, lmaxi, lmaxj, iatom, iatom)
720 drij(i) = rij(i) + ddr
722 llm, lmaxi, lmaxj, iatom, iatom)
724 dsblock = dsblock - dsblockm
725 dsblock = dsblock/(2.0_dp*ddr)
727 clist%acint(1:n1, 1:n2, i) = -dsblock(1:n1, 1:n2)
729 DEALLOCATE (dsblock, dsblockm)
736 CALL timestop(handle)
738 END SUBROUTINE dftb_dsint_list
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.
elemental logical function, public is_hydrogen(atomic_kind)
Determines if the atomic_kind is HYDROGEN.
Holds information on atomic properties.
subroutine, public atprop_array_init(atarray, natom)
...
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.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
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_get_block_p(matrix, row, col, block, found, row_size, col_size)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_add(matrix_a, matrix_b, alpha_scalar, beta_scalar)
...
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
subroutine, public ewald_env_get(ewald_env, ewald_type, alpha, eps_pol, epsilon, gmax, ns_max, o_spline, group, para_env, poisson_section, precs, rcut, do_multipoles, max_multipole, do_ipol, max_ipol_iter, interaction_cutoffs, cell_hmat)
Purpose: Get the EWALD environment.
Calculation of Ewald contributions in DFTB.
subroutine, public tb_ewald_overlap(gmcharge, mcharge, alpha, n_list, virial, use_virial)
...
subroutine, public tb_spme_evaluate(ewald_env, ewald_pw, particle_set, box, gmcharge, mcharge, calculate_forces, virial, use_virial)
...
Calculation of the incomplete Gamma function F_n(t) for multi-center integrals over Cartesian Gaussia...
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.
Definition of mathematical constants and functions.
real(kind=dp), parameter, public oorootpi
real(kind=dp), parameter, public pi
real(kind=dp), dimension(0:maxfac), parameter, public fac
Interface to the message passing library MPI.
Define the data structure for the particle information.
functions related to the poisson solver on regular grids
integer, parameter, public do_ewald_pme
integer, parameter, public do_ewald_ewald
integer, parameter, public do_ewald_none
integer, parameter, public do_ewald_spme
Calculation of QMMM Coulomb contributions in TB.
subroutine, public build_tb_coulomb_qmqm(qs_env, ks_matrix, rho, mcharge, energy, calculate_forces, just_energy)
...
Calculation of DFTB3 Terms.
subroutine, public build_dftb3_diagonal(qs_env, ks_matrix, rho, mcharge, energy, xgamma, zeff, sap_int, calculate_forces, just_energy)
...
Calculation of Coulomb contributions in DFTB.
subroutine, public build_dftb_coulomb(qs_env, ks_matrix, rho, mcharge, energy, calculate_forces, just_energy)
...
real(dp) function, public gamma_rab_sr(r, ga, gb, hb_para)
Computes the short-range gamma parameter from exact Coulomb interaction of normalized exp(-a*r) charg...
Definition of the DFTB parameter types.
Working with the DFTB parameter types.
subroutine, public compute_block_sk(block, smatij, smatji, rij, ngrd, ngrdcut, dgrd, llm, lmaxi, lmaxj, irow, iatom)
...
subroutine, public get_dftb_atom_param(dftb_parameter, name, typ, defined, z, zeff, natorb, lmax, skself, occupation, eta, energy, cutoff, xi, di, rcdisp, dudq)
...
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zatom, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
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)
...
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public qs_rho_get(rho_struct, rho_ao, rho_ao_im, rho_ao_kp, rho_ao_im_kp, rho_r, drho_r, rho_g, drho_g, tau_r, tau_g, rho_r_valid, drho_r_valid, rho_g_valid, drho_g_valid, tau_r_valid, tau_g_valid, tot_rho_r, tot_rho_g, rho_r_sccs, soft_valid, complex_rho_ao)
returns info about the density described by this object. If some representation is not available an e...
General overlap type integrals containers.
subroutine, public release_sap_int(sap_int)
...
pure subroutine, public virial_pair_force(pv_virial, f0, force, rab)
Computes the contribution to the stress tensor from two-body pair-wise forces.
Provides all information about an atomic kind.
type for the atomic properties
Type defining parameters related to the simulation cell.
structure to store local (to a processor) ordered lists of integers.
to build arrays of pointers
Contains information about kpoints.
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.
keeps the density in various representations, keeping track of which ones are valid.