(git:374b731)
Loading...
Searching...
No Matches
rt_propagation_methods.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 Routines for propagating the orbitals
10!> \author Florian Schiffmann (02.09)
11! **************************************************************************************************
13 USE bibliography, ONLY: kolafa2004,&
14 kuhne2007,&
15 cite_reference
16 USE cell_types, ONLY: cell_type
19 USE cp_cfm_types, ONLY: cp_cfm_create,&
34 USE cp_fm_types, ONLY: cp_fm_create,&
43 USE dbcsr_api, ONLY: &
44 dbcsr_add, dbcsr_copy, dbcsr_create, dbcsr_deallocate_matrix, dbcsr_desymmetrize, &
45 dbcsr_filter, dbcsr_frobenius_norm, dbcsr_get_block_p, dbcsr_init_p, &
46 dbcsr_iterator_blocks_left, dbcsr_iterator_next_block, dbcsr_iterator_start, &
47 dbcsr_iterator_stop, dbcsr_iterator_type, dbcsr_multiply, dbcsr_p_type, dbcsr_release, &
48 dbcsr_scale, dbcsr_set, dbcsr_transposed, dbcsr_type, dbcsr_type_antisymmetric
50 USE input_constants, ONLY: do_arnoldi,&
51 do_bch,&
52 do_em,&
53 do_pade,&
56 USE kinds, ONLY: dp
58 USE mathlib, ONLY: binomial
65 USE qs_ks_types, ONLY: set_ks_env
73 USE rt_propagation_types, ONLY: get_rtp,&
80#include "../base/base_uses.f90"
81
82 IMPLICIT NONE
83
84 PRIVATE
85
86 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'rt_propagation_methods'
87
88 PUBLIC :: propagation_step, &
90 calc_sinvh, &
92
93CONTAINS
94
95! **************************************************************************************************
96!> \brief performs a single propagation step a(t+Dt)=U(t+Dt,t)*a(0)
97!> and calculates the new exponential
98!> \param qs_env ...
99!> \param rtp ...
100!> \param rtp_control ...
101!> \author Florian Schiffmann (02.09)
102! **************************************************************************************************
103
104 SUBROUTINE propagation_step(qs_env, rtp, rtp_control)
105
106 TYPE(qs_environment_type), POINTER :: qs_env
107 TYPE(rt_prop_type), POINTER :: rtp
108 TYPE(rtp_control_type), POINTER :: rtp_control
109
110 CHARACTER(len=*), PARAMETER :: routinen = 'propagation_step'
111
112 INTEGER :: aspc_order, handle, i, im, re, unit_nr
113 TYPE(cell_type), POINTER :: cell
114 TYPE(cp_fm_type), DIMENSION(:), POINTER :: delta_mos, mos_new
115 TYPE(cp_logger_type), POINTER :: logger
116 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: delta_p, h_last_iter, ks_mix, ks_mix_im, &
117 matrix_ks, matrix_ks_im, matrix_s, &
118 rho_new
119 TYPE(dft_control_type), POINTER :: dft_control
120
121 CALL timeset(routinen, handle)
122
123 logger => cp_get_default_logger()
124 IF (logger%para_env%is_source()) THEN
125 unit_nr = cp_logger_get_default_unit_nr(logger, local=.true.)
126 ELSE
127 unit_nr = -1
128 END IF
129
130 NULLIFY (cell, delta_p, rho_new, delta_mos, mos_new)
131 NULLIFY (ks_mix, ks_mix_im)
132 ! get everything needed and set some values
133 CALL get_qs_env(qs_env, cell=cell, matrix_s=matrix_s, dft_control=dft_control)
134
135 IF (rtp%iter == 1) THEN
136 CALL qs_energies_init(qs_env, .false.)
137 !the above recalculates matrix_s, but matrix not changed if ions are fixed
138 IF (rtp_control%fixed_ions) CALL set_ks_env(qs_env%ks_env, s_mstruct_changed=.false.)
139
140 ! add additional terms to matrix_h and matrix_h_im in the case of applied electric field,
141 ! either in the lengh or velocity gauge.
142 ! should be called after qs_energies_init and before qs_ks_update_qs_env
143 IF (dft_control%apply_efield_field) THEN
144 IF (any(cell%perd(1:3) /= 0)) &
145 cpabort("Length gauge (efield) and periodicity are not compatible")
147 ELSE IF (rtp_control%velocity_gauge) THEN
148 IF (dft_control%apply_vector_potential) &
149 CALL update_vector_potential(qs_env, dft_control)
150 CALL velocity_gauge_ks_matrix(qs_env, subtract_nl_term=.false.)
151 END IF
152
153 CALL get_qs_env(qs_env, matrix_s=matrix_s)
154 IF (.NOT. rtp_control%fixed_ions) THEN
155 CALL s_matrices_create(matrix_s, rtp)
156 END IF
157 rtp%delta_iter = 100.0_dp
158 rtp%mixing_factor = 1.0_dp
159 rtp%mixing = .false.
160 aspc_order = rtp_control%aspc_order
161 CALL aspc_extrapolate(rtp, matrix_s, aspc_order)
162 IF (rtp%linear_scaling) THEN
163 CALL calc_update_rho_sparse(qs_env)
164 ELSE
165 CALL calc_update_rho(qs_env)
166 END IF
167 CALL qs_ks_update_qs_env(qs_env, calculate_forces=.false.)
168 END IF
169 IF (.NOT. rtp_control%fixed_ions) THEN
170 CALL calc_s_derivs(qs_env)
171 END IF
172 rtp%converged = .false.
173
174 IF (rtp%linear_scaling) THEN
175 ! keep temporary copy of the starting density matrix to check for convergence
176 CALL get_rtp(rtp=rtp, rho_new=rho_new)
177 NULLIFY (delta_p)
178 CALL dbcsr_allocate_matrix_set(delta_p, SIZE(rho_new))
179 DO i = 1, SIZE(rho_new)
180 CALL dbcsr_init_p(delta_p(i)%matrix)
181 CALL dbcsr_create(delta_p(i)%matrix, template=rho_new(i)%matrix)
182 CALL dbcsr_copy(delta_p(i)%matrix, rho_new(i)%matrix)
183 END DO
184 ELSE
185 ! keep temporary copy of the starting mos to check for convergence
186 CALL get_rtp(rtp=rtp, mos_new=mos_new)
187 ALLOCATE (delta_mos(SIZE(mos_new)))
188 DO i = 1, SIZE(mos_new)
189 CALL cp_fm_create(delta_mos(i), &
190 matrix_struct=mos_new(i)%matrix_struct, &
191 name="delta_mos"//trim(adjustl(cp_to_string(i))))
192 CALL cp_fm_to_fm(mos_new(i), delta_mos(i))
193 END DO
194 END IF
195
196 CALL get_qs_env(qs_env, &
197 matrix_ks=matrix_ks, &
198 matrix_ks_im=matrix_ks_im)
199
200 CALL get_rtp(rtp=rtp, h_last_iter=h_last_iter)
201 IF (rtp%mixing) THEN
202 IF (unit_nr > 0) THEN
203 WRITE (unit_nr, '(t3,a,2f16.8)') "Mixing the Hamiltonians to improve robustness, mixing factor: ", rtp%mixing_factor
204 END IF
205 CALL dbcsr_allocate_matrix_set(ks_mix, SIZE(matrix_ks))
206 CALL dbcsr_allocate_matrix_set(ks_mix_im, SIZE(matrix_ks))
207 DO i = 1, SIZE(matrix_ks)
208 CALL dbcsr_init_p(ks_mix(i)%matrix)
209 CALL dbcsr_create(ks_mix(i)%matrix, template=matrix_ks(1)%matrix)
210 CALL dbcsr_init_p(ks_mix_im(i)%matrix)
211 CALL dbcsr_create(ks_mix_im(i)%matrix, template=matrix_ks(1)%matrix, matrix_type=dbcsr_type_antisymmetric)
212 END DO
213 DO i = 1, SIZE(matrix_ks)
214 re = 2*i - 1
215 im = 2*i
216 CALL dbcsr_add(ks_mix(i)%matrix, matrix_ks(i)%matrix, 0.0_dp, rtp%mixing_factor)
217 CALL dbcsr_add(ks_mix(i)%matrix, h_last_iter(re)%matrix, 1.0_dp, 1.0_dp - rtp%mixing_factor)
218 IF (rtp%propagate_complex_ks) THEN
219 CALL dbcsr_add(ks_mix_im(i)%matrix, matrix_ks_im(i)%matrix, 0.0_dp, rtp%mixing_factor)
220 CALL dbcsr_add(ks_mix_im(i)%matrix, h_last_iter(im)%matrix, 1.0_dp, 1.0_dp - rtp%mixing_factor)
221 END IF
222 END DO
223 CALL calc_sinvh(rtp, ks_mix, ks_mix_im, rtp_control)
224 DO i = 1, SIZE(matrix_ks)
225 re = 2*i - 1
226 im = 2*i
227 CALL dbcsr_copy(h_last_iter(re)%matrix, ks_mix(i)%matrix)
228 IF (rtp%propagate_complex_ks) THEN
229 CALL dbcsr_copy(h_last_iter(im)%matrix, ks_mix_im(i)%matrix)
230 END IF
231 END DO
232 CALL dbcsr_deallocate_matrix_set(ks_mix)
233 CALL dbcsr_deallocate_matrix_set(ks_mix_im)
234 ELSE
235 CALL calc_sinvh(rtp, matrix_ks, matrix_ks_im, rtp_control)
236 DO i = 1, SIZE(matrix_ks)
237 re = 2*i - 1
238 im = 2*i
239 CALL dbcsr_copy(h_last_iter(re)%matrix, matrix_ks(i)%matrix)
240 IF (rtp%propagate_complex_ks) THEN
241 CALL dbcsr_copy(h_last_iter(im)%matrix, matrix_ks_im(i)%matrix)
242 END IF
243 END DO
244 END IF
245
246 CALL compute_propagator_matrix(rtp, rtp_control%propagator)
247
248 SELECT CASE (rtp_control%mat_exp)
249 CASE (do_pade, do_taylor)
250 IF (rtp%linear_scaling) THEN
251 CALL propagate_exp_density(rtp, rtp_control)
252 CALL calc_update_rho_sparse(qs_env)
253 ELSE
254 CALL propagate_exp(rtp, rtp_control)
255 CALL calc_update_rho(qs_env)
256 END IF
257 CASE (do_arnoldi)
258 CALL propagate_arnoldi(rtp, rtp_control)
259 CALL calc_update_rho(qs_env)
260 CASE (do_bch)
261 CALL propagate_bch(rtp, rtp_control)
262 CALL calc_update_rho_sparse(qs_env)
263 END SELECT
264 CALL step_finalize(qs_env, rtp_control, delta_mos, delta_p)
265 IF (rtp%linear_scaling) THEN
266 CALL dbcsr_deallocate_matrix_set(delta_p)
267 ELSE
268 CALL cp_fm_release(delta_mos)
269 END IF
270
271 CALL timestop(handle)
272
273 END SUBROUTINE propagation_step
274
275! **************************************************************************************************
276!> \brief Performs all the stuff to finish the step:
277!> convergence checks
278!> copying stuff into right place for the next step
279!> updating the history for extrapolation
280!> \param qs_env ...
281!> \param rtp_control ...
282!> \param delta_mos ...
283!> \param delta_P ...
284!> \author Florian Schiffmann (02.09)
285! **************************************************************************************************
286
287 SUBROUTINE step_finalize(qs_env, rtp_control, delta_mos, delta_P)
288 TYPE(qs_environment_type), POINTER :: qs_env
289 TYPE(rtp_control_type), POINTER :: rtp_control
290 TYPE(cp_fm_type), DIMENSION(:), POINTER :: delta_mos
291 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: delta_p
292
293 CHARACTER(len=*), PARAMETER :: routinen = 'step_finalize'
294
295 INTEGER :: handle, i, ihist
296 TYPE(cp_fm_type), DIMENSION(:), POINTER :: mos_new, mos_old
297 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: exp_h_new, exp_h_old, matrix_ks, &
298 matrix_ks_im, rho_new, rho_old, s_mat
299 TYPE(qs_energy_type), POINTER :: energy
300 TYPE(rt_prop_type), POINTER :: rtp
301
302 CALL timeset(routinen, handle)
303
304 CALL get_qs_env(qs_env=qs_env, rtp=rtp, matrix_s=s_mat, &
305 matrix_ks=matrix_ks, matrix_ks_im=matrix_ks_im, energy=energy)
306 CALL get_rtp(rtp=rtp, exp_h_old=exp_h_old, exp_h_new=exp_h_new)
307
308 IF (rtp_control%sc_check_start .LT. rtp%iter) THEN
309 rtp%delta_iter_old = rtp%delta_iter
310 IF (rtp%linear_scaling) THEN
311 CALL rt_convergence_density(rtp, delta_p, rtp%delta_iter)
312 ELSE
313 CALL rt_convergence(rtp, s_mat(1)%matrix, delta_mos, rtp%delta_iter)
314 END IF
315 rtp%converged = (rtp%delta_iter .LT. rtp_control%eps_ener)
316 !Apply mixing if scf loop is not converging
317
318 !It would be better to redo the current step with mixixng,
319 !but currently the decision is made to use mixing from the next step on
320 IF (rtp_control%sc_check_start .LT. rtp%iter + 1) THEN
321 IF (rtp%delta_iter/rtp%delta_iter_old > 0.9) THEN
322 rtp%mixing_factor = max(rtp%mixing_factor/2.0_dp, 0.125_dp)
323 rtp%mixing = .true.
324 END IF
325 END IF
326 END IF
327
328 IF (rtp%converged) THEN
329 IF (rtp%linear_scaling) THEN
330 CALL get_rtp(rtp=rtp, rho_old=rho_old, rho_new=rho_new)
331 CALL purify_mcweeny_complex_nonorth(rho_new, s_mat, rtp%filter_eps, rtp%filter_eps_small, &
332 rtp_control%mcweeny_max_iter, rtp_control%mcweeny_eps)
333 IF (rtp_control%mcweeny_max_iter > 0) CALL calc_update_rho_sparse(qs_env)
334 CALL report_density_occupation(rtp%filter_eps, rho_new)
335 DO i = 1, SIZE(rho_new)
336 CALL dbcsr_copy(rho_old(i)%matrix, rho_new(i)%matrix)
337 END DO
338 ELSE
339 CALL get_rtp(rtp=rtp, mos_old=mos_old, mos_new=mos_new)
340 DO i = 1, SIZE(mos_new)
341 CALL cp_fm_to_fm(mos_new(i), mos_old(i))
342 END DO
343 END IF
344 IF (rtp_control%propagator == do_em) CALL calc_sinvh(rtp, matrix_ks, matrix_ks_im, rtp_control)
345 DO i = 1, SIZE(exp_h_new)
346 CALL dbcsr_copy(exp_h_old(i)%matrix, exp_h_new(i)%matrix)
347 END DO
348 ihist = mod(rtp%istep, rtp_control%aspc_order) + 1
349 IF (rtp_control%fixed_ions) THEN
350 CALL put_data_to_history(rtp, rho=rho_new, mos=mos_new, ihist=ihist)
351 ELSE
352 CALL put_data_to_history(rtp, rho=rho_new, mos=mos_new, s_mat=s_mat, ihist=ihist)
353 END IF
354 END IF
355
356 rtp%energy_new = energy%total
357
358 CALL timestop(handle)
359
360 END SUBROUTINE step_finalize
361
362! **************************************************************************************************
363!> \brief computes the propagator matrix for EM/ETRS, RTP/EMD
364!> \param rtp ...
365!> \param propagator ...
366!> \author Florian Schiffmann (02.09)
367! **************************************************************************************************
368
369 SUBROUTINE compute_propagator_matrix(rtp, propagator)
370 TYPE(rt_prop_type), POINTER :: rtp
371 INTEGER :: propagator
372
373 CHARACTER(len=*), PARAMETER :: routinen = 'compute_propagator_matrix'
374
375 INTEGER :: handle, i
376 REAL(kind=dp) :: dt, prefac
377 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: exp_h_new, exp_h_old, propagator_matrix
378
379 CALL timeset(routinen, handle)
380 CALL get_rtp(rtp=rtp, exp_h_new=exp_h_new, exp_h_old=exp_h_old, &
381 propagator_matrix=propagator_matrix, dt=dt)
382
383 prefac = -0.5_dp*dt
384
385 DO i = 1, SIZE(exp_h_new)
386 CALL dbcsr_add(propagator_matrix(i)%matrix, exp_h_new(i)%matrix, 0.0_dp, prefac)
387 IF (propagator == do_em) &
388 CALL dbcsr_add(propagator_matrix(i)%matrix, exp_h_old(i)%matrix, 1.0_dp, prefac)
389 END DO
390
391 CALL timestop(handle)
392
393 END SUBROUTINE compute_propagator_matrix
394
395! **************************************************************************************************
396!> \brief computes S_inv*H, if needed Sinv*B and S_inv*H_imag and store these quantities to the
397!> \brief exp_H for the real and imag part (for RTP and EMD)
398!> \param rtp ...
399!> \param matrix_ks ...
400!> \param matrix_ks_im ...
401!> \param rtp_control ...
402!> \author Florian Schiffmann (02.09)
403! **************************************************************************************************
404
405 SUBROUTINE calc_sinvh(rtp, matrix_ks, matrix_ks_im, rtp_control)
406 TYPE(rt_prop_type), POINTER :: rtp
407 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_ks, matrix_ks_im
408 TYPE(rtp_control_type), POINTER :: rtp_control
409
410 CHARACTER(len=*), PARAMETER :: routinen = 'calc_SinvH'
411 REAL(kind=dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
412
413 INTEGER :: handle, im, ispin, re
414 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: exp_h, sinvb, sinvh, sinvh_imag
415 TYPE(dbcsr_type) :: matrix_ks_nosym
416 TYPE(dbcsr_type), POINTER :: b_mat, s_inv
417
418 CALL timeset(routinen, handle)
419 CALL get_rtp(rtp=rtp, s_inv=s_inv, exp_h_new=exp_h)
420 DO ispin = 1, SIZE(matrix_ks)
421 re = ispin*2 - 1
422 im = ispin*2
423 CALL dbcsr_set(exp_h(re)%matrix, zero)
424 CALL dbcsr_set(exp_h(im)%matrix, zero)
425 END DO
426 CALL dbcsr_create(matrix_ks_nosym, template=matrix_ks(1)%matrix, matrix_type="N")
427
428 ! Real part of S_inv x H -> imag part of exp_H
429 DO ispin = 1, SIZE(matrix_ks)
430 re = ispin*2 - 1
431 im = ispin*2
432 CALL dbcsr_desymmetrize(matrix_ks(ispin)%matrix, matrix_ks_nosym)
433 CALL dbcsr_multiply("N", "N", one, s_inv, matrix_ks_nosym, zero, exp_h(im)%matrix, &
434 filter_eps=rtp%filter_eps)
435 IF (.NOT. rtp_control%fixed_ions) THEN
436 CALL get_rtp(rtp=rtp, sinvh=sinvh)
437 CALL dbcsr_copy(sinvh(ispin)%matrix, exp_h(im)%matrix)
438 END IF
439 END DO
440
441 ! Imag part of S_inv x H -> real part of exp_H
442 IF (rtp%propagate_complex_ks) THEN
443 DO ispin = 1, SIZE(matrix_ks)
444 re = ispin*2 - 1
445 im = ispin*2
446 CALL dbcsr_set(matrix_ks_nosym, 0.0_dp)
447 CALL dbcsr_desymmetrize(matrix_ks_im(ispin)%matrix, matrix_ks_nosym)
448 ! - SinvH_imag is added to exp_H(re)%matrix
449 CALL dbcsr_multiply("N", "N", -one, s_inv, matrix_ks_nosym, zero, exp_h(re)%matrix, &
450 filter_eps=rtp%filter_eps)
451 IF (.NOT. rtp_control%fixed_ions) THEN
452 CALL get_rtp(rtp=rtp, sinvh_imag=sinvh_imag)
453 ! -SinvH_imag is saved
454 CALL dbcsr_copy(sinvh_imag(ispin)%matrix, exp_h(re)%matrix)
455 END IF
456 END DO
457 END IF
458 ! EMD case: the real part of exp_H should be updated with B
459 IF (.NOT. rtp_control%fixed_ions) THEN
460 CALL get_rtp(rtp=rtp, b_mat=b_mat, sinvb=sinvb)
461 CALL dbcsr_set(matrix_ks_nosym, 0.0_dp)
462 CALL dbcsr_multiply("N", "N", one, s_inv, b_mat, zero, matrix_ks_nosym, filter_eps=rtp%filter_eps)
463 DO ispin = 1, SIZE(matrix_ks)
464 re = ispin*2 - 1
465 im = ispin*2
466 ! + SinvB is added to exp_H(re)%matrix
467 CALL dbcsr_add(exp_h(re)%matrix, matrix_ks_nosym, 1.0_dp, 1.0_dp)
468 ! + SinvB is saved
469 CALL dbcsr_copy(sinvb(ispin)%matrix, matrix_ks_nosym)
470 END DO
471 END IF
472 ! Otherwise no real part for exp_H
473
474 CALL dbcsr_release(matrix_ks_nosym)
475 CALL timestop(handle)
476
477 END SUBROUTINE calc_sinvh
478
479! **************************************************************************************************
480!> \brief calculates the needed overlap-like matrices
481!> depending on the way the exponential is calculated, only S^-1 is needed
482!> \param s_mat ...
483!> \param rtp ...
484!> \author Florian Schiffmann (02.09)
485! **************************************************************************************************
486
487 SUBROUTINE s_matrices_create(s_mat, rtp)
488
489 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: s_mat
490 TYPE(rt_prop_type), POINTER :: rtp
491
492 CHARACTER(len=*), PARAMETER :: routinen = 's_matrices_create'
493 REAL(kind=dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
494
495 INTEGER :: handle
496 TYPE(dbcsr_type), POINTER :: s_half, s_inv, s_minus_half
497
498 CALL timeset(routinen, handle)
499
500 CALL get_rtp(rtp=rtp, s_inv=s_inv)
501
502 IF (rtp%linear_scaling) THEN
503 CALL get_rtp(rtp=rtp, s_half=s_half, s_minus_half=s_minus_half)
504 CALL matrix_sqrt_newton_schulz(s_half, s_minus_half, s_mat(1)%matrix, rtp%filter_eps, &
505 rtp%newton_schulz_order, rtp%lanzcos_threshold, rtp%lanzcos_max_iter)
506 CALL dbcsr_multiply("N", "N", one, s_minus_half, s_minus_half, zero, s_inv, &
507 filter_eps=rtp%filter_eps)
508 ELSE
509 CALL dbcsr_copy(s_inv, s_mat(1)%matrix)
510 CALL cp_dbcsr_cholesky_decompose(s_inv, para_env=rtp%ao_ao_fmstruct%para_env, &
511 blacs_env=rtp%ao_ao_fmstruct%context)
512 CALL cp_dbcsr_cholesky_invert(s_inv, para_env=rtp%ao_ao_fmstruct%para_env, &
513 blacs_env=rtp%ao_ao_fmstruct%context, upper_to_full=.true.)
514 END IF
515
516 CALL timestop(handle)
517 END SUBROUTINE s_matrices_create
518
519! **************************************************************************************************
520!> \brief Calculates the frobenius norm of a complex matrix represented by two real matrices
521!> \param frob_norm ...
522!> \param mat_re ...
523!> \param mat_im ...
524!> \author Samuel Andermatt (04.14)
525! **************************************************************************************************
526
527 SUBROUTINE complex_frobenius_norm(frob_norm, mat_re, mat_im)
528
529 REAL(kind=dp), INTENT(out) :: frob_norm
530 TYPE(dbcsr_type), POINTER :: mat_re, mat_im
531
532 CHARACTER(len=*), PARAMETER :: routinen = 'complex_frobenius_norm'
533 REAL(kind=dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
534
535 INTEGER :: col_atom, handle, row_atom
536 LOGICAL :: found
537 REAL(dp), DIMENSION(:), POINTER :: block_values, block_values2
538 TYPE(dbcsr_iterator_type) :: iter
539 TYPE(dbcsr_type), POINTER :: tmp
540
541 CALL timeset(routinen, handle)
542
543 NULLIFY (tmp)
544 ALLOCATE (tmp)
545 CALL dbcsr_create(tmp, template=mat_re)
546 !make sure the tmp has the same sparsity pattern as the real and the complex part combined
547 CALL dbcsr_add(tmp, mat_re, zero, one)
548 CALL dbcsr_add(tmp, mat_im, zero, one)
549 CALL dbcsr_set(tmp, zero)
550 !calculate the hadamard product
551 CALL dbcsr_iterator_start(iter, tmp)
552 DO WHILE (dbcsr_iterator_blocks_left(iter))
553 CALL dbcsr_iterator_next_block(iter, row_atom, col_atom, block_values)
554 CALL dbcsr_get_block_p(mat_re, row_atom, col_atom, block_values2, found=found)
555 IF (found) THEN
556 block_values = block_values2*block_values2
557 END IF
558 CALL dbcsr_get_block_p(mat_im, row_atom, col_atom, block_values2, found=found)
559 IF (found) THEN
560 block_values = block_values + block_values2*block_values2
561 END IF
562 block_values = sqrt(block_values)
563 END DO
564 CALL dbcsr_iterator_stop(iter)
565 frob_norm = dbcsr_frobenius_norm(tmp)
566
567 CALL dbcsr_deallocate_matrix(tmp)
568
569 CALL timestop(handle)
570
571 END SUBROUTINE complex_frobenius_norm
572
573! **************************************************************************************************
574!> \brief Does McWeeny for complex matrices in the non-orthogonal basis
575!> \param P ...
576!> \param s_mat ...
577!> \param eps ...
578!> \param eps_small ...
579!> \param max_iter ...
580!> \param threshold ...
581!> \author Samuel Andermatt (04.14)
582! **************************************************************************************************
583
584 SUBROUTINE purify_mcweeny_complex_nonorth(P, s_mat, eps, eps_small, max_iter, threshold)
585
586 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: p, s_mat
587 REAL(kind=dp), INTENT(in) :: eps, eps_small
588 INTEGER, INTENT(in) :: max_iter
589 REAL(kind=dp), INTENT(in) :: threshold
590
591 CHARACTER(len=*), PARAMETER :: routinen = 'purify_mcweeny_complex_nonorth'
592 REAL(kind=dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
593
594 INTEGER :: handle, i, im, imax, ispin, re, unit_nr
595 REAL(kind=dp) :: frob_norm
596 TYPE(cp_logger_type), POINTER :: logger
597 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: ps, psp, tmp
598
599 CALL timeset(routinen, handle)
600
601 logger => cp_get_default_logger()
602 IF (logger%para_env%is_source()) THEN
603 unit_nr = cp_logger_get_default_unit_nr(logger, local=.true.)
604 ELSE
605 unit_nr = -1
606 END IF
607
608 NULLIFY (tmp, ps, psp)
609 CALL dbcsr_allocate_matrix_set(tmp, SIZE(p))
610 CALL dbcsr_allocate_matrix_set(psp, SIZE(p))
611 CALL dbcsr_allocate_matrix_set(ps, SIZE(p))
612 DO i = 1, SIZE(p)
613 CALL dbcsr_init_p(ps(i)%matrix)
614 CALL dbcsr_create(ps(i)%matrix, template=p(1)%matrix)
615 CALL dbcsr_init_p(psp(i)%matrix)
616 CALL dbcsr_create(psp(i)%matrix, template=p(1)%matrix)
617 CALL dbcsr_init_p(tmp(i)%matrix)
618 CALL dbcsr_create(tmp(i)%matrix, template=p(1)%matrix)
619 END DO
620 IF (SIZE(p) == 2) THEN
621 CALL dbcsr_scale(p(1)%matrix, one/2)
622 CALL dbcsr_scale(p(2)%matrix, one/2)
623 END IF
624 DO ispin = 1, SIZE(p)/2
625 re = 2*ispin - 1
626 im = 2*ispin
627 imax = max(max_iter, 1) !if max_iter is 0 then only the deviation from idempotency needs to be calculated
628 DO i = 1, imax
629 CALL dbcsr_multiply("N", "N", one, p(re)%matrix, s_mat(1)%matrix, &
630 zero, ps(re)%matrix, filter_eps=eps_small)
631 CALL dbcsr_multiply("N", "N", one, p(im)%matrix, s_mat(1)%matrix, &
632 zero, ps(im)%matrix, filter_eps=eps_small)
633 CALL cp_complex_dbcsr_gemm_3("N", "N", one, ps(re)%matrix, ps(im)%matrix, &
634 p(re)%matrix, p(im)%matrix, zero, psp(re)%matrix, psp(im)%matrix, &
635 filter_eps=eps_small)
636 CALL dbcsr_copy(tmp(re)%matrix, psp(re)%matrix)
637 CALL dbcsr_copy(tmp(im)%matrix, psp(im)%matrix)
638 CALL dbcsr_add(tmp(re)%matrix, p(re)%matrix, 1.0_dp, -1.0_dp)
639 CALL dbcsr_add(tmp(im)%matrix, p(im)%matrix, 1.0_dp, -1.0_dp)
640 CALL complex_frobenius_norm(frob_norm, tmp(re)%matrix, tmp(im)%matrix)
641 IF (unit_nr .GT. 0) WRITE (unit_nr, '(t3,a,2f16.8)') "Deviation from idempotency: ", frob_norm
642 IF (frob_norm .GT. threshold .AND. max_iter > 0) THEN
643 CALL dbcsr_copy(p(re)%matrix, psp(re)%matrix)
644 CALL dbcsr_copy(p(im)%matrix, psp(im)%matrix)
645 CALL cp_complex_dbcsr_gemm_3("N", "N", -2.0_dp, ps(re)%matrix, ps(im)%matrix, &
646 psp(re)%matrix, psp(im)%matrix, 3.0_dp, p(re)%matrix, p(im)%matrix, &
647 filter_eps=eps_small)
648 CALL dbcsr_filter(p(re)%matrix, eps)
649 CALL dbcsr_filter(p(im)%matrix, eps)
650 !make sure P is exactly hermitian
651 CALL dbcsr_transposed(tmp(re)%matrix, p(re)%matrix)
652 CALL dbcsr_add(p(re)%matrix, tmp(re)%matrix, one/2, one/2)
653 CALL dbcsr_transposed(tmp(im)%matrix, p(im)%matrix)
654 CALL dbcsr_add(p(im)%matrix, tmp(im)%matrix, one/2, -one/2)
655 ELSE
656 EXIT
657 END IF
658 END DO
659 !make sure P is hermitian
660 CALL dbcsr_transposed(tmp(re)%matrix, p(re)%matrix)
661 CALL dbcsr_add(p(re)%matrix, tmp(re)%matrix, one/2, one/2)
662 CALL dbcsr_transposed(tmp(im)%matrix, p(im)%matrix)
663 CALL dbcsr_add(p(im)%matrix, tmp(im)%matrix, one/2, -one/2)
664 END DO
665 IF (SIZE(p) == 2) THEN
666 CALL dbcsr_scale(p(1)%matrix, one*2)
667 CALL dbcsr_scale(p(2)%matrix, one*2)
668 END IF
672
673 CALL timestop(handle)
674
675 END SUBROUTINE purify_mcweeny_complex_nonorth
676
677! **************************************************************************************************
678!> \brief ...
679!> \param rtp ...
680!> \param matrix_s ...
681!> \param aspc_order ...
682! **************************************************************************************************
683 SUBROUTINE aspc_extrapolate(rtp, matrix_s, aspc_order)
684 TYPE(rt_prop_type), POINTER :: rtp
685 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s
686 INTEGER, INTENT(in) :: aspc_order
687
688 CHARACTER(len=*), PARAMETER :: routinen = 'aspc_extrapolate'
689 COMPLEX(KIND=dp), PARAMETER :: cone = (1.0_dp, 0.0_dp), &
690 czero = (0.0_dp, 0.0_dp)
691 REAL(kind=dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
692
693 INTEGER :: handle, i, iaspc, icol_local, ihist, &
694 imat, k, kdbl, n, naspc, ncol_local, &
695 nmat
696 REAL(kind=dp) :: alpha
697 TYPE(cp_cfm_type) :: cfm_tmp, cfm_tmp1, csc
698 TYPE(cp_fm_struct_type), POINTER :: matrix_struct, matrix_struct_new
699 TYPE(cp_fm_type) :: fm_tmp, fm_tmp1, fm_tmp2
700 TYPE(cp_fm_type), DIMENSION(:), POINTER :: mos_new
701 TYPE(cp_fm_type), DIMENSION(:, :), POINTER :: mo_hist
702 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rho_new, s_hist
703 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: rho_hist
704
705 NULLIFY (rho_hist)
706 CALL timeset(routinen, handle)
707 CALL cite_reference(kolafa2004)
708 CALL cite_reference(kuhne2007)
709
710 IF (rtp%linear_scaling) THEN
711 CALL get_rtp(rtp=rtp, rho_new=rho_new)
712 ELSE
713 CALL get_rtp(rtp=rtp, mos_new=mos_new)
714 END IF
715
716 naspc = min(rtp%istep, aspc_order)
717 IF (rtp%linear_scaling) THEN
718 nmat = SIZE(rho_new)
719 rho_hist => rtp%history%rho_history
720 DO imat = 1, nmat
721 DO iaspc = 1, naspc
722 alpha = (-1.0_dp)**(iaspc + 1)*real(iaspc, kind=dp)* &
723 binomial(2*naspc, naspc - iaspc)/binomial(2*naspc - 2, naspc - 1)
724 ihist = mod(rtp%istep - iaspc, aspc_order) + 1
725 IF (iaspc == 1) THEN
726 CALL dbcsr_add(rho_new(imat)%matrix, rho_hist(imat, ihist)%matrix, zero, alpha)
727 ELSE
728 CALL dbcsr_add(rho_new(imat)%matrix, rho_hist(imat, ihist)%matrix, one, alpha)
729 END IF
730 END DO
731 END DO
732 ELSE
733 mo_hist => rtp%history%mo_history
734 nmat = SIZE(mos_new)
735 DO imat = 1, nmat
736 DO iaspc = 1, naspc
737 alpha = (-1.0_dp)**(iaspc + 1)*real(iaspc, kind=dp)* &
738 binomial(2*naspc, naspc - iaspc)/binomial(2*naspc - 2, naspc - 1)
739 ihist = mod(rtp%istep - iaspc, aspc_order) + 1
740 IF (iaspc == 1) THEN
741 CALL cp_fm_scale_and_add(zero, mos_new(imat), alpha, mo_hist(imat, ihist))
742 ELSE
743 CALL cp_fm_scale_and_add(one, mos_new(imat), alpha, mo_hist(imat, ihist))
744 END IF
745 END DO
746 END DO
747
748 mo_hist => rtp%history%mo_history
749 s_hist => rtp%history%s_history
750 DO i = 1, SIZE(mos_new)/2
751 NULLIFY (matrix_struct, matrix_struct_new)
752
753 CALL cp_fm_struct_double(matrix_struct, &
754 mos_new(2*i)%matrix_struct, &
755 mos_new(2*i)%matrix_struct%context, &
756 .true., .false.)
757
758 CALL cp_fm_create(fm_tmp, matrix_struct)
759 CALL cp_fm_create(fm_tmp1, matrix_struct)
760 CALL cp_fm_create(fm_tmp2, mos_new(2*i)%matrix_struct)
761 CALL cp_cfm_create(cfm_tmp, mos_new(2*i)%matrix_struct)
762 CALL cp_cfm_create(cfm_tmp1, mos_new(2*i)%matrix_struct)
763
764 CALL cp_fm_get_info(fm_tmp, ncol_global=kdbl)
765
766 CALL cp_fm_get_info(mos_new(2*i), &
767 nrow_global=n, &
768 ncol_global=k, &
769 ncol_local=ncol_local)
770
771 CALL cp_fm_struct_create(matrix_struct_new, &
772 template_fmstruct=mos_new(2*i)%matrix_struct, &
773 nrow_global=k, &
774 ncol_global=k)
775 CALL cp_cfm_create(csc, matrix_struct_new)
776
777 CALL cp_fm_struct_release(matrix_struct_new)
778 CALL cp_fm_struct_release(matrix_struct)
779
780 ! first the most recent
781
782! reorthogonalize vectors
783 DO icol_local = 1, ncol_local
784 fm_tmp%local_data(:, icol_local) = mos_new(2*i - 1)%local_data(:, icol_local)
785 fm_tmp%local_data(:, icol_local + ncol_local) = mos_new(2*i)%local_data(:, icol_local)
786 END DO
787
788 CALL cp_dbcsr_sm_fm_multiply(matrix_s(1)%matrix, fm_tmp, fm_tmp1, kdbl)
789
790 DO icol_local = 1, ncol_local
791 cfm_tmp%local_data(:, icol_local) = cmplx(fm_tmp1%local_data(:, icol_local), &
792 fm_tmp1%local_data(:, icol_local + ncol_local), dp)
793 cfm_tmp1%local_data(:, icol_local) = cmplx(mos_new(2*i - 1)%local_data(:, icol_local), &
794 mos_new(2*i)%local_data(:, icol_local), dp)
795 END DO
796 CALL parallel_gemm('C', 'N', k, k, n, cone, cfm_tmp1, cfm_tmp, czero, csc)
798 CALL cp_cfm_triangular_multiply(csc, cfm_tmp1, n_cols=k, side='R', invert_tr=.true.)
799 DO icol_local = 1, ncol_local
800 mos_new(2*i - 1)%local_data(:, icol_local) = real(cfm_tmp1%local_data(:, icol_local), dp)
801 mos_new(2*i)%local_data(:, icol_local) = aimag(cfm_tmp1%local_data(:, icol_local))
802 END DO
803
804! deallocate work matrices
805 CALL cp_cfm_release(csc)
806 CALL cp_fm_release(fm_tmp)
807 CALL cp_fm_release(fm_tmp1)
808 CALL cp_fm_release(fm_tmp2)
809 CALL cp_cfm_release(cfm_tmp)
810 CALL cp_cfm_release(cfm_tmp1)
811 END DO
812
813 END IF
814
815 CALL timestop(handle)
816
817 END SUBROUTINE aspc_extrapolate
818
819! **************************************************************************************************
820!> \brief ...
821!> \param rtp ...
822!> \param mos ...
823!> \param rho ...
824!> \param s_mat ...
825!> \param ihist ...
826! **************************************************************************************************
827 SUBROUTINE put_data_to_history(rtp, mos, rho, s_mat, ihist)
828 TYPE(rt_prop_type), POINTER :: rtp
829 TYPE(cp_fm_type), DIMENSION(:), POINTER :: mos
830 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rho
831 TYPE(dbcsr_p_type), DIMENSION(:), OPTIONAL, &
832 POINTER :: s_mat
833 INTEGER :: ihist
834
835 INTEGER :: i
836
837 IF (rtp%linear_scaling) THEN
838 DO i = 1, SIZE(rho)
839 CALL dbcsr_copy(rtp%history%rho_history(i, ihist)%matrix, rho(i)%matrix)
840 END DO
841 ELSE
842 DO i = 1, SIZE(mos)
843 CALL cp_fm_to_fm(mos(i), rtp%history%mo_history(i, ihist))
844 END DO
845 IF (PRESENT(s_mat)) THEN
846 IF (ASSOCIATED(rtp%history%s_history(ihist)%matrix)) THEN ! the sparsity might be different
847 ! (future struct:check)
848 CALL dbcsr_deallocate_matrix(rtp%history%s_history(ihist)%matrix)
849 END IF
850 ALLOCATE (rtp%history%s_history(ihist)%matrix)
851 CALL dbcsr_copy(rtp%history%s_history(ihist)%matrix, s_mat(1)%matrix)
852 END IF
853 END IF
854
855 END SUBROUTINE put_data_to_history
856
857END MODULE rt_propagation_methods
static int imax(int x, int y)
Returns the larger of two given integer (missing from the C standard)
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public kuhne2007
integer, save, public kolafa2004
Handles all functions related to the CELL.
Definition cell_types.F:15
Basic linear algebra operations for complex full matrices.
subroutine, public cp_cfm_triangular_multiply(triangular_matrix, matrix_b, side, transa_tr, invert_tr, uplo_tr, unit_diag_tr, n_rows, n_cols, alpha)
Multiplies in place by a triangular matrix: matrix_b = alpha op(triangular_matrix) matrix_b or (if si...
subroutine, public cp_cfm_cholesky_decompose(matrix, n, info_out)
Used to replace a symmetric positive definite matrix M with its Cholesky decomposition U: M = U^T * U...
Represents a complex full matrix distributed on many processors.
subroutine, public cp_cfm_create(matrix, matrix_struct, name)
Creates a new full matrix with the given structure.
subroutine, public cp_cfm_release(matrix)
Releases a full matrix.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
Interface to (sca)lapack for the Cholesky based procedures.
subroutine, public cp_dbcsr_cholesky_decompose(matrix, n, para_env, blacs_env)
used to replace a symmetric positive def. matrix M with its cholesky decomposition U: M = U^T * U,...
subroutine, public cp_dbcsr_cholesky_invert(matrix, n, para_env, blacs_env, upper_to_full)
used to replace the cholesky decomposition by the inverse
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
basic linear algebra operations for full matrices
subroutine, public cp_fm_scale_and_add(alpha, matrix_a, beta, matrix_b)
calc A <- alpha*A + beta*B optimized for alpha == 1.0 (just add beta*B) and beta == 0....
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_double(fmstruct, struct, context, col, row)
creates a struct with twice the number of blocks on each core. If matrix A has to be multiplied with ...
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_create(matrix, matrix_struct, name, use_sp)
creates a new full matrix with the given structure
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
all routins needed for a nonperiodic electric field
subroutine, public efield_potential_lengh_gauge(qs_env)
Replace the original implementation of the electric-electronic interaction in the length gauge....
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_bch
integer, parameter, public do_pade
integer, parameter, public do_taylor
integer, parameter, public do_em
integer, parameter, public do_arnoldi
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)
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
Routines for calculating a complex matrix exponential with dbcsr matrices. Based on the code in matri...
subroutine, public cp_complex_dbcsr_gemm_3(transa, transb, alpha, a_re, a_im, b_re, b_im, beta, c_re, c_im, filter_eps)
Convenience function. Computes the matrix multiplications needed for the multiplication of complex sp...
Collection of simple mathematical functions and subroutines.
Definition mathlib.F:15
elemental real(kind=dp) function, public binomial(n, k)
The binomial coefficient n over k for 0 <= k <= n is calculated, otherwise zero is returned.
Definition mathlib.F:205
basic linear algebra operations for full matrixes
Utility subroutine for qs energy calculation.
subroutine, public qs_energies_init(qs_env, calc_forces)
Refactoring of qs_energies_scf. Driver routine for the initial setup and calculations for a qs energy...
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.
routines that build the Kohn-Sham matrix (i.e calculate the coulomb and xc parts
subroutine, public qs_ks_update_qs_env(qs_env, calculate_forces, just_energy, print_active)
updates the Kohn Sham matrix of the given qs_env (facility method)
subroutine, public set_ks_env(ks_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, complex_ks, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, kinetic, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_ks_im_kp, vppl, rho_core, rho_nlcc, rho_nlcc_g, vee, neighbor_list_id, kpoints, sab_orb, sab_all, sac_ae, sac_ppl, sac_lri, sap_ppnl, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_nonbond, sab_vdw, sab_scp, sab_almo, sab_kp, sab_kp_nosym, task_list, task_list_soft, subsys, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env)
...
Routines for calculating a complex matrix exponential.
subroutine, public propagate_exp(rtp, rtp_control)
performs propagations if explicit matrix exponentials are used ETRS: exp(i*H(t+dt)*dt/2)*exp(i*H(t)*d...
subroutine, public propagate_bch(rtp, rtp_control)
Propagation using the Baker-Campbell-Hausdorff expansion, currently only works for rtp.
subroutine, public propagate_exp_density(rtp, rtp_control)
Propagation of the density matrix instead of the atomic orbitals via a matrix exponential.
subroutine, public propagate_arnoldi(rtp, rtp_control)
computes U_prop*MOs using arnoldi subspace algorithm
Routines for propagating the orbitals.
subroutine, public calc_sinvh(rtp, matrix_ks, matrix_ks_im, rtp_control)
computes S_inv*H, if needed Sinv*B and S_inv*H_imag and store these quantities to the
subroutine, public propagation_step(qs_env, rtp, rtp_control)
performs a single propagation step a(t+Dt)=U(t+Dt,t)*a(0) and calculates the new exponential
subroutine, public put_data_to_history(rtp, mos, rho, s_mat, ihist)
...
subroutine, public s_matrices_create(s_mat, rtp)
calculates the needed overlap-like matrices depending on the way the exponential is calculated,...
Routine for the real time propagation output.
subroutine, public report_density_occupation(filter_eps, rho)
Reports the sparsity pattern of the complex density matrix.
subroutine, public rt_convergence(rtp, matrix_s, delta_mos, delta_eps)
computes the convergence criterion for RTP and EMD
subroutine, public rt_convergence_density(rtp, delta_p, delta_eps)
computes the convergence criterion for RTP and EMD based on the density matrix
Types and set_get for real time propagation depending on runtype and diagonalization method different...
subroutine, public get_rtp(rtp, exp_h_old, exp_h_new, h_last_iter, rho_old, rho_next, rho_new, mos, mos_new, mos_old, mos_next, s_inv, s_half, s_minus_half, b_mat, c_mat, propagator_matrix, mixing, mixing_factor, s_der, dt, nsteps, sinvh, sinvh_imag, sinvb, admm_mos)
...
Routines needed for EMD.
subroutine, public calc_update_rho_sparse(qs_env)
Copies the density matrix back into the qs_envrhorho_ao.
subroutine, public calc_s_derivs(qs_env)
Calculates dS/dR respectily the velocity weighted derivatves only needed for ehrenfest MD.
subroutine, public calc_update_rho(qs_env)
calculates the density from the complex MOs and passes the density to qs_env.
Routines to perform the RTP in the velocity gauge.
subroutine, public velocity_gauge_ks_matrix(qs_env, subtract_nl_term)
...
subroutine, public update_vector_potential(qs_env, dft_control)
Update the vector potential in the case where a time-dependant electric field is apply.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
Represent a complex full matrix.
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...