(git:ed6f26b)
Loading...
Searching...
No Matches
xtb_coulomb.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Calculation of Coulomb contributions in xTB
10!> \author JGH
11! **************************************************************************************************
13 USE ai_contraction, ONLY: block_add,&
15 USE ai_overlap, ONLY: overlap_ab
22 USE cell_types, ONLY: cell_type,&
23 get_cell,&
24 pbc
27 USE cp_dbcsr_api, ONLY: dbcsr_add,&
41 USE kinds, ONLY: dp
42 USE kpoint_types, ONLY: get_kpoint_info,&
44 USE mathconstants, ONLY: oorootpi,&
45 pi
47 USE orbital_pointers, ONLY: ncoset
61 USE qs_kind_types, ONLY: get_qs_kind,&
69 USE qs_rho_types, ONLY: qs_rho_get,&
71 USE sap_kind_types, ONLY: clist_type,&
75 USE virial_types, ONLY: virial_type
76 USE xtb_types, ONLY: get_xtb_atom_param,&
78#include "./base/base_uses.f90"
79
80 IMPLICIT NONE
81
82 PRIVATE
83
84 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'xtb_coulomb'
85
87
88CONTAINS
89
90! **************************************************************************************************
91!> \brief ...
92!> \param qs_env ...
93!> \param ks_matrix ...
94!> \param rho ...
95!> \param charges ...
96!> \param mcharge ...
97!> \param energy ...
98!> \param calculate_forces ...
99!> \param just_energy ...
100! **************************************************************************************************
101 SUBROUTINE build_xtb_coulomb(qs_env, ks_matrix, rho, charges, mcharge, energy, &
102 calculate_forces, just_energy)
103
104 TYPE(qs_environment_type), POINTER :: qs_env
105 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: ks_matrix
106 TYPE(qs_rho_type), POINTER :: rho
107 REAL(dp), DIMENSION(:, :), INTENT(in) :: charges
108 REAL(dp), DIMENSION(:), INTENT(in) :: mcharge
109 TYPE(qs_energy_type), POINTER :: energy
110 LOGICAL, INTENT(in) :: calculate_forces, just_energy
111
112 CHARACTER(len=*), PARAMETER :: routinen = 'build_xtb_coulomb'
113
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, &
116 nkind, nmat, za, zb
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, &
123 found, use_virial
124 REAL(kind=dp) :: alpha, deth, dr, ecsr, etaa, etab, f1, &
125 f2, fi, gmij, kg, rcut, rcuta, rcutb, &
126 zeff
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
135 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
136 TYPE(atprop_type), POINTER :: atprop
137 TYPE(cell_type), POINTER :: cell
138 TYPE(dbcsr_iterator_type) :: iter
139 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_p, matrix_s
140 TYPE(dft_control_type), POINTER :: dft_control
141 TYPE(distribution_1d_type), POINTER :: local_particles
142 TYPE(ewald_environment_type), POINTER :: ewald_env
143 TYPE(ewald_pw_type), POINTER :: ewald_pw
144 TYPE(kpoint_type), POINTER :: kpoints
145 TYPE(mp_para_env_type), POINTER :: para_env
147 DIMENSION(:), POINTER :: nl_iterator
148 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
149 POINTER :: n_list
150 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
151 TYPE(qs_force_type), DIMENSION(:), POINTER :: force
152 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
153 TYPE(sap_int_type), DIMENSION(:), POINTER :: sap_int
154 TYPE(virial_type), POINTER :: virial
155 TYPE(xtb_atom_type), POINTER :: xtb_atom_a, xtb_atom_b, xtb_kind
156 TYPE(xtb_control_type), POINTER :: xtb_control
157
158 CALL timeset(routinen, handle)
159
160 NULLIFY (matrix_p, matrix_s, virial, atprop, dft_control)
161
162 CALL get_qs_env(qs_env, &
163 qs_kind_set=qs_kind_set, &
164 particle_set=particle_set, &
165 cell=cell, &
166 virial=virial, &
167 atprop=atprop, &
168 dft_control=dft_control)
169
170 xtb_control => dft_control%qs_control%xtb_control
171
172 use_virial = .false.
173 IF (calculate_forces) THEN
174 use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
175 END IF
176
177 do_gamma_stress = .false.
178 IF (.NOT. just_energy .AND. use_virial) THEN
179 IF (dft_control%nimages == 1) do_gamma_stress = .true.
180 END IF
181
182 IF (atprop%energy) THEN
183 CALL get_qs_env(qs_env=qs_env, particle_set=particle_set)
184 natom = SIZE(particle_set)
185 CALL atprop_array_init(atprop%atecoul, natom)
186 END IF
187
188 IF (calculate_forces) THEN
189 nmat = 4
190 ELSE
191 nmat = 1
192 END IF
193
194 CALL get_qs_env(qs_env, nkind=nkind, natom=natom)
195 ALLOCATE (gchrg(natom, 5, nmat))
196 gchrg = 0._dp
197 ALLOCATE (gmcharge(natom, nmat))
198 gmcharge = 0._dp
199
200 ! short range contribution (gamma)
201 ! loop over all atom pairs (sab_xtbe)
202 kg = xtb_control%kg
203 NULLIFY (n_list)
204 CALL get_qs_env(qs_env=qs_env, sab_xtbe=n_list)
205 CALL neighbor_list_iterator_create(nl_iterator, n_list)
206 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
207 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
208 iatom=iatom, jatom=jatom, r=rij, cell=cellind)
209 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
210 CALL get_xtb_atom_param(xtb_atom_a, defined=defined, natorb=natorb_a)
211 IF (.NOT. defined .OR. natorb_a < 1) cycle
212 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
213 CALL get_xtb_atom_param(xtb_atom_b, defined=defined, natorb=natorb_b)
214 IF (.NOT. defined .OR. natorb_b < 1) cycle
215 ! atomic parameters
216 CALL get_xtb_atom_param(xtb_atom_a, eta=etaa, lmax=lmaxa, kappa=kappaa, rcut=rcuta)
217 CALL get_xtb_atom_param(xtb_atom_b, eta=etab, lmax=lmaxb, kappa=kappab, rcut=rcutb)
218 ! gamma matrix
219 ni = lmaxa + 1
220 nj = lmaxb + 1
221 ALLOCATE (gammab(ni, nj))
222 rcut = rcuta + rcutb
223 dr = sqrt(sum(rij(:)**2))
224 CALL gamma_rab_sr(gammab, dr, ni, kappaa, etaa, nj, kappab, etab, kg, rcut)
225 gchrg(iatom, 1:ni, 1) = gchrg(iatom, 1:ni, 1) + matmul(gammab, charges(jatom, 1:nj))
226 IF (iatom /= jatom) THEN
227 gchrg(jatom, 1:nj, 1) = gchrg(jatom, 1:nj, 1) + matmul(charges(iatom, 1:ni), gammab)
228 END IF
229 IF (calculate_forces) THEN
230 IF (dr > 1.e-6_dp) THEN
231 CALL dgamma_rab_sr(gammab, dr, ni, kappaa, etaa, nj, kappab, etab, kg, rcut)
232 DO i = 1, 3
233 gchrg(iatom, 1:ni, i + 1) = gchrg(iatom, 1:ni, i + 1) &
234 + matmul(gammab, charges(jatom, 1:nj))*rij(i)/dr
235 IF (iatom /= jatom) THEN
236 gchrg(jatom, 1:nj, i + 1) = gchrg(jatom, 1:nj, i + 1) &
237 - matmul(charges(iatom, 1:ni), gammab)*rij(i)/dr
238 END IF
239 END DO
240 IF (use_virial) THEN
241 gcint(1:ni) = matmul(gammab, charges(jatom, 1:nj))
242 DO i = 1, 3
243 fij(i) = -sum(charges(iatom, 1:ni)*gcint(1:ni))*rij(i)/dr
244 END DO
245 fi = 1.0_dp
246 IF (iatom == jatom) fi = 0.5_dp
247 CALL virial_pair_force(virial%pv_virial, fi, fij, rij)
248 END IF
249 END IF
250 END IF
251 DEALLOCATE (gammab)
252 END DO
253 CALL neighbor_list_iterator_release(nl_iterator)
254
255 ! 1/R contribution
256
257 IF (xtb_control%coulomb_lr) THEN
258 do_ewald = xtb_control%do_ewald
259 IF (do_ewald) THEN
260 ! Ewald sum
261 NULLIFY (ewald_env, ewald_pw)
262 CALL get_qs_env(qs_env=qs_env, &
263 ewald_env=ewald_env, ewald_pw=ewald_pw)
264 CALL get_cell(cell=cell, periodic=periodic, deth=deth)
265 CALL ewald_env_get(ewald_env, alpha=alpha, ewald_type=ewald_type)
266 CALL get_qs_env(qs_env=qs_env, sab_tbe=n_list)
267 CALL tb_ewald_overlap(gmcharge, mcharge, alpha, n_list, virial, use_virial)
268 SELECT CASE (ewald_type)
269 CASE DEFAULT
270 cpabort("Invalid Ewald type")
271 CASE (do_ewald_none)
272 cpabort("Not allowed with xTB/DFTB")
273 CASE (do_ewald_ewald)
274 cpabort("Standard Ewald not implemented in xTB/DFTB")
275 CASE (do_ewald_pme)
276 cpabort("PME not implemented in xTB/DFTB")
277 CASE (do_ewald_spme)
278 CALL tb_spme_evaluate(ewald_env, ewald_pw, particle_set, cell, &
279 gmcharge, mcharge, calculate_forces, virial, use_virial)
280 END SELECT
281 ELSE
282 ! direct sum
283 CALL get_qs_env(qs_env=qs_env, &
284 local_particles=local_particles)
285 DO ikind = 1, SIZE(local_particles%n_el)
286 DO ia = 1, local_particles%n_el(ikind)
287 iatom = local_particles%list(ikind)%array(ia)
288 DO jatom = 1, iatom - 1
289 rij = particle_set(iatom)%r - particle_set(jatom)%r
290 rij = pbc(rij, cell)
291 dr = sqrt(sum(rij(:)**2))
292 IF (dr > 1.e-6_dp) THEN
293 gmcharge(iatom, 1) = gmcharge(iatom, 1) + mcharge(jatom)/dr
294 gmcharge(jatom, 1) = gmcharge(jatom, 1) + mcharge(iatom)/dr
295 DO i = 2, nmat
296 gmcharge(iatom, i) = gmcharge(iatom, i) + rij(i - 1)*mcharge(jatom)/dr**3
297 gmcharge(jatom, i) = gmcharge(jatom, i) - rij(i - 1)*mcharge(iatom)/dr**3
298 END DO
299 END IF
300 END DO
301 END DO
302 END DO
303 cpassert(.NOT. use_virial)
304 END IF
305 END IF
306
307 ! global sum of gamma*p arrays
308 CALL get_qs_env(qs_env=qs_env, &
309 atomic_kind_set=atomic_kind_set, &
310 force=force, para_env=para_env)
311 CALL para_env%sum(gmcharge(:, 1))
312 CALL para_env%sum(gchrg(:, :, 1))
313
314 IF (xtb_control%coulomb_lr) THEN
315 IF (do_ewald) THEN
316 ! add self charge interaction and background charge contribution
317 gmcharge(:, 1) = gmcharge(:, 1) - 2._dp*alpha*oorootpi*mcharge(:)
318 IF (any(periodic(:) == 1)) THEN
319 gmcharge(:, 1) = gmcharge(:, 1) - pi/alpha**2/deth
320 END IF
321 END IF
322 END IF
323
324 ! energy
325 CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, &
326 kind_of=kind_of, &
327 atom_of_kind=atom_of_kind)
328 ecsr = 0.0_dp
329 DO iatom = 1, natom
330 ikind = kind_of(iatom)
331 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_kind)
332 CALL get_xtb_atom_param(xtb_kind, lmax=ni)
333 ni = ni + 1
334 ecsr = ecsr + sum(charges(iatom, 1:ni)*gchrg(iatom, 1:ni, 1))
335 END DO
336
337 energy%hartree = energy%hartree + 0.5_dp*ecsr
338 energy%hartree = energy%hartree + 0.5_dp*sum(mcharge(:)*gmcharge(:, 1))
339
340 IF (atprop%energy) THEN
341 CALL get_qs_env(qs_env=qs_env, local_particles=local_particles)
342 DO ikind = 1, SIZE(local_particles%n_el)
343 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_kind)
344 CALL get_xtb_atom_param(xtb_kind, lmax=ni, occupation=occ)
345 ni = ni + 1
346 zeff = sum(real(occ, kind=dp))
347 DO ia = 1, local_particles%n_el(ikind)
348 iatom = local_particles%list(ikind)%array(ia)
349 atprop%atecoul(iatom) = atprop%atecoul(iatom) + &
350 0.5_dp*sum(real(occ(1:ni), kind=dp)*gchrg(iatom, 1:ni, 1))
351 atprop%atecoul(iatom) = atprop%atecoul(iatom) + &
352 0.5_dp*zeff*gmcharge(iatom, 1)
353 END DO
354 END DO
355 END IF
356
357 IF (calculate_forces) THEN
358 DO iatom = 1, natom
359 ikind = kind_of(iatom)
360 atom_i = atom_of_kind(iatom)
361 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_kind)
362 CALL get_xtb_atom_param(xtb_kind, lmax=ni)
363 ! short range
364 ni = ni + 1
365 DO i = 1, 3
366 fij(i) = sum(charges(iatom, 1:ni)*gchrg(iatom, 1:ni, i + 1))
367 END DO
368 force(ikind)%rho_elec(1, atom_i) = force(ikind)%rho_elec(1, atom_i) - fij(1)
369 force(ikind)%rho_elec(2, atom_i) = force(ikind)%rho_elec(2, atom_i) - fij(2)
370 force(ikind)%rho_elec(3, atom_i) = force(ikind)%rho_elec(3, atom_i) - fij(3)
371 ! long range
372 DO i = 1, 3
373 fij(i) = gmcharge(iatom, i + 1)*mcharge(iatom)
374 END DO
375 force(ikind)%rho_elec(1, atom_i) = force(ikind)%rho_elec(1, atom_i) - fij(1)
376 force(ikind)%rho_elec(2, atom_i) = force(ikind)%rho_elec(2, atom_i) - fij(2)
377 force(ikind)%rho_elec(3, atom_i) = force(ikind)%rho_elec(3, atom_i) - fij(3)
378 END DO
379 END IF
380
381 IF (.NOT. just_energy) THEN
382 CALL get_qs_env(qs_env=qs_env, matrix_s_kp=matrix_s)
383 CALL qs_rho_get(rho, rho_ao_kp=matrix_p)
384
385 nimg = dft_control%nimages
386 NULLIFY (cell_to_index)
387 IF (nimg > 1) THEN
388 NULLIFY (kpoints)
389 CALL get_qs_env(qs_env=qs_env, kpoints=kpoints)
390 CALL get_kpoint_info(kpoint=kpoints, cell_to_index=cell_to_index)
391 END IF
392
393 IF (calculate_forces .AND. SIZE(matrix_p, 1) == 2) THEN
394 DO img = 1, nimg
395 CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
396 alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
397 END DO
398 END IF
399
400 NULLIFY (sap_int)
401 IF (do_gamma_stress) THEN
402 ! derivative overlap integral (non collapsed)
403 CALL xtb_dsint_list(qs_env, sap_int)
404 END IF
405
406 IF (nimg == 1) THEN
407 ! no k-points; all matrices have been transformed to periodic bsf
408 CALL dbcsr_iterator_start(iter, matrix_s(1, 1)%matrix)
409 DO WHILE (dbcsr_iterator_blocks_left(iter))
410 CALL dbcsr_iterator_next_block(iter, irow, icol, sblock)
411 ikind = kind_of(irow)
412 jkind = kind_of(icol)
413
414 ! atomic parameters
415 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
416 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
417 CALL get_xtb_atom_param(xtb_atom_a, z=za, lao=laoa)
418 CALL get_xtb_atom_param(xtb_atom_b, z=zb, lao=laob)
419
420 ni = SIZE(sblock, 1)
421 nj = SIZE(sblock, 2)
422 ALLOCATE (gcij(ni, nj))
423 DO i = 1, ni
424 DO j = 1, nj
425 la = laoa(i) + 1
426 lb = laob(j) + 1
427 gcij(i, j) = 0.5_dp*(gchrg(irow, la, 1) + gchrg(icol, lb, 1))
428 END DO
429 END DO
430 gmij = 0.5_dp*(gmcharge(irow, 1) + gmcharge(icol, 1))
431 DO is = 1, SIZE(ks_matrix, 1)
432 NULLIFY (ksblock)
433 CALL dbcsr_get_block_p(matrix=ks_matrix(is, 1)%matrix, &
434 row=irow, col=icol, block=ksblock, found=found)
435 cpassert(found)
436 ksblock = ksblock - gcij*sblock
437 ksblock = ksblock - gmij*sblock
438 END DO
439 IF (calculate_forces) THEN
440 atom_i = atom_of_kind(irow)
441 atom_j = atom_of_kind(icol)
442 NULLIFY (pblock)
443 CALL dbcsr_get_block_p(matrix=matrix_p(1, 1)%matrix, &
444 row=irow, col=icol, block=pblock, found=found)
445 cpassert(found)
446 DO i = 1, 3
447 NULLIFY (dsblock)
448 CALL dbcsr_get_block_p(matrix=matrix_s(1 + i, 1)%matrix, &
449 row=irow, col=icol, block=dsblock, found=found)
450 cpassert(found)
451 fij(i) = 0.0_dp
452 ! short range
453 fi = -2.0_dp*sum(pblock*dsblock*gcij)
454 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
455 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
456 fij(i) = fij(i) + fi
457 ! long range
458 fi = -2.0_dp*gmij*sum(pblock*dsblock)
459 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
460 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
461 fij(i) = fij(i) + fi
462 END DO
463 END IF
464 DEALLOCATE (gcij)
465 END DO
466 CALL dbcsr_iterator_stop(iter)
467 ! stress tensor (needs recalculation of overlap integrals)
468 IF (do_gamma_stress) THEN
469 DO ikind = 1, nkind
470 DO jkind = 1, nkind
471 iac = ikind + nkind*(jkind - 1)
472 IF (.NOT. ASSOCIATED(sap_int(iac)%alist)) cycle
473 ! atomic parameters
474 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
475 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
476 CALL get_xtb_atom_param(xtb_atom_a, lao=laoa, natorb=ni)
477 CALL get_xtb_atom_param(xtb_atom_b, lao=laob, natorb=nj)
478 DO ia = 1, sap_int(iac)%nalist
479 IF (.NOT. ASSOCIATED(sap_int(iac)%alist(ia)%clist)) cycle
480 iatom = sap_int(iac)%alist(ia)%aatom
481 DO ic = 1, sap_int(iac)%alist(ia)%nclist
482 jatom = sap_int(iac)%alist(ia)%clist(ic)%catom
483 rij = sap_int(iac)%alist(ia)%clist(ic)%rac
484 dr = sqrt(sum(rij(:)**2))
485 IF (dr > 1.e-6_dp) THEN
486 dsint => sap_int(iac)%alist(ia)%clist(ic)%acint
487 ALLOCATE (gcij(ni, nj))
488 DO i = 1, ni
489 DO j = 1, nj
490 la = laoa(i) + 1
491 lb = laob(j) + 1
492 gcij(i, j) = 0.5_dp*(gchrg(iatom, la, 1) + gchrg(jatom, lb, 1))
493 END DO
494 END DO
495 gmij = 0.5_dp*(gmcharge(iatom, 1) + gmcharge(jatom, 1))
496 icol = max(iatom, jatom)
497 irow = min(iatom, jatom)
498 NULLIFY (pblock)
499 CALL dbcsr_get_block_p(matrix=matrix_p(1, 1)%matrix, &
500 row=irow, col=icol, block=pblock, found=found)
501 cpassert(found)
502 fij = 0.0_dp
503 DO i = 1, 3
504 ! short/long range
505 IF (irow == iatom) THEN
506 f1 = -2.0_dp*sum(pblock*dsint(:, :, i)*gcij)
507 f2 = -2.0_dp*gmij*sum(pblock*dsint(:, :, i))
508 ELSE
509 f1 = -2.0_dp*sum(transpose(pblock)*dsint(:, :, i)*gcij)
510 f2 = -2.0_dp*gmij*sum(transpose(pblock)*dsint(:, :, i))
511 END IF
512 fij(i) = f1 + f2
513 END DO
514 DEALLOCATE (gcij)
515 fi = 1.0_dp
516 IF (iatom == jatom) fi = 0.5_dp
517 CALL virial_pair_force(virial%pv_virial, fi, fij, rij)
518 END IF
519 END DO
520 END DO
521 END DO
522 END DO
523 END IF
524 ELSE
525 NULLIFY (n_list)
526 CALL get_qs_env(qs_env=qs_env, sab_orb=n_list)
527 CALL neighbor_list_iterator_create(nl_iterator, n_list)
528 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
529 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
530 iatom=iatom, jatom=jatom, r=rij, cell=cellind)
531
532 icol = max(iatom, jatom)
533 irow = min(iatom, jatom)
534
535 ic = cell_to_index(cellind(1), cellind(2), cellind(3))
536 cpassert(ic > 0)
537
538 NULLIFY (sblock)
539 CALL dbcsr_get_block_p(matrix=matrix_s(1, ic)%matrix, &
540 row=irow, col=icol, block=sblock, found=found)
541 cpassert(found)
542
543 ! atomic parameters
544 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
545 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
546 CALL get_xtb_atom_param(xtb_atom_a, z=za, lao=laoa)
547 CALL get_xtb_atom_param(xtb_atom_b, z=zb, lao=laob)
548
549 ni = SIZE(sblock, 1)
550 nj = SIZE(sblock, 2)
551 ALLOCATE (gcij(ni, nj))
552 DO i = 1, ni
553 DO j = 1, nj
554 IF (irow == iatom) THEN
555 la = laoa(i) + 1
556 lb = laob(j) + 1
557 gcij(i, j) = 0.5_dp*(gchrg(iatom, la, 1) + gchrg(jatom, lb, 1))
558 ELSE
559 la = laoa(j) + 1
560 lb = laob(i) + 1
561 gcij(i, j) = 0.5_dp*(gchrg(iatom, la, 1) + gchrg(jatom, lb, 1))
562 END IF
563 END DO
564 END DO
565 gmij = 0.5_dp*(gmcharge(iatom, 1) + gmcharge(jatom, 1))
566 DO is = 1, SIZE(ks_matrix, 1)
567 NULLIFY (ksblock)
568 CALL dbcsr_get_block_p(matrix=ks_matrix(is, ic)%matrix, &
569 row=irow, col=icol, block=ksblock, found=found)
570 cpassert(found)
571 ksblock = ksblock - gcij*sblock
572 ksblock = ksblock - gmij*sblock
573 END DO
574
575 IF (calculate_forces) THEN
576 atom_i = atom_of_kind(iatom)
577 atom_j = atom_of_kind(jatom)
578 IF (irow /= iatom) THEN
579 gmij = -gmij
580 gcij = -gcij
581 END IF
582 NULLIFY (pblock)
583 CALL dbcsr_get_block_p(matrix=matrix_p(1, ic)%matrix, &
584 row=irow, col=icol, block=pblock, found=found)
585 cpassert(found)
586 DO i = 1, 3
587 NULLIFY (dsblock)
588 CALL dbcsr_get_block_p(matrix=matrix_s(1 + i, ic)%matrix, &
589 row=irow, col=icol, block=dsblock, found=found)
590 cpassert(found)
591 fij(i) = 0.0_dp
592 ! short range
593 fi = -2.0_dp*sum(pblock*dsblock*gcij)
594 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
595 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
596 fij(i) = fij(i) + fi
597 ! long range
598 fi = -2.0_dp*gmij*sum(pblock*dsblock)
599 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
600 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
601 fij(i) = fij(i) + fi
602 END DO
603 IF (use_virial) THEN
604 dr = sqrt(sum(rij(:)**2))
605 IF (dr > 1.e-6_dp) THEN
606 fi = 1.0_dp
607 IF (iatom == jatom) fi = 0.5_dp
608 CALL virial_pair_force(virial%pv_virial, fi, fij, rij)
609 END IF
610 END IF
611 END IF
612 DEALLOCATE (gcij)
613
614 END DO
615 CALL neighbor_list_iterator_release(nl_iterator)
616 END IF
617
618 IF (calculate_forces .AND. SIZE(matrix_p, 1) == 2) THEN
619 DO img = 1, nimg
620 CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
621 alpha_scalar=1.0_dp, beta_scalar=-1.0_dp)
622 END DO
623 END IF
624 END IF
625
626 IF (xtb_control%tb3_interaction) THEN
627 CALL get_qs_env(qs_env, nkind=nkind)
628 ALLOCATE (zeffk(nkind), xgamma(nkind))
629 DO ikind = 1, nkind
630 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_kind)
631 CALL get_xtb_atom_param(xtb_kind, xgamma=xgamma(ikind), zeff=zeffk(ikind))
632 END DO
633 ! Diagonal 3rd order correction (DFTB3)
634 CALL build_dftb3_diagonal(qs_env, ks_matrix, rho, mcharge, energy, xgamma, zeffk, &
635 sap_int, calculate_forces, just_energy)
636 DEALLOCATE (zeffk, xgamma)
637 END IF
638
639 ! QMMM
640 IF (qs_env%qmmm .AND. qs_env%qmmm_periodic) THEN
641 CALL build_tb_coulomb_qmqm(qs_env, ks_matrix, rho, mcharge, energy, &
642 calculate_forces, just_energy)
643 END IF
644
645 IF (do_gamma_stress) THEN
646 CALL release_sap_int(sap_int)
647 END IF
648
649 CALL timestop(handle)
650
651 END SUBROUTINE build_xtb_coulomb
652
653! **************************************************************************************************
654!> \brief Computes the short-range gamma parameter from
655!> Nataga-Mishimoto-Ohno-Klopman formula for xTB
656!> WARNING: The xTB function (gamma - 1/r) has still an l-dependent longrange (1/r^3)
657!> behaviour. We use a cutoff function to smoothly remove this part.
658!> However, this will change energies and effect final results.
659!>
660!> \param gmat ...
661!> \param rab ...
662!> \param nla ...
663!> \param kappaa ...
664!> \param etaa ...
665!> \param nlb ...
666!> \param kappab ...
667!> \param etab ...
668!> \param kg ...
669!> \param rcut ...
670!> \par History
671!> 10.2018 JGH
672!> \version 1.1
673! **************************************************************************************************
674 SUBROUTINE gamma_rab_sr(gmat, rab, nla, kappaa, etaa, nlb, kappab, etab, kg, rcut)
675 REAL(dp), DIMENSION(:, :), INTENT(INOUT) :: gmat
676 REAL(dp), INTENT(IN) :: rab
677 INTEGER, INTENT(IN) :: nla
678 REAL(dp), DIMENSION(:), INTENT(IN) :: kappaa
679 REAL(dp), INTENT(IN) :: etaa
680 INTEGER, INTENT(IN) :: nlb
681 REAL(dp), DIMENSION(:), INTENT(IN) :: kappab
682 REAL(dp), INTENT(IN) :: etab, kg, rcut
683
684 REAL(kind=dp), PARAMETER :: rsmooth = 1.0_dp
685
686 INTEGER :: i, j
687 REAL(kind=dp) :: fcut, r, rk, x
688 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: eta
689
690 ALLOCATE (eta(nla, nlb))
691 eta = 0.0_dp
692
693 DO j = 1, nlb
694 DO i = 1, nla
695 eta(i, j) = 1._dp/(etaa*(1._dp + kappaa(i))) + 1._dp/(etab*(1._dp + kappab(j)))
696 eta(i, j) = 2._dp/eta(i, j)
697 END DO
698 END DO
699
700 gmat = 0.0_dp
701 IF (rab < 1.e-6_dp) THEN
702 ! on site terms
703 gmat(:, :) = eta(:, :)
704 ELSEIF (rab > rcut) THEN
705 ! do nothing
706 ELSE
707 rk = rab**kg
708 eta = eta**(-kg)
709 IF (rab < rcut - rsmooth) THEN
710 fcut = 1.0_dp
711 ELSE
712 r = rab - (rcut - rsmooth)
713 x = r/rsmooth
714 fcut = -6._dp*x**5 + 15._dp*x**4 - 10._dp*x**3 + 1._dp
715 END IF
716 gmat(:, :) = fcut*(1._dp/(rk + eta(:, :)))**(1._dp/kg) - fcut/rab
717 END IF
718
719 DEALLOCATE (eta)
720
721 END SUBROUTINE gamma_rab_sr
722
723! **************************************************************************************************
724!> \brief Computes the derivative of the short-range gamma parameter from
725!> Nataga-Mishimoto-Ohno-Klopman formula for xTB
726!> WARNING: The xTB function (gamma - 1/r) has still an l-dependent longrange (1/r^3)
727!> behaviour. We use a cutoff function to smoothly remove this part.
728!> However, this will change energies and effect final results.
729!>
730!> \param dgmat ...
731!> \param rab ...
732!> \param nla ...
733!> \param kappaa ...
734!> \param etaa ...
735!> \param nlb ...
736!> \param kappab ...
737!> \param etab ...
738!> \param kg ...
739!> \param rcut ...
740!> \par History
741!> 10.2018 JGH
742!> \version 1.1
743! **************************************************************************************************
744 SUBROUTINE dgamma_rab_sr(dgmat, rab, nla, kappaa, etaa, nlb, kappab, etab, kg, rcut)
745 REAL(dp), DIMENSION(:, :), INTENT(INOUT) :: dgmat
746 REAL(dp), INTENT(IN) :: rab
747 INTEGER, INTENT(IN) :: nla
748 REAL(dp), DIMENSION(:), INTENT(IN) :: kappaa
749 REAL(dp), INTENT(IN) :: etaa
750 INTEGER, INTENT(IN) :: nlb
751 REAL(dp), DIMENSION(:), INTENT(IN) :: kappab
752 REAL(dp), INTENT(IN) :: etab, kg, rcut
753
754 REAL(kind=dp), PARAMETER :: rsmooth = 1.0_dp
755
756 INTEGER :: i, j
757 REAL(kind=dp) :: dfcut, fcut, r, rk, x
758 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: eta
759
760 ALLOCATE (eta(nla, nlb))
761
762 DO j = 1, nlb
763 DO i = 1, nla
764 eta(i, j) = 1._dp/(etaa*(1._dp + kappaa(i))) + 1._dp/(etab*(1._dp + kappab(j)))
765 eta(i, j) = 2._dp/eta(i, j)
766 END DO
767 END DO
768
769 IF (rab < 1.e-6) THEN
770 ! on site terms
771 dgmat(:, :) = 0.0_dp
772 ELSEIF (rab > rcut) THEN
773 dgmat(:, :) = 0.0_dp
774 ELSE
775 eta = eta**(-kg)
776 rk = rab**kg
777 IF (rab < rcut - rsmooth) THEN
778 fcut = 1.0_dp
779 dfcut = 0.0_dp
780 ELSE
781 r = rab - (rcut - rsmooth)
782 x = r/rsmooth
783 fcut = -6._dp*x**5 + 15._dp*x**4 - 10._dp*x**3 + 1._dp
784 dfcut = -30._dp*x**4 + 60._dp*x**3 - 30._dp*x**2
785 dfcut = dfcut/rsmooth
786 END IF
787 dgmat(:, :) = dfcut*(1._dp/(rk + eta(:, :)))**(1._dp/kg)
788 dgmat(:, :) = dgmat(:, :) - dfcut/rab + fcut/rab**2
789 dgmat(:, :) = dgmat(:, :) - fcut/(rk + eta(:, :))*(1._dp/(rk + eta(:, :)))**(1._dp/kg)*rk/rab
790 END IF
791
792 DEALLOCATE (eta)
793
794 END SUBROUTINE dgamma_rab_sr
795
796! **************************************************************************************************
797!> \brief ...
798!> \param qs_env ...
799!> \param sap_int ...
800! **************************************************************************************************
801 SUBROUTINE xtb_dsint_list(qs_env, sap_int)
802
803 TYPE(qs_environment_type), POINTER :: qs_env
804 TYPE(sap_int_type), DIMENSION(:), POINTER :: sap_int
805
806 CHARACTER(LEN=*), PARAMETER :: routinen = 'xtb_dsint_list'
807
808 INTEGER :: handle, i, iac, iatom, ikind, ilist, iset, jatom, jkind, jneighbor, jset, ldsab, &
809 n1, n2, natorb_a, natorb_b, ncoa, ncob, nkind, nlist, nneighbor, nseta, nsetb, sgfa, sgfb
810 INTEGER, DIMENSION(3) :: cell
811 INTEGER, DIMENSION(:), POINTER :: la_max, la_min, lb_max, lb_min, npgfa, &
812 npgfb, nsgfa, nsgfb
813 INTEGER, DIMENSION(:, :), POINTER :: first_sgfa, first_sgfb
814 LOGICAL :: defined
815 REAL(kind=dp) :: dr
816 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: owork
817 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: oint, sint
818 REAL(kind=dp), DIMENSION(3) :: rij
819 REAL(kind=dp), DIMENSION(:), POINTER :: set_radius_a, set_radius_b
820 REAL(kind=dp), DIMENSION(:, :), POINTER :: rpgfa, rpgfb, scon_a, scon_b, zeta, zetb
821 TYPE(clist_type), POINTER :: clist
822 TYPE(dft_control_type), POINTER :: dft_control
823 TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER :: basis_set_list
824 TYPE(gto_basis_set_type), POINTER :: basis_set_a, basis_set_b
826 DIMENSION(:), POINTER :: nl_iterator
827 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
828 POINTER :: sab_orb
829 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
830 TYPE(xtb_atom_type), POINTER :: xtb_atom_a, xtb_atom_b
831
832 CALL timeset(routinen, handle)
833
834 CALL get_qs_env(qs_env=qs_env, nkind=nkind)
835 cpassert(.NOT. ASSOCIATED(sap_int))
836 ALLOCATE (sap_int(nkind*nkind))
837 DO i = 1, nkind*nkind
838 NULLIFY (sap_int(i)%alist, sap_int(i)%asort, sap_int(i)%aindex)
839 sap_int(i)%nalist = 0
840 END DO
841
842 CALL get_qs_env(qs_env=qs_env, &
843 qs_kind_set=qs_kind_set, &
844 dft_control=dft_control, &
845 sab_orb=sab_orb)
846
847 ! set up basis set lists
848 ALLOCATE (basis_set_list(nkind))
849 CALL basis_set_list_setup(basis_set_list, "ORB", qs_kind_set)
850
851 ! loop over all atom pairs with a non-zero overlap (sab_orb)
852 CALL neighbor_list_iterator_create(nl_iterator, sab_orb)
853 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
854 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, iatom=iatom, &
855 jatom=jatom, nlist=nlist, ilist=ilist, nnode=nneighbor, &
856 inode=jneighbor, cell=cell, r=rij)
857 iac = ikind + nkind*(jkind - 1)
858 !
859 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
860 CALL get_xtb_atom_param(xtb_atom_a, defined=defined, natorb=natorb_a)
861 IF (.NOT. defined .OR. natorb_a < 1) cycle
862 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
863 CALL get_xtb_atom_param(xtb_atom_b, defined=defined, natorb=natorb_b)
864 IF (.NOT. defined .OR. natorb_b < 1) cycle
865
866 dr = sqrt(sum(rij(:)**2))
867
868 ! integral list
869 IF (.NOT. ASSOCIATED(sap_int(iac)%alist)) THEN
870 sap_int(iac)%a_kind = ikind
871 sap_int(iac)%p_kind = jkind
872 sap_int(iac)%nalist = nlist
873 ALLOCATE (sap_int(iac)%alist(nlist))
874 DO i = 1, nlist
875 NULLIFY (sap_int(iac)%alist(i)%clist)
876 sap_int(iac)%alist(i)%aatom = 0
877 sap_int(iac)%alist(i)%nclist = 0
878 END DO
879 END IF
880 IF (.NOT. ASSOCIATED(sap_int(iac)%alist(ilist)%clist)) THEN
881 sap_int(iac)%alist(ilist)%aatom = iatom
882 sap_int(iac)%alist(ilist)%nclist = nneighbor
883 ALLOCATE (sap_int(iac)%alist(ilist)%clist(nneighbor))
884 DO i = 1, nneighbor
885 sap_int(iac)%alist(ilist)%clist(i)%catom = 0
886 END DO
887 END IF
888 clist => sap_int(iac)%alist(ilist)%clist(jneighbor)
889 clist%catom = jatom
890 clist%cell = cell
891 clist%rac = rij
892 ALLOCATE (clist%acint(natorb_a, natorb_b, 3))
893 NULLIFY (clist%achint)
894 clist%acint = 0._dp
895 clist%nsgf_cnt = 0
896 NULLIFY (clist%sgf_list)
897
898 ! overlap
899 basis_set_a => basis_set_list(ikind)%gto_basis_set
900 IF (.NOT. ASSOCIATED(basis_set_a)) cycle
901 basis_set_b => basis_set_list(jkind)%gto_basis_set
902 IF (.NOT. ASSOCIATED(basis_set_b)) cycle
903 ! basis ikind
904 first_sgfa => basis_set_a%first_sgf
905 la_max => basis_set_a%lmax
906 la_min => basis_set_a%lmin
907 npgfa => basis_set_a%npgf
908 nseta = basis_set_a%nset
909 nsgfa => basis_set_a%nsgf_set
910 rpgfa => basis_set_a%pgf_radius
911 set_radius_a => basis_set_a%set_radius
912 scon_a => basis_set_a%scon
913 zeta => basis_set_a%zet
914 ! basis jkind
915 first_sgfb => basis_set_b%first_sgf
916 lb_max => basis_set_b%lmax
917 lb_min => basis_set_b%lmin
918 npgfb => basis_set_b%npgf
919 nsetb = basis_set_b%nset
920 nsgfb => basis_set_b%nsgf_set
921 rpgfb => basis_set_b%pgf_radius
922 set_radius_b => basis_set_b%set_radius
923 scon_b => basis_set_b%scon
924 zetb => basis_set_b%zet
925
926 ldsab = get_memory_usage(qs_kind_set, "ORB", "ORB")
927 ALLOCATE (oint(ldsab, ldsab, 4), owork(ldsab, ldsab))
928 ALLOCATE (sint(natorb_a, natorb_b, 4))
929 sint = 0.0_dp
930
931 DO iset = 1, nseta
932 ncoa = npgfa(iset)*ncoset(la_max(iset))
933 n1 = npgfa(iset)*(ncoset(la_max(iset)) - ncoset(la_min(iset) - 1))
934 sgfa = first_sgfa(1, iset)
935 DO jset = 1, nsetb
936 IF (set_radius_a(iset) + set_radius_b(jset) < dr) cycle
937 ncob = npgfb(jset)*ncoset(lb_max(jset))
938 n2 = npgfb(jset)*(ncoset(lb_max(jset)) - ncoset(lb_min(jset) - 1))
939 sgfb = first_sgfb(1, jset)
940 CALL overlap_ab(la_max(iset), la_min(iset), npgfa(iset), rpgfa(:, iset), zeta(:, iset), &
941 lb_max(jset), lb_min(jset), npgfb(jset), rpgfb(:, jset), zetb(:, jset), &
942 rij, sab=oint(:, :, 1), dab=oint(:, :, 2:4))
943 ! Contraction
944 DO i = 1, 4
945 CALL contraction(oint(:, :, i), owork, ca=scon_a(:, sgfa:), na=n1, ma=nsgfa(iset), &
946 cb=scon_b(:, sgfb:), nb=n2, mb=nsgfb(jset), fscale=1.0_dp, trans=.false.)
947 CALL block_add("IN", owork, nsgfa(iset), nsgfb(jset), sint(:, :, i), &
948 sgfa, sgfb, trans=.false.)
949 END DO
950 END DO
951 END DO
952 ! update dS/dR matrix
953 clist%acint(1:natorb_a, 1:natorb_b, 1:3) = sint(1:natorb_a, 1:natorb_b, 2:4)
954
955 DEALLOCATE (oint, owork, sint)
956
957 END DO
958 CALL neighbor_list_iterator_release(nl_iterator)
959
960 DEALLOCATE (basis_set_list)
961
962 CALL timestop(handle)
963
964 END SUBROUTINE xtb_dsint_list
965
966END MODULE xtb_coulomb
967
Set of routines to: Contract integrals over primitive Gaussians Decontract (density) matrices Trace m...
Calculation of the overlap integrals over Cartesian Gaussian-type functions.
Definition ai_overlap.F:18
subroutine, public overlap_ab(la_max, la_min, npgfa, rpgfa, zeta, lb_max, lb_min, npgfb, rpgfb, zetb, rab, sab, dab, ddab)
Calculation of the two-center overlap integrals [a|b] over Cartesian Gaussian-type functions....
Definition ai_overlap.F:680
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.
Holds information on atomic properties.
subroutine, public atprop_array_init(atarray, natom)
...
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, h, h_inv, symmetry_id, tag)
Get informations about a simulation cell.
Definition cell_types.F:195
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)
...
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
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
Interface to the message passing library MPI.
Provides Cartesian and spherical orbital pointers and indices.
integer, dimension(:), allocatable, public ncoset
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)
...
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.
Some utility functions for the calculation of integrals.
subroutine, public basis_set_list_setup(basis_set_list, basis_type, qs_kind_set)
Set up an easy accessible list of the basis sets for all kinds.
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.
Calculation of Coulomb contributions in xTB.
Definition xtb_coulomb.F:12
subroutine, public gamma_rab_sr(gmat, rab, nla, kappaa, etaa, nlb, kappab, etab, kg, rcut)
Computes the short-range gamma parameter from Nataga-Mishimoto-Ohno-Klopman formula for xTB WARNING: ...
subroutine, public dgamma_rab_sr(dgmat, rab, nla, kappaa, etaa, nlb, kappab, etab, kg, rcut)
Computes the derivative of the short-range gamma parameter from Nataga-Mishimoto-Ohno-Klopman formula...
subroutine, public build_xtb_coulomb(qs_env, ks_matrix, rho, charges, mcharge, energy, calculate_forces, just_energy)
...
subroutine, public xtb_dsint_list(qs_env, sap_int)
...
Definition of the xTB parameter types.
Definition xtb_types.F:20
subroutine, public get_xtb_atom_param(xtb_parameter, symbol, aname, typ, defined, z, zeff, natorb, lmax, nao, lao, rcut, rcov, kx, eta, xgamma, alpha, zneff, nshell, nval, lval, kpoly, kappa, hen, zeta, xi, kappa0, alpg, occupation, electronegativity, chmax, en, kqat2, kcn, kq)
...
Definition xtb_types.F:199
Provides all information about an atomic kind.
type for the atomic properties
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
structure to store local (to a processor) ordered lists of integers.
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.