(git:e7e05ae)
rt_propagation_output.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 Routine for the real time propagation output.
10 !> \author Florian Schiffmann (02.09)
11 ! **************************************************************************************************
12 
14  USE atomic_kind_types, ONLY: atomic_kind_type
15  USE cp_control_types, ONLY: dft_control_type
23  cp_fm_struct_type
24  USE cp_fm_types, ONLY: cp_fm_create,&
26  cp_fm_release,&
28  cp_fm_type
32  cp_logger_type
34  cp_p_file,&
39  USE dbcsr_api, ONLY: &
40  dbcsr_add, dbcsr_binary_write, dbcsr_checksum, dbcsr_copy, dbcsr_create, &
41  dbcsr_deallocate_matrix, dbcsr_desymmetrize, dbcsr_filter, dbcsr_get_info, &
42  dbcsr_get_occupation, dbcsr_init_p, dbcsr_iterator_blocks_left, dbcsr_iterator_next_block, &
43  dbcsr_iterator_start, dbcsr_iterator_stop, dbcsr_iterator_type, dbcsr_p_type, dbcsr_scale, &
44  dbcsr_set, dbcsr_type
45  USE efield_utils, ONLY: make_field
46  USE input_constants, ONLY: ehrenfest,&
50  section_vals_type
51  USE kahan_sum, ONLY: accurate_sum
52  USE kinds, ONLY: default_path_length,&
53  dp
54  USE machine, ONLY: m_flush
55  USE message_passing, ONLY: mp_comm_type
56  USE parallel_gemm_api, ONLY: parallel_gemm
57  USE particle_list_types, ONLY: particle_list_type
58  USE particle_types, ONLY: particle_type
59  USE physcon, ONLY: femtoseconds
60  USE pw_env_types, ONLY: pw_env_get,&
61  pw_env_type
62  USE pw_methods, ONLY: pw_zero
63  USE pw_pool_types, ONLY: pw_pool_type
64  USE pw_types, ONLY: pw_c1d_gs_type,&
65  pw_r3d_rs_type
66  USE qs_energy_types, ONLY: qs_energy_type
67  USE qs_environment_types, ONLY: get_qs_env,&
68  qs_environment_type
69  USE qs_kind_types, ONLY: get_qs_kind_set,&
70  qs_kind_type
72  USE qs_linres_types, ONLY: current_env_type
74  USE qs_rho_types, ONLY: qs_rho_get,&
75  qs_rho_type
80  USE qs_scf_types, ONLY: qs_scf_env_type
81  USE qs_subsys_types, ONLY: qs_subsys_get,&
82  qs_subsys_type
84  USE rt_propagation_types, ONLY: get_rtp,&
85  rt_prop_type
89 #include "../base/base_uses.f90"
90 
91  IMPLICIT NONE
92 
93  PRIVATE
94 
95  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'rt_propagation_output'
96 
97  PUBLIC :: rt_prop_output, &
101 
102 CONTAINS
103 
104 ! **************************************************************************************************
105 !> \brief ...
106 !> \param qs_env ...
107 !> \param run_type ...
108 !> \param delta_iter ...
109 !> \param used_time ...
110 ! **************************************************************************************************
111  SUBROUTINE rt_prop_output(qs_env, run_type, delta_iter, used_time)
112  TYPE(qs_environment_type), POINTER :: qs_env
113  INTEGER, INTENT(in) :: run_type
114  REAL(dp), INTENT(in), OPTIONAL :: delta_iter, used_time
115 
116  INTEGER :: n_electrons, n_proj, nspin, output_unit, &
117  spin
118  REAL(dp) :: orthonormality, tot_rho_r
119  REAL(kind=dp), DIMENSION(:), POINTER :: qs_tot_rho_r
120  TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
121  TYPE(cp_fm_type), DIMENSION(:), POINTER :: mos_new
122  TYPE(cp_logger_type), POINTER :: logger
123  TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s, p_im, rho_new
124  TYPE(dft_control_type), POINTER :: dft_control
125  TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
126  TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
127  TYPE(qs_rho_type), POINTER :: rho
128  TYPE(rt_prop_type), POINTER :: rtp
129  TYPE(section_vals_type), POINTER :: dft_section, input, rtp_section
130 
131  NULLIFY (logger, dft_control)
132 
133  logger => cp_get_default_logger()
134  CALL get_qs_env(qs_env, &
135  rtp=rtp, &
136  matrix_s=matrix_s, &
137  input=input, &
138  rho=rho, &
139  particle_set=particle_set, &
140  atomic_kind_set=atomic_kind_set, &
141  qs_kind_set=qs_kind_set, &
142  dft_control=dft_control)
143 
144  rtp_section => section_vals_get_subs_vals(input, "DFT%REAL_TIME_PROPAGATION")
145 
146  CALL get_qs_kind_set(qs_kind_set, nelectron=n_electrons)
147  n_electrons = n_electrons - dft_control%charge
148 
149  CALL qs_rho_get(rho_struct=rho, tot_rho_r=qs_tot_rho_r)
150 
151  tot_rho_r = accurate_sum(qs_tot_rho_r)
152 
153  output_unit = cp_print_key_unit_nr(logger, rtp_section, "PRINT%PROGRAM_RUN_INFO", &
154  extension=".scfLog")
155 
156  IF (output_unit > 0) THEN
157  WRITE (output_unit, fmt="(/,(T3,A,T40,I5))") &
158  "Information at iteration step:", rtp%iter
159  WRITE (unit=output_unit, fmt="((T3,A,T41,2F20.10))") &
160  "Total electronic density (r-space): ", &
161  tot_rho_r, &
162  tot_rho_r + &
163  REAL(n_electrons, dp)
164  WRITE (unit=output_unit, fmt="((T3,A,T59,F22.14))") &
165  "Total energy:", rtp%energy_new
166  IF (run_type == ehrenfest) &
167  WRITE (unit=output_unit, fmt="((T3,A,T61,F20.14))") &
168  "Energy difference to previous iteration step:", rtp%energy_new - rtp%energy_old
169  IF (run_type == real_time_propagation) &
170  WRITE (unit=output_unit, fmt="((T3,A,T61,F20.14))") &
171  "Energy difference to initial state:", rtp%energy_new - rtp%energy_old
172  IF (PRESENT(delta_iter)) &
173  WRITE (unit=output_unit, fmt="((T3,A,T61,E20.6))") &
174  "Convergence:", delta_iter
175  IF (rtp%converged) THEN
176  IF (run_type == real_time_propagation) &
177  WRITE (unit=output_unit, fmt="((T3,A,T61,F12.2))") &
178  "Time needed for propagation:", used_time
179  WRITE (unit=output_unit, fmt="(/,(T3,A,3X,F16.14))") &
180  "CONVERGENCE REACHED", rtp%energy_new - rtp%energy_old
181  END IF
182  END IF
183 
184  IF (rtp%converged) THEN
185  IF (.NOT. rtp%linear_scaling) THEN
186  CALL get_rtp(rtp=rtp, mos_new=mos_new)
187  CALL rt_calculate_orthonormality(orthonormality, &
188  mos_new, matrix_s(1)%matrix)
189  IF (output_unit > 0) &
190  WRITE (output_unit, fmt="(/,(T3,A,T60,F20.10))") &
191  "Max deviation from orthonormalization:", orthonormality
192  END IF
193  END IF
194 
195  IF (output_unit > 0) &
196  CALL m_flush(output_unit)
197  CALL cp_print_key_finished_output(output_unit, logger, rtp_section, &
198  "PRINT%PROGRAM_RUN_INFO")
199 
200  IF (rtp%converged) THEN
201  dft_section => section_vals_get_subs_vals(input, "DFT")
202  IF (btest(cp_print_key_should_output(logger%iter_info, &
203  dft_section, "REAL_TIME_PROPAGATION%PRINT%FIELD"), cp_p_file)) &
204  CALL print_field_applied(qs_env, dft_section)
205  CALL make_moment(qs_env)
206  IF (btest(cp_print_key_should_output(logger%iter_info, &
207  dft_section, "REAL_TIME_PROPAGATION%PRINT%E_CONSTITUENTS"), cp_p_file)) THEN
208  CALL print_rtp_energy_components(qs_env, dft_section)
209  END IF
210  IF (.NOT. dft_control%qs_control%dftb) THEN
211  CALL write_available_results(qs_env=qs_env, rtp=rtp)
212  END IF
213  IF (rtp%linear_scaling) THEN
214  CALL get_rtp(rtp=rtp, rho_new=rho_new)
215  IF (btest(cp_print_key_should_output(logger%iter_info, &
216  dft_section, "REAL_TIME_PROPAGATION%PRINT%RESTART"), cp_p_file)) THEN
217  CALL write_rt_p_to_restart(rho_new, .false.)
218  END IF
219  IF (btest(cp_print_key_should_output(logger%iter_info, &
220  dft_section, "REAL_TIME_PROPAGATION%PRINT%RESTART_HISTORY"), cp_p_file)) THEN
221  CALL write_rt_p_to_restart(rho_new, .true.)
222  END IF
223  IF (.NOT. dft_control%qs_control%dftb) THEN
224  !Not sure if these things could also work with dftb or not
225  IF (btest(cp_print_key_should_output(logger%iter_info, &
226  dft_section, "REAL_TIME_PROPAGATION%PRINT%CURRENT"), cp_p_file)) THEN
227  DO spin = 1, SIZE(rho_new)/2
228  CALL rt_current(qs_env, rho_new(2*spin)%matrix, dft_section, spin, SIZE(rho_new)/2)
229  END DO
230  END IF
231  END IF
232  ELSE
233  CALL get_rtp(rtp=rtp, mos_new=mos_new)
234  IF (.NOT. dft_control%qs_control%dftb) THEN
235  IF (btest(cp_print_key_should_output(logger%iter_info, &
236  dft_section, "REAL_TIME_PROPAGATION%PRINT%CURRENT"), cp_p_file)) THEN
237  NULLIFY (p_im)
238  nspin = SIZE(mos_new)/2
239  CALL dbcsr_allocate_matrix_set(p_im, nspin)
240  DO spin = 1, nspin
241  CALL dbcsr_init_p(p_im(spin)%matrix)
242  CALL dbcsr_create(p_im(spin)%matrix, template=matrix_s(1)%matrix, matrix_type="N")
243  END DO
244  CALL calculate_p_imaginary(qs_env, rtp, p_im)
245  DO spin = 1, nspin
246  CALL rt_current(qs_env, p_im(spin)%matrix, dft_section, spin, nspin)
247  END DO
248  CALL dbcsr_deallocate_matrix_set(p_im)
249  END IF
250  IF (dft_control%rtp_control%is_proj_mo) THEN
251  DO n_proj = 1, SIZE(dft_control%rtp_control%proj_mo_list)
252  CALL compute_and_write_proj_mo(qs_env, mos_new, &
253  dft_control%rtp_control%proj_mo_list(n_proj)%proj_mo, n_proj)
254  END DO
255  END IF
256  END IF
257  CALL write_rt_mos_to_restart(qs_env%mos, mos_new, particle_set, &
258  dft_section, qs_kind_set)
259  END IF
260  END IF
261 
262  rtp%energy_old = rtp%energy_new
263 
264  IF (.NOT. rtp%converged .AND. rtp%iter >= dft_control%rtp_control%max_iter) &
265  CALL cp_abort(__location__, "EMD did not converge, either increase MAX_ITER "// &
266  "or use a smaller TIMESTEP")
267 
268  END SUBROUTINE rt_prop_output
269 
270 ! **************************************************************************************************
271 !> \brief computes the effective orthonormality of a set of mos given an s-matrix
272 !> orthonormality is the max deviation from unity of the C^T S C
273 !> \param orthonormality ...
274 !> \param mos_new ...
275 !> \param matrix_s ...
276 !> \author Florian Schiffmann (02.09)
277 ! **************************************************************************************************
278  SUBROUTINE rt_calculate_orthonormality(orthonormality, mos_new, matrix_s)
279  REAL(kind=dp), INTENT(out) :: orthonormality
280  TYPE(cp_fm_type), DIMENSION(:), POINTER :: mos_new
281  TYPE(dbcsr_type), OPTIONAL, POINTER :: matrix_s
282 
283  CHARACTER(len=*), PARAMETER :: routinen = 'rt_calculate_orthonormality'
284 
285  INTEGER :: handle, i, im, ispin, j, k, n, &
286  ncol_local, nrow_local, nspin, re
287  INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
288  REAL(kind=dp) :: alpha, max_alpha, max_beta
289  TYPE(cp_fm_struct_type), POINTER :: tmp_fm_struct
290  TYPE(cp_fm_type) :: overlap_re, svec_im, svec_re
291 
292  NULLIFY (tmp_fm_struct)
293 
294  CALL timeset(routinen, handle)
295 
296  nspin = SIZE(mos_new)/2
297  max_alpha = 0.0_dp
298  max_beta = 0.0_dp
299  DO ispin = 1, nspin
300  re = ispin*2 - 1
301  im = ispin*2
302  ! get S*C
303  CALL cp_fm_create(svec_re, mos_new(im)%matrix_struct)
304  CALL cp_fm_create(svec_im, mos_new(im)%matrix_struct)
305  CALL cp_fm_get_info(mos_new(im), &
306  nrow_global=n, ncol_global=k)
307  CALL cp_dbcsr_sm_fm_multiply(matrix_s, mos_new(re), &
308  svec_re, k)
309  CALL cp_dbcsr_sm_fm_multiply(matrix_s, mos_new(im), &
310  svec_im, k)
311 
312  ! get C^T (S*C)
313  CALL cp_fm_struct_create(tmp_fm_struct, nrow_global=k, ncol_global=k, &
314  para_env=mos_new(re)%matrix_struct%para_env, &
315  context=mos_new(re)%matrix_struct%context)
316  CALL cp_fm_create(overlap_re, tmp_fm_struct)
317 
318  CALL cp_fm_struct_release(tmp_fm_struct)
319 
320  CALL parallel_gemm('T', 'N', k, k, n, 1.0_dp, mos_new(re), &
321  svec_re, 0.0_dp, overlap_re)
322  CALL parallel_gemm('T', 'N', k, k, n, 1.0_dp, mos_new(im), &
323  svec_im, 1.0_dp, overlap_re)
324 
325  CALL cp_fm_release(svec_re)
326  CALL cp_fm_release(svec_im)
327 
328  CALL cp_fm_get_info(overlap_re, nrow_local=nrow_local, ncol_local=ncol_local, &
329  row_indices=row_indices, col_indices=col_indices)
330  DO i = 1, nrow_local
331  DO j = 1, ncol_local
332  alpha = overlap_re%local_data(i, j)
333  IF (row_indices(i) .EQ. col_indices(j)) alpha = alpha - 1.0_dp
334  max_alpha = max(max_alpha, abs(alpha))
335  END DO
336  END DO
337  CALL cp_fm_release(overlap_re)
338  END DO
339  CALL mos_new(1)%matrix_struct%para_env%max(max_alpha)
340  CALL mos_new(1)%matrix_struct%para_env%max(max_beta)
341  orthonormality = max_alpha
342 
343  CALL timestop(handle)
344 
345  END SUBROUTINE rt_calculate_orthonormality
346 
347 ! **************************************************************************************************
348 !> \brief computes the convergence criterion for RTP and EMD
349 !> \param rtp ...
350 !> \param matrix_s Overlap matrix without the derivatives
351 !> \param delta_mos ...
352 !> \param delta_eps ...
353 !> \author Florian Schiffmann (02.09)
354 ! **************************************************************************************************
355 
356  SUBROUTINE rt_convergence(rtp, matrix_s, delta_mos, delta_eps)
357  TYPE(rt_prop_type), POINTER :: rtp
358  TYPE(dbcsr_type), POINTER :: matrix_s
359  TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: delta_mos
360  REAL(dp), INTENT(out) :: delta_eps
361 
362  CHARACTER(len=*), PARAMETER :: routinen = 'rt_convergence'
363  REAL(kind=dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
364 
365  INTEGER :: handle, i, icol, im, ispin, j, lcol, &
366  lrow, nao, newdim, nmo, nspin, re
367  LOGICAL :: double_col, double_row
368  REAL(kind=dp) :: alpha, max_alpha
369  TYPE(cp_fm_struct_type), POINTER :: newstruct, newstruct1, tmp_fm_struct
370  TYPE(cp_fm_type) :: work, work1, work2
371  TYPE(cp_fm_type), DIMENSION(:), POINTER :: mos_new
372 
373  NULLIFY (tmp_fm_struct)
374 
375  CALL timeset(routinen, handle)
376 
377  CALL get_rtp(rtp=rtp, mos_new=mos_new)
378 
379  nspin = SIZE(delta_mos)/2
380  max_alpha = 0.0_dp
381 
382  DO i = 1, SIZE(mos_new)
383  CALL cp_fm_scale_and_add(-one, delta_mos(i), one, mos_new(i))
384  END DO
385 
386  DO ispin = 1, nspin
387  re = ispin*2 - 1
388  im = ispin*2
389 
390  double_col = .true.
391  double_row = .false.
392  CALL cp_fm_struct_double(newstruct, &
393  delta_mos(re)%matrix_struct, &
394  delta_mos(re)%matrix_struct%context, &
395  double_col, &
396  double_row)
397 
398  CALL cp_fm_create(work, matrix_struct=newstruct)
399  CALL cp_fm_create(work1, matrix_struct=newstruct)
400 
401  CALL cp_fm_get_info(delta_mos(re), ncol_local=lcol, ncol_global=nmo, &
402  nrow_global=nao)
403  CALL cp_fm_get_info(work, ncol_global=newdim)
404 
405  CALL cp_fm_set_all(work, zero, zero)
406 
407  DO icol = 1, lcol
408  work%local_data(:, icol) = delta_mos(re)%local_data(:, icol)
409  work%local_data(:, icol + lcol) = delta_mos(im)%local_data(:, icol)
410  END DO
411 
412  CALL cp_dbcsr_sm_fm_multiply(matrix_s, work, work1, ncol=newdim)
413 
414  CALL cp_fm_release(work)
415 
416  CALL cp_fm_struct_create(tmp_fm_struct, nrow_global=nmo, ncol_global=nmo, &
417  para_env=delta_mos(re)%matrix_struct%para_env, &
418  context=delta_mos(re)%matrix_struct%context)
419  CALL cp_fm_struct_double(newstruct1, &
420  tmp_fm_struct, &
421  delta_mos(re)%matrix_struct%context, &
422  double_col, &
423  double_row)
424 
425  CALL cp_fm_create(work, matrix_struct=newstruct1)
426  CALL cp_fm_create(work2, matrix_struct=newstruct1)
427 
428  CALL parallel_gemm("T", "N", nmo, newdim, nao, one, delta_mos(re), &
429  work1, zero, work)
430 
431  CALL parallel_gemm("T", "N", nmo, newdim, nao, one, delta_mos(im), &
432  work1, zero, work2)
433 
434  CALL cp_fm_get_info(work, nrow_local=lrow)
435  DO i = 1, lrow
436  DO j = 1, lcol
437  alpha = sqrt((work%local_data(i, j) + work2%local_data(i, j + lcol))**2 + &
438  (work%local_data(i, j + lcol) - work2%local_data(i, j))**2)
439  max_alpha = max(max_alpha, abs(alpha))
440  END DO
441  END DO
442 
443  CALL cp_fm_release(work)
444  CALL cp_fm_release(work1)
445  CALL cp_fm_release(work2)
446  CALL cp_fm_struct_release(tmp_fm_struct)
447  CALL cp_fm_struct_release(newstruct)
448  CALL cp_fm_struct_release(newstruct1)
449 
450  END DO
451 
452  CALL delta_mos(1)%matrix_struct%para_env%max(max_alpha)
453  delta_eps = sqrt(max_alpha)
454 
455  CALL timestop(handle)
456 
457  END SUBROUTINE rt_convergence
458 
459 ! **************************************************************************************************
460 !> \brief computes the convergence criterion for RTP and EMD based on the density matrix
461 !> \param rtp ...
462 !> \param delta_P ...
463 !> \param delta_eps ...
464 !> \author Samuel Andermatt (02.14)
465 ! **************************************************************************************************
466 
467  SUBROUTINE rt_convergence_density(rtp, delta_P, delta_eps)
468 
469  TYPE(rt_prop_type), POINTER :: rtp
470  TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: delta_p
471  REAL(dp), INTENT(out) :: delta_eps
472 
473  CHARACTER(len=*), PARAMETER :: routinen = 'rt_convergence_density'
474  REAL(kind=dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
475 
476  INTEGER :: col_atom, group_handle, handle, i, &
477  ispin, row_atom
478  REAL(dp) :: alpha, max_alpha
479  REAL(dp), DIMENSION(:), POINTER :: block_values
480  TYPE(dbcsr_iterator_type) :: iter
481  TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rho_new
482  TYPE(dbcsr_type), POINTER :: tmp
483  TYPE(mp_comm_type) :: group
484 
485  CALL timeset(routinen, handle)
486 
487  CALL get_rtp(rtp=rtp, rho_new=rho_new)
488 
489  DO i = 1, SIZE(rho_new)
490  CALL dbcsr_add(delta_p(i)%matrix, rho_new(i)%matrix, one, -one)
491  END DO
492  !get the maximum value of delta_P
493  DO i = 1, SIZE(delta_p)
494  !square all entries of both matrices
495  CALL dbcsr_iterator_start(iter, delta_p(i)%matrix)
496  DO WHILE (dbcsr_iterator_blocks_left(iter))
497  CALL dbcsr_iterator_next_block(iter, row_atom, col_atom, block_values)
498  block_values = block_values*block_values
499  END DO
500  CALL dbcsr_iterator_stop(iter)
501  END DO
502  NULLIFY (tmp)
503  ALLOCATE (tmp)
504  CALL dbcsr_create(tmp, template=delta_p(1)%matrix, matrix_type="N")
505  DO ispin = 1, SIZE(delta_p)/2
506  CALL dbcsr_desymmetrize(delta_p(2*ispin - 1)%matrix, tmp)
507  CALL dbcsr_add(delta_p(2*ispin)%matrix, tmp, one, one)
508  END DO
509  !the absolute values are now in the even entries of delta_P
510  max_alpha = zero
511  DO ispin = 1, SIZE(delta_p)/2
512  CALL dbcsr_iterator_start(iter, delta_p(2*ispin)%matrix)
513  DO WHILE (dbcsr_iterator_blocks_left(iter))
514  CALL dbcsr_iterator_next_block(iter, row_atom, col_atom, block_values)
515  alpha = maxval(block_values)
516  IF (alpha > max_alpha) max_alpha = alpha
517  END DO
518  CALL dbcsr_iterator_stop(iter)
519  END DO
520  CALL dbcsr_get_info(delta_p(1)%matrix, group=group_handle)
521  CALL group%set_handle(group_handle)
522  CALL group%max(max_alpha)
523  delta_eps = sqrt(max_alpha)
524  CALL dbcsr_deallocate_matrix(tmp)
525  CALL timestop(handle)
526 
527  END SUBROUTINE rt_convergence_density
528 
529 ! **************************************************************************************************
530 !> \brief interface to qs_moments. Does only work for nonperiodic dipole
531 !> \param qs_env ...
532 !> \author Florian Schiffmann (02.09)
533 ! **************************************************************************************************
534 
535  SUBROUTINE make_moment(qs_env)
536 
537  TYPE(qs_environment_type), POINTER :: qs_env
538 
539  CHARACTER(len=*), PARAMETER :: routinen = 'make_moment'
540 
541  INTEGER :: handle, output_unit
542  TYPE(cp_logger_type), POINTER :: logger
543  TYPE(dft_control_type), POINTER :: dft_control
544 
545  CALL timeset(routinen, handle)
546 
547  NULLIFY (dft_control)
548 
549  logger => cp_get_default_logger()
550  output_unit = cp_logger_get_default_io_unit(logger)
551  CALL get_qs_env(qs_env, dft_control=dft_control)
552  IF (dft_control%qs_control%dftb) THEN
553  CALL scf_post_calculation_tb(qs_env, "DFTB", .false.)
554  ELSE IF (dft_control%qs_control%xtb) THEN
555  CALL scf_post_calculation_tb(qs_env, "xTB", .false.)
556  ELSE
557  CALL qs_scf_post_moments(qs_env%input, logger, qs_env, output_unit)
558  END IF
559  CALL timestop(handle)
560 
561  END SUBROUTINE make_moment
562 
563 ! **************************************************************************************************
564 !> \brief Reports the sparsity pattern of the complex density matrix
565 !> \param filter_eps ...
566 !> \param rho ...
567 !> \author Samuel Andermatt (09.14)
568 ! **************************************************************************************************
569 
570  SUBROUTINE report_density_occupation(filter_eps, rho)
571 
572  REAL(kind=dp) :: filter_eps
573  TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rho
574 
575  CHARACTER(len=*), PARAMETER :: routinen = 'report_density_occupation'
576 
577  INTEGER :: handle, i, im, ispin, re, unit_nr
578  REAL(kind=dp) :: eps, occ
579  TYPE(cp_logger_type), POINTER :: logger
580  TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: tmp
581 
582  CALL timeset(routinen, handle)
583 
584  logger => cp_get_default_logger()
585  unit_nr = cp_logger_get_default_io_unit(logger)
586  NULLIFY (tmp)
587  CALL dbcsr_allocate_matrix_set(tmp, SIZE(rho))
588  DO i = 1, SIZE(rho)
589  CALL dbcsr_init_p(tmp(i)%matrix)
590  CALL dbcsr_create(tmp(i)%matrix, template=rho(i)%matrix)
591  CALL dbcsr_copy(tmp(i)%matrix, rho(i)%matrix)
592  END DO
593  DO ispin = 1, SIZE(rho)/2
594  re = 2*ispin - 1
595  im = 2*ispin
596  eps = max(filter_eps, 1.0e-11_dp)
597  DO WHILE (eps < 1.1_dp)
598  CALL dbcsr_filter(tmp(re)%matrix, eps)
599  occ = dbcsr_get_occupation(tmp(re)%matrix)
600  IF (unit_nr > 0) WRITE (unit_nr, fmt="((T3,A,I1,A,F15.12,A,T61,F20.10))") "Occupation of rho spin ", &
601  ispin, " eps ", eps, " real: ", occ
602  eps = eps*10
603  END DO
604  eps = max(filter_eps, 1.0e-11_dp)
605  DO WHILE (eps < 1.1_dp)
606  CALL dbcsr_filter(tmp(im)%matrix, eps)
607  occ = dbcsr_get_occupation(tmp(im)%matrix)
608  IF (unit_nr > 0) WRITE (unit_nr, fmt="((T3,A,I1,A,F15.12,A,T61,F20.10))") "Occupation of rho spin ", &
609  ispin, " eps ", eps, " imag: ", occ
610  eps = eps*10.0_dp
611  END DO
612  END DO
614  CALL timestop(handle)
615 
616  END SUBROUTINE report_density_occupation
617 
618 ! **************************************************************************************************
619 !> \brief Writes the density matrix and the atomic positions to a restart file
620 !> \param rho_new ...
621 !> \param history ...
622 !> \author Samuel Andermatt (09.14)
623 ! **************************************************************************************************
624 
625  SUBROUTINE write_rt_p_to_restart(rho_new, history)
626 
627  TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rho_new
628  LOGICAL :: history
629 
630  CHARACTER(LEN=*), PARAMETER :: routinen = 'write_rt_p_to_restart'
631 
632  CHARACTER(LEN=default_path_length) :: file_name, project_name
633  INTEGER :: handle, im, ispin, re, unit_nr
634  REAL(kind=dp) :: cs_pos
635  TYPE(cp_logger_type), POINTER :: logger
636 
637  CALL timeset(routinen, handle)
638  logger => cp_get_default_logger()
639  IF (logger%para_env%is_source()) THEN
640  unit_nr = cp_logger_get_default_unit_nr(logger, local=.true.)
641  ELSE
642  unit_nr = -1
643  END IF
644 
645  project_name = logger%iter_info%project_name
646  DO ispin = 1, SIZE(rho_new)/2
647  re = 2*ispin - 1
648  im = 2*ispin
649  IF (history) THEN
650  WRITE (file_name, '(A,I0,A)') &
651  trim(project_name)//"_LS_DM_SPIN_RE", ispin, "_"//trim(cp_iter_string(logger%iter_info))//"_RESTART.dm"
652  ELSE
653  WRITE (file_name, '(A,I0,A)') trim(project_name)//"_LS_DM_SPIN_RE", ispin, "_RESTART.dm"
654  END IF
655  cs_pos = dbcsr_checksum(rho_new(re)%matrix, pos=.true.)
656  IF (unit_nr > 0) THEN
657  WRITE (unit_nr, '(T2,A,E20.8)') "Writing restart DM "//trim(file_name)//" with checksum: ", cs_pos
658  END IF
659  CALL dbcsr_binary_write(rho_new(re)%matrix, file_name)
660  IF (history) THEN
661  WRITE (file_name, '(A,I0,A)') &
662  trim(project_name)//"_LS_DM_SPIN_IM", ispin, "_"//trim(cp_iter_string(logger%iter_info))//"_RESTART.dm"
663  ELSE
664  WRITE (file_name, '(A,I0,A)') trim(project_name)//"_LS_DM_SPIN_IM", ispin, "_RESTART.dm"
665  END IF
666  cs_pos = dbcsr_checksum(rho_new(im)%matrix, pos=.true.)
667  IF (unit_nr > 0) THEN
668  WRITE (unit_nr, '(T2,A,E20.8)') "Writing restart DM "//trim(file_name)//" with checksum: ", cs_pos
669  END IF
670  CALL dbcsr_binary_write(rho_new(im)%matrix, file_name)
671  END DO
672 
673  CALL timestop(handle)
674 
675  END SUBROUTINE write_rt_p_to_restart
676 
677 ! **************************************************************************************************
678 !> \brief Collocation of the current and printing of it in a cube file
679 !> \param qs_env ...
680 !> \param P_im ...
681 !> \param dft_section ...
682 !> \param spin ...
683 !> \param nspin ...
684 !> \author Samuel Andermatt (06.15)
685 ! **************************************************************************************************
686  SUBROUTINE rt_current(qs_env, P_im, dft_section, spin, nspin)
687  TYPE(qs_environment_type), POINTER :: qs_env
688  TYPE(dbcsr_type), POINTER :: p_im
689  TYPE(section_vals_type), POINTER :: dft_section
690  INTEGER :: spin, nspin
691 
692  CHARACTER(len=*), PARAMETER :: routinen = 'rt_current'
693 
694  CHARACTER(len=1) :: char_spin
695  CHARACTER(len=14) :: ext
696  CHARACTER(len=2) :: sdir
697  INTEGER :: dir, handle, print_unit
698  INTEGER, DIMENSION(:), POINTER :: stride(:)
699  LOGICAL :: mpi_io
700  TYPE(cp_logger_type), POINTER :: logger
701  TYPE(current_env_type) :: current_env
702  TYPE(dbcsr_type), POINTER :: tmp, zero
703  TYPE(particle_list_type), POINTER :: particles
704  TYPE(pw_c1d_gs_type) :: gs
705  TYPE(pw_env_type), POINTER :: pw_env
706  TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
707  TYPE(pw_r3d_rs_type) :: rs
708  TYPE(qs_subsys_type), POINTER :: subsys
709 
710  CALL timeset(routinen, handle)
711 
712  logger => cp_get_default_logger()
713  CALL get_qs_env(qs_env=qs_env, subsys=subsys, pw_env=pw_env)
714  CALL qs_subsys_get(subsys, particles=particles)
715  CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
716 
717  NULLIFY (zero, tmp)
718  ALLOCATE (zero, tmp)
719  CALL dbcsr_create(zero, template=p_im)
720  CALL dbcsr_copy(zero, p_im)
721  CALL dbcsr_set(zero, 0.0_dp)
722  CALL dbcsr_create(tmp, template=p_im)
723  CALL dbcsr_copy(tmp, p_im)
724  IF (nspin == 1) THEN
725  CALL dbcsr_scale(tmp, 0.5_dp)
726  END IF
727  current_env%gauge = -1
728  current_env%gauge_init = .false.
729  CALL auxbas_pw_pool%create_pw(rs)
730  CALL auxbas_pw_pool%create_pw(gs)
731 
732  NULLIFY (stride)
733  ALLOCATE (stride(3))
734 
735  DO dir = 1, 3
736 
737  CALL pw_zero(rs)
738  CALL pw_zero(gs)
739 
740  CALL calculate_jrho_resp(zero, tmp, zero, zero, dir, dir, rs, gs, qs_env, current_env, retain_rsgrid=.true.)
741 
742  stride = section_get_ivals(dft_section, "REAL_TIME_PROPAGATION%PRINT%CURRENT%STRIDE")
743 
744  IF (dir == 1) THEN
745  sdir = "-x"
746  ELSEIF (dir == 2) THEN
747  sdir = "-y"
748  ELSE
749  sdir = "-z"
750  END IF
751  WRITE (char_spin, "(I1)") spin
752 
753  ext = "-SPIN-"//char_spin//sdir//".cube"
754  mpi_io = .true.
755  print_unit = cp_print_key_unit_nr(logger, dft_section, "REAL_TIME_PROPAGATION%PRINT%CURRENT", &
756  extension=ext, file_status="REPLACE", file_action="WRITE", &
757  log_filename=.false., mpi_io=mpi_io)
758 
759  CALL cp_pw_to_cube(rs, print_unit, "EMD current", particles=particles, stride=stride, &
760  mpi_io=mpi_io)
761 
762  CALL cp_print_key_finished_output(print_unit, logger, dft_section, "REAL_TIME_PROPAGATION%PRINT%CURRENT", &
763  mpi_io=mpi_io)
764 
765  END DO
766 
767  CALL auxbas_pw_pool%give_back_pw(rs)
768  CALL auxbas_pw_pool%give_back_pw(gs)
769 
770  CALL dbcsr_deallocate_matrix(zero)
771  CALL dbcsr_deallocate_matrix(tmp)
772 
773  DEALLOCATE (stride)
774 
775  CALL timestop(handle)
776 
777  END SUBROUTINE rt_current
778 
779 ! **************************************************************************************************
780 !> \brief Interface routine to trigger writing of results available from normal
781 !> SCF. Can write MO-dependent and MO free results (needed for call from
782 !> the linear scaling code)
783 !> Update: trigger also some of prints for time-dependent runs
784 !> \param qs_env ...
785 !> \param rtp ...
786 !> \par History
787 !> 2022-11 Update [Guillaume Le Breton]
788 ! **************************************************************************************************
789  SUBROUTINE write_available_results(qs_env, rtp)
790  TYPE(qs_environment_type), POINTER :: qs_env
791  TYPE(rt_prop_type), POINTER :: rtp
792 
793  CHARACTER(len=*), PARAMETER :: routinen = 'write_available_results'
794 
795  INTEGER :: handle
796  TYPE(qs_scf_env_type), POINTER :: scf_env
797 
798  CALL timeset(routinen, handle)
799 
800  CALL get_qs_env(qs_env, scf_env=scf_env)
801  IF (rtp%linear_scaling) THEN
802  CALL write_mo_free_results(qs_env)
803  ELSE
804  CALL write_mo_free_results(qs_env)
805  CALL write_mo_dependent_results(qs_env, scf_env)
806  ! Time-dependent MO print
807  CALL write_rtp_mos_to_output_unit(qs_env, rtp)
808  CALL write_rtp_mo_cubes(qs_env, rtp)
809  END IF
810 
811  CALL timestop(handle)
812 
813  END SUBROUTINE write_available_results
814 
815 ! **************************************************************************************************
816 !> \brief Print the field applied to the system. Either the electric
817 !> field or the vector potential depending on the gauge used
818 !> \param qs_env ...
819 !> \param dft_section ...
820 !> \par History
821 !> 2023-01 Created [Guillaume Le Breton]
822 ! **************************************************************************************************
823  SUBROUTINE print_field_applied(qs_env, dft_section)
824  TYPE(qs_environment_type), POINTER :: qs_env
825  TYPE(section_vals_type), POINTER :: dft_section
826 
827  CHARACTER(LEN=3), DIMENSION(3) :: rlab
828  CHARACTER(LEN=default_path_length) :: filename
829  INTEGER :: i, i_step, output_unit, unit_nr
830  LOGICAL :: new_file
831  REAL(kind=dp) :: field(3)
832  TYPE(cp_logger_type), POINTER :: logger
833  TYPE(dft_control_type), POINTER :: dft_control
834  TYPE(rt_prop_type), POINTER :: rtp
835 
836  NULLIFY (dft_control)
837 
838  logger => cp_get_default_logger()
839  output_unit = cp_logger_get_default_io_unit(logger)
840 
841  CALL get_qs_env(qs_env, dft_control=dft_control, rtp=rtp)
842 
843  i_step = rtp%istep
844 
845  unit_nr = cp_print_key_unit_nr(logger, dft_section, &
846  "REAL_TIME_PROPAGATION%PRINT%FIELD", extension=".dat", is_new_file=new_file)
847 
848  IF (output_unit > 0) THEN
849  rlab = [CHARACTER(LEN=3) :: "X", "Y", "Z"]
850  IF (unit_nr /= output_unit) THEN
851  INQUIRE (unit=unit_nr, name=filename)
852  WRITE (unit=output_unit, fmt="(/,T2,A,2(/,T3,A),/)") &
853  "FIELD", "The field applied is written to the file:", &
854  trim(filename)
855  ELSE
856  WRITE (unit=output_unit, fmt="(/,T2,A)") "FIELD APPLIED [a.u.]"
857  WRITE (unit=output_unit, fmt="(T5,3(A,A,E16.8,1X))") &
858  (trim(rlab(i)), "=", dft_control%rtp_control%field(i), i=1, 3)
859  END IF
860 
861  IF (new_file) THEN
862  IF (dft_control%apply_efield_field) THEN
863  WRITE (unit=unit_nr, fmt='("#",5X,A,8X,A,3(6X,A))') "Step Nr.", "Time[fs]", " Field X", " Field Y", " Field Z"
864  ELSE IF (dft_control%apply_vector_potential) THEN
865  WRITE (unit=unit_nr, fmt='("#",5X,A,8X,A,6(6X,A))') "Step Nr.", "Time[fs]", " Field X", " Field Y", " Field Z", &
866  "Vec. Pot. X", "Vec. Pot. Y", "Vec. Pot. Z"
867  END IF
868  END IF
869 
870  IF (dft_control%apply_efield_field) THEN
871  CALL make_field(dft_control, field, qs_env%sim_step, qs_env%sim_time)
872  WRITE (unit=unit_nr, fmt="(I10,F20.6,3(E16.8,1X))") i_step, i_step*rtp%dt*femtoseconds, field(1), field(2), field(3)
873  ELSE IF (dft_control%apply_vector_potential) THEN
874  WRITE (unit=unit_nr, fmt="(I10,F20.6,6(E16.8,1X))") i_step, i_step*rtp%dt*femtoseconds, &
875  dft_control%rtp_control%field(1), dft_control%rtp_control%field(2), dft_control%rtp_control%field(3), &
876  dft_control%rtp_control%vec_pot(1), dft_control%rtp_control%vec_pot(2), dft_control%rtp_control%vec_pot(3)
877  END IF
878 
879  END IF
880 
881  CALL cp_print_key_finished_output(unit_nr, logger, dft_section, &
882  "REAL_TIME_PROPAGATION%PRINT%FIELD")
883 
884  END SUBROUTINE print_field_applied
885 
886 ! **************************************************************************************************
887 !> \brief Print the components of the total energy used in an RTP calculation
888 !> \param qs_env ...
889 !> \param dft_section ...
890 !> \par History
891 !> 2024-02 Created [ANB]
892 ! **************************************************************************************************
893  SUBROUTINE print_rtp_energy_components(qs_env, dft_section)
894  TYPE(qs_environment_type), POINTER :: qs_env
895  TYPE(section_vals_type), POINTER :: dft_section
896 
897  CHARACTER(LEN=default_path_length) :: filename
898  INTEGER :: i_step, output_unit, unit_nr
899  LOGICAL :: new_file
900  TYPE(cp_logger_type), POINTER :: logger
901  TYPE(dft_control_type), POINTER :: dft_control
902  TYPE(qs_energy_type), POINTER :: energy
903  TYPE(rt_prop_type), POINTER :: rtp
904 
905  NULLIFY (dft_control, energy, rtp)
906 
907  logger => cp_get_default_logger()
908  output_unit = cp_logger_get_default_io_unit(logger)
909 
910  CALL get_qs_env(qs_env, dft_control=dft_control, rtp=rtp, energy=energy)
911  i_step = rtp%istep
912 
913  unit_nr = cp_print_key_unit_nr(logger, dft_section, &
914  "REAL_TIME_PROPAGATION%PRINT%E_CONSTITUENTS", extension=".ener", &
915  file_action="WRITE", is_new_file=new_file)
916 
917  IF (output_unit > 0) THEN
918  IF (unit_nr /= output_unit) THEN
919  INQUIRE (unit=unit_nr, name=filename)
920  WRITE (unit=output_unit, fmt="(/,T2,A,2(/,T3,A),/)") &
921  "ENERGY_CONSTITUENTS", "Total Energy constituents written to file:", &
922  trim(filename)
923  ELSE
924  WRITE (unit=output_unit, fmt="(/,T2,A)") "ENERGY_CONSTITUENTS"
925  END IF
926 
927  IF (new_file) THEN
928  ! NOTE that these are not all terms contributing to the total energy for RTP, only a selection of those
929  ! most significant / impactful. Therefore the printed components likely will not add up to the total energy.
930  WRITE (unit=unit_nr, fmt='("#",5X,A,8X,A,10(6X,A))') "Step Nr.", "Time[fs]", &
931  "Total ener.[a.u.]", "core[a.u.] ", " overlap [a.u.]", "hartree[a.u.]", " exc. [a.u.] ", &
932  " hartree 1c[a.u.]", "exc 1c[a.u.] ", "exc admm[a.u.]", "exc 1c admm[a.u.]", "efield LG"
933 
934  END IF
935  WRITE (unit=unit_nr, fmt="(I10,F20.6,10(F20.9))") &
936  i_step, i_step*rtp%dt*femtoseconds, energy%total, energy%core, energy%core_overlap, energy%hartree, energy%exc, &
937  energy%hartree_1c, energy%exc1, energy%exc_aux_fit, energy%exc1_aux_fit, energy%efield_core
938 
939  END IF
940 
941  CALL cp_print_key_finished_output(unit_nr, logger, dft_section, &
942  "REAL_TIME_PROPAGATION%PRINT%E_CONSTITUENTS")
943 
944  END SUBROUTINE print_rtp_energy_components
945 
946 END MODULE rt_propagation_output
Define the atomic kind types and their sub types.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
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
Definition: cp_fm_struct.F:14
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
Definition: cp_fm_struct.F:125
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 ...
Definition: cp_fm_struct.F:536
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
Definition: cp_fm_struct.F:320
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
Definition: cp_fm_types.F:1016
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
Definition: cp_fm_types.F:535
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp)
creates a new full matrix with the given structure
Definition: cp_fm_types.F:167
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
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
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
character(len=default_string_length) function, public cp_iter_string(iter_info, print_key, for_file)
returns the iteration string, a string that is useful to create unique filenames (once you trim it)
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
integer, parameter, public cp_p_file
integer function, public cp_print_key_should_output(iteration_info, basis_section, print_key_path, used_print_key, first_time)
returns what should be done with the given property if btest(res,cp_p_store) then the property should...
A wrapper around pw_to_cube() which accepts particle_list_type.
subroutine, public cp_pw_to_cube(pw, unit_nr, title, particles, stride, zero_tails, silent, mpi_io)
...
all routins needed for a nonperiodic electric field
Definition: efield_utils.F:12
subroutine, public make_field(dft_control, field, sim_step, sim_time)
computes the amplitude of the efield within a given envelop
Definition: efield_utils.F:118
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public ehrenfest
integer, parameter, public real_time_propagation
objects that represent the structure of input sections and the data contained in an input section
integer function, dimension(:), pointer, public section_get_ivals(section_vals, keyword_name)
...
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
sums arrays of real/complex numbers with much reduced round-off as compared to a naive implementation...
Definition: kahan_sum.F:29
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
integer, parameter, public default_path_length
Definition: kinds.F:58
Machine interface based on Fortran 2003 and POSIX.
Definition: machine.F:17
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition: machine.F:106
Interface to the message passing library MPI.
basic linear algebra operations for full matrixes
represent a simple array based list of the given type
Define the data structure for the particle information.
Definition of physical constants:
Definition: physcon.F:68
real(kind=dp), parameter, public femtoseconds
Definition: physcon.F:153
container for various plainwaves related things
Definition: pw_env_types.F:14
subroutine, public pw_env_get(pw_env, pw_pools, cube_info, gridlevel_info, auxbas_pw_pool, auxbas_grid, auxbas_rs_desc, auxbas_rs_grid, rs_descs, rs_grids, xc_pw_pool, vdw_pw_pool, poisson_env, interp_section)
returns the various attributes of the pw env
Definition: pw_env_types.F:113
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Definition: pw_pool_types.F:24
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_RI_aux_kp, matrix_s, matrix_s_RI_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, WannierCentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, rhs)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
Definition: qs_kind_types.F:23
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)
Get attributes of an atomic kind set.
given the response wavefunctions obtained by the application of the (rxp), p, and ((dk-dl)xp) operato...
subroutine, public calculate_jrho_resp(mat_d0, mat_jp, mat_jp_rii, mat_jp_riii, iB, idir, current_rs, current_gs, qs_env, current_env, soft_valid, retain_rsgrid)
Calculation of the idir component of the response current density in the presence of a constant magne...
Type definitiona for linear response calculations.
Definition and initialisation of the mo data type.
Definition: qs_mo_io.F:21
subroutine, public write_rt_mos_to_restart(mo_array, rt_mos, particle_set, dft_section, qs_kind_set)
...
Definition: qs_mo_io.F:233
superstucture that hold various representations of the density and keeps track of which ones are vali...
Definition: qs_rho_types.F:18
subroutine, public qs_rho_get(rho_struct, rho_ao, rho_ao_im, rho_ao_kp, rho_ao_im_kp, rho_r, drho_r, rho_g, drho_g, tau_r, tau_g, rho_r_valid, drho_r_valid, rho_g_valid, drho_g_valid, tau_r_valid, tau_g_valid, tot_rho_r, tot_rho_g, rho_r_sccs, soft_valid, complex_rho_ao)
returns info about the density described by this object. If some representation is not available an e...
Definition: qs_rho_types.F:229
Does all kind of post scf calculations for GPW/GAPW.
subroutine, public write_mo_free_results(qs_env)
Write QS results always available (if switched on through the print_keys) Can be called from ls_scf.
subroutine, public qs_scf_post_moments(input, logger, qs_env, output_unit)
Computes and prints electric moments.
subroutine, public write_mo_dependent_results(qs_env, scf_env)
Write QS results available if MO's are present (if switched on through the print_keys) Writes only MO...
Does all kind of post scf calculations for DFTB.
subroutine, public scf_post_calculation_tb(qs_env, tb_type, no_mos)
collects possible post - scf calculations and prints info / computes properties.
module that contains the definitions of the scf types
Definition: qs_scf_types.F:14
types that represent a quickstep subsys
subroutine, public qs_subsys_get(subsys, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell, cell_ref, use_ref_cell, energy, force, qs_kind_set, cp_subsys, nelectron_total, nelectron_spin)
...
Function related to MO projection in RTP calculations.
subroutine, public compute_and_write_proj_mo(qs_env, mos_new, proj_mo, n_proj)
Compute the projection of the current MO coefficients on reference ones and write the results.
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_prop_output(qs_env, run_type, delta_iter, used_time)
...
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 write_rtp_mos_to_output_unit(qs_env, rtp)
...
subroutine, public write_rtp_mo_cubes(qs_env, rtp)
Write the time dependent amplitude of the MOs in real grid. Very close to qs_scf_post_gpw/qs_scf_post...
subroutine, public calculate_p_imaginary(qs_env, rtp, matrix_p_im, keep_sparsity)
...