(git:374b731)
Loading...
Searching...
No Matches
xtb_ehess_force.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 forces for Coulomb contributions in response xTB
10!> \author JGH
11! **************************************************************************************************
15 USE cell_types, ONLY: cell_type,&
16 get_cell,&
17 pbc
23 USE dbcsr_api, ONLY: &
24 dbcsr_add, dbcsr_get_block_p, dbcsr_iterator_blocks_left, dbcsr_iterator_next_block, &
25 dbcsr_iterator_start, dbcsr_iterator_stop, dbcsr_iterator_type, dbcsr_p_type, dbcsr_type
32 USE kinds, ONLY: dp
33 USE mathconstants, ONLY: oorootpi,&
34 pi
45 USE qs_kind_types, ONLY: get_qs_kind,&
53 USE qs_rho_types, ONLY: qs_rho_type
54 USE virial_types, ONLY: virial_type
55 USE xtb_coulomb, ONLY: dgamma_rab_sr,&
57 USE xtb_types, ONLY: get_xtb_atom_param,&
59#include "./base/base_uses.f90"
60
61 IMPLICIT NONE
62
63 PRIVATE
64
65 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'xtb_ehess_force'
66
67 PUBLIC :: calc_xtb_ehess_force
68
69! **************************************************************************************************
70
71CONTAINS
72
73! **************************************************************************************************
74!> \brief ...
75!> \param qs_env ...
76!> \param matrix_p0 ...
77!> \param matrix_p1 ...
78!> \param charges0 ...
79!> \param mcharge0 ...
80!> \param charges1 ...
81!> \param mcharge1 ...
82!> \param debug_forces ...
83! **************************************************************************************************
84 SUBROUTINE calc_xtb_ehess_force(qs_env, matrix_p0, matrix_p1, charges0, mcharge0, &
85 charges1, mcharge1, debug_forces)
86
87 TYPE(qs_environment_type), POINTER :: qs_env
88 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_p0, matrix_p1
89 REAL(kind=dp), DIMENSION(:, :), INTENT(in) :: charges0
90 REAL(kind=dp), DIMENSION(:), INTENT(in) :: mcharge0
91 REAL(kind=dp), DIMENSION(:, :), INTENT(in) :: charges1
92 REAL(kind=dp), DIMENSION(:), INTENT(in) :: mcharge1
93 LOGICAL, INTENT(IN) :: debug_forces
94
95 CHARACTER(len=*), PARAMETER :: routinen = 'calc_xtb_ehess_force'
96
97 INTEGER :: atom_i, atom_j, blk, ewald_type, handle, i, ia, iatom, icol, ikind, iounit, irow, &
98 j, jatom, jkind, la, lb, lmaxa, lmaxb, natom, natorb_a, natorb_b, ni, nimg, nj, nkind, &
99 nmat, za, zb
100 INTEGER, ALLOCATABLE, DIMENSION(:) :: atom_of_kind, kind_of
101 INTEGER, DIMENSION(25) :: laoa, laob
102 INTEGER, DIMENSION(3) :: cellind, periodic
103 LOGICAL :: calculate_forces, defined, do_ewald, &
104 found, just_energy, use_virial
105 REAL(kind=dp) :: alpha, deth, dr, etaa, etab, fi, gmij0, &
106 gmij1, kg, rcut, rcuta, rcutb
107 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: xgamma
108 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: gammab, gcij0, gcij1, gmcharge0, &
109 gmcharge1
110 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: gchrg0, gchrg1
111 REAL(kind=dp), DIMENSION(3) :: fij, fodeb, rij
112 REAL(kind=dp), DIMENSION(5) :: kappaa, kappab
113 REAL(kind=dp), DIMENSION(:, :), POINTER :: dsblock, pblock0, pblock1, sblock
114 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
115 TYPE(cell_type), POINTER :: cell
116 TYPE(cp_logger_type), POINTER :: logger
117 TYPE(dbcsr_iterator_type) :: iter
118 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_s
119 TYPE(dft_control_type), POINTER :: dft_control
120 TYPE(distribution_1d_type), POINTER :: local_particles
121 TYPE(ewald_environment_type), POINTER :: ewald_env
122 TYPE(ewald_pw_type), POINTER :: ewald_pw
123 TYPE(mp_para_env_type), POINTER :: para_env
125 DIMENSION(:), POINTER :: nl_iterator
126 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
127 POINTER :: n_list
128 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
129 TYPE(qs_energy_type), POINTER :: energy
130 TYPE(qs_force_type), DIMENSION(:), POINTER :: force
131 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
132 TYPE(qs_rho_type), POINTER :: rho
133 TYPE(virial_type), POINTER :: virial
134 TYPE(xtb_atom_type), POINTER :: xtb_atom_a, xtb_atom_b, xtb_kind
135 TYPE(xtb_control_type), POINTER :: xtb_control
136
137 CALL timeset(routinen, handle)
138
139 logger => cp_get_default_logger()
140 IF (logger%para_env%is_source()) THEN
141 iounit = cp_logger_get_default_unit_nr(logger, local=.true.)
142 ELSE
143 iounit = -1
144 END IF
145
146 cpassert(ASSOCIATED(matrix_p1))
147
148 CALL get_qs_env(qs_env, &
149 qs_kind_set=qs_kind_set, &
150 particle_set=particle_set, &
151 cell=cell, &
152 rho=rho, &
153 energy=energy, &
154 virial=virial, &
155 dft_control=dft_control)
156
157 xtb_control => dft_control%qs_control%xtb_control
158
159 calculate_forces = .true.
160 just_energy = .false.
161 use_virial = .false.
162 nmat = 4
163 nimg = dft_control%nimages
164 IF (nimg > 1) THEN
165 cpabort('xTB-sTDA forces for k-points not available')
166 END IF
167
168 CALL get_qs_env(qs_env, nkind=nkind, natom=natom)
169 ALLOCATE (gchrg0(natom, 5, nmat))
170 gchrg0 = 0._dp
171 ALLOCATE (gmcharge0(natom, nmat))
172 gmcharge0 = 0._dp
173 ALLOCATE (gchrg1(natom, 5, nmat))
174 gchrg1 = 0._dp
175 ALLOCATE (gmcharge1(natom, nmat))
176 gmcharge1 = 0._dp
177
178 ! short range contribution (gamma)
179 ! loop over all atom pairs (sab_xtbe)
180 kg = xtb_control%kg
181 NULLIFY (n_list)
182 IF (xtb_control%old_coulomb_damping) THEN
183 CALL get_qs_env(qs_env=qs_env, sab_orb=n_list)
184 ELSE
185 CALL get_qs_env(qs_env=qs_env, sab_xtbe=n_list)
186 END IF
187 CALL neighbor_list_iterator_create(nl_iterator, n_list)
188 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
189 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
190 iatom=iatom, jatom=jatom, r=rij, cell=cellind)
191 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
192 CALL get_xtb_atom_param(xtb_atom_a, defined=defined, natorb=natorb_a)
193 IF (.NOT. defined .OR. natorb_a < 1) cycle
194 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
195 CALL get_xtb_atom_param(xtb_atom_b, defined=defined, natorb=natorb_b)
196 IF (.NOT. defined .OR. natorb_b < 1) cycle
197 ! atomic parameters
198 CALL get_xtb_atom_param(xtb_atom_a, eta=etaa, lmax=lmaxa, kappa=kappaa, rcut=rcuta)
199 CALL get_xtb_atom_param(xtb_atom_b, eta=etab, lmax=lmaxb, kappa=kappab, rcut=rcutb)
200 ! gamma matrix
201 ni = lmaxa + 1
202 nj = lmaxb + 1
203 ALLOCATE (gammab(ni, nj))
204 rcut = rcuta + rcutb
205 dr = sqrt(sum(rij(:)**2))
206 CALL gamma_rab_sr(gammab, dr, ni, kappaa, etaa, nj, kappab, etab, kg, rcut)
207 gchrg0(iatom, 1:ni, 1) = gchrg0(iatom, 1:ni, 1) + matmul(gammab, charges0(jatom, 1:nj))
208 gchrg1(iatom, 1:ni, 1) = gchrg1(iatom, 1:ni, 1) + matmul(gammab, charges1(jatom, 1:nj))
209 IF (iatom /= jatom) THEN
210 gchrg0(jatom, 1:nj, 1) = gchrg0(jatom, 1:nj, 1) + matmul(charges0(iatom, 1:ni), gammab)
211 gchrg1(jatom, 1:nj, 1) = gchrg1(jatom, 1:nj, 1) + matmul(charges1(iatom, 1:ni), gammab)
212 END IF
213 IF (dr > 1.e-6_dp) THEN
214 CALL dgamma_rab_sr(gammab, dr, ni, kappaa, etaa, nj, kappab, etab, kg, rcut)
215 DO i = 1, 3
216 gchrg0(iatom, 1:ni, i + 1) = gchrg0(iatom, 1:ni, i + 1) &
217 + matmul(gammab, charges0(jatom, 1:nj))*rij(i)/dr
218 gchrg1(iatom, 1:ni, i + 1) = gchrg1(iatom, 1:ni, i + 1) &
219 + matmul(gammab, charges1(jatom, 1:nj))*rij(i)/dr
220 IF (iatom /= jatom) THEN
221 gchrg0(jatom, 1:nj, i + 1) = gchrg0(jatom, 1:nj, i + 1) &
222 - matmul(charges0(iatom, 1:ni), gammab)*rij(i)/dr
223 gchrg1(jatom, 1:nj, i + 1) = gchrg1(jatom, 1:nj, i + 1) &
224 - matmul(charges1(iatom, 1:ni), gammab)*rij(i)/dr
225 END IF
226 END DO
227 END IF
228 DEALLOCATE (gammab)
229 END DO
230 CALL neighbor_list_iterator_release(nl_iterator)
231
232 ! 1/R contribution
233
234 IF (xtb_control%coulomb_lr) THEN
235 do_ewald = xtb_control%do_ewald
236 IF (do_ewald) THEN
237 ! Ewald sum
238 NULLIFY (ewald_env, ewald_pw)
239 CALL get_qs_env(qs_env=qs_env, &
240 ewald_env=ewald_env, ewald_pw=ewald_pw)
241 CALL get_cell(cell=cell, periodic=periodic, deth=deth)
242 CALL ewald_env_get(ewald_env, alpha=alpha, ewald_type=ewald_type)
243 CALL get_qs_env(qs_env=qs_env, sab_tbe=n_list)
244 CALL tb_ewald_overlap(gmcharge0, mcharge0, alpha, n_list, virial, use_virial)
245 CALL tb_ewald_overlap(gmcharge1, mcharge1, alpha, n_list, virial, use_virial)
246 SELECT CASE (ewald_type)
247 CASE DEFAULT
248 cpabort("Invalid Ewald type")
249 CASE (do_ewald_none)
250 cpabort("Not allowed with DFTB")
251 CASE (do_ewald_ewald)
252 cpabort("Standard Ewald not implemented in DFTB")
253 CASE (do_ewald_pme)
254 cpabort("PME not implemented in DFTB")
255 CASE (do_ewald_spme)
256 CALL tb_spme_zforce(ewald_env, ewald_pw, particle_set, cell, gmcharge0, mcharge0)
257 CALL tb_spme_zforce(ewald_env, ewald_pw, particle_set, cell, gmcharge1, mcharge1)
258 END SELECT
259 ELSE
260 ! direct sum
261 CALL get_qs_env(qs_env=qs_env, local_particles=local_particles)
262 DO ikind = 1, SIZE(local_particles%n_el)
263 DO ia = 1, local_particles%n_el(ikind)
264 iatom = local_particles%list(ikind)%array(ia)
265 DO jatom = 1, iatom - 1
266 rij = particle_set(iatom)%r - particle_set(jatom)%r
267 rij = pbc(rij, cell)
268 dr = sqrt(sum(rij(:)**2))
269 IF (dr > 1.e-6_dp) THEN
270 gmcharge0(iatom, 1) = gmcharge0(iatom, 1) + mcharge0(jatom)/dr
271 gmcharge0(jatom, 1) = gmcharge0(jatom, 1) + mcharge0(iatom)/dr
272 gmcharge1(iatom, 1) = gmcharge1(iatom, 1) + mcharge1(jatom)/dr
273 gmcharge1(jatom, 1) = gmcharge1(jatom, 1) + mcharge1(iatom)/dr
274 DO i = 2, nmat
275 gmcharge0(iatom, i) = gmcharge0(iatom, i) + rij(i - 1)*mcharge0(jatom)/dr**3
276 gmcharge0(jatom, i) = gmcharge0(jatom, i) - rij(i - 1)*mcharge0(iatom)/dr**3
277 gmcharge1(iatom, i) = gmcharge1(iatom, i) + rij(i - 1)*mcharge1(jatom)/dr**3
278 gmcharge1(jatom, i) = gmcharge1(jatom, i) - rij(i - 1)*mcharge1(iatom)/dr**3
279 END DO
280 END IF
281 END DO
282 END DO
283 END DO
284 cpassert(.NOT. use_virial)
285 END IF
286 END IF
287
288 ! global sum of gamma*p arrays
289 CALL get_qs_env(qs_env=qs_env, &
290 atomic_kind_set=atomic_kind_set, &
291 force=force, para_env=para_env)
292 CALL para_env%sum(gmcharge0(:, 1))
293 CALL para_env%sum(gchrg0(:, :, 1))
294 CALL para_env%sum(gmcharge1(:, 1))
295 CALL para_env%sum(gchrg1(:, :, 1))
296
297 IF (xtb_control%coulomb_lr) THEN
298 IF (do_ewald) THEN
299 ! add self charge interaction and background charge contribution
300 gmcharge0(:, 1) = gmcharge0(:, 1) - 2._dp*alpha*oorootpi*mcharge0(:)
301 IF (any(periodic(:) == 1)) THEN
302 gmcharge0(:, 1) = gmcharge0(:, 1) - pi/alpha**2/deth
303 END IF
304 gmcharge1(:, 1) = gmcharge1(:, 1) - 2._dp*alpha*oorootpi*mcharge1(:)
305 IF (any(periodic(:) == 1)) THEN
306 gmcharge1(:, 1) = gmcharge1(:, 1) - pi/alpha**2/deth
307 END IF
308 END IF
309 END IF
310
311 CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, &
312 kind_of=kind_of, &
313 atom_of_kind=atom_of_kind)
314
315 IF (debug_forces) fodeb(1:3) = force(1)%rho_elec(1:3, 1)
316 DO iatom = 1, natom
317 ikind = kind_of(iatom)
318 atom_i = atom_of_kind(iatom)
319 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_kind)
320 CALL get_xtb_atom_param(xtb_kind, lmax=ni)
321 ni = ni + 1
322 ! short range
323 fij = 0.0_dp
324 DO i = 1, 3
325 fij(i) = sum(charges0(iatom, 1:ni)*gchrg1(iatom, 1:ni, i + 1)) + &
326 sum(charges1(iatom, 1:ni)*gchrg0(iatom, 1:ni, i + 1))
327 END DO
328 force(ikind)%rho_elec(1, atom_i) = force(ikind)%rho_elec(1, atom_i) - fij(1)
329 force(ikind)%rho_elec(2, atom_i) = force(ikind)%rho_elec(2, atom_i) - fij(2)
330 force(ikind)%rho_elec(3, atom_i) = force(ikind)%rho_elec(3, atom_i) - fij(3)
331 ! long range
332 fij = 0.0_dp
333 DO i = 1, 3
334 fij(i) = gmcharge1(iatom, i + 1)*mcharge0(iatom) + &
335 gmcharge0(iatom, i + 1)*mcharge1(iatom)
336 END DO
337 force(ikind)%rho_elec(1, atom_i) = force(ikind)%rho_elec(1, atom_i) - fij(1)
338 force(ikind)%rho_elec(2, atom_i) = force(ikind)%rho_elec(2, atom_i) - fij(2)
339 force(ikind)%rho_elec(3, atom_i) = force(ikind)%rho_elec(3, atom_i) - fij(3)
340 END DO
341 IF (debug_forces) THEN
342 fodeb(1:3) = force(1)%rho_elec(1:3, 1) - fodeb(1:3)
343 CALL para_env%sum(fodeb)
344 IF (iounit > 0) WRITE (iounit, "(T3,A,T33,3F16.8)") "DEBUG:: P*dH[Pz] ", fodeb
345 END IF
346
347 CALL get_qs_env(qs_env=qs_env, matrix_s_kp=matrix_s)
348
349 IF (SIZE(matrix_p0) == 2) THEN
350 CALL dbcsr_add(matrix_p0(1)%matrix, matrix_p0(2)%matrix, &
351 alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
352 CALL dbcsr_add(matrix_p1(1)%matrix, matrix_p1(2)%matrix, &
353 alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
354 END IF
355
356 ! no k-points; all matrices have been transformed to periodic bsf
357 IF (debug_forces) fodeb(1:3) = force(1)%rho_elec(1:3, 1)
358 CALL dbcsr_iterator_start(iter, matrix_s(1, 1)%matrix)
359 DO WHILE (dbcsr_iterator_blocks_left(iter))
360 CALL dbcsr_iterator_next_block(iter, irow, icol, sblock, blk)
361 ikind = kind_of(irow)
362 jkind = kind_of(icol)
363
364 ! atomic parameters
365 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
366 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
367 CALL get_xtb_atom_param(xtb_atom_a, z=za, lao=laoa)
368 CALL get_xtb_atom_param(xtb_atom_b, z=zb, lao=laob)
369
370 ni = SIZE(sblock, 1)
371 nj = SIZE(sblock, 2)
372 ALLOCATE (gcij0(ni, nj))
373 ALLOCATE (gcij1(ni, nj))
374 DO i = 1, ni
375 DO j = 1, nj
376 la = laoa(i) + 1
377 lb = laob(j) + 1
378 gcij0(i, j) = 0.5_dp*(gchrg0(irow, la, 1) + gchrg0(icol, lb, 1))
379 gcij1(i, j) = 0.5_dp*(gchrg1(irow, la, 1) + gchrg1(icol, lb, 1))
380 END DO
381 END DO
382 gmij0 = 0.5_dp*(gmcharge0(irow, 1) + gmcharge0(icol, 1))
383 gmij1 = 0.5_dp*(gmcharge1(irow, 1) + gmcharge1(icol, 1))
384 atom_i = atom_of_kind(irow)
385 atom_j = atom_of_kind(icol)
386 NULLIFY (pblock0)
387 CALL dbcsr_get_block_p(matrix=matrix_p0(1)%matrix, &
388 row=irow, col=icol, block=pblock0, found=found)
389 cpassert(found)
390 NULLIFY (pblock1)
391 CALL dbcsr_get_block_p(matrix=matrix_p1(1)%matrix, &
392 row=irow, col=icol, block=pblock1, found=found)
393 cpassert(found)
394 DO i = 1, 3
395 NULLIFY (dsblock)
396 CALL dbcsr_get_block_p(matrix=matrix_s(1 + i, 1)%matrix, &
397 row=irow, col=icol, block=dsblock, found=found)
398 cpassert(found)
399 ! short range
400 fi = -2.0_dp*sum(pblock0*dsblock*gcij1) - 2.0_dp*sum(pblock1*dsblock*gcij0)
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 ! long range
404 fi = -2.0_dp*gmij1*sum(pblock0*dsblock) - 2.0_dp*gmij0*sum(pblock1*dsblock)
405 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
406 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
407 END DO
408 DEALLOCATE (gcij0, gcij1)
409 END DO
410 CALL dbcsr_iterator_stop(iter)
411 IF (debug_forces) THEN
412 fodeb(1:3) = force(1)%rho_elec(1:3, 1) - fodeb(1:3)
413 CALL para_env%sum(fodeb)
414 IF (iounit > 0) WRITE (iounit, "(T3,A,T33,3F16.8)") "DEBUG:: Pz*H[P]*dS ", fodeb
415 END IF
416
417 IF (xtb_control%tb3_interaction) THEN
418 CALL get_qs_env(qs_env, nkind=nkind)
419 ALLOCATE (xgamma(nkind))
420 DO ikind = 1, nkind
421 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_kind)
422 CALL get_xtb_atom_param(xtb_kind, xgamma=xgamma(ikind))
423 END DO
424 ! Diagonal 3rd order correction (DFTB3)
425 IF (debug_forces) fodeb(1:3) = force(1)%rho_elec(1:3, 1)
426 CALL dftb3_diagonal_hessian_force(qs_env, mcharge0, mcharge1, &
427 matrix_p0(1)%matrix, matrix_p1(1)%matrix, xgamma)
428 IF (debug_forces) THEN
429 fodeb(1:3) = force(1)%rho_elec(1:3, 1) - fodeb(1:3)
430 CALL para_env%sum(fodeb)
431 IF (iounit > 0) WRITE (iounit, "(T3,A,T33,3F16.8)") "DEBUG:: Pz*H3[P] ", fodeb
432 END IF
433 DEALLOCATE (xgamma)
434 END IF
435
436 IF (SIZE(matrix_p0) == 2) THEN
437 CALL dbcsr_add(matrix_p0(1)%matrix, matrix_p0(2)%matrix, &
438 alpha_scalar=1.0_dp, beta_scalar=-1.0_dp)
439 CALL dbcsr_add(matrix_p1(1)%matrix, matrix_p1(2)%matrix, &
440 alpha_scalar=1.0_dp, beta_scalar=-1.0_dp)
441 END IF
442
443 ! QMMM
444 IF (qs_env%qmmm .AND. qs_env%qmmm_periodic) THEN
445 cpabort("Not Available")
446 END IF
447
448 DEALLOCATE (gmcharge0, gchrg0, gmcharge1, gchrg1)
449
450 CALL timestop(handle)
451
452 END SUBROUTINE calc_xtb_ehess_force
453
454! **************************************************************************************************
455!> \brief ...
456!> \param qs_env ...
457!> \param mcharge0 ...
458!> \param mcharge1 ...
459!> \param matrixp0 ...
460!> \param matrixp1 ...
461!> \param xgamma ...
462! **************************************************************************************************
463 SUBROUTINE dftb3_diagonal_hessian_force(qs_env, mcharge0, mcharge1, &
464 matrixp0, matrixp1, xgamma)
465
466 TYPE(qs_environment_type), POINTER :: qs_env
467 REAL(dp), DIMENSION(:) :: mcharge0, mcharge1
468 TYPE(dbcsr_type), POINTER :: matrixp0, matrixp1
469 REAL(dp), DIMENSION(:) :: xgamma
470
471 CHARACTER(len=*), PARAMETER :: routinen = 'dftb3_diagonal_hessian_force'
472
473 INTEGER :: atom_i, atom_j, blk, handle, i, icol, &
474 ikind, irow, jkind
475 INTEGER, ALLOCATABLE, DIMENSION(:) :: atom_of_kind, kind_of
476 LOGICAL :: found
477 REAL(kind=dp) :: fi, gmijp, gmijq, ui, uj
478 REAL(kind=dp), DIMENSION(:, :), POINTER :: dsblock, p0block, p1block, sblock
479 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
480 TYPE(dbcsr_iterator_type) :: iter
481 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s
482 TYPE(qs_force_type), DIMENSION(:), POINTER :: force
483
484 CALL timeset(routinen, handle)
485 CALL get_qs_env(qs_env=qs_env, matrix_s=matrix_s)
486 CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set)
487 CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, &
488 kind_of=kind_of, atom_of_kind=atom_of_kind)
489 CALL get_qs_env(qs_env=qs_env, force=force)
490 ! no k-points; all matrices have been transformed to periodic bsf
491 CALL dbcsr_iterator_start(iter, matrix_s(1)%matrix)
492 DO WHILE (dbcsr_iterator_blocks_left(iter))
493 CALL dbcsr_iterator_next_block(iter, irow, icol, sblock, blk)
494 ikind = kind_of(irow)
495 atom_i = atom_of_kind(irow)
496 ui = xgamma(ikind)
497 jkind = kind_of(icol)
498 atom_j = atom_of_kind(icol)
499 uj = xgamma(jkind)
500 !
501 gmijp = ui*mcharge0(irow)*mcharge1(irow) + uj*mcharge0(icol)*mcharge1(icol)
502 gmijq = 0.5_dp*ui*mcharge0(irow)**2 + 0.5_dp*uj*mcharge0(icol)**2
503 !
504 NULLIFY (p0block)
505 CALL dbcsr_get_block_p(matrix=matrixp0, &
506 row=irow, col=icol, block=p0block, found=found)
507 cpassert(found)
508 NULLIFY (p1block)
509 CALL dbcsr_get_block_p(matrix=matrixp1, &
510 row=irow, col=icol, block=p1block, found=found)
511 cpassert(found)
512 DO i = 1, 3
513 NULLIFY (dsblock)
514 CALL dbcsr_get_block_p(matrix=matrix_s(1 + i)%matrix, &
515 row=irow, col=icol, block=dsblock, found=found)
516 cpassert(found)
517 fi = gmijp*sum(p0block*dsblock) + gmijq*sum(p1block*dsblock)
518 force(ikind)%rho_elec(i, atom_i) = force(ikind)%rho_elec(i, atom_i) + fi
519 force(jkind)%rho_elec(i, atom_j) = force(jkind)%rho_elec(i, atom_j) - fi
520 END DO
521 END DO
522 CALL dbcsr_iterator_stop(iter)
523
524 CALL timestop(handle)
525
526 END SUBROUTINE dftb3_diagonal_hessian_force
527
528END MODULE xtb_ehess_force
529
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.
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...
various routines to log and control the output. The idea is that decisions about where to log should ...
recursive integer function, public cp_logger_get_default_unit_nr(logger, local, skip_not_ionode)
asks the default unit number of the given logger. try to use cp_logger_get_unit_nr
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
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_spme_zforce(ewald_env, ewald_pw, particle_set, box, gmcharge, mcharge)
...
subroutine, public tb_ewald_overlap(gmcharge, mcharge, alpha, n_list, virial, use_virial, atprop)
...
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
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.
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
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...
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...
Calculation of forces for Coulomb contributions in response xTB.
subroutine, public calc_xtb_ehess_force(qs_env, matrix_p0, matrix_p1, charges0, mcharge0, charges1, mcharge1, debug_forces)
...
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, occupation, electronegativity, chmax)
...
Definition xtb_types.F:175
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
type of a logger, at the moment it contains just a print level starting at which level it should be l...
structure to store local (to a processor) ordered lists of integers.
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.