Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2026 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 : ! **************************************************************************************************
12 : MODULE qs_dftb_coulomb
13 : USE atomic_kind_types, ONLY: atomic_kind_type,&
14 : get_atomic_kind_set,&
15 : is_hydrogen
16 : USE atprop_types, ONLY: atprop_array_init,&
17 : atprop_type
18 : USE cell_types, ONLY: cell_type,&
19 : get_cell,&
20 : pbc
21 : USE cp_control_types, ONLY: dft_control_type,&
22 : dftb_control_type
23 : USE cp_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
31 : USE distribution_1d_types, ONLY: distribution_1d_type
32 : USE ewald_environment_types, ONLY: ewald_env_get,&
33 : ewald_environment_type
34 : USE ewald_methods_tb, ONLY: tb_ewald_overlap,&
35 : tb_spme_evaluate
36 : USE ewald_pw_types, ONLY: ewald_pw_type
37 : USE kinds, ONLY: dp
38 : USE kpoint_types, ONLY: get_kpoint_info,&
39 : kpoint_type
40 : USE mathconstants, ONLY: oorootpi,&
41 : pi
42 : USE message_passing, ONLY: mp_para_env_type
43 : USE particle_types, ONLY: particle_type
44 : USE pw_poisson_types, ONLY: do_ewald_ewald,&
45 : do_ewald_none,&
46 : do_ewald_pme,&
47 : do_ewald_spme
48 : USE qmmm_tb_coulomb, ONLY: build_tb_coulomb_qmqm
49 : USE qs_dftb3_methods, ONLY: build_dftb3_diagonal
50 : USE qs_dftb_types, ONLY: qs_dftb_atom_type,&
51 : qs_dftb_pairpot_type
52 : USE qs_dftb_utils, ONLY: compute_block_sk,&
53 : get_dftb_atom_param
54 : USE qs_energy_types, ONLY: qs_energy_type
55 : USE qs_environment_types, ONLY: get_qs_env,&
56 : qs_environment_type
57 : USE qs_force_types, ONLY: qs_force_type
58 : USE qs_kind_types, ONLY: get_qs_kind,&
59 : qs_kind_type
60 : USE qs_neighbor_list_types, ONLY: get_iterator_info,&
61 : neighbor_list_iterate,&
62 : neighbor_list_iterator_create,&
63 : neighbor_list_iterator_p_type,&
64 : neighbor_list_iterator_release,&
65 : neighbor_list_set_p_type
66 : USE qs_rho_types, ONLY: qs_rho_get,&
67 : qs_rho_type
68 : USE sap_kind_types, ONLY: clist_type,&
69 : release_sap_int,&
70 : sap_int_type
71 : USE virial_methods, ONLY: virial_pair_force
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 : REAL(KIND=dp), PARAMETER, PRIVATE :: dftb_fd_deriv_step = 1.0e-3_dp
86 :
87 : PUBLIC :: build_dftb_coulomb, gamma_rab_sr
88 :
89 : CONTAINS
90 :
91 : ! **************************************************************************************************
92 : !> \brief ...
93 : !> \param qs_env ...
94 : !> \param ks_matrix ...
95 : !> \param rho ...
96 : !> \param mcharge ...
97 : !> \param energy ...
98 : !> \param calculate_forces ...
99 : !> \param just_energy ...
100 : ! **************************************************************************************************
101 16708 : SUBROUTINE build_dftb_coulomb(qs_env, ks_matrix, rho, 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(:) :: mcharge
108 : TYPE(qs_energy_type), POINTER :: energy
109 : LOGICAL, INTENT(in) :: calculate_forces, just_energy
110 :
111 : CHARACTER(len=*), PARAMETER :: routineN = 'build_dftb_coulomb'
112 :
113 : INTEGER :: atom_i, atom_j, ewald_type, handle, i, ia, iac, iatom, ic, icol, ikind, img, &
114 : irow, is, jatom, jkind, llm, lmaxi, lmaxj, n1, n2, natom, natorb_a, natorb_b, ngrd, &
115 : ngrdcut, nimg, nkind, nmat
116 16708 : INTEGER, ALLOCATABLE, DIMENSION(:) :: atom_of_kind, kind_of
117 : INTEGER, DIMENSION(3) :: cellind, periodic
118 16708 : INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
119 : LOGICAL :: defined, defined_a, defined_b, do_ewald, &
120 : do_gamma_stress, found, hb_sr_damp, &
121 : use_virial
122 : REAL(KIND=dp) :: alpha, background, ddr, deth, dgam, &
123 : dgrd, dr, drm, drp, fi, ga, gb, gmat, &
124 : gmij, gmij_virial, hb_para, zeff
125 16708 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: xgamma, zeffk
126 : REAL(KIND=dp), DIMENSION(0:3) :: eta_a, eta_b
127 : REAL(KIND=dp), DIMENSION(3) :: drij, fij, rij
128 16708 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: dsblock, dsblock_vir, dsblockm_vir, &
129 16708 : gmcharge, ksblock, pblock, sblock, &
130 16708 : smatij, smatji
131 16708 : REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: dsint
132 16708 : TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
133 : TYPE(atprop_type), POINTER :: atprop
134 : TYPE(cell_type), POINTER :: cell
135 : TYPE(dbcsr_iterator_type) :: iter
136 16708 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_p, matrix_s
137 : TYPE(dft_control_type), POINTER :: dft_control
138 : TYPE(distribution_1d_type), POINTER :: local_particles
139 : TYPE(ewald_environment_type), POINTER :: ewald_env
140 : TYPE(ewald_pw_type), POINTER :: ewald_pw
141 : TYPE(kpoint_type), POINTER :: kpoints
142 : TYPE(mp_para_env_type), POINTER :: para_env
143 : TYPE(neighbor_list_iterator_p_type), &
144 16708 : DIMENSION(:), POINTER :: nl_iterator
145 : TYPE(neighbor_list_set_p_type), DIMENSION(:), &
146 16708 : POINTER :: n_list
147 16708 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
148 : TYPE(qs_dftb_atom_type), POINTER :: dftb_kind, dftb_kind_a, dftb_kind_b
149 : TYPE(qs_dftb_pairpot_type), DIMENSION(:, :), &
150 16708 : POINTER :: dftb_potential
151 : TYPE(qs_dftb_pairpot_type), POINTER :: dftb_param_ij, dftb_param_ji
152 16708 : TYPE(qs_force_type), DIMENSION(:), POINTER :: force
153 16708 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
154 16708 : TYPE(sap_int_type), DIMENSION(:), POINTER :: sap_int
155 : TYPE(virial_type), POINTER :: virial
156 :
157 16708 : CALL timeset(routineN, handle)
158 :
159 16708 : NULLIFY (matrix_p, matrix_s, virial, atprop, dft_control)
160 :
161 16708 : use_virial = .FALSE.
162 :
163 16708 : IF (calculate_forces) THEN
164 : nmat = 4
165 : ELSE
166 16076 : nmat = 1
167 : END IF
168 :
169 16708 : CALL get_qs_env(qs_env, nkind=nkind, natom=natom)
170 66832 : ALLOCATE (gmcharge(natom, nmat))
171 261114 : gmcharge = 0._dp
172 :
173 : CALL get_qs_env(qs_env, &
174 : particle_set=particle_set, &
175 : cell=cell, &
176 : virial=virial, &
177 : atprop=atprop, &
178 : dft_control=dft_control, &
179 : atomic_kind_set=atomic_kind_set, &
180 : qs_kind_set=qs_kind_set, &
181 16708 : force=force, para_env=para_env)
182 :
183 16708 : IF (calculate_forces) THEN
184 1148 : use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
185 : END IF
186 :
187 16708 : NULLIFY (dftb_potential)
188 16708 : CALL get_qs_env(qs_env=qs_env, dftb_potential=dftb_potential)
189 16708 : NULLIFY (n_list)
190 16708 : CALL get_qs_env(qs_env=qs_env, sab_orb=n_list)
191 16708 : CALL neighbor_list_iterator_create(nl_iterator, n_list)
192 2807522 : DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
193 : CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
194 2790814 : iatom=iatom, jatom=jatom, r=rij, cell=cellind)
195 :
196 2790814 : CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind_a)
197 : CALL get_dftb_atom_param(dftb_kind_a, &
198 2790814 : defined=defined, eta=eta_a, natorb=natorb_a)
199 2790814 : IF (.NOT. defined .OR. natorb_a < 1) CYCLE
200 2790814 : CALL get_qs_kind(qs_kind_set(jkind), dftb_parameter=dftb_kind_b)
201 : CALL get_dftb_atom_param(dftb_kind_b, &
202 2790814 : defined=defined, eta=eta_b, natorb=natorb_b)
203 2790814 : IF (.NOT. defined .OR. natorb_b < 1) CYCLE
204 :
205 : ! gamma matrix
206 2790814 : hb_sr_damp = dft_control%qs_control%dftb_control%hb_sr_damp
207 2790814 : IF (hb_sr_damp) THEN
208 : ! short range correction enabled only when iatom XOR jatom are hydrogens
209 : hb_sr_damp = is_hydrogen(particle_set(iatom)%atomic_kind) .NEQV. &
210 1089438 : is_hydrogen(particle_set(jatom)%atomic_kind)
211 : END IF
212 1089438 : IF (hb_sr_damp) THEN
213 435556 : hb_para = dft_control%qs_control%dftb_control%hb_sr_para
214 : ELSE
215 2355258 : hb_para = 0.0_dp
216 : END IF
217 2790814 : ga = eta_a(0)
218 2790814 : gb = eta_b(0)
219 11163256 : dr = SQRT(SUM(rij(:)**2))
220 2790814 : gmat = gamma_rab_sr(dr, ga, gb, hb_para)
221 2790814 : gmcharge(jatom, 1) = gmcharge(jatom, 1) + gmat*mcharge(iatom)
222 2790814 : IF (iatom /= jatom) THEN
223 2530487 : gmcharge(iatom, 1) = gmcharge(iatom, 1) + gmat*mcharge(jatom)
224 : END IF
225 2807522 : IF (calculate_forces .AND. (iatom /= jatom .OR. dr > 0.001_dp)) THEN
226 254964 : ddr = dftb_fd_deriv_step*dftb_potential(ikind, jkind)%dgrd
227 254964 : drp = dr + ddr
228 254964 : drm = dr - ddr
229 254964 : dgam = 0.5_dp*(gamma_rab_sr(drp, ga, gb, hb_para) - gamma_rab_sr(drm, ga, gb, hb_para))/ddr
230 1019856 : DO i = 1, 3
231 764892 : gmcharge(jatom, i + 1) = gmcharge(jatom, i + 1) - dgam*mcharge(iatom)*rij(i)/dr
232 764892 : IF (dr > 0.001_dp) THEN
233 764892 : gmcharge(iatom, i + 1) = gmcharge(iatom, i + 1) + dgam*mcharge(jatom)*rij(i)/dr
234 : END IF
235 1019856 : IF (use_virial) THEN
236 512889 : fij(i) = -mcharge(iatom)*mcharge(jatom)*dgam*rij(i)/dr
237 : END IF
238 : END DO
239 254964 : IF (use_virial) THEN
240 170963 : fi = 1.0_dp
241 170963 : IF (iatom == jatom) fi = 0.5_dp
242 170963 : CALL virial_pair_force(virial%pv_virial, fi, fij, rij)
243 : END IF
244 : END IF
245 : END DO
246 16708 : CALL neighbor_list_iterator_release(nl_iterator)
247 :
248 16708 : IF (atprop%energy) THEN
249 466 : CALL get_qs_env(qs_env=qs_env, particle_set=particle_set)
250 466 : natom = SIZE(particle_set)
251 466 : CALL atprop_array_init(atprop%atecoul, natom)
252 : END IF
253 :
254 : ! 1/R contribution
255 16708 : do_ewald = dft_control%qs_control%dftb_control%do_ewald
256 16708 : IF (do_ewald) THEN
257 : ! Ewald sum
258 9482 : NULLIFY (ewald_env, ewald_pw)
259 : CALL get_qs_env(qs_env=qs_env, &
260 9482 : ewald_env=ewald_env, ewald_pw=ewald_pw)
261 9482 : CALL get_cell(cell=cell, periodic=periodic, deth=deth)
262 9482 : CALL ewald_env_get(ewald_env, alpha=alpha, ewald_type=ewald_type)
263 9482 : CALL get_qs_env(qs_env=qs_env, sab_tbe=n_list)
264 9482 : CALL tb_ewald_overlap(gmcharge, mcharge, alpha, n_list, virial, use_virial)
265 0 : SELECT CASE (ewald_type)
266 : CASE DEFAULT
267 0 : CPABORT("Invalid Ewald type")
268 : CASE (do_ewald_none)
269 0 : CPABORT("Not allowed with DFTB")
270 : CASE (do_ewald_ewald)
271 0 : CPABORT("Standard Ewald not implemented in DFTB")
272 : CASE (do_ewald_pme)
273 0 : CPABORT("PME not implemented in DFTB")
274 : CASE (do_ewald_spme)
275 : CALL tb_spme_evaluate(ewald_env, ewald_pw, particle_set, cell, &
276 9482 : gmcharge, mcharge, calculate_forces, virial, use_virial)
277 : END SELECT
278 : ELSE
279 : ! direct sum
280 : CALL get_qs_env(qs_env=qs_env, &
281 7226 : local_particles=local_particles)
282 24240 : DO ikind = 1, SIZE(local_particles%n_el)
283 36355 : DO ia = 1, local_particles%n_el(ikind)
284 12115 : iatom = local_particles%list(ikind)%array(ia)
285 43801 : DO jatom = 1, iatom - 1
286 58688 : rij = particle_set(iatom)%r - particle_set(jatom)%r
287 58688 : rij = pbc(rij, cell)
288 58688 : dr = SQRT(SUM(rij(:)**2))
289 14672 : gmcharge(iatom, 1) = gmcharge(iatom, 1) + mcharge(jatom)/dr
290 14672 : gmcharge(jatom, 1) = gmcharge(jatom, 1) + mcharge(iatom)/dr
291 27726 : DO i = 2, nmat
292 939 : gmcharge(iatom, i) = gmcharge(iatom, i) + rij(i - 1)*mcharge(jatom)/dr**3
293 15611 : 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 7226 : CPASSERT(.NOT. use_virial)
299 : END IF
300 :
301 382536 : CALL para_env%sum(gmcharge(:, 1))
302 :
303 16708 : background = 0.0_dp
304 16708 : IF (do_ewald) THEN
305 : ! add self charge interaction and background charge contribution
306 168166 : gmcharge(:, 1) = gmcharge(:, 1) - 2._dp*alpha*oorootpi*mcharge(:)
307 10304 : IF (ANY(periodic(:) == 1)) THEN
308 9208 : background = pi/alpha**2/deth
309 : END IF
310 : END IF
311 :
312 199622 : energy%hartree = energy%hartree + 0.5_dp*SUM(mcharge(:)*(gmcharge(:, 1) - background))
313 16708 : IF (atprop%energy) THEN
314 : CALL get_qs_env(qs_env=qs_env, &
315 466 : local_particles=local_particles)
316 1526 : DO ikind = 1, SIZE(local_particles%n_el)
317 1060 : CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind)
318 1060 : CALL get_dftb_atom_param(dftb_kind, zeff=zeff)
319 18006 : DO ia = 1, local_particles%n_el(ikind)
320 16480 : iatom = local_particles%list(ikind)%array(ia)
321 : atprop%atecoul(iatom) = atprop%atecoul(iatom) + &
322 17540 : 0.5_dp*zeff*gmcharge(iatom, 1)
323 : END DO
324 : END DO
325 : END IF
326 :
327 16708 : IF (calculate_forces) THEN
328 : CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, &
329 : kind_of=kind_of, &
330 632 : atom_of_kind=atom_of_kind)
331 :
332 14928 : gmcharge(:, 2) = gmcharge(:, 2)*mcharge(:)
333 14928 : gmcharge(:, 3) = gmcharge(:, 3)*mcharge(:)
334 14928 : gmcharge(:, 4) = gmcharge(:, 4)*mcharge(:)
335 14928 : DO iatom = 1, natom
336 14296 : ikind = kind_of(iatom)
337 14296 : atom_i = atom_of_kind(iatom)
338 14296 : force(ikind)%rho_elec(1, atom_i) = force(ikind)%rho_elec(1, atom_i) - gmcharge(iatom, 2)
339 14296 : force(ikind)%rho_elec(2, atom_i) = force(ikind)%rho_elec(2, atom_i) - gmcharge(iatom, 3)
340 14928 : force(ikind)%rho_elec(3, atom_i) = force(ikind)%rho_elec(3, atom_i) - gmcharge(iatom, 4)
341 : END DO
342 : END IF
343 :
344 16708 : do_gamma_stress = .FALSE.
345 16708 : IF (.NOT. just_energy .AND. use_virial) THEN
346 116 : IF (dft_control%nimages == 1) do_gamma_stress = .TRUE.
347 : END IF
348 :
349 16708 : IF (.NOT. just_energy) THEN
350 16432 : CALL get_qs_env(qs_env=qs_env, matrix_s_kp=matrix_s)
351 16432 : CALL qs_rho_get(rho, rho_ao_kp=matrix_p)
352 :
353 16432 : nimg = dft_control%nimages
354 16432 : NULLIFY (cell_to_index)
355 16432 : IF (nimg > 1) THEN
356 4442 : NULLIFY (kpoints)
357 4442 : CALL get_qs_env(qs_env=qs_env, kpoints=kpoints)
358 4442 : CALL get_kpoint_info(kpoint=kpoints, cell_to_index=cell_to_index)
359 : END IF
360 :
361 16432 : IF (calculate_forces .AND. SIZE(matrix_p, 1) == 2) THEN
362 328 : DO img = 1, nimg
363 : CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
364 328 : alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
365 : END DO
366 : END IF
367 :
368 16432 : NULLIFY (sap_int)
369 16432 : IF (do_gamma_stress) THEN
370 : ! derivative overlap integral (non collapsed)
371 90 : CALL dftb_dsint_list(qs_env, sap_int)
372 : END IF
373 :
374 16432 : IF (nimg == 1) THEN
375 : ! no k-points; all matrices have been transformed to periodic bsf
376 11990 : CALL dbcsr_iterator_start(iter, matrix_s(1, 1)%matrix)
377 1750722 : DO WHILE (dbcsr_iterator_blocks_left(iter))
378 1738732 : CALL dbcsr_iterator_next_block(iter, irow, icol, sblock)
379 1738732 : gmij = 0.5_dp*(gmcharge(irow, 1) + gmcharge(icol, 1))
380 3478734 : DO is = 1, SIZE(ks_matrix, 1)
381 1740002 : NULLIFY (ksblock)
382 : CALL dbcsr_get_block_p(matrix=ks_matrix(is, 1)%matrix, &
383 1740002 : row=irow, col=icol, block=ksblock, found=found)
384 1740002 : CPASSERT(found)
385 26229190 : ksblock = ksblock - gmij*sblock
386 : END DO
387 1750722 : IF (calculate_forces) THEN
388 230509 : ikind = kind_of(irow)
389 230509 : atom_i = atom_of_kind(irow)
390 230509 : jkind = kind_of(icol)
391 230509 : atom_j = atom_of_kind(icol)
392 230509 : NULLIFY (pblock)
393 : CALL dbcsr_get_block_p(matrix=matrix_p(1, 1)%matrix, &
394 230509 : row=irow, col=icol, block=pblock, found=found)
395 230509 : CPASSERT(found)
396 922036 : DO i = 1, 3
397 691527 : NULLIFY (dsblock)
398 : CALL dbcsr_get_block_p(matrix=matrix_s(1 + i, 1)%matrix, &
399 691527 : row=irow, col=icol, block=dsblock, found=found)
400 691527 : CPASSERT(found)
401 4854888 : fi = -2.0_dp*gmij*SUM(pblock*dsblock)
402 691527 : force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
403 691527 : force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
404 1613563 : fij(i) = fi
405 : END DO
406 : END IF
407 : END DO
408 11990 : CALL dbcsr_iterator_stop(iter)
409 : !
410 11990 : IF (do_gamma_stress) THEN
411 270 : DO ikind = 1, nkind
412 630 : DO jkind = 1, nkind
413 360 : iac = ikind + nkind*(jkind - 1)
414 360 : IF (.NOT. ASSOCIATED(sap_int(iac)%alist)) CYCLE
415 8610 : DO ia = 1, sap_int(iac)%nalist
416 8082 : IF (.NOT. ASSOCIATED(sap_int(iac)%alist(ia)%clist)) CYCLE
417 8079 : iatom = sap_int(iac)%alist(ia)%aatom
418 177675 : DO ic = 1, sap_int(iac)%alist(ia)%nclist
419 169236 : jatom = sap_int(iac)%alist(ia)%clist(ic)%catom
420 676944 : rij = sap_int(iac)%alist(ia)%clist(ic)%rac
421 676944 : dr = SQRT(SUM(rij(:)**2))
422 177318 : IF (dr > 1.e-6_dp) THEN
423 165195 : dsint => sap_int(iac)%alist(ia)%clist(ic)%acint
424 165195 : gmij = 0.5_dp*(gmcharge(iatom, 1) + gmcharge(jatom, 1))
425 165195 : icol = MAX(iatom, jatom)
426 165195 : irow = MIN(iatom, jatom)
427 165195 : NULLIFY (pblock)
428 : CALL dbcsr_get_block_p(matrix=matrix_p(1, 1)%matrix, &
429 165195 : row=irow, col=icol, block=pblock, found=found)
430 165195 : CPASSERT(found)
431 660780 : DO i = 1, 3
432 660780 : IF (irow == iatom) THEN
433 1722726 : fij(i) = -2.0_dp*gmij*SUM(pblock*dsint(:, :, i))
434 : ELSE
435 1731411 : fij(i) = -2.0_dp*gmij*SUM(TRANSPOSE(pblock)*dsint(:, :, i))
436 : END IF
437 : END DO
438 165195 : fi = 1.0_dp
439 165195 : IF (iatom == jatom) fi = 0.5_dp
440 165195 : CALL virial_pair_force(virial%pv_virial, fi, fij, rij)
441 : END IF
442 : END DO
443 : END DO
444 : END DO
445 : END DO
446 : END IF
447 : ELSE
448 4442 : NULLIFY (n_list)
449 4442 : CALL get_qs_env(qs_env=qs_env, sab_orb=n_list)
450 4442 : CALL neighbor_list_iterator_create(nl_iterator, n_list)
451 550446 : DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
452 : CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
453 546004 : iatom=iatom, jatom=jatom, r=rij, cell=cellind)
454 :
455 546004 : icol = MAX(iatom, jatom)
456 546004 : irow = MIN(iatom, jatom)
457 546004 : ic = cell_to_index(cellind(1), cellind(2), cellind(3))
458 546004 : CPASSERT(ic > 0)
459 :
460 546004 : gmij = 0.5_dp*(gmcharge(iatom, 1) + gmcharge(jatom, 1))
461 :
462 546004 : NULLIFY (sblock)
463 : CALL dbcsr_get_block_p(matrix=matrix_s(1, ic)%matrix, &
464 546004 : row=irow, col=icol, block=sblock, found=found)
465 546004 : CPASSERT(found)
466 1092008 : DO is = 1, SIZE(ks_matrix, 1)
467 546004 : NULLIFY (ksblock)
468 : CALL dbcsr_get_block_p(matrix=ks_matrix(is, ic)%matrix, &
469 546004 : row=irow, col=icol, block=ksblock, found=found)
470 546004 : CPASSERT(found)
471 41390038 : ksblock = ksblock - gmij*sblock
472 : END DO
473 :
474 550446 : IF (calculate_forces) THEN
475 5962 : ikind = kind_of(iatom)
476 5962 : atom_i = atom_of_kind(iatom)
477 5962 : jkind = kind_of(jatom)
478 5962 : atom_j = atom_of_kind(jatom)
479 5962 : IF (irow == jatom) gmij = -gmij
480 5962 : NULLIFY (pblock)
481 : CALL dbcsr_get_block_p(matrix=matrix_p(1, ic)%matrix, &
482 5962 : row=irow, col=icol, block=pblock, found=found)
483 5962 : CPASSERT(found)
484 23848 : DO i = 1, 3
485 17886 : NULLIFY (dsblock)
486 : CALL dbcsr_get_block_p(matrix=matrix_s(1 + i, ic)%matrix, &
487 17886 : row=irow, col=icol, block=dsblock, found=found)
488 17886 : CPASSERT(found)
489 387837 : fi = -2.0_dp*gmij*SUM(pblock*dsblock)
490 17886 : force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
491 17886 : force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
492 41734 : fij(i) = fi
493 : END DO
494 5962 : IF (use_virial) THEN
495 : ! The image-resolved overlap derivative blocks are matrix oriented. Rebuild the
496 : ! atom-oriented derivative for the virial contraction, matching the Gamma path.
497 5868 : fij = 0.0_dp
498 5868 : CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind_a)
499 : CALL get_dftb_atom_param(dftb_kind_a, &
500 5868 : defined=defined_a, lmax=lmaxi, natorb=natorb_a)
501 5868 : CALL get_qs_kind(qs_kind_set(jkind), dftb_parameter=dftb_kind_b)
502 : CALL get_dftb_atom_param(dftb_kind_b, &
503 5868 : defined=defined_b, lmax=lmaxj, natorb=natorb_b)
504 5868 : IF (defined_a .AND. defined_b .AND. natorb_a > 0 .AND. natorb_b > 0) THEN
505 5868 : dftb_param_ij => dftb_potential(ikind, jkind)
506 5868 : dftb_param_ji => dftb_potential(jkind, ikind)
507 5868 : ngrd = dftb_param_ij%ngrd
508 5868 : ngrdcut = dftb_param_ij%ngrdcut
509 5868 : dgrd = dftb_param_ij%dgrd
510 5868 : ddr = dgrd*dftb_fd_deriv_step
511 5868 : CPASSERT(dftb_param_ij%llm == dftb_param_ji%llm)
512 5868 : llm = dftb_param_ij%llm
513 5868 : smatij => dftb_param_ij%smat
514 5868 : smatji => dftb_param_ji%smat
515 23472 : dr = SQRT(SUM(rij(:)**2))
516 5868 : IF (NINT(dr/dgrd) <= ngrdcut .AND. dr > 1.e-6_dp) THEN
517 5757 : n1 = natorb_a
518 5757 : n2 = natorb_b
519 34542 : ALLOCATE (dsblock_vir(n1, n2), dsblockm_vir(n1, n2))
520 5757 : gmij_virial = 0.5_dp*(gmcharge(iatom, 1) + gmcharge(jatom, 1))
521 23028 : DO i = 1, 3
522 379353 : dsblock_vir = 0._dp
523 379353 : dsblockm_vir = 0._dp
524 17271 : drij = rij
525 17271 : drij(i) = rij(i) - ddr
526 : CALL compute_block_sk(dsblockm_vir, smatij, smatji, drij, ngrd, ngrdcut, dgrd, &
527 17271 : llm, lmaxi, lmaxj, iatom, iatom)
528 17271 : drij(i) = rij(i) + ddr
529 : CALL compute_block_sk(dsblock_vir, smatij, smatji, drij, ngrd, ngrdcut, dgrd, &
530 17271 : llm, lmaxi, lmaxj, iatom, iatom)
531 741435 : dsblock_vir = (dsblock_vir - dsblockm_vir)/(2.0_dp*ddr)
532 23028 : IF (irow == iatom) THEN
533 208056 : fij(i) = 2.0_dp*gmij_virial*SUM(pblock*dsblock_vir)
534 : ELSE
535 171297 : fij(i) = 2.0_dp*gmij_virial*SUM(TRANSPOSE(pblock)*dsblock_vir)
536 : END IF
537 : END DO
538 5757 : DEALLOCATE (dsblock_vir, dsblockm_vir)
539 5757 : fi = 1.0_dp
540 5757 : IF (iatom == jatom) fi = 0.5_dp
541 5757 : CALL virial_pair_force(virial%pv_virial, fi, fij, rij)
542 : END IF
543 : END IF
544 : END IF
545 : END IF
546 : END DO
547 4442 : CALL neighbor_list_iterator_release(nl_iterator)
548 : END IF
549 :
550 16432 : IF (calculate_forces .AND. SIZE(matrix_p, 1) == 2) THEN
551 328 : DO img = 1, nimg
552 : CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
553 328 : alpha_scalar=1.0_dp, beta_scalar=-1.0_dp)
554 : END DO
555 : END IF
556 : END IF
557 :
558 16708 : IF (dft_control%qs_control%dftb_control%dftb3_diagonal) THEN
559 4634 : CALL get_qs_env(qs_env, nkind=nkind)
560 18536 : ALLOCATE (zeffk(nkind), xgamma(nkind))
561 13802 : DO ikind = 1, nkind
562 9168 : CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind)
563 13802 : CALL get_dftb_atom_param(dftb_kind, dudq=xgamma(ikind), zeff=zeffk(ikind))
564 : END DO
565 : ! Diagonal 3rd order correction (DFTB3)
566 : CALL build_dftb3_diagonal(qs_env, ks_matrix, rho, mcharge, energy, xgamma, zeffk, &
567 4634 : sap_int, calculate_forces, just_energy)
568 4634 : DEALLOCATE (zeffk, xgamma)
569 : END IF
570 :
571 : ! QMMM
572 16708 : IF (qs_env%qmmm .AND. qs_env%qmmm_periodic) THEN
573 : CALL build_tb_coulomb_qmqm(qs_env, ks_matrix, rho, mcharge, energy, &
574 1626 : calculate_forces, just_energy)
575 : END IF
576 :
577 16708 : IF (do_gamma_stress) THEN
578 90 : CALL release_sap_int(sap_int)
579 : END IF
580 :
581 16708 : DEALLOCATE (gmcharge)
582 :
583 16708 : CALL timestop(handle)
584 :
585 33416 : END SUBROUTINE build_dftb_coulomb
586 :
587 : ! **************************************************************************************************
588 : !> \brief Computes the short-range gamma parameter from exact Coulomb
589 : !> interaction of normalized exp(-a*r) charge distribution - 1/r
590 : !> \param r ...
591 : !> \param ga ...
592 : !> \param gb ...
593 : !> \param hb_para ...
594 : !> \return ...
595 : !> \par Literature
596 : !> Elstner et al, PRB 58 (1998) 7260
597 : !> \par History
598 : !> 10.2008 Axel Kohlmeyer - adding sr_damp
599 : !> 08.2014 JGH - adding flexible exponent for damping
600 : !> \version 1.1
601 : ! **************************************************************************************************
602 3330838 : FUNCTION gamma_rab_sr(r, ga, gb, hb_para) RESULT(gamma)
603 : REAL(dp), INTENT(in) :: r, ga, gb, hb_para
604 : REAL(dp) :: gamma
605 :
606 : REAL(dp) :: a, b, fac, g_sum
607 :
608 3330838 : gamma = 0.0_dp
609 3330838 : a = 3.2_dp*ga ! 3.2 = 16/5 in Eq. 18 and ff.
610 3330838 : b = 3.2_dp*gb
611 3330838 : g_sum = a + b
612 3330838 : IF (g_sum < tol_gamma) RETURN ! hardness screening
613 3330838 : IF (r < rtiny) THEN ! This is for short distances but non-onsite terms
614 : ! This gives also correct diagonal elements (a=b, r=0)
615 91457 : gamma = 0.5_dp*(a*b/g_sum + (a*b)**2/g_sum**3)
616 91457 : RETURN
617 : END IF
618 : !
619 : ! distinguish two cases: Gamma's are very close, e.g. for the same atom type,
620 : ! and Gamma's are different
621 : !
622 3239381 : IF (ABS(a - b) < rtiny) THEN
623 1905313 : fac = 1.6_dp*r*a*b/g_sum*(1.0_dp + a*b/g_sum**2)
624 1905313 : gamma = -(48.0_dp + 33._dp*fac + (9.0_dp + fac)*fac**2)*EXP(-fac)/(48._dp*r)
625 : ELSE
626 : gamma = -EXP(-a*r)*(0.5_dp*a*b**4/(a**2 - b**2)**2 - &
627 : (b**6 - 3._dp*a**2*b**4)/(r*(a**2 - b**2)**3)) - & ! a-> b
628 : EXP(-b*r)*(0.5_dp*b*a**4/(b**2 - a**2)**2 - &
629 1334068 : (a**6 - 3._dp*b**2*a**4)/(r*(b**2 - a**2)**3)) ! b-> a
630 : END IF
631 : !
632 : ! damping function for better short range hydrogen bonds.
633 : ! functional form from Hu H. et al., J. Phys. Chem. A 2007, 111, 5685-5691
634 : ! according to Elstner M, Theor. Chem. Acc. 2006, 116, 316-325,
635 : ! this should only be applied to a-b pairs involving hydrogen.
636 3239381 : IF (hb_para > 0.0_dp) THEN
637 493380 : gamma = gamma*EXP(-(0.5_dp*(ga + gb))**hb_para*r*r)
638 : END IF
639 : END FUNCTION gamma_rab_sr
640 :
641 : ! **************************************************************************************************
642 : !> \brief ...
643 : !> \param qs_env ...
644 : !> \param sap_int ...
645 : ! **************************************************************************************************
646 90 : SUBROUTINE dftb_dsint_list(qs_env, sap_int)
647 :
648 : TYPE(qs_environment_type), POINTER :: qs_env
649 : TYPE(sap_int_type), DIMENSION(:), POINTER :: sap_int
650 :
651 : CHARACTER(LEN=*), PARAMETER :: routineN = 'dftb_dsint_list'
652 :
653 : INTEGER :: handle, i, iac, iatom, ikind, ilist, jatom, jkind, jneighbor, llm, lmaxi, lmaxj, &
654 : n1, n2, natorb_a, natorb_b, ngrd, ngrdcut, nkind, nlist, nneighbor
655 : INTEGER, DIMENSION(3) :: cell
656 : LOGICAL :: defined
657 : REAL(KIND=dp) :: ddr, dgrd, dr
658 : REAL(KIND=dp), DIMENSION(3) :: drij, rij
659 90 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: dsblock, dsblockm, smatij, smatji
660 : TYPE(clist_type), POINTER :: clist
661 : TYPE(dft_control_type), POINTER :: dft_control
662 : TYPE(dftb_control_type), POINTER :: dftb_control
663 : TYPE(neighbor_list_iterator_p_type), &
664 90 : DIMENSION(:), POINTER :: nl_iterator
665 : TYPE(neighbor_list_set_p_type), DIMENSION(:), &
666 90 : POINTER :: sab_orb
667 : TYPE(qs_dftb_atom_type), POINTER :: dftb_kind_a, dftb_kind_b
668 : TYPE(qs_dftb_pairpot_type), DIMENSION(:, :), &
669 90 : POINTER :: dftb_potential
670 : TYPE(qs_dftb_pairpot_type), POINTER :: dftb_param_ij, dftb_param_ji
671 90 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
672 :
673 90 : CALL timeset(routineN, handle)
674 :
675 90 : CALL get_qs_env(qs_env=qs_env, nkind=nkind)
676 90 : CPASSERT(.NOT. ASSOCIATED(sap_int))
677 630 : ALLOCATE (sap_int(nkind*nkind))
678 450 : DO i = 1, nkind*nkind
679 360 : NULLIFY (sap_int(i)%alist, sap_int(i)%asort, sap_int(i)%aindex)
680 450 : sap_int(i)%nalist = 0
681 : END DO
682 :
683 : CALL get_qs_env(qs_env=qs_env, &
684 : qs_kind_set=qs_kind_set, &
685 : dft_control=dft_control, &
686 90 : sab_orb=sab_orb)
687 :
688 90 : dftb_control => dft_control%qs_control%dftb_control
689 :
690 90 : NULLIFY (dftb_potential)
691 : CALL get_qs_env(qs_env=qs_env, &
692 90 : dftb_potential=dftb_potential)
693 :
694 : ! loop over all atom pairs with a non-zero overlap (sab_orb)
695 90 : CALL neighbor_list_iterator_create(nl_iterator, sab_orb)
696 169326 : DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
697 : CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, iatom=iatom, &
698 : jatom=jatom, nlist=nlist, ilist=ilist, nnode=nneighbor, &
699 169236 : inode=jneighbor, cell=cell, r=rij)
700 169236 : iac = ikind + nkind*(jkind - 1)
701 : !
702 169236 : CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_kind_a)
703 : CALL get_dftb_atom_param(dftb_kind_a, &
704 169236 : defined=defined, lmax=lmaxi, natorb=natorb_a)
705 169236 : IF (.NOT. defined .OR. natorb_a < 1) CYCLE
706 169236 : CALL get_qs_kind(qs_kind_set(jkind), dftb_parameter=dftb_kind_b)
707 : CALL get_dftb_atom_param(dftb_kind_b, &
708 169236 : defined=defined, lmax=lmaxj, natorb=natorb_b)
709 169236 : IF (.NOT. defined .OR. natorb_b < 1) CYCLE
710 :
711 676944 : dr = SQRT(SUM(rij(:)**2))
712 :
713 : ! integral list
714 169236 : IF (.NOT. ASSOCIATED(sap_int(iac)%alist)) THEN
715 348 : sap_int(iac)%a_kind = ikind
716 348 : sap_int(iac)%p_kind = jkind
717 348 : sap_int(iac)%nalist = nlist
718 9126 : ALLOCATE (sap_int(iac)%alist(nlist))
719 8430 : DO i = 1, nlist
720 8082 : NULLIFY (sap_int(iac)%alist(i)%clist)
721 8082 : sap_int(iac)%alist(i)%aatom = 0
722 8430 : sap_int(iac)%alist(i)%nclist = 0
723 : END DO
724 : END IF
725 169236 : IF (.NOT. ASSOCIATED(sap_int(iac)%alist(ilist)%clist)) THEN
726 8079 : sap_int(iac)%alist(ilist)%aatom = iatom
727 8079 : sap_int(iac)%alist(ilist)%nclist = nneighbor
728 241947 : ALLOCATE (sap_int(iac)%alist(ilist)%clist(nneighbor))
729 177315 : DO i = 1, nneighbor
730 177315 : sap_int(iac)%alist(ilist)%clist(i)%catom = 0
731 : END DO
732 : END IF
733 169236 : clist => sap_int(iac)%alist(ilist)%clist(jneighbor)
734 169236 : clist%catom = jatom
735 676944 : clist%cell = cell
736 676944 : clist%rac = rij
737 846180 : ALLOCATE (clist%acint(natorb_a, natorb_b, 3))
738 169236 : NULLIFY (clist%achint)
739 3732480 : clist%acint = 0._dp
740 169236 : clist%nsgf_cnt = 0
741 169236 : NULLIFY (clist%sgf_list)
742 :
743 : ! retrieve information on S matrix
744 169236 : dftb_param_ij => dftb_potential(ikind, jkind)
745 169236 : dftb_param_ji => dftb_potential(jkind, ikind)
746 : ! assume table size and type is symmetric
747 169236 : ngrd = dftb_param_ij%ngrd
748 169236 : ngrdcut = dftb_param_ij%ngrdcut
749 169236 : dgrd = dftb_param_ij%dgrd
750 169236 : ddr = dgrd*dftb_fd_deriv_step
751 169236 : CPASSERT(dftb_param_ij%llm == dftb_param_ji%llm)
752 169236 : llm = dftb_param_ij%llm
753 169236 : smatij => dftb_param_ij%smat
754 169236 : smatji => dftb_param_ji%smat
755 :
756 507798 : IF (NINT(dr/dgrd) <= ngrdcut) THEN
757 168862 : IF (iatom == jatom .AND. dr < 0.001_dp) THEN
758 : ! diagonal block
759 : ELSE
760 : ! off-diagonal block
761 164821 : n1 = natorb_a
762 164821 : n2 = natorb_b
763 988926 : ALLOCATE (dsblock(n1, n2), dsblockm(n1, n2))
764 659284 : DO i = 1, 3
765 3445497 : dsblock = 0._dp
766 3445497 : dsblockm = 0.0_dp
767 494463 : drij = rij
768 :
769 494463 : drij(i) = rij(i) - ddr
770 : CALL compute_block_sk(dsblockm, smatij, smatji, drij, ngrd, ngrdcut, dgrd, &
771 494463 : llm, lmaxi, lmaxj, iatom, iatom)
772 :
773 494463 : drij(i) = rij(i) + ddr
774 : CALL compute_block_sk(dsblock, smatij, smatji, drij, ngrd, ngrdcut, dgrd, &
775 494463 : llm, lmaxi, lmaxj, iatom, iatom)
776 :
777 6396531 : dsblock = dsblock - dsblockm
778 3445497 : dsblock = dsblock/(2.0_dp*ddr)
779 :
780 6561352 : clist%acint(1:n1, 1:n2, i) = -dsblock(1:n1, 1:n2)
781 : END DO
782 164821 : DEALLOCATE (dsblock, dsblockm)
783 : END IF
784 : END IF
785 :
786 : END DO
787 90 : CALL neighbor_list_iterator_release(nl_iterator)
788 :
789 90 : CALL timestop(handle)
790 :
791 90 : END SUBROUTINE dftb_dsint_list
792 :
793 : END MODULE qs_dftb_coulomb
|