(git:374b731)
Loading...
Searching...
No Matches
qs_dftb_coulomb.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2024 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 DFTB
10!> \author JGH
11! **************************************************************************************************
18 USE cell_types, ONLY: cell_type,&
19 get_cell,&
20 pbc
23 USE dbcsr_api, ONLY: dbcsr_add,&
24 dbcsr_get_block_p,&
25 dbcsr_iterator_blocks_left,&
26 dbcsr_iterator_next_block,&
27 dbcsr_iterator_start,&
28 dbcsr_iterator_stop,&
29 dbcsr_iterator_type,&
30 dbcsr_p_type
37 USE kinds, ONLY: dp
38 USE kpoint_types, ONLY: get_kpoint_info,&
40 USE mathconstants, ONLY: oorootpi,&
41 pi
58 USE qs_kind_types, ONLY: get_qs_kind,&
66 USE qs_rho_types, ONLY: qs_rho_get,&
68 USE sap_kind_types, ONLY: clist_type,&
72 USE virial_types, ONLY: virial_type
73#include "./base/base_uses.f90"
74
75 IMPLICIT NONE
76
77 PRIVATE
78
79 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_dftb_coulomb'
80
81 ! screening for gamma function
82 REAL(dp), PARAMETER :: tol_gamma = 1.e-4_dp
83 ! small real number
84 REAL(dp), PARAMETER :: rtiny = 1.e-10_dp
85
87
88CONTAINS
89
90! **************************************************************************************************
91!> \brief ...
92!> \param qs_env ...
93!> \param ks_matrix ...
94!> \param rho ...
95!> \param mcharge ...
96!> \param energy ...
97!> \param calculate_forces ...
98!> \param just_energy ...
99! **************************************************************************************************
100 SUBROUTINE build_dftb_coulomb(qs_env, ks_matrix, rho, mcharge, energy, &
101 calculate_forces, just_energy)
102
103 TYPE(qs_environment_type), POINTER :: qs_env
104 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: ks_matrix
105 TYPE(qs_rho_type), POINTER :: rho
106 REAL(dp), DIMENSION(:) :: mcharge
107 TYPE(qs_energy_type), POINTER :: energy
108 LOGICAL, INTENT(in) :: calculate_forces, just_energy
109
110 CHARACTER(len=*), PARAMETER :: routinen = 'build_dftb_coulomb'
111
112 INTEGER :: atom_i, atom_j, blk, 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, &
125 sblock
126 REAL(kind=dp), DIMENSION(:, :, :), POINTER :: dsint
127 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
128 TYPE(atprop_type), POINTER :: atprop
129 TYPE(cell_type), POINTER :: cell
130 TYPE(dbcsr_iterator_type) :: iter
131 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_p, matrix_s
132 TYPE(dft_control_type), POINTER :: dft_control
133 TYPE(distribution_1d_type), POINTER :: local_particles
134 TYPE(ewald_environment_type), POINTER :: ewald_env
135 TYPE(ewald_pw_type), POINTER :: ewald_pw
136 TYPE(kpoint_type), POINTER :: kpoints
137 TYPE(mp_para_env_type), POINTER :: para_env
139 DIMENSION(:), POINTER :: nl_iterator
140 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
141 POINTER :: n_list
142 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
143 TYPE(qs_dftb_atom_type), POINTER :: dftb_kind, dftb_kind_a, dftb_kind_b
144 TYPE(qs_dftb_pairpot_type), DIMENSION(:, :), &
145 POINTER :: dftb_potential
146 TYPE(qs_force_type), DIMENSION(:), POINTER :: force
147 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
148 TYPE(sap_int_type), DIMENSION(:), POINTER :: sap_int
149 TYPE(virial_type), POINTER :: virial
150
151 CALL timeset(routinen, handle)
152
153 NULLIFY (matrix_p, matrix_s, virial, atprop, dft_control)
154
155 use_virial = .false.
156
157 IF (calculate_forces) THEN
158 nmat = 4
159 ELSE
160 nmat = 1
161 END IF
162
163 CALL get_qs_env(qs_env, nkind=nkind, natom=natom)
164 ALLOCATE (gmcharge(natom, nmat))
165 gmcharge = 0._dp
166
167 CALL get_qs_env(qs_env, &
168 particle_set=particle_set, &
169 cell=cell, &
170 virial=virial, &
171 atprop=atprop, &
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)
176
177 IF (calculate_forces) THEN
178 use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
179 END IF
180
181 NULLIFY (dftb_potential)
182 CALL get_qs_env(qs_env=qs_env, dftb_potential=dftb_potential)
183 NULLIFY (n_list)
184 CALL get_qs_env(qs_env=qs_env, sab_orb=n_list)
185 CALL neighbor_list_iterator_create(nl_iterator, n_list)
186 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
187 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
188 iatom=iatom, jatom=jatom, r=rij, cell=cellind)
189
190 CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind_a)
191 CALL get_dftb_atom_param(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)
195 CALL get_dftb_atom_param(dftb_kind_b, &
196 defined=defined, eta=eta_b, natorb=natorb_b)
197 IF (.NOT. defined .OR. natorb_b < 1) cycle
198
199 ! gamma matrix
200 hb_sr_damp = dft_control%qs_control%dftb_control%hb_sr_damp
201 IF (hb_sr_damp) THEN
202 ! short range correction enabled only when iatom XOR jatom are hydrogens
203 hb_sr_damp = is_hydrogen(particle_set(iatom)%atomic_kind) .NEQV. &
204 is_hydrogen(particle_set(jatom)%atomic_kind)
205 END IF
206 IF (hb_sr_damp) THEN
207 hb_para = dft_control%qs_control%dftb_control%hb_sr_para
208 ELSE
209 hb_para = 0.0_dp
210 END IF
211 ga = eta_a(0)
212 gb = eta_b(0)
213 dr = sqrt(sum(rij(:)**2))
214 gmat = gamma_rab_sr(dr, ga, gb, hb_para)
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)
218 END IF
219 IF (calculate_forces .AND. (iatom /= jatom .OR. dr > 0.001_dp)) THEN
220 ddr = 0.1_dp*dftb_potential(ikind, jkind)%dgrd
221 drp = dr + ddr
222 drm = dr - ddr
223 dgam = 0.5_dp*(gamma_rab_sr(drp, ga, gb, hb_para) - gamma_rab_sr(drm, ga, gb, hb_para))/ddr
224 DO i = 1, 3
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
228 END IF
229 IF (use_virial) THEN
230 fij(i) = -mcharge(iatom)*mcharge(jatom)*dgam*rij(i)/dr
231 END IF
232 END DO
233 IF (use_virial) THEN
234 fi = 1.0_dp
235 IF (iatom == jatom) fi = 0.5_dp
236 CALL virial_pair_force(virial%pv_virial, fi, fij, rij)
237 IF (atprop%stress) THEN
238 CALL virial_pair_force(atprop%atstress(:, :, iatom), fi*0.5_dp, fij, rij)
239 CALL virial_pair_force(atprop%atstress(:, :, jatom), fi*0.5_dp, fij, rij)
240 END IF
241 END IF
242 END IF
243 END DO
244 CALL neighbor_list_iterator_release(nl_iterator)
245
246 IF (atprop%energy) THEN
247 CALL get_qs_env(qs_env=qs_env, particle_set=particle_set)
248 natom = SIZE(particle_set)
249 CALL atprop_array_init(atprop%atecoul, natom)
250 END IF
251
252 ! 1/R contribution
253 do_ewald = dft_control%qs_control%dftb_control%do_ewald
254 IF (do_ewald) THEN
255 ! Ewald sum
256 NULLIFY (ewald_env, ewald_pw)
257 CALL get_qs_env(qs_env=qs_env, &
258 ewald_env=ewald_env, ewald_pw=ewald_pw)
259 CALL get_cell(cell=cell, periodic=periodic, deth=deth)
260 CALL ewald_env_get(ewald_env, alpha=alpha, ewald_type=ewald_type)
261 CALL get_qs_env(qs_env=qs_env, sab_tbe=n_list)
262 CALL tb_ewald_overlap(gmcharge, mcharge, alpha, n_list, &
263 virial, use_virial, atprop=atprop)
264 SELECT CASE (ewald_type)
265 CASE DEFAULT
266 cpabort("Invalid Ewald type")
267 CASE (do_ewald_none)
268 cpabort("Not allowed with DFTB")
269 CASE (do_ewald_ewald)
270 cpabort("Standard Ewald not implemented in DFTB")
271 CASE (do_ewald_pme)
272 cpabort("PME not implemented in DFTB")
273 CASE (do_ewald_spme)
274 CALL tb_spme_evaluate(ewald_env, ewald_pw, particle_set, cell, &
275 gmcharge, mcharge, calculate_forces, virial, &
276 use_virial, atprop=atprop)
277 END SELECT
278 ELSE
279 ! direct sum
280 CALL get_qs_env(qs_env=qs_env, &
281 local_particles=local_particles)
282 DO ikind = 1, SIZE(local_particles%n_el)
283 DO ia = 1, local_particles%n_el(ikind)
284 iatom = local_particles%list(ikind)%array(ia)
285 DO jatom = 1, iatom - 1
286 rij = particle_set(iatom)%r - particle_set(jatom)%r
287 rij = pbc(rij, cell)
288 dr = sqrt(sum(rij(:)**2))
289 gmcharge(iatom, 1) = gmcharge(iatom, 1) + mcharge(jatom)/dr
290 gmcharge(jatom, 1) = gmcharge(jatom, 1) + mcharge(iatom)/dr
291 DO i = 2, nmat
292 gmcharge(iatom, i) = gmcharge(iatom, i) + rij(i - 1)*mcharge(jatom)/dr**3
293 gmcharge(jatom, i) = gmcharge(jatom, i) - rij(i - 1)*mcharge(iatom)/dr**3
294 END DO
295 END DO
296 END DO
297 END DO
298 cpassert(.NOT. use_virial)
299 END IF
300
301 CALL para_env%sum(gmcharge(:, 1))
302
303 IF (do_ewald) THEN
304 ! add self charge interaction and background charge contribution
305 gmcharge(:, 1) = gmcharge(:, 1) - 2._dp*alpha*oorootpi*mcharge(:)
306 IF (any(periodic(:) == 1)) THEN
307 gmcharge(:, 1) = gmcharge(:, 1) - pi/alpha**2/deth
308 END IF
309 END IF
310
311 energy%hartree = energy%hartree + 0.5_dp*sum(mcharge(:)*gmcharge(:, 1))
312 IF (atprop%energy) THEN
313 CALL get_qs_env(qs_env=qs_env, &
314 local_particles=local_particles)
315 DO ikind = 1, SIZE(local_particles%n_el)
316 CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind)
317 CALL get_dftb_atom_param(dftb_kind, zeff=zeff)
318 DO ia = 1, local_particles%n_el(ikind)
319 iatom = local_particles%list(ikind)%array(ia)
320 atprop%atecoul(iatom) = atprop%atecoul(iatom) + &
321 0.5_dp*zeff*gmcharge(iatom, 1)
322 END DO
323 END DO
324 END IF
325
326 IF (calculate_forces) THEN
327 CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, &
328 kind_of=kind_of, &
329 atom_of_kind=atom_of_kind)
330
331 gmcharge(:, 2) = gmcharge(:, 2)*mcharge(:)
332 gmcharge(:, 3) = gmcharge(:, 3)*mcharge(:)
333 gmcharge(:, 4) = gmcharge(:, 4)*mcharge(:)
334 DO iatom = 1, natom
335 ikind = kind_of(iatom)
336 atom_i = atom_of_kind(iatom)
337 force(ikind)%rho_elec(1, atom_i) = force(ikind)%rho_elec(1, atom_i) - gmcharge(iatom, 2)
338 force(ikind)%rho_elec(2, atom_i) = force(ikind)%rho_elec(2, atom_i) - gmcharge(iatom, 3)
339 force(ikind)%rho_elec(3, atom_i) = force(ikind)%rho_elec(3, atom_i) - gmcharge(iatom, 4)
340 END DO
341 END IF
342
343 do_gamma_stress = .false.
344 IF (.NOT. just_energy .AND. use_virial) THEN
345 IF (dft_control%nimages == 1) do_gamma_stress = .true.
346 END IF
347
348 IF (.NOT. just_energy) THEN
349 CALL get_qs_env(qs_env=qs_env, matrix_s_kp=matrix_s)
350 CALL qs_rho_get(rho, rho_ao_kp=matrix_p)
351
352 nimg = dft_control%nimages
353 NULLIFY (cell_to_index)
354 IF (nimg > 1) THEN
355 NULLIFY (kpoints)
356 CALL get_qs_env(qs_env=qs_env, kpoints=kpoints)
357 CALL get_kpoint_info(kpoint=kpoints, cell_to_index=cell_to_index)
358 END IF
359
360 IF (calculate_forces .AND. SIZE(matrix_p, 1) == 2) THEN
361 DO img = 1, nimg
362 CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
363 alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
364 END DO
365 END IF
366
367 NULLIFY (sap_int)
368 IF (do_gamma_stress) THEN
369 ! derivative overlap integral (non collapsed)
370 CALL dftb_dsint_list(qs_env, sap_int)
371 END IF
372
373 IF (nimg == 1) THEN
374 ! no k-points; all matrices have been transformed to periodic bsf
375 CALL dbcsr_iterator_start(iter, matrix_s(1, 1)%matrix)
376 DO WHILE (dbcsr_iterator_blocks_left(iter))
377 CALL dbcsr_iterator_next_block(iter, irow, icol, sblock, blk)
378 gmij = 0.5_dp*(gmcharge(irow, 1) + gmcharge(icol, 1))
379 DO is = 1, SIZE(ks_matrix, 1)
380 NULLIFY (ksblock)
381 CALL dbcsr_get_block_p(matrix=ks_matrix(is, 1)%matrix, &
382 row=irow, col=icol, block=ksblock, found=found)
383 cpassert(found)
384 ksblock = ksblock - gmij*sblock
385 END DO
386 IF (calculate_forces) THEN
387 ikind = kind_of(irow)
388 atom_i = atom_of_kind(irow)
389 jkind = kind_of(icol)
390 atom_j = atom_of_kind(icol)
391 NULLIFY (pblock)
392 CALL dbcsr_get_block_p(matrix=matrix_p(1, 1)%matrix, &
393 row=irow, col=icol, block=pblock, found=found)
394 cpassert(found)
395 DO i = 1, 3
396 NULLIFY (dsblock)
397 CALL dbcsr_get_block_p(matrix=matrix_s(1 + i, 1)%matrix, &
398 row=irow, col=icol, block=dsblock, found=found)
399 cpassert(found)
400 fi = -2.0_dp*gmij*sum(pblock*dsblock)
401 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
402 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
403 fij(i) = fi
404 END DO
405 END IF
406 END DO
407 CALL dbcsr_iterator_stop(iter)
408 !
409 IF (do_gamma_stress) THEN
410 DO ikind = 1, nkind
411 DO jkind = 1, nkind
412 iac = ikind + nkind*(jkind - 1)
413 IF (.NOT. ASSOCIATED(sap_int(iac)%alist)) cycle
414 DO ia = 1, sap_int(iac)%nalist
415 IF (.NOT. ASSOCIATED(sap_int(iac)%alist(ia)%clist)) cycle
416 iatom = sap_int(iac)%alist(ia)%aatom
417 DO ic = 1, sap_int(iac)%alist(ia)%nclist
418 jatom = sap_int(iac)%alist(ia)%clist(ic)%catom
419 rij = sap_int(iac)%alist(ia)%clist(ic)%rac
420 dr = sqrt(sum(rij(:)**2))
421 IF (dr > 1.e-6_dp) THEN
422 dsint => sap_int(iac)%alist(ia)%clist(ic)%acint
423 gmij = 0.5_dp*(gmcharge(iatom, 1) + gmcharge(jatom, 1))
424 icol = max(iatom, jatom)
425 irow = min(iatom, jatom)
426 NULLIFY (pblock)
427 CALL dbcsr_get_block_p(matrix=matrix_p(1, 1)%matrix, &
428 row=irow, col=icol, block=pblock, found=found)
429 cpassert(found)
430 DO i = 1, 3
431 IF (irow == iatom) THEN
432 fij(i) = -2.0_dp*gmij*sum(pblock*dsint(:, :, i))
433 ELSE
434 fij(i) = -2.0_dp*gmij*sum(transpose(pblock)*dsint(:, :, i))
435 END IF
436 END DO
437 fi = 1.0_dp
438 IF (iatom == jatom) fi = 0.5_dp
439 CALL virial_pair_force(virial%pv_virial, fi, fij, rij)
440 IF (atprop%stress) THEN
441 CALL virial_pair_force(atprop%atstress(:, :, iatom), 0.5_dp*fi, fij, rij)
442 CALL virial_pair_force(atprop%atstress(:, :, jatom), 0.5_dp*fi, fij, rij)
443 END IF
444 END IF
445 END DO
446 END DO
447 END DO
448 END DO
449 END IF
450 ELSE
451 NULLIFY (n_list)
452 CALL get_qs_env(qs_env=qs_env, sab_orb=n_list)
453 CALL neighbor_list_iterator_create(nl_iterator, n_list)
454 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
455 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
456 iatom=iatom, jatom=jatom, r=rij, cell=cellind)
457
458 icol = max(iatom, jatom)
459 irow = min(iatom, jatom)
460 ic = cell_to_index(cellind(1), cellind(2), cellind(3))
461 cpassert(ic > 0)
462
463 gmij = 0.5_dp*(gmcharge(iatom, 1) + gmcharge(jatom, 1))
464
465 NULLIFY (sblock)
466 CALL dbcsr_get_block_p(matrix=matrix_s(1, ic)%matrix, &
467 row=irow, col=icol, block=sblock, found=found)
468 cpassert(found)
469 DO is = 1, SIZE(ks_matrix, 1)
470 NULLIFY (ksblock)
471 CALL dbcsr_get_block_p(matrix=ks_matrix(is, ic)%matrix, &
472 row=irow, col=icol, block=ksblock, found=found)
473 cpassert(found)
474 ksblock = ksblock - gmij*sblock
475 END DO
476
477 IF (calculate_forces) THEN
478 ikind = kind_of(iatom)
479 atom_i = atom_of_kind(iatom)
480 jkind = kind_of(jatom)
481 atom_j = atom_of_kind(jatom)
482 IF (irow == jatom) gmij = -gmij
483 NULLIFY (pblock)
484 CALL dbcsr_get_block_p(matrix=matrix_p(1, ic)%matrix, &
485 row=irow, col=icol, block=pblock, found=found)
486 cpassert(found)
487 DO i = 1, 3
488 NULLIFY (dsblock)
489 CALL dbcsr_get_block_p(matrix=matrix_s(1 + i, ic)%matrix, &
490 row=irow, col=icol, block=dsblock, found=found)
491 cpassert(found)
492 fi = -2.0_dp*gmij*sum(pblock*dsblock)
493 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
494 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
495 fij(i) = fi
496 END DO
497 IF (use_virial) THEN
498 fi = 1.0_dp
499 IF (iatom == jatom) fi = 0.5_dp
500 CALL virial_pair_force(virial%pv_virial, fi, fij, rij)
501 IF (atprop%stress) THEN
502 CALL virial_pair_force(atprop%atstress(:, :, iatom), fi*0.5_dp, fij, rij)
503 CALL virial_pair_force(atprop%atstress(:, :, jatom), fi*0.5_dp, fij, rij)
504 END IF
505 END IF
506 END IF
507 END DO
508 CALL neighbor_list_iterator_release(nl_iterator)
509 END IF
510
511 IF (calculate_forces .AND. SIZE(matrix_p, 1) == 2) THEN
512 DO img = 1, nimg
513 CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
514 alpha_scalar=1.0_dp, beta_scalar=-1.0_dp)
515 END DO
516 END IF
517 END IF
518
519 IF (dft_control%qs_control%dftb_control%dftb3_diagonal) THEN
520 CALL get_qs_env(qs_env, nkind=nkind)
521 ALLOCATE (zeffk(nkind), xgamma(nkind))
522 DO ikind = 1, nkind
523 CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind)
524 CALL get_dftb_atom_param(dftb_kind, dudq=xgamma(ikind), zeff=zeffk(ikind))
525 END DO
526 ! Diagonal 3rd order correction (DFTB3)
527 CALL build_dftb3_diagonal(qs_env, ks_matrix, rho, mcharge, energy, xgamma, zeffk, &
528 sap_int, calculate_forces, just_energy)
529 DEALLOCATE (zeffk, xgamma)
530 END IF
531
532 ! QMMM
533 IF (qs_env%qmmm .AND. qs_env%qmmm_periodic) THEN
534 CALL build_tb_coulomb_qmqm(qs_env, ks_matrix, rho, mcharge, energy, &
535 calculate_forces, just_energy)
536 END IF
537
538 IF (do_gamma_stress) THEN
539 CALL release_sap_int(sap_int)
540 END IF
541
542 DEALLOCATE (gmcharge)
543
544 CALL timestop(handle)
545
546 END SUBROUTINE build_dftb_coulomb
547
548! **************************************************************************************************
549!> \brief Computes the short-range gamma parameter from exact Coulomb
550!> interaction of normalized exp(-a*r) charge distribution - 1/r
551!> \param r ...
552!> \param ga ...
553!> \param gb ...
554!> \param hb_para ...
555!> \return ...
556!> \par Literature
557!> Elstner et al, PRB 58 (1998) 7260
558!> \par History
559!> 10.2008 Axel Kohlmeyer - adding sr_damp
560!> 08.2014 JGH - adding flexible exponent for damping
561!> \version 1.1
562! **************************************************************************************************
563 FUNCTION gamma_rab_sr(r, ga, gb, hb_para) RESULT(gamma)
564 REAL(dp), INTENT(in) :: r, ga, gb, hb_para
565 REAL(dp) :: gamma
566
567 REAL(dp) :: a, b, fac, g_sum
568
569 gamma = 0.0_dp
570 a = 3.2_dp*ga ! 3.2 = 16/5 in Eq. 18 and ff.
571 b = 3.2_dp*gb
572 g_sum = a + b
573 IF (g_sum < tol_gamma) RETURN ! hardness screening
574 IF (r < rtiny) THEN ! This is for short distances but non-onsite terms
575 ! This gives also correct diagonal elements (a=b, r=0)
576 gamma = 0.5_dp*(a*b/g_sum + (a*b)**2/g_sum**3)
577 RETURN
578 END IF
579 !
580 ! distinguish two cases: Gamma's are very close, e.g. for the same atom type,
581 ! and Gamma's are different
582 !
583 IF (abs(a - b) < rtiny) THEN
584 fac = 1.6_dp*r*a*b/g_sum*(1.0_dp + a*b/g_sum**2)
585 gamma = -(48.0_dp + 33._dp*fac + (9.0_dp + fac)*fac**2)*exp(-fac)/(48._dp*r)
586 ELSE
587 gamma = -exp(-a*r)*(0.5_dp*a*b**4/(a**2 - b**2)**2 - &
588 (b**6 - 3._dp*a**2*b**4)/(r*(a**2 - b**2)**3)) - & ! a-> b
589 exp(-b*r)*(0.5_dp*b*a**4/(b**2 - a**2)**2 - &
590 (a**6 - 3._dp*b**2*a**4)/(r*(b**2 - a**2)**3)) ! b-> a
591 END IF
592 !
593 ! damping function for better short range hydrogen bonds.
594 ! functional form from Hu H. et al., J. Phys. Chem. A 2007, 111, 5685-5691
595 ! according to Elstner M, Theor. Chem. Acc. 2006, 116, 316-325,
596 ! this should only be applied to a-b pairs involving hydrogen.
597 IF (hb_para > 0.0_dp) THEN
598 gamma = gamma*exp(-(0.5_dp*(ga + gb))**hb_para*r*r)
599 END IF
600 END FUNCTION gamma_rab_sr
601
602! **************************************************************************************************
603!> \brief ...
604!> \param qs_env ...
605!> \param sap_int ...
606! **************************************************************************************************
607 SUBROUTINE dftb_dsint_list(qs_env, sap_int)
608
609 TYPE(qs_environment_type), POINTER :: qs_env
610 TYPE(sap_int_type), DIMENSION(:), POINTER :: sap_int
611
612 CHARACTER(LEN=*), PARAMETER :: routinen = 'dftb_dsint_list'
613
614 INTEGER :: handle, i, iac, iatom, ikind, ilist, jatom, jkind, jneighbor, llm, lmaxi, lmaxj, &
615 n1, n2, natorb_a, natorb_b, ngrd, ngrdcut, nkind, nlist, nneighbor
616 INTEGER, DIMENSION(3) :: cell
617 LOGICAL :: defined
618 REAL(kind=dp) :: ddr, dgrd, dr
619 REAL(kind=dp), DIMENSION(3) :: drij, rij
620 REAL(kind=dp), DIMENSION(:, :), POINTER :: dsblock, dsblockm, smatij, smatji
621 TYPE(clist_type), POINTER :: clist
622 TYPE(dft_control_type), POINTER :: dft_control
623 TYPE(dftb_control_type), POINTER :: dftb_control
625 DIMENSION(:), POINTER :: nl_iterator
626 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
627 POINTER :: sab_orb
628 TYPE(qs_dftb_atom_type), POINTER :: dftb_kind_a, dftb_kind_b
629 TYPE(qs_dftb_pairpot_type), DIMENSION(:, :), &
630 POINTER :: dftb_potential
631 TYPE(qs_dftb_pairpot_type), POINTER :: dftb_param_ij, dftb_param_ji
632 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
633
634 CALL timeset(routinen, handle)
635
636 CALL get_qs_env(qs_env=qs_env, nkind=nkind)
637 cpassert(.NOT. ASSOCIATED(sap_int))
638 ALLOCATE (sap_int(nkind*nkind))
639 DO i = 1, nkind*nkind
640 NULLIFY (sap_int(i)%alist, sap_int(i)%asort, sap_int(i)%aindex)
641 sap_int(i)%nalist = 0
642 END DO
643
644 CALL get_qs_env(qs_env=qs_env, &
645 qs_kind_set=qs_kind_set, &
646 dft_control=dft_control, &
647 sab_orb=sab_orb)
648
649 dftb_control => dft_control%qs_control%dftb_control
650
651 NULLIFY (dftb_potential)
652 CALL get_qs_env(qs_env=qs_env, &
653 dftb_potential=dftb_potential)
654
655 ! loop over all atom pairs with a non-zero overlap (sab_orb)
656 CALL neighbor_list_iterator_create(nl_iterator, sab_orb)
657 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
658 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, iatom=iatom, &
659 jatom=jatom, nlist=nlist, ilist=ilist, nnode=nneighbor, &
660 inode=jneighbor, cell=cell, r=rij)
661 iac = ikind + nkind*(jkind - 1)
662 !
663 CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind_a)
664 CALL get_dftb_atom_param(dftb_kind_a, &
665 defined=defined, lmax=lmaxi, natorb=natorb_a)
666 IF (.NOT. defined .OR. natorb_a < 1) cycle
667 CALL get_qs_kind(qs_kind_set(jkind), dftb_parameter=dftb_kind_b)
668 CALL get_dftb_atom_param(dftb_kind_b, &
669 defined=defined, lmax=lmaxj, natorb=natorb_b)
670 IF (.NOT. defined .OR. natorb_b < 1) cycle
671
672 dr = sqrt(sum(rij(:)**2))
673
674 ! integral list
675 IF (.NOT. ASSOCIATED(sap_int(iac)%alist)) THEN
676 sap_int(iac)%a_kind = ikind
677 sap_int(iac)%p_kind = jkind
678 sap_int(iac)%nalist = nlist
679 ALLOCATE (sap_int(iac)%alist(nlist))
680 DO i = 1, nlist
681 NULLIFY (sap_int(iac)%alist(i)%clist)
682 sap_int(iac)%alist(i)%aatom = 0
683 sap_int(iac)%alist(i)%nclist = 0
684 END DO
685 END IF
686 IF (.NOT. ASSOCIATED(sap_int(iac)%alist(ilist)%clist)) THEN
687 sap_int(iac)%alist(ilist)%aatom = iatom
688 sap_int(iac)%alist(ilist)%nclist = nneighbor
689 ALLOCATE (sap_int(iac)%alist(ilist)%clist(nneighbor))
690 DO i = 1, nneighbor
691 sap_int(iac)%alist(ilist)%clist(i)%catom = 0
692 END DO
693 END IF
694 clist => sap_int(iac)%alist(ilist)%clist(jneighbor)
695 clist%catom = jatom
696 clist%cell = cell
697 clist%rac = rij
698 ALLOCATE (clist%acint(natorb_a, natorb_b, 3))
699 NULLIFY (clist%achint)
700 clist%acint = 0._dp
701 clist%nsgf_cnt = 0
702 NULLIFY (clist%sgf_list)
703
704 ! retrieve information on S matrix
705 dftb_param_ij => dftb_potential(ikind, jkind)
706 dftb_param_ji => dftb_potential(jkind, ikind)
707 ! assume table size and type is symmetric
708 ngrd = dftb_param_ij%ngrd
709 ngrdcut = dftb_param_ij%ngrdcut
710 dgrd = dftb_param_ij%dgrd
711 ddr = dgrd*0.1_dp
712 cpassert(dftb_param_ij%llm == dftb_param_ji%llm)
713 llm = dftb_param_ij%llm
714 smatij => dftb_param_ij%smat
715 smatji => dftb_param_ji%smat
716
717 IF (nint(dr/dgrd) <= ngrdcut) THEN
718 IF (iatom == jatom .AND. dr < 0.001_dp) THEN
719 ! diagonal block
720 ELSE
721 ! off-diagonal block
722 n1 = natorb_a
723 n2 = natorb_b
724 ALLOCATE (dsblock(n1, n2), dsblockm(n1, n2))
725 DO i = 1, 3
726 dsblock = 0._dp
727 dsblockm = 0.0_dp
728 drij = rij
729
730 drij(i) = rij(i) - ddr
731 CALL compute_block_sk(dsblockm, smatij, smatji, drij, ngrd, ngrdcut, dgrd, &
732 llm, lmaxi, lmaxj, iatom, iatom)
733
734 drij(i) = rij(i) + ddr
735 CALL compute_block_sk(dsblock, smatij, smatji, drij, ngrd, ngrdcut, dgrd, &
736 llm, lmaxi, lmaxj, iatom, iatom)
737
738 dsblock = dsblock - dsblockm
739 dsblock = dsblock/(2.0_dp*ddr)
740
741 clist%acint(1:n1, 1:n2, i) = -dsblock(1:n1, 1:n2)
742 END DO
743 DEALLOCATE (dsblock, dsblockm)
744 END IF
745 END IF
746
747 END DO
748 CALL neighbor_list_iterator_release(nl_iterator)
749
750 CALL timestop(handle)
751
752 END SUBROUTINE dftb_dsint_list
753
754END MODULE qs_dftb_coulomb
755
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.
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...
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, atprop)
...
subroutine, public tb_spme_evaluate(ewald_env, ewald_pw, particle_set, box, gmcharge, mcharge, calculate_forces, virial, use_virial, atprop)
...
Calculation of the incomplete Gamma function F_n(t) for multi-center integrals over Cartesian Gaussia...
Definition gamma.F:15
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
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_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, 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, 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, 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_r3d_rs_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_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.
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.