(git:06f838d)
Loading...
Searching...
No Matches
qs_tddfpt2_stda_utils.F
Go to the documentation of this file.
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!> \brief Simplified Tamm Dancoff approach (sTDA).
9! **************************************************************************************************
11
14 USE cell_types, ONLY: cell_type,&
15 pbc
18 USE cp_dbcsr_api, ONLY: &
22 dbcsr_type, dbcsr_type_antisymmetric, dbcsr_type_no_symmetry, dbcsr_type_symmetric
37 USE cp_fm_types, ONLY: cp_fm_create,&
60 USE kinds, ONLY: dp
61 USE mathconstants, ONLY: oorootpi
79 USE util, ONLY: get_limit
80 USE virial_types, ONLY: virial_type
81#include "./base/base_uses.f90"
82
83 IMPLICIT NONE
84
85 PRIVATE
86
87 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_tddfpt2_stda_utils'
88
91
92CONTAINS
93
94! **************************************************************************************************
95!> \brief Calculate sTDA matrices
96!> \param qs_env ...
97!> \param stda_kernel ...
98!> \param sub_env ...
99!> \param work ...
100!> \param tddfpt_control ...
101! **************************************************************************************************
102 SUBROUTINE stda_init_matrices(qs_env, stda_kernel, sub_env, work, tddfpt_control)
103
104 TYPE(qs_environment_type), POINTER :: qs_env
105 TYPE(stda_env_type) :: stda_kernel
106 TYPE(tddfpt_subgroup_env_type), INTENT(in) :: sub_env
107 TYPE(tddfpt_work_matrices) :: work
108 TYPE(tddfpt2_control_type), POINTER :: tddfpt_control
109
110 CHARACTER(len=*), PARAMETER :: routinen = 'stda_init_matrices'
111
112 INTEGER :: handle
113 LOGICAL :: do_coulomb
114 TYPE(cell_type), POINTER :: cell, cell_ref
115 TYPE(ewald_environment_type), POINTER :: ewald_env
116 TYPE(ewald_pw_type), POINTER :: ewald_pw
117 TYPE(section_vals_type), POINTER :: ewald_section, poisson_section, &
118 print_section
119
120 CALL timeset(routinen, handle)
121
122 do_coulomb = .NOT. tddfpt_control%rks_triplets
123 IF (do_coulomb) THEN
124 ! calculate exchange gamma matrix
125 CALL setup_gamma(qs_env, stda_kernel, sub_env, work%gamma_exchange)
126 END IF
127
128 ! calculate S_half and Lowdin MO coefficients
129 CALL get_lowdin_mo_coefficients(qs_env, sub_env, work)
130
131 ! initialize Ewald for sTDA
132 IF (tddfpt_control%stda_control%do_ewald) THEN
133 NULLIFY (ewald_env, ewald_pw)
134 ALLOCATE (ewald_env)
135 CALL ewald_env_create(ewald_env, sub_env%para_env)
136 poisson_section => section_vals_get_subs_vals(qs_env%input, "DFT%POISSON")
137 CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
138 ewald_section => section_vals_get_subs_vals(poisson_section, "EWALD")
139 print_section => section_vals_get_subs_vals(qs_env%input, "PRINT%GRID_INFORMATION")
140 CALL get_qs_env(qs_env, cell=cell, cell_ref=cell_ref)
141 CALL read_ewald_section_tb(ewald_env, ewald_section, cell_ref%hmat, &
142 cell_periodic=cell%perd)
143 ALLOCATE (ewald_pw)
144 CALL ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section=print_section)
145 work%ewald_env => ewald_env
146 work%ewald_pw => ewald_pw
147 END IF
148
149 CALL timestop(handle)
150
151 END SUBROUTINE stda_init_matrices
152! **************************************************************************************************
153!> \brief Calculate sTDA exchange-type contributions
154!> \param qs_env ...
155!> \param stda_env ...
156!> \param sub_env ...
157!> \param gamma_matrix sTDA exchange-type contributions
158!> \param ndim ...
159!> \note Note the specific sTDA notation exchange-type integrals (ia|jb) refer to Coulomb interaction
160! **************************************************************************************************
161 SUBROUTINE setup_gamma(qs_env, stda_env, sub_env, gamma_matrix, ndim)
162
163 TYPE(qs_environment_type), POINTER :: qs_env
164 TYPE(stda_env_type) :: stda_env
165 TYPE(tddfpt_subgroup_env_type), INTENT(in) :: sub_env
166 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: gamma_matrix
167 INTEGER, INTENT(IN), OPTIONAL :: ndim
168
169 CHARACTER(len=*), PARAMETER :: routinen = 'setup_gamma'
170 REAL(kind=dp), PARAMETER :: rsmooth = 1.0_dp
171
172 INTEGER :: handle, i, iatom, icol, ikind, imat, &
173 irow, jatom, jkind, natom, nmat
174 INTEGER, DIMENSION(:), POINTER :: row_blk_sizes
175 LOGICAL :: found
176 REAL(kind=dp) :: dfcut, dgb, dr, eta, fcut, r, rcut, &
177 rcuta, rcutb, x
178 REAL(kind=dp), DIMENSION(3) :: rij
179 REAL(kind=dp), DIMENSION(:, :), POINTER :: dgblock, gblock
180 TYPE(dbcsr_distribution_type), POINTER :: dbcsr_dist
182 DIMENSION(:), POINTER :: nl_iterator
183 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
184 POINTER :: n_list
185
186 CALL timeset(routinen, handle)
187
188 CALL get_qs_env(qs_env=qs_env, natom=natom)
189 dbcsr_dist => sub_env%dbcsr_dist
190 ! Using the overlap list here can have a considerable effect on the number of
191 ! terms calculated. This makes gamma also dependent on EPS_DEFAULT -> Overlap
192 n_list => sub_env%sab_orb
193
194 IF (PRESENT(ndim)) THEN
195 nmat = ndim
196 ELSE
197 nmat = 1
198 END IF
199 cpassert(nmat == 1 .OR. nmat == 4)
200 cpassert(.NOT. ASSOCIATED(gamma_matrix))
201 CALL dbcsr_allocate_matrix_set(gamma_matrix, nmat)
202
203 ALLOCATE (row_blk_sizes(natom))
204 row_blk_sizes(1:natom) = 1
205 DO imat = 1, nmat
206 ALLOCATE (gamma_matrix(imat)%matrix)
207 END DO
208
209 CALL dbcsr_create(gamma_matrix(1)%matrix, name="gamma", dist=dbcsr_dist, &
210 matrix_type=dbcsr_type_symmetric, row_blk_size=row_blk_sizes, &
211 col_blk_size=row_blk_sizes)
212 DO imat = 2, nmat
213 CALL dbcsr_create(gamma_matrix(imat)%matrix, name="dgamma", dist=dbcsr_dist, &
214 matrix_type=dbcsr_type_antisymmetric, row_blk_size=row_blk_sizes, &
215 col_blk_size=row_blk_sizes)
216 END DO
217
218 DEALLOCATE (row_blk_sizes)
219
220 ! setup the matrices using the neighbor list
221 DO imat = 1, nmat
222 CALL cp_dbcsr_alloc_block_from_nbl(gamma_matrix(imat)%matrix, n_list)
223 CALL dbcsr_set(gamma_matrix(imat)%matrix, 0.0_dp)
224 END DO
225
226 NULLIFY (nl_iterator)
227 CALL neighbor_list_iterator_create(nl_iterator, n_list)
228 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
229 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
230 iatom=iatom, jatom=jatom, r=rij)
231
232 dr = sqrt(sum(rij(:)**2)) ! interatomic distance
233
234 eta = (stda_env%kind_param_set(ikind)%kind_param%hardness_param + &
235 stda_env%kind_param_set(jkind)%kind_param%hardness_param)/2.0_dp
236
237 icol = max(iatom, jatom)
238 irow = min(iatom, jatom)
239
240 NULLIFY (gblock)
241 CALL dbcsr_get_block_p(matrix=gamma_matrix(1)%matrix, &
242 row=irow, col=icol, block=gblock, found=found)
243 cpassert(found)
244
245 ! get rcuta and rcutb
246 rcuta = stda_env%kind_param_set(ikind)%kind_param%rcut
247 rcutb = stda_env%kind_param_set(jkind)%kind_param%rcut
248 rcut = rcuta + rcutb
249
250 !> Computes the short-range gamma parameter from
251 !> Nataga-Mishimoto-Ohno-Klopman formula equivalently as it is done for xTB
252 IF (dr < 1.e-6) THEN
253 ! on site terms
254 gblock(:, :) = gblock(:, :) + eta
255 ELSEIF (dr > rcut) THEN
256 ! do nothing
257 ELSE
258 IF (dr < rcut - rsmooth) THEN
259 fcut = 1.0_dp
260 ELSE
261 r = dr - (rcut - rsmooth)
262 x = r/rsmooth
263 fcut = -6._dp*x**5 + 15._dp*x**4 - 10._dp*x**3 + 1._dp
264 END IF
265 gblock(:, :) = gblock(:, :) + &
266 fcut*(1._dp/(dr**(stda_env%alpha_param) + eta**(-stda_env%alpha_param))) &
267 **(1._dp/stda_env%alpha_param) - fcut/dr
268 END IF
269
270 IF (nmat > 1) THEN
271 !> Computes the short-range gamma parameter from
272 !> Nataga-Mishimoto-Ohno-Klopman formula equivalently as it is done for xTB
273 !> Derivatives
274 IF (dr < 1.e-6 .OR. dr > rcut) THEN
275 ! on site terms or beyond cutoff
276 dgb = 0.0_dp
277 ELSE
278 IF (dr < rcut - rsmooth) THEN
279 fcut = 1.0_dp
280 dfcut = 0.0_dp
281 ELSE
282 r = dr - (rcut - rsmooth)
283 x = r/rsmooth
284 fcut = -6._dp*x**5 + 15._dp*x**4 - 10._dp*x**3 + 1._dp
285 dfcut = -30._dp*x**4 + 60._dp*x**3 - 30._dp*x**2
286 dfcut = dfcut/rsmooth
287 END IF
288 dgb = dfcut*(1._dp/(dr**(stda_env%alpha_param) + eta**(-stda_env%alpha_param))) &
289 **(1._dp/stda_env%alpha_param)
290 dgb = dgb - dfcut/dr + fcut/dr**2
291 dgb = dgb - fcut*(1._dp/(dr**(stda_env%alpha_param) + eta**(-stda_env%alpha_param))) &
292 **(1._dp/stda_env%alpha_param + 1._dp)*dr**(stda_env%alpha_param - 1._dp)
293 END IF
294 DO imat = 2, nmat
295 NULLIFY (dgblock)
296 CALL dbcsr_get_block_p(matrix=gamma_matrix(imat)%matrix, &
297 row=irow, col=icol, block=dgblock, found=found)
298 IF (found) THEN
299 IF (dr > 1.e-6) THEN
300 i = imat - 1
301 IF (irow == iatom) THEN
302 dgblock(:, :) = dgblock(:, :) + dgb*rij(i)/dr
303 ELSE
304 dgblock(:, :) = dgblock(:, :) - dgb*rij(i)/dr
305 END IF
306 END IF
307 END IF
308 END DO
309 END IF
310
311 END DO
312
313 CALL neighbor_list_iterator_release(nl_iterator)
314
315 DO imat = 1, nmat
316 CALL dbcsr_finalize(gamma_matrix(imat)%matrix)
317 END DO
318
319 CALL timestop(handle)
320
321 END SUBROUTINE setup_gamma
322
323! **************************************************************************************************
324!> \brief Calculate Lowdin MO coefficients
325!> \param qs_env ...
326!> \param sub_env ...
327!> \param work ...
328! **************************************************************************************************
329 SUBROUTINE get_lowdin_mo_coefficients(qs_env, sub_env, work)
330
331 TYPE(qs_environment_type), POINTER :: qs_env
332 TYPE(tddfpt_subgroup_env_type), INTENT(in) :: sub_env
333 TYPE(tddfpt_work_matrices) :: work
334
335 CHARACTER(len=*), PARAMETER :: routinen = 'get_lowdin_mo_coefficients'
336
337 INTEGER :: handle, i, iounit, ispin, j, &
338 max_iter_lanczos, nactive, ndep, nsgf, &
339 nspins, order_lanczos
340 LOGICAL :: converged
341 REAL(kind=dp) :: eps_lanczos, sij, threshold
342 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: slam
343 REAL(kind=dp), CONTIGUOUS, DIMENSION(:, :), &
344 POINTER :: local_data
345 TYPE(cp_fm_struct_type), POINTER :: fmstruct
346 TYPE(cp_fm_type) :: fm_s_half, fm_work1
347 TYPE(cp_logger_type), POINTER :: logger
348 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrixkp_s
349 TYPE(dbcsr_type) :: sm_hinv
350 TYPE(dbcsr_type), POINTER :: sm_h, sm_s
351 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
352 TYPE(scf_control_type), POINTER :: scf_control
353
354 CALL timeset(routinen, handle)
355
356 NULLIFY (logger) !get output_unit
357 logger => cp_get_default_logger()
358 iounit = cp_logger_get_default_io_unit(logger)
359
360 ! Calculate S^1/2 matrix
361 IF (iounit > 0) THEN
362 WRITE (iounit, "(1X,A)") "", &
363 "-------------------------------------------------------------------------------", &
364 "- Create Matrix SQRT(S) -", &
365 "-------------------------------------------------------------------------------"
366 END IF
367
368 IF (sub_env%is_split) THEN
369 cpabort('SPLIT')
370 ELSE
371 CALL get_qs_env(qs_env=qs_env, matrix_s_kp=matrixkp_s)
372 cpassert(ASSOCIATED(matrixkp_s))
373 cpwarn_if(SIZE(matrixkp_s, 2) > 1, "not implemented for k-points.")
374 sm_s => matrixkp_s(1, 1)%matrix
375 END IF
376 sm_h => work%shalf
377
378 CALL dbcsr_create(sm_hinv, template=sm_s)
379 CALL dbcsr_add_on_diag(sm_h, 1.0_dp)
380 threshold = 1.0e-8_dp
381 order_lanczos = 3
382 eps_lanczos = 1.0e-4_dp
383 max_iter_lanczos = 40
384 CALL matrix_sqrt_newton_schulz(sm_h, sm_hinv, sm_s, &
385 threshold, order_lanczos, eps_lanczos, max_iter_lanczos, &
386 converged=converged)
387 CALL dbcsr_release(sm_hinv)
388 !
389 NULLIFY (qs_kind_set)
390 CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set)
391 ! Get the total number of contracted spherical Gaussian basis functions
392 CALL get_qs_kind_set(qs_kind_set, nsgf=nsgf)
393 !
394 IF (.NOT. converged) THEN
395 IF (iounit > 0) THEN
396 WRITE (iounit, "(T3,A)") "STDA| Newton-Schulz iteration did not converge"
397 WRITE (iounit, "(T3,A)") "STDA| Calculate SQRT(S) from diagonalization"
398 END IF
399 CALL get_qs_env(qs_env=qs_env, scf_control=scf_control)
400 ! Provide full size work matrices
401 CALL cp_fm_struct_create(fmstruct=fmstruct, &
402 para_env=sub_env%para_env, &
403 context=sub_env%blacs_env, &
404 nrow_global=nsgf, &
405 ncol_global=nsgf)
406 CALL cp_fm_create(matrix=fm_s_half, matrix_struct=fmstruct, name="S^(1/2) MATRIX")
407 CALL cp_fm_create(matrix=fm_work1, matrix_struct=fmstruct, name="TMP MATRIX")
408 CALL cp_fm_struct_release(fmstruct=fmstruct)
409 CALL copy_dbcsr_to_fm(sm_s, fm_s_half)
410 CALL cp_fm_power(fm_s_half, fm_work1, 0.5_dp, scf_control%eps_eigval, ndep)
411 IF (ndep /= 0) &
412 CALL cp_warn(__location__, &
413 "Overlap matrix exhibits linear dependencies. At least some "// &
414 "eigenvalues have been quenched.")
415 CALL copy_fm_to_dbcsr(fm_s_half, sm_h)
416 CALL cp_fm_release(fm_s_half)
417 CALL cp_fm_release(fm_work1)
418 IF (iounit > 0) WRITE (iounit, *)
419 END IF
420
421 nspins = SIZE(sub_env%mos_occ)
422
423 DO ispin = 1, nspins
424 CALL cp_fm_get_info(work%ctransformed(ispin), ncol_global=nactive)
425 CALL cp_dbcsr_sm_fm_multiply(work%shalf, sub_env%mos_active(ispin), &
426 work%ctransformed(ispin), nactive, alpha=1.0_dp, beta=0.0_dp)
427 END DO
428
429 ! for Lowdin forces
430 CALL cp_fm_create(matrix=fm_work1, matrix_struct=work%S_eigenvectors%matrix_struct, name="TMP MATRIX")
431 CALL copy_dbcsr_to_fm(sm_s, fm_work1)
432 CALL choose_eigv_solver(fm_work1, work%S_eigenvectors, work%S_eigenvalues)
433 CALL cp_fm_release(fm_work1)
434 !
435 ALLOCATE (slam(nsgf, 1))
436 DO i = 1, nsgf
437 IF (work%S_eigenvalues(i) > 0._dp) THEN
438 slam(i, 1) = sqrt(work%S_eigenvalues(i))
439 ELSE
440 cpabort("S matrix not positive definit")
441 END IF
442 END DO
443 DO i = 1, nsgf
444 CALL cp_fm_set_submatrix(work%slambda, slam, 1, i, nsgf, 1, 1.0_dp, 0.0_dp)
445 END DO
446 DO i = 1, nsgf
447 CALL cp_fm_set_submatrix(work%slambda, slam, i, 1, 1, nsgf, 1.0_dp, 1.0_dp, .true.)
448 END DO
449 CALL cp_fm_get_info(work%slambda, local_data=local_data)
450 DO i = 1, SIZE(local_data, 2)
451 DO j = 1, SIZE(local_data, 1)
452 sij = local_data(j, i)
453 IF (sij > 0.0_dp) sij = 1.0_dp/sij
454 local_data(j, i) = sij
455 END DO
456 END DO
457 DEALLOCATE (slam)
458
459 CALL timestop(handle)
460
461 END SUBROUTINE get_lowdin_mo_coefficients
462
463! **************************************************************************************************
464!> \brief Calculate Lowdin transformed Davidson trial vector X
465!> shalf (dbcsr), xvec, xt (fm) are defined in the same sub_env
466!> \param shalf ...
467!> \param xvec ...
468!> \param xt ...
469! **************************************************************************************************
470 SUBROUTINE get_lowdin_x(shalf, xvec, xt)
471
472 TYPE(dbcsr_type), INTENT(IN) :: shalf
473 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: xvec
474 TYPE(cp_fm_type), DIMENSION(:), INTENT(INOUT) :: xt
475
476 CHARACTER(len=*), PARAMETER :: routinen = 'get_lowdin_x'
477
478 INTEGER :: handle, ispin, nactive, nspins
479
480 CALL timeset(routinen, handle)
481
482 nspins = SIZE(xvec)
483
484 ! Build Lowdin transformed tilde(X)= S^1/2 X for each spin
485 DO ispin = 1, nspins
486 CALL cp_fm_get_info(xt(ispin), ncol_global=nactive)
487 CALL cp_dbcsr_sm_fm_multiply(shalf, xvec(ispin), &
488 xt(ispin), nactive, alpha=1.0_dp, beta=0.0_dp)
489 END DO
490
491 CALL timestop(handle)
492
493 END SUBROUTINE get_lowdin_x
494
495! **************************************************************************************************
496!> \brief ...Calculate the sTDA kernel contribution by contracting the Lowdin MO coefficients --
497!> transition charges with the Coulomb-type or exchange-type integrals
498!> \param qs_env ...
499!> \param stda_control ...
500!> \param stda_env ...
501!> \param sub_env ...
502!> \param work ...
503!> \param is_rks_triplets ...
504!> \param X ...
505!> \param res ... vector AX with A being the sTDA matrix and X the Davidson trial vector of the
506!> eigenvalue problem A*X = omega*X
507! **************************************************************************************************
508 SUBROUTINE stda_calculate_kernel(qs_env, stda_control, stda_env, sub_env, &
509 work, is_rks_triplets, X, res)
510
511 TYPE(qs_environment_type), POINTER :: qs_env
512 TYPE(stda_control_type) :: stda_control
513 TYPE(stda_env_type) :: stda_env
514 TYPE(tddfpt_subgroup_env_type) :: sub_env
515 TYPE(tddfpt_work_matrices) :: work
516 LOGICAL, INTENT(IN) :: is_rks_triplets
517 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: x
518 TYPE(cp_fm_type), DIMENSION(:), INTENT(INOUT) :: res
519
520 CHARACTER(len=*), PARAMETER :: routinen = 'stda_calculate_kernel'
521
522 INTEGER :: ewald_type, handle, ia, iatom, ikind, &
523 is, ispin, jatom, jkind, jspin, natom, &
524 nsgf, nspins
525 INTEGER, ALLOCATABLE, DIMENSION(:) :: first_sgf, kind_of, last_sgf
526 INTEGER, DIMENSION(2) :: nactive, nlim
527 LOGICAL :: calculate_forces, do_coulomb, do_ewald, &
528 do_exchange, use_virial
529 REAL(kind=dp) :: alpha, bp, dr, eta, gabr, hfx, rbeta, &
530 spinfac
531 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: tcharge, tv
532 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: gtcharge
533 REAL(kind=dp), DIMENSION(3) :: rij
534 REAL(kind=dp), DIMENSION(:, :), POINTER :: gab, pblock
535 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
536 TYPE(cell_type), POINTER :: cell
537 TYPE(cp_fm_struct_type), POINTER :: fmstruct, fmstructjspin
538 TYPE(cp_fm_type) :: cvec, cvecjspin
539 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:) :: xtransformed
540 TYPE(cp_fm_type), POINTER :: ct, ctjspin
541 TYPE(dbcsr_iterator_type) :: iter
542 TYPE(dbcsr_type) :: pdens
543 TYPE(dbcsr_type), POINTER :: tempmat
544 TYPE(ewald_environment_type), POINTER :: ewald_env
545 TYPE(ewald_pw_type), POINTER :: ewald_pw
546 TYPE(mp_para_env_type), POINTER :: para_env
547 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
548 POINTER :: n_list
549 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
550 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
551 TYPE(virial_type), POINTER :: virial
552
553 CALL timeset(routinen, handle)
554
555 nactive(:) = stda_env%nactive(:)
556 nspins = SIZE(x)
557 spinfac = 2.0_dp
558 IF (nspins == 2) spinfac = 1.0_dp
559
560 IF (nspins == 1 .AND. is_rks_triplets) THEN
561 do_coulomb = .false.
562 ELSE
563 do_coulomb = .true.
564 END IF
565 do_ewald = stda_control%do_ewald
566 do_exchange = stda_control%do_exchange
567
568 para_env => sub_env%para_env
569
570 CALL get_qs_env(qs_env, natom=natom, cell=cell, &
571 particle_set=particle_set, qs_kind_set=qs_kind_set)
572 ALLOCATE (first_sgf(natom))
573 ALLOCATE (last_sgf(natom))
574 CALL get_particle_set(particle_set, qs_kind_set, first_sgf=first_sgf, last_sgf=last_sgf)
575
576 ! calculate Loewdin transformed Davidson trial vector tilde(X)=S^1/2*X
577 ! and tilde(tilde(X))=S^1/2_A*tilde(X)_A
578 ALLOCATE (xtransformed(nspins))
579 DO ispin = 1, nspins
580 NULLIFY (fmstruct)
581 ct => work%ctransformed(ispin)
582 CALL cp_fm_get_info(ct, matrix_struct=fmstruct)
583 CALL cp_fm_create(matrix=xtransformed(ispin), matrix_struct=fmstruct, name="XTRANSFORMED")
584 END DO
585 CALL get_lowdin_x(work%shalf, x, xtransformed)
586
587 ALLOCATE (tcharge(natom), gtcharge(natom, 1))
588
589 DO ispin = 1, nspins
590 CALL cp_fm_set_all(res(ispin), 0.0_dp)
591 END DO
592
593 DO ispin = 1, nspins
594 ct => work%ctransformed(ispin)
595 CALL cp_fm_get_info(ct, matrix_struct=fmstruct, nrow_global=nsgf)
596 ALLOCATE (tv(nsgf))
597 CALL cp_fm_create(cvec, fmstruct)
598 !
599 ! *** Coulomb contribution
600 !
601 IF (do_coulomb) THEN
602 tcharge(:) = 0.0_dp
603 DO jspin = 1, nspins
604 ctjspin => work%ctransformed(jspin)
605 CALL cp_fm_get_info(ctjspin, matrix_struct=fmstructjspin)
606 CALL cp_fm_get_info(ctjspin, matrix_struct=fmstructjspin, nrow_global=nsgf)
607 CALL cp_fm_create(cvecjspin, fmstructjspin)
608 ! CV(mu,j) = CT(mu,j)*XT(mu,j)
609 CALL cp_fm_schur_product(ctjspin, xtransformed(jspin), cvecjspin)
610 ! TV(mu) = SUM_j CV(mu,j)
611 CALL cp_fm_vectorssum(cvecjspin, tv, "R")
612 ! contract charges
613 ! TC(a) = SUM_(mu of a) TV(mu)
614 DO ia = 1, natom
615 DO is = first_sgf(ia), last_sgf(ia)
616 tcharge(ia) = tcharge(ia) + tv(is)
617 END DO
618 END DO
619 CALL cp_fm_release(cvecjspin)
620 END DO !jspin
621 ! Apply tcharge*gab -> gtcharge
622 ! gT(b) = SUM_a g(a,b)*TC(a)
623 ! gab = work%gamma_exchange(1)%matrix
624 gtcharge = 0.0_dp
625 ! short range contribution
626 tempmat => work%gamma_exchange(1)%matrix
627 CALL dbcsr_iterator_start(iter, tempmat)
628 DO WHILE (dbcsr_iterator_blocks_left(iter))
629 CALL dbcsr_iterator_next_block(iter, iatom, jatom, gab)
630 gtcharge(iatom, 1) = gtcharge(iatom, 1) + gab(1, 1)*tcharge(jatom)
631 IF (iatom /= jatom) THEN
632 gtcharge(jatom, 1) = gtcharge(jatom, 1) + gab(1, 1)*tcharge(iatom)
633 END IF
634 END DO
635 CALL dbcsr_iterator_stop(iter)
636 ! Ewald long range contribution
637 IF (do_ewald) THEN
638 ewald_env => work%ewald_env
639 ewald_pw => work%ewald_pw
640 CALL ewald_env_get(ewald_env, alpha=alpha, ewald_type=ewald_type)
641 CALL get_qs_env(qs_env=qs_env, virial=virial)
642 use_virial = .false.
643 calculate_forces = .false.
644 n_list => sub_env%sab_orb
645 CALL tb_ewald_overlap(gtcharge, tcharge, alpha, n_list, virial, use_virial)
646 CALL tb_spme_evaluate(ewald_env, ewald_pw, particle_set, cell, &
647 gtcharge, tcharge, calculate_forces, virial, use_virial)
648 ! add self charge interaction contribution
649 IF (para_env%is_source()) THEN
650 gtcharge(:, 1) = gtcharge(:, 1) - 2._dp*alpha*oorootpi*tcharge(:)
651 END IF
652 ELSE
653 nlim = get_limit(natom, para_env%num_pe, para_env%mepos)
654 DO iatom = nlim(1), nlim(2)
655 DO jatom = 1, iatom - 1
656 rij = particle_set(iatom)%r - particle_set(jatom)%r
657 rij = pbc(rij, cell)
658 dr = sqrt(sum(rij(:)**2))
659 IF (dr > 1.e-6_dp) THEN
660 gtcharge(iatom, 1) = gtcharge(iatom, 1) + tcharge(jatom)/dr
661 gtcharge(jatom, 1) = gtcharge(jatom, 1) + tcharge(iatom)/dr
662 END IF
663 END DO
664 END DO
665 END IF
666 CALL para_env%sum(gtcharge)
667 ! expand charges
668 ! TV(mu) = TC(a of mu)
669 tv(1:nsgf) = 0.0_dp
670 DO ia = 1, natom
671 DO is = first_sgf(ia), last_sgf(ia)
672 tv(is) = gtcharge(ia, 1)
673 END DO
674 END DO
675 ! CV(mu,i) = TV(mu)*CV(mu,i)
676 ct => work%ctransformed(ispin)
677 CALL cp_fm_to_fm(ct, cvec)
678 CALL cp_fm_row_scale(cvec, tv)
679 ! rho(nu,i) = rho(nu,i) + Shalf(nu,mu)*CV(mu,i)
680 CALL cp_dbcsr_sm_fm_multiply(work%shalf, cvec, res(ispin), nactive(ispin), spinfac, 1.0_dp)
681 END IF
682 !
683 ! *** Exchange contribution
684 !
685 IF (do_exchange) THEN ! option to explicitly switch off exchange
686 ! (exchange contributes also if hfx_fraction=0)
687 CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set)
688 CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, kind_of=kind_of)
689 !
690 tempmat => work%shalf
691 CALL dbcsr_create(pdens, template=tempmat, matrix_type=dbcsr_type_no_symmetry)
692 ! P(nu,mu) = SUM_j XT(nu,j)*CT(mu,j)
693 ct => work%ctransformed(ispin)
694 CALL dbcsr_set(pdens, 0.0_dp)
695 CALL cp_dbcsr_plus_fm_fm_t(pdens, xtransformed(ispin), ct, nactive(ispin), &
696 1.0_dp, keep_sparsity=.false.)
697 CALL dbcsr_filter(pdens, stda_env%eps_td_filter)
698 ! Apply PP*gab -> PP; gab = gamma_coulomb
699 ! P(nu,mu) = P(nu,mu)*g(a of nu,b of mu)
700 bp = stda_env%beta_param
701 hfx = stda_env%hfx_fraction
702 CALL dbcsr_iterator_start(iter, pdens)
703 DO WHILE (dbcsr_iterator_blocks_left(iter))
704 CALL dbcsr_iterator_next_block(iter, iatom, jatom, pblock)
705 rij = particle_set(iatom)%r - particle_set(jatom)%r
706 rij = pbc(rij, cell)
707 dr = sqrt(sum(rij(:)**2))
708 ikind = kind_of(iatom)
709 jkind = kind_of(jatom)
710 eta = (stda_env%kind_param_set(ikind)%kind_param%hardness_param + &
711 stda_env%kind_param_set(jkind)%kind_param%hardness_param)/2.0_dp
712 rbeta = dr**bp
713 IF (hfx > 0.0_dp) THEN
714 gabr = (1._dp/(rbeta + (hfx*eta)**(-bp)))**(1._dp/bp)
715 ELSE
716 IF (dr < 1.e-6) THEN
717 gabr = 0.0_dp
718 ELSE
719 gabr = 1._dp/dr
720 END IF
721 END IF
722 pblock = gabr*pblock
723 END DO
724 CALL dbcsr_iterator_stop(iter)
725 ! CV(mu,i) = P(nu,mu)*CT(mu,i)
726 CALL cp_dbcsr_sm_fm_multiply(pdens, ct, cvec, nactive(ispin), 1.0_dp, 0.0_dp)
727 ! rho(nu,i) = rho(nu,i) + ShalfP(nu,mu)*CV(mu,i)
728 CALL cp_dbcsr_sm_fm_multiply(work%shalf, cvec, res(ispin), nactive(ispin), -1.0_dp, 1.0_dp)
729 !
730 CALL dbcsr_release(pdens)
731 DEALLOCATE (kind_of)
732 END IF
733 !
734 CALL cp_fm_release(cvec)
735 DEALLOCATE (tv)
736 END DO
737
738 CALL cp_fm_release(xtransformed)
739 DEALLOCATE (tcharge, gtcharge)
740 DEALLOCATE (first_sgf, last_sgf)
741
742 CALL timestop(handle)
743
744 END SUBROUTINE stda_calculate_kernel
745
746! **************************************************************************************************
747
748END MODULE qs_tddfpt2_stda_utils
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
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
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_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset, transposed)
...
subroutine, public dbcsr_filter(matrix, eps)
...
subroutine, public dbcsr_finalize(matrix)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_release(matrix)
...
subroutine, public dbcsr_add_on_diag(matrix, alpha)
Adds the given scalar to the diagonal of the matrix. Reserves any missing diagonal blocks.
DBCSR operations in CP2K.
subroutine, public cp_dbcsr_sm_fm_multiply(matrix, fm_in, fm_out, ncol, alpha, beta)
multiply a dbcsr with a fm matrix
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
subroutine, public cp_dbcsr_plus_fm_fm_t(sparse_matrix, matrix_v, matrix_g, ncol, alpha, keep_sparsity, symmetry_mode)
performs the multiplication sparse_matrix+dense_mat*dens_mat^T if matrix_g is not explicitly given,...
subroutine, public copy_fm_to_dbcsr(fm, matrix, keep_sparsity)
Copy a BLACS matrix to a dbcsr matrix.
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_row_scale(matrixa, scaling)
scales row i of matrix a with scaling(i)
subroutine, public cp_fm_schur_product(matrix_a, matrix_b, matrix_c)
computes the schur product of two matrices c_ij = a_ij * b_ij
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
Definition cp_fm_diag.F:17
subroutine, public cp_fm_power(matrix, work, exponent, threshold, n_dependent, verbose, eigvals)
...
subroutine, public choose_eigv_solver(matrix, eigenvectors, eigenvalues, info)
Choose the Eigensolver depending on which library is available ELPA seems to be unstable for small sy...
Definition cp_fm_diag.F:245
represent the structure of a full matrix
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
subroutine, public cp_fm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
subroutine, public cp_fm_vectorssum(matrix, sum_array, dir)
summing up all the elements along the matrix's i-th index or
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
subroutine, public cp_fm_set_submatrix(fm, new_values, start_row, start_col, n_rows, n_cols, alpha, beta, transpose)
sets a submatrix of a full matrix fm(start_row:start_row+n_rows,start_col:start_col+n_cols) = alpha*o...
subroutine, public cp_fm_set_all(matrix, alpha, beta)
set all elements of a matrix to the same value, and optionally the diagonal to a different one
various routines to log and control the output. The idea is that decisions about where to log should ...
integer function, public cp_logger_get_default_io_unit(logger)
returns the unit nr for the ionode (-1 on all other processors) skips as well checks if the procs cal...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
subroutine, public ewald_env_set(ewald_env, ewald_type, alpha, epsilon, eps_pol, gmax, ns_max, precs, o_spline, para_env, poisson_section, interaction_cutoffs, cell_hmat)
Purpose: Set the EWALD environment.
subroutine, public ewald_env_create(ewald_env, para_env)
allocates and intitializes a ewald_env
subroutine, public read_ewald_section_tb(ewald_env, ewald_section, hmat, silent, pset, cell_periodic)
Purpose: read the EWALD section for TB methods.
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)
...
subroutine, public ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section)
creates the structure ewald_pw_type
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
Routines useful for iterative matrix calculations.
subroutine, public matrix_sqrt_newton_schulz(matrix_sqrt, matrix_sqrt_inv, matrix, threshold, order, eps_lanczos, max_iter_lanczos, symmetrize, converged, iounit)
compute the sqrt of a matrix via the sign function and the corresponding Newton-Schulz iterations the...
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
Interface to the message passing library MPI.
Define methods related to particle_type.
subroutine, public get_particle_set(particle_set, qs_kind_set, first_sgf, last_sgf, nsgf, nmao, basis, ncgf)
Get the components of a particle set.
Define the data structure for the particle information.
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, mimic, 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, sab_cneo, 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, xcint_weights, 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, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_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, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind_set(qs_kind_set, all_potential_present, tnadd_potential_present, gth_potential_present, sgp_potential_present, paw_atom_present, dft_plus_u_atom_present, maxcgf, maxsgf, maxco, maxco_proj, maxgtops, maxlgto, maxlprj, maxnset, maxsgf_set, ncgf, npgf, nset, nsgf, nshell, maxpol, maxlppl, maxlppnl, maxppnl, nelectron, maxder, max_ngrid_rad, max_sph_harm, maxg_iso_not0, lmax_rho0, basis_rcut, basis_type, total_zeff_corr, npgf_seg, cneo_potential_present, nkind_q, natom_q)
Get attributes of an atomic kind set.
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)
...
Simplified Tamm Dancoff approach (sTDA).
Simplified Tamm Dancoff approach (sTDA).
subroutine, public get_lowdin_mo_coefficients(qs_env, sub_env, work)
Calculate Lowdin MO coefficients.
subroutine, public get_lowdin_x(shalf, xvec, xt)
Calculate Lowdin transformed Davidson trial vector X shalf (dbcsr), xvec, xt (fm) are defined in the ...
subroutine, public setup_gamma(qs_env, stda_env, sub_env, gamma_matrix, ndim)
Calculate sTDA exchange-type contributions.
subroutine, public stda_init_matrices(qs_env, stda_kernel, sub_env, work, tddfpt_control)
Calculate sTDA matrices.
subroutine, public stda_calculate_kernel(qs_env, stda_control, stda_env, sub_env, work, is_rks_triplets, x, res)
...Calculate the sTDA kernel contribution by contracting the Lowdin MO coefficients – transition char...
parameters that control an scf iteration
All kind of helpful little routines.
Definition util.F:14
pure integer function, dimension(2), public get_limit(m, n, me)
divide m entries into n parts, return size of part me
Definition util.F:333
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:60
keeps the information about the structure of a full matrix
represent a full matrix
type of a logger, at the moment it contains just a print level starting at which level it should be l...
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.
Set of temporary ("work") matrices.