(git:ed6f26b)
Loading...
Searching...
No Matches
xtb_ehess.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Calculation of Coulomb Hessian contributions in xTB
10!> \author JGH
11! **************************************************************************************************
15 USE cell_types, ONLY: cell_type,&
16 get_cell,&
17 pbc
33 USE kinds, ONLY: dp
34 USE mathconstants, ONLY: oorootpi,&
35 pi
44 USE qs_kind_types, ONLY: get_qs_kind,&
52 USE virial_types, ONLY: virial_type
53 USE xtb_coulomb, ONLY: gamma_rab_sr
54 USE xtb_types, ONLY: get_xtb_atom_param,&
56#include "./base/base_uses.f90"
57
58 IMPLICIT NONE
59
60 PRIVATE
61
62 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'xtb_ehess'
63
64 PUBLIC :: xtb_coulomb_hessian
65
66CONTAINS
67
68! **************************************************************************************************
69!> \brief ...
70!> \param qs_env ...
71!> \param ks_matrix ...
72!> \param charges1 ...
73!> \param mcharge1 ...
74!> \param mcharge ...
75! **************************************************************************************************
76 SUBROUTINE xtb_coulomb_hessian(qs_env, ks_matrix, charges1, mcharge1, mcharge)
77
78 TYPE(qs_environment_type), POINTER :: qs_env
79 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: ks_matrix
80 REAL(dp), DIMENSION(:, :) :: charges1
81 REAL(dp), DIMENSION(:) :: mcharge1, mcharge
82
83 CHARACTER(len=*), PARAMETER :: routinen = 'xtb_coulomb_hessian'
84
85 INTEGER :: ewald_type, handle, i, ia, iatom, icol, ikind, irow, is, j, jatom, jkind, la, lb, &
86 lmaxa, lmaxb, natom, natorb_a, natorb_b, ni, nj, nkind, nmat, za, zb
87 INTEGER, ALLOCATABLE, DIMENSION(:) :: kind_of
88 INTEGER, DIMENSION(25) :: laoa, laob
89 INTEGER, DIMENSION(3) :: cellind, periodic
90 LOGICAL :: defined, do_ewald, found
91 REAL(kind=dp) :: alpha, deth, dr, etaa, etab, gmij, kg, &
92 rcut, rcuta, rcutb
93 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: xgamma
94 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: gammab, gcij, gmcharge
95 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: gchrg
96 REAL(kind=dp), DIMENSION(3) :: rij
97 REAL(kind=dp), DIMENSION(5) :: kappaa, kappab
98 REAL(kind=dp), DIMENSION(:, :), POINTER :: ksblock, sblock
99 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
100 TYPE(cell_type), POINTER :: cell
101 TYPE(dbcsr_iterator_type) :: iter
102 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_s
103 TYPE(dft_control_type), POINTER :: dft_control
104 TYPE(distribution_1d_type), POINTER :: local_particles
105 TYPE(ewald_environment_type), POINTER :: ewald_env
106 TYPE(ewald_pw_type), POINTER :: ewald_pw
107 TYPE(mp_para_env_type), POINTER :: para_env
109 DIMENSION(:), POINTER :: nl_iterator
110 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
111 POINTER :: n_list
112 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
113 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
114 TYPE(virial_type), POINTER :: virial
115 TYPE(xtb_atom_type), POINTER :: xtb_atom_a, xtb_atom_b, xtb_kind
116 TYPE(xtb_control_type), POINTER :: xtb_control
117
118 CALL timeset(routinen, handle)
119
120 CALL get_qs_env(qs_env, &
121 matrix_s_kp=matrix_s, &
122 qs_kind_set=qs_kind_set, &
123 particle_set=particle_set, &
124 cell=cell, &
125 dft_control=dft_control)
126
127 xtb_control => dft_control%qs_control%xtb_control
128
129 IF (dft_control%nimages /= 1) THEN
130 cpabort("No kpoints allowed in xTB response calculation")
131 END IF
132
133 CALL get_qs_env(qs_env, nkind=nkind, natom=natom)
134 nmat = 1
135 ALLOCATE (gchrg(natom, 5, nmat))
136 gchrg = 0._dp
137 ALLOCATE (gmcharge(natom, nmat))
138 gmcharge = 0._dp
139
140 ! short range contribution (gamma)
141 ! loop over all atom pairs (sab_xtbe)
142 kg = xtb_control%kg
143 NULLIFY (n_list)
144 CALL get_qs_env(qs_env=qs_env, sab_xtbe=n_list)
145 CALL neighbor_list_iterator_create(nl_iterator, n_list)
146 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
147 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
148 iatom=iatom, jatom=jatom, r=rij, cell=cellind)
149 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
150 CALL get_xtb_atom_param(xtb_atom_a, defined=defined, natorb=natorb_a)
151 IF (.NOT. defined .OR. natorb_a < 1) cycle
152 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
153 CALL get_xtb_atom_param(xtb_atom_b, defined=defined, natorb=natorb_b)
154 IF (.NOT. defined .OR. natorb_b < 1) cycle
155 ! atomic parameters
156 CALL get_xtb_atom_param(xtb_atom_a, eta=etaa, lmax=lmaxa, kappa=kappaa, rcut=rcuta)
157 CALL get_xtb_atom_param(xtb_atom_b, eta=etab, lmax=lmaxb, kappa=kappab, rcut=rcutb)
158 ! gamma matrix
159 ni = lmaxa + 1
160 nj = lmaxb + 1
161 ALLOCATE (gammab(ni, nj))
162 rcut = rcuta + rcutb
163 dr = sqrt(sum(rij(:)**2))
164 CALL gamma_rab_sr(gammab, dr, ni, kappaa, etaa, nj, kappab, etab, kg, rcut)
165 gchrg(iatom, 1:ni, 1) = gchrg(iatom, 1:ni, 1) + matmul(gammab, charges1(jatom, 1:nj))
166 IF (iatom /= jatom) THEN
167 gchrg(jatom, 1:nj, 1) = gchrg(jatom, 1:nj, 1) + matmul(charges1(iatom, 1:ni), gammab)
168 END IF
169 DEALLOCATE (gammab)
170 END DO
171 CALL neighbor_list_iterator_release(nl_iterator)
172
173 ! 1/R contribution
174
175 IF (xtb_control%coulomb_lr) THEN
176 do_ewald = xtb_control%do_ewald
177 IF (do_ewald) THEN
178 ! Ewald sum
179 NULLIFY (ewald_env, ewald_pw)
180 NULLIFY (virial)
181 CALL get_qs_env(qs_env=qs_env, &
182 ewald_env=ewald_env, ewald_pw=ewald_pw)
183 CALL get_cell(cell=cell, periodic=periodic, deth=deth)
184 CALL ewald_env_get(ewald_env, alpha=alpha, ewald_type=ewald_type)
185 CALL get_qs_env(qs_env=qs_env, sab_tbe=n_list)
186 CALL tb_ewald_overlap(gmcharge, mcharge1, alpha, n_list, virial, .false.)
187 SELECT CASE (ewald_type)
188 CASE DEFAULT
189 cpabort("Invalid Ewald type")
190 CASE (do_ewald_none)
191 cpabort("Not allowed with DFTB")
192 CASE (do_ewald_ewald)
193 cpabort("Standard Ewald not implemented in DFTB")
194 CASE (do_ewald_pme)
195 cpabort("PME not implemented in DFTB")
196 CASE (do_ewald_spme)
197 CALL tb_spme_evaluate(ewald_env, ewald_pw, particle_set, cell, &
198 gmcharge, mcharge1, .false., virial, .false.)
199 END SELECT
200 ELSE
201 ! direct sum
202 CALL get_qs_env(qs_env=qs_env, &
203 local_particles=local_particles)
204 DO ikind = 1, SIZE(local_particles%n_el)
205 DO ia = 1, local_particles%n_el(ikind)
206 iatom = local_particles%list(ikind)%array(ia)
207 DO jatom = 1, iatom - 1
208 rij = particle_set(iatom)%r - particle_set(jatom)%r
209 rij = pbc(rij, cell)
210 dr = sqrt(sum(rij(:)**2))
211 IF (dr > 1.e-6_dp) THEN
212 gmcharge(iatom, 1) = gmcharge(iatom, 1) + mcharge1(jatom)/dr
213 gmcharge(jatom, 1) = gmcharge(jatom, 1) + mcharge1(iatom)/dr
214 END IF
215 END DO
216 END DO
217 END DO
218 END IF
219 END IF
220
221 ! global sum of gamma*p arrays
222 CALL get_qs_env(qs_env=qs_env, para_env=para_env)
223 CALL para_env%sum(gmcharge(:, 1))
224 CALL para_env%sum(gchrg(:, :, 1))
225
226 IF (xtb_control%coulomb_lr) THEN
227 IF (do_ewald) THEN
228 ! add self charge interaction and background charge contribution
229 gmcharge(:, 1) = gmcharge(:, 1) - 2._dp*alpha*oorootpi*mcharge1(:)
230 IF (any(periodic(:) == 1)) THEN
231 gmcharge(:, 1) = gmcharge(:, 1) - pi/alpha**2/deth
232 END IF
233 END IF
234 END IF
235
236 CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set)
237 CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, kind_of=kind_of)
238
239 ! no k-points; all matrices have been transformed to periodic bsf
240 CALL dbcsr_iterator_start(iter, matrix_s(1, 1)%matrix)
241 DO WHILE (dbcsr_iterator_blocks_left(iter))
242 CALL dbcsr_iterator_next_block(iter, irow, icol, sblock)
243 ikind = kind_of(irow)
244 jkind = kind_of(icol)
245
246 ! atomic parameters
247 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom_a)
248 CALL get_qs_kind(qs_kind_set(jkind), xtb_parameter=xtb_atom_b)
249 CALL get_xtb_atom_param(xtb_atom_a, z=za, lao=laoa)
250 CALL get_xtb_atom_param(xtb_atom_b, z=zb, lao=laob)
251
252 ni = SIZE(sblock, 1)
253 nj = SIZE(sblock, 2)
254 ALLOCATE (gcij(ni, nj))
255 DO i = 1, ni
256 DO j = 1, nj
257 la = laoa(i) + 1
258 lb = laob(j) + 1
259 gcij(i, j) = gchrg(irow, la, 1) + gchrg(icol, lb, 1)
260 END DO
261 END DO
262 gmij = gmcharge(irow, 1) + gmcharge(icol, 1)
263 DO is = 1, SIZE(ks_matrix)
264 NULLIFY (ksblock)
265 CALL dbcsr_get_block_p(matrix=ks_matrix(is)%matrix, &
266 row=irow, col=icol, block=ksblock, found=found)
267 cpassert(found)
268 ksblock = ksblock - gcij*sblock
269 ksblock = ksblock - gmij*sblock
270 END DO
271 DEALLOCATE (gcij)
272 END DO
273 CALL dbcsr_iterator_stop(iter)
274
275 IF (xtb_control%tb3_interaction) THEN
276 CALL get_qs_env(qs_env, nkind=nkind)
277 ALLOCATE (xgamma(nkind))
278 DO ikind = 1, nkind
279 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_kind)
280 CALL get_xtb_atom_param(xtb_kind, xgamma=xgamma(ikind))
281 END DO
282 ! Diagonal 3rd order correction (DFTB3)
283 CALL dftb3_diagonal_hessian(qs_env, ks_matrix, mcharge, mcharge1, xgamma)
284 DEALLOCATE (xgamma)
285 END IF
286
287 IF (qs_env%qmmm .AND. qs_env%qmmm_periodic) THEN
288 cpabort("QMMM not available in xTB response calculations")
289 END IF
290
291 DEALLOCATE (gmcharge, gchrg)
292
293 CALL timestop(handle)
294
295 END SUBROUTINE xtb_coulomb_hessian
296
297! **************************************************************************************************
298!> \brief ...
299!> \param qs_env ...
300!> \param ks_matrix ...
301!> \param mcharge ...
302!> \param mcharge1 ...
303!> \param xgamma ...
304! **************************************************************************************************
305 SUBROUTINE dftb3_diagonal_hessian(qs_env, ks_matrix, mcharge, mcharge1, xgamma)
306
307 TYPE(qs_environment_type), POINTER :: qs_env
308 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: ks_matrix
309 REAL(dp), DIMENSION(:) :: mcharge, mcharge1, xgamma
310
311 CHARACTER(len=*), PARAMETER :: routinen = 'dftb3_diagonal_hessian'
312
313 INTEGER :: handle, icol, ikind, irow, is, jkind
314 INTEGER, ALLOCATABLE, DIMENSION(:) :: kind_of
315 LOGICAL :: found
316 REAL(kind=dp) :: gmij, ui, uj
317 REAL(kind=dp), DIMENSION(:, :), POINTER :: ksblock, sblock
318 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
319 TYPE(dbcsr_iterator_type) :: iter
320 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_s
321 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
322
323 CALL timeset(routinen, handle)
324
325 CALL get_qs_env(qs_env=qs_env, matrix_s_kp=matrix_s)
326 CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set, qs_kind_set=qs_kind_set)
327 CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, kind_of=kind_of)
328 ! no k-points; all matrices have been transformed to periodic bsf
329 CALL dbcsr_iterator_start(iter, matrix_s(1, 1)%matrix)
330 DO WHILE (dbcsr_iterator_blocks_left(iter))
331 CALL dbcsr_iterator_next_block(iter, irow, icol, sblock)
332 ikind = kind_of(irow)
333 ui = xgamma(ikind)
334 jkind = kind_of(icol)
335 uj = xgamma(jkind)
336 gmij = ui*mcharge(irow)*mcharge1(irow) + uj*mcharge(icol)*mcharge1(icol)
337 DO is = 1, SIZE(ks_matrix)
338 NULLIFY (ksblock)
339 CALL dbcsr_get_block_p(matrix=ks_matrix(is)%matrix, &
340 row=irow, col=icol, block=ksblock, found=found)
341 cpassert(found)
342 ksblock = ksblock + gmij*sblock
343 END DO
344 END DO
345 CALL dbcsr_iterator_stop(iter)
346
347 CALL timestop(handle)
348
349 END SUBROUTINE dftb3_diagonal_hessian
350
351END MODULE xtb_ehess
352
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...
subroutine, public dbcsr_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset)
...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_get_block_p(matrix, row, col, block, found, row_size, col_size)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
subroutine, public ewald_env_get(ewald_env, ewald_type, alpha, eps_pol, epsilon, gmax, ns_max, o_spline, group, para_env, poisson_section, precs, rcut, do_multipoles, max_multipole, do_ipol, max_ipol_iter, interaction_cutoffs, cell_hmat)
Purpose: Get the EWALD environment.
Calculation of Ewald contributions in DFTB.
subroutine, public tb_ewald_overlap(gmcharge, mcharge, alpha, n_list, virial, use_virial)
...
subroutine, public tb_spme_evaluate(ewald_env, ewald_pw, particle_set, box, gmcharge, mcharge, calculate_forces, virial, use_virial)
...
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
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_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zatom, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
Define the neighbor list data types and the corresponding functionality.
subroutine, public neighbor_list_iterator_create(iterator_set, nl, search, nthread)
Neighbor list iterator functions.
subroutine, public neighbor_list_iterator_release(iterator_set)
...
integer function, public neighbor_list_iterate(iterator_set, mepos)
...
subroutine, public get_iterator_info(iterator_set, mepos, ikind, jkind, nkind, ilist, nlist, inode, nnode, iatom, jatom, r, cell)
...
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: ...
Calculation of Coulomb Hessian contributions in xTB.
Definition xtb_ehess.F:12
subroutine, public xtb_coulomb_hessian(qs_env, ks_matrix, charges1, mcharge1, mcharge)
...
Definition xtb_ehess.F:77
Definition of the xTB parameter types.
Definition xtb_types.F:20
subroutine, public get_xtb_atom_param(xtb_parameter, symbol, aname, typ, defined, z, zeff, natorb, lmax, nao, lao, rcut, rcov, kx, eta, xgamma, alpha, zneff, nshell, nval, lval, kpoly, kappa, hen, zeta, xi, kappa0, alpg, occupation, electronegativity, chmax, en, kqat2, kcn, kq)
...
Definition xtb_types.F:199
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
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.