(git:f56c6e3)
Loading...
Searching...
No Matches
rt_propagation_utils.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Routines needed for EMD
10!> \author Florian Schiffmann (02.09)
11! **************************************************************************************************
12
15 USE cell_types, ONLY: cell_type
19 USE cp_dbcsr_api, ONLY: &
27 USE cp_files, ONLY: close_file,&
31 USE cp_fm_types, ONLY: cp_fm_create,&
41 USE cp_output_handling, ONLY: cp_p_file,&
47 USE efield_utils, ONLY: make_field
57 USE kinds, ONLY: default_path_length,&
59 dp
60 USE mathconstants, ONLY: zero
63 USE orbital_pointers, ONLY: ncoset
66 USE physcon, ONLY: femtoseconds
67 USE pw_env_types, ONLY: pw_env_get,&
69 USE pw_methods, ONLY: pw_multiply,&
71 USE pw_pool_types, ONLY: pw_pool_p_type,&
73 USE pw_types, ONLY: pw_c1d_gs_type,&
81 USE qs_ks_types, ONLY: qs_ks_did_change,&
86 USE qs_mo_types, ONLY: allocate_mo_set,&
94 USE qs_rho_types, ONLY: qs_rho_get,&
97 USE qs_scf_wfn_mix, ONLY: wfn_mix
100 USE rt_propagation_types, ONLY: get_rtp,&
102#include "../base/base_uses.f90"
103
104 IMPLICIT NONE
105 PRIVATE
106
107 PUBLIC :: get_restart_wfn, &
115 read_moments, &
117
118 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'rt_propagation_utils'
119
120CONTAINS
121
122! **************************************************************************************************
123!> \brief Calculates dS/dR respectily the velocity weighted derivatves
124!> only needed for ehrenfest MD.
125!>
126!> \param qs_env the qs environment
127!> \par History
128!> 02.2009 created [Manuel Guidon]
129!> 02.2014 switched to dbcsr matrices [Samuel Andermatt]
130!> \author Florian Schiffmann
131! **************************************************************************************************
132 SUBROUTINE calc_s_derivs(qs_env)
133 TYPE(qs_environment_type), POINTER :: qs_env
134
135 CHARACTER(LEN=*), PARAMETER :: routinen = 'calc_S_derivs'
136 REAL(kind=dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
137
138 INTEGER :: col_atom, handle, i, j, m, maxder, n, &
139 nder, row_atom
140 INTEGER, DIMENSION(6, 2) :: c_map_mat
141 LOGICAL :: return_s_derivatives
142 REAL(dp), DIMENSION(:, :), POINTER :: block_values
143 TYPE(dbcsr_iterator_type) :: iter
144 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: c_mat, s_der, s_derivs
145 TYPE(dbcsr_type), POINTER :: b_mat, tmp_mat, tmp_mat2
146 TYPE(dft_control_type), POINTER :: dft_control
147 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
148 POINTER :: sab_orb
149 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
150 TYPE(qs_ks_env_type), POINTER :: ks_env
151 TYPE(rt_prop_type), POINTER :: rtp
152
153 CALL timeset(routinen, handle)
154
155 return_s_derivatives = .true.
156
157 NULLIFY (particle_set)
158 NULLIFY (rtp)
159 NULLIFY (s_derivs)
160 NULLIFY (dft_control)
161 NULLIFY (ks_env)
162
163 CALL get_qs_env(qs_env=qs_env, &
164 rtp=rtp, &
165 particle_set=particle_set, &
166 sab_orb=sab_orb, &
167 dft_control=dft_control, &
168 ks_env=ks_env)
169
170 CALL get_rtp(rtp=rtp, b_mat=b_mat, c_mat=c_mat, s_der=s_der)
171
172 nder = 2
173 maxder = ncoset(nder)
174
175 NULLIFY (tmp_mat)
176 ALLOCATE (tmp_mat)
177 CALL dbcsr_create(tmp_mat, template=s_der(1)%matrix, matrix_type="N")
178
179 IF (rtp%iter < 2) THEN
180 ! calculate the overlap derivative matrices
181 IF (dft_control%qs_control%dftb) THEN
182 CALL build_dftb_overlap(qs_env, nder, s_derivs)
183 ELSE
184 CALL build_overlap_matrix(ks_env, nderivative=nder, matrix_s=s_derivs, &
185 basis_type_a="ORB", basis_type_b="ORB", sab_nl=sab_orb)
186 END IF
187
188 NULLIFY (tmp_mat2)
189 ALLOCATE (tmp_mat2)
190 CALL dbcsr_create(tmp_mat2, template=s_der(1)%matrix, matrix_type="S")
191 DO m = 1, 9
192 CALL dbcsr_copy(tmp_mat2, s_derivs(m + 1)%matrix)
193 CALL dbcsr_desymmetrize(tmp_mat2, s_der(m)%matrix)
194 CALL dbcsr_scale(s_der(m)%matrix, -one)
195 CALL dbcsr_filter(s_der(m)%matrix, rtp%filter_eps)
196 !The diagonal should be zero
197 CALL dbcsr_iterator_start(iter, s_der(m)%matrix)
198 DO WHILE (dbcsr_iterator_blocks_left(iter))
199 CALL dbcsr_iterator_next_block(iter, row_atom, col_atom, block_values)
200 IF (row_atom == col_atom) block_values = 0.0_dp
201 END DO
202 CALL dbcsr_iterator_stop(iter)
203 END DO
204 CALL dbcsr_deallocate_matrix_set(s_derivs)
205 CALL dbcsr_deallocate_matrix(tmp_mat2)
206 END IF
207
208 !calculate scalar product v(Rb)*<alpha|d/dRb beta> (B_mat), and store the first derivatives
209
210 CALL dbcsr_set(b_mat, zero)
211 DO m = 1, 3
212 CALL dbcsr_copy(tmp_mat, s_der(m)%matrix)
213 CALL dbcsr_iterator_start(iter, tmp_mat)
214 DO WHILE (dbcsr_iterator_blocks_left(iter))
215 CALL dbcsr_iterator_next_block(iter, row_atom, col_atom, block_values)
216 IF (row_atom == col_atom) block_values = 0.0_dp
217 block_values = block_values*particle_set(col_atom)%v(m)
218 END DO
219 CALL dbcsr_iterator_stop(iter)
220 CALL dbcsr_add(b_mat, tmp_mat, one, one)
221 END DO
222 CALL dbcsr_filter(b_mat, rtp%filter_eps)
223 !calculate C matrix: v(Rb)*<d/dRa alpha| d/dRb beta>
224
225 c_map_mat = 0
226 n = 0
227 DO j = 1, 3
228 DO m = j, 3
229 n = n + 1
230 c_map_mat(n, 1) = j
231 IF (m == j) cycle
232 c_map_mat(n, 2) = m
233 END DO
234 END DO
235
236 DO i = 1, 3
237 CALL dbcsr_set(c_mat(i)%matrix, zero)
238 END DO
239 DO m = 1, 6
240 CALL dbcsr_copy(tmp_mat, s_der(m + 3)%matrix)
241 DO j = 1, 2
242 IF (c_map_mat(m, j) == 0) cycle
243 CALL dbcsr_add(c_mat(c_map_mat(m, j))%matrix, tmp_mat, one, one)
244 END DO
245 END DO
246
247 DO m = 1, 3
248 CALL dbcsr_iterator_start(iter, c_mat(m)%matrix)
249 DO WHILE (dbcsr_iterator_blocks_left(iter))
250 CALL dbcsr_iterator_next_block(iter, row_atom, col_atom, block_values)
251 block_values = block_values*particle_set(row_atom)%v(m)
252 END DO
253 CALL dbcsr_iterator_stop(iter)
254 CALL dbcsr_filter(c_mat(m)%matrix, rtp%filter_eps)
255 END DO
256
257 CALL dbcsr_deallocate_matrix(tmp_mat)
258 CALL timestop(handle)
259 END SUBROUTINE calc_s_derivs
260
261! **************************************************************************************************
262!> \brief reads the restart file. At the moment only SCF (means only real)
263!> \param qs_env ...
264!> \author Florian Schiffmann (02.09)
265! **************************************************************************************************
266
267 SUBROUTINE get_restart_wfn(qs_env)
268 TYPE(qs_environment_type), POINTER :: qs_env
269
270 CHARACTER(LEN=default_path_length) :: file_name, project_name
271 INTEGER :: i, id_nr, im, ispin, ncol, nspin, &
272 output_unit, re, unit_nr
273 REAL(kind=dp) :: alpha, cs_pos
274 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
275 TYPE(cp_fm_type) :: mos_occ
276 TYPE(cp_fm_type), DIMENSION(:), POINTER :: mos_old
277 TYPE(cp_logger_type), POINTER :: logger
278 TYPE(dbcsr_distribution_type) :: dist
279 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: p_rmpv, rho_new, rho_old
280 TYPE(dft_control_type), POINTER :: dft_control
281 TYPE(mo_set_type), DIMENSION(:), POINTER :: mo_array
282 TYPE(mp_para_env_type), POINTER :: para_env
283 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
284 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
285 TYPE(qs_rho_type), POINTER :: rho_struct
286 TYPE(rt_prop_type), POINTER :: rtp
287 TYPE(section_vals_type), POINTER :: dft_section, input
288
289 NULLIFY (atomic_kind_set, qs_kind_set, mo_array, particle_set, rho_struct, para_env)
290
291 CALL get_qs_env(qs_env, &
292 qs_kind_set=qs_kind_set, &
293 atomic_kind_set=atomic_kind_set, &
294 particle_set=particle_set, &
295 mos=mo_array, &
296 input=input, &
297 rtp=rtp, &
298 dft_control=dft_control, &
299 rho=rho_struct, &
300 para_env=para_env)
301 logger => cp_get_default_logger()
302 output_unit = cp_logger_get_default_io_unit(logger)
303
304 IF (logger%para_env%is_source()) THEN
305 unit_nr = cp_logger_get_default_unit_nr(logger, local=.true.)
306 ELSE
307 unit_nr = -1
308 END IF
309
310 id_nr = 0
311 nspin = SIZE(mo_array)
312 CALL qs_rho_get(rho_struct, rho_ao=p_rmpv)
313 dft_section => section_vals_get_subs_vals(input, "DFT")
314 SELECT CASE (dft_control%rtp_control%initial_wfn)
315 CASE (use_restart_wfn)
316 CALL read_mo_set_from_restart(mo_array, qs_kind_set, particle_set, para_env, &
317 id_nr=id_nr, multiplicity=dft_control%multiplicity, &
318 dft_section=dft_section)
319 CALL set_uniform_occupation_mo_array(mo_array, nspin)
320
321 IF (dft_control%rtp_control%apply_wfn_mix_init_restart) &
322 CALL wfn_mix(mo_array, particle_set, dft_section, qs_kind_set, para_env, output_unit, &
323 for_rtp=.true.)
324
325 DO ispin = 1, nspin
326 CALL calculate_density_matrix(mo_array(ispin), p_rmpv(ispin)%matrix)
327 END DO
328 IF (rtp%linear_scaling) THEN
329 CALL get_rtp(rtp=rtp, rho_old=rho_old, rho_new=rho_new)
330 DO ispin = 1, nspin
331 re = 2*ispin - 1
332 im = 2*ispin
333 CALL cp_fm_get_info(mo_array(ispin)%mo_coeff, ncol_global=ncol)
334 CALL cp_fm_create(mos_occ, &
335 matrix_struct=mo_array(ispin)%mo_coeff%matrix_struct, &
336 name="mos_occ")
337 CALL cp_fm_to_fm(mo_array(ispin)%mo_coeff, mos_occ)
338 IF (mo_array(ispin)%uniform_occupation) THEN
339 alpha = 3.0_dp - real(nspin, dp)
340 CALL cp_fm_column_scale(mos_occ, mo_array(ispin)%occupation_numbers/alpha)
341 CALL cp_dbcsr_plus_fm_fm_t(sparse_matrix=rho_old(re)%matrix, &
342 matrix_v=mos_occ, &
343 ncol=ncol, &
344 alpha=alpha, keep_sparsity=.false.)
345 ELSE
346 alpha = 1.0_dp
347 CALL cp_fm_column_scale(mos_occ, mo_array(ispin)%occupation_numbers/alpha)
348 CALL cp_dbcsr_plus_fm_fm_t(sparse_matrix=rho_old(re)%matrix, &
349 matrix_v=mo_array(ispin)%mo_coeff, &
350 matrix_g=mos_occ, &
351 ncol=ncol, &
352 alpha=alpha, keep_sparsity=.false.)
353 END IF
354 CALL dbcsr_filter(rho_old(re)%matrix, rtp%filter_eps)
355 CALL dbcsr_copy(rho_new(re)%matrix, rho_old(re)%matrix)
356 CALL cp_fm_release(mos_occ)
357 END DO
358 CALL calc_update_rho_sparse(qs_env)
359 ELSE
360 CALL get_rtp(rtp=rtp, mos_old=mos_old)
361 DO i = 1, SIZE(qs_env%mos)
362 CALL cp_fm_to_fm(mo_array(i)%mo_coeff, mos_old(2*i - 1))
363 CALL cp_fm_set_all(mos_old(2*i), zero, zero)
364 END DO
365 END IF
366 CASE (use_rt_restart)
367 IF (rtp%linear_scaling) THEN
368 CALL get_rtp(rtp=rtp, rho_old=rho_old, rho_new=rho_new)
369 project_name = logger%iter_info%project_name
370 DO ispin = 1, nspin
371 re = 2*ispin - 1
372 im = 2*ispin
373 WRITE (file_name, '(A,I0,A)') trim(project_name)//"_LS_DM_SPIN_RE", ispin, "_RESTART.dm"
374 CALL dbcsr_get_info(rho_old(re)%matrix, distribution=dist)
375 CALL dbcsr_binary_read(file_name, distribution=dist, matrix_new=rho_old(re)%matrix)
376 cs_pos = dbcsr_checksum(rho_old(re)%matrix, pos=.true.)
377 IF (unit_nr > 0) THEN
378 WRITE (unit_nr, '(T2,A,E20.8)') "Read restart DM "//trim(file_name)//" with checksum: ", cs_pos
379 END IF
380 WRITE (file_name, '(A,I0,A)') trim(project_name)//"_LS_DM_SPIN_IM", ispin, "_RESTART.dm"
381 CALL dbcsr_get_info(rho_old(im)%matrix, distribution=dist)
382 CALL dbcsr_binary_read(file_name, distribution=dist, matrix_new=rho_old(im)%matrix)
383 cs_pos = dbcsr_checksum(rho_old(im)%matrix, pos=.true.)
384 IF (unit_nr > 0) THEN
385 WRITE (unit_nr, '(T2,A,E20.8)') "Read restart DM "//trim(file_name)//" with checksum: ", cs_pos
386 END IF
387 END DO
388 DO i = 1, SIZE(rho_new)
389 CALL dbcsr_copy(rho_new(i)%matrix, rho_old(i)%matrix)
390 END DO
391 CALL calc_update_rho_sparse(qs_env)
392 ELSE
393 CALL get_rtp(rtp=rtp, mos_old=mos_old)
394 CALL read_rt_mos_from_restart(mo_array, mos_old, qs_kind_set, particle_set, para_env, &
395 id_nr, dft_control%multiplicity, dft_section)
396 CALL set_uniform_occupation_mo_array(mo_array, nspin)
397 DO ispin = 1, nspin
398 CALL calculate_density_matrix(mo_array(ispin), &
399 p_rmpv(ispin)%matrix)
400 END DO
401 END IF
402 END SELECT
403
404 END SUBROUTINE get_restart_wfn
405
406! **************************************************************************************************
407!> \brief Set mo_array(ispin)%uniform_occupation after a restart
408!> \param mo_array ...
409!> \param nspin ...
410!> \author Guillaume Le Breton (03.23)
411! **************************************************************************************************
412
413 SUBROUTINE set_uniform_occupation_mo_array(mo_array, nspin)
414
415 TYPE(mo_set_type), DIMENSION(:), POINTER :: mo_array
416 INTEGER :: nspin
417
418 INTEGER :: ispin, mo
419 LOGICAL :: is_uniform
420
421 DO ispin = 1, nspin
422 is_uniform = .true.
423 DO mo = 1, mo_array(ispin)%nmo
424 IF (mo_array(ispin)%occupation_numbers(mo) /= 0.0 .AND. &
425 mo_array(ispin)%occupation_numbers(mo) /= 1.0 .AND. &
426 mo_array(ispin)%occupation_numbers(mo) /= 2.0) &
427 is_uniform = .false.
428 END DO
429 mo_array(ispin)%uniform_occupation = is_uniform
430 END DO
431
432 END SUBROUTINE set_uniform_occupation_mo_array
433
434! **************************************************************************************************
435!> \brief calculates the density from the complex MOs and passes the density to
436!> qs_env.
437!> \param qs_env ...
438!> \author Florian Schiffmann (02.09)
439! **************************************************************************************************
440
441 SUBROUTINE calc_update_rho(qs_env)
442
443 TYPE(qs_environment_type), POINTER :: qs_env
444
445 CHARACTER(len=*), PARAMETER :: routinen = 'calc_update_rho'
446 REAL(kind=dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
447
448 INTEGER :: handle, i, im, ncol, re
449 REAL(kind=dp) :: alpha
450 TYPE(cp_fm_type) :: mos_occ
451 TYPE(cp_fm_type), DIMENSION(:), POINTER :: mos_new
452 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rho_ao, rho_ao_im
453 TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
454 TYPE(qs_ks_env_type), POINTER :: ks_env
455 TYPE(qs_rho_type), POINTER :: rho
456 TYPE(rt_prop_type), POINTER :: rtp
457
458 CALL timeset(routinen, handle)
459
460 NULLIFY (rho, ks_env, mos_new, rtp)
461 CALL get_qs_env(qs_env, &
462 ks_env=ks_env, &
463 rho=rho, &
464 rtp=rtp, &
465 mos=mos)
466 CALL get_rtp(rtp=rtp, mos_new=mos_new)
467 CALL qs_rho_get(rho_struct=rho, rho_ao=rho_ao)
468 DO i = 1, SIZE(mos_new)/2
469 re = 2*i - 1; im = 2*i
470 CALL dbcsr_set(rho_ao(i)%matrix, zero)
471 CALL cp_fm_get_info(mos_new(re), ncol_global=ncol)
472 CALL cp_fm_create(mos_occ, &
473 matrix_struct=mos(i)%mo_coeff%matrix_struct, &
474 name="mos_occ")
475 CALL cp_fm_to_fm(mos_new(re), mos_occ)
476 IF (mos(i)%uniform_occupation) THEN
477 alpha = 3*one - real(SIZE(mos_new)/2, dp)
478 CALL cp_fm_column_scale(mos_occ, mos(i)%occupation_numbers/alpha)
479 CALL cp_dbcsr_plus_fm_fm_t(sparse_matrix=rho_ao(i)%matrix, &
480 matrix_v=mos_occ, &
481 ncol=ncol, &
482 alpha=alpha)
483 ELSE
484 alpha = 1.0_dp
485 CALL cp_fm_column_scale(mos_occ, mos(i)%occupation_numbers/alpha)
486 CALL cp_dbcsr_plus_fm_fm_t(sparse_matrix=rho_ao(i)%matrix, &
487 matrix_v=mos_new(re), &
488 matrix_g=mos_occ, &
489 ncol=ncol, &
490 alpha=alpha)
491 END IF
492
493 ! It is actually complex conjugate but i*i=-1 therefore it must be added
494 CALL cp_fm_to_fm(mos_new(im), mos_occ)
495 IF (mos(i)%uniform_occupation) THEN
496 alpha = 3*one - real(SIZE(mos_new)/2, dp)
497 CALL cp_fm_column_scale(mos_occ, mos(i)%occupation_numbers/alpha)
498 CALL cp_dbcsr_plus_fm_fm_t(sparse_matrix=rho_ao(i)%matrix, &
499 matrix_v=mos_occ, &
500 ncol=ncol, &
501 alpha=alpha)
502 ELSE
503 alpha = 1.0_dp
504 CALL cp_fm_column_scale(mos_occ, mos(i)%occupation_numbers/alpha)
505 CALL cp_dbcsr_plus_fm_fm_t(sparse_matrix=rho_ao(i)%matrix, &
506 matrix_v=mos_new(im), &
507 matrix_g=mos_occ, &
508 ncol=ncol, &
509 alpha=alpha)
510 END IF
511 CALL cp_fm_release(mos_occ)
512 END DO
513 CALL qs_rho_update_rho(rho, qs_env)
514
515 IF (rtp%track_imag_density) THEN
516 CALL qs_rho_get(rho_struct=rho, rho_ao_im=rho_ao_im)
517 CALL calculate_p_imaginary(qs_env, rtp, rho_ao_im, keep_sparsity=.true.)
518 CALL qs_rho_set(rho, rho_ao_im=rho_ao_im)
519 END IF
520
521 CALL qs_ks_did_change(ks_env, rho_changed=.true.)
522
523 CALL timestop(handle)
524
525 END SUBROUTINE calc_update_rho
526
527! **************************************************************************************************
528!> \brief Copies the density matrix back into the qs_env%rho%rho_ao
529!> \param qs_env ...
530!> \author Samuel Andermatt (3.14)
531! **************************************************************************************************
532
533 SUBROUTINE calc_update_rho_sparse(qs_env)
534
535 TYPE(qs_environment_type), POINTER :: qs_env
536
537 CHARACTER(len=*), PARAMETER :: routinen = 'calc_update_rho_sparse'
538 REAL(kind=dp), PARAMETER :: zero = 0.0_dp
539
540 INTEGER :: handle, ispin
541 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rho_ao, rho_ao_im, rho_new
542 TYPE(dft_control_type), POINTER :: dft_control
543 TYPE(qs_ks_env_type), POINTER :: ks_env
544 TYPE(qs_rho_type), POINTER :: rho
545 TYPE(rt_prop_type), POINTER :: rtp
546 TYPE(rtp_control_type), POINTER :: rtp_control
547
548 NULLIFY (rho, ks_env, rtp, dft_control)
549 CALL timeset(routinen, handle)
550 CALL get_qs_env(qs_env, &
551 ks_env=ks_env, &
552 rho=rho, &
553 rtp=rtp, &
554 dft_control=dft_control)
555 rtp_control => dft_control%rtp_control
556 CALL get_rtp(rtp=rtp, rho_new=rho_new)
557 CALL qs_rho_get(rho_struct=rho, rho_ao=rho_ao)
558 IF (rtp%track_imag_density) CALL qs_rho_get(rho_struct=rho, rho_ao_im=rho_ao_im)
559 DO ispin = 1, SIZE(rho_ao)
560 CALL dbcsr_set(rho_ao(ispin)%matrix, zero)
561 CALL dbcsr_copy(rho_ao(ispin)%matrix, rho_new(ispin*2 - 1)%matrix, keep_sparsity=.true.)
562 IF (rtp%track_imag_density) THEN
563 CALL dbcsr_copy(rho_ao_im(ispin)%matrix, rho_new(ispin*2)%matrix, keep_sparsity=.true.)
564 END IF
565 END DO
566
567 CALL qs_rho_update_rho(rho, qs_env)
568 CALL qs_ks_did_change(ks_env, rho_changed=.true.)
569
570 CALL timestop(handle)
571
572 END SUBROUTINE calc_update_rho_sparse
573
574! **************************************************************************************************
575!> \brief ...
576!> \param qs_env ...
577!> \param rtp ...
578!> \param matrix_p_im ...
579!> \param keep_sparsity ...
580! **************************************************************************************************
581 SUBROUTINE calculate_p_imaginary(qs_env, rtp, matrix_p_im, keep_sparsity)
582 TYPE(qs_environment_type), POINTER :: qs_env
583 TYPE(rt_prop_type), POINTER :: rtp
584 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_p_im
585 LOGICAL, OPTIONAL :: keep_sparsity
586
587 INTEGER :: i, im, ncol, re
588 LOGICAL :: my_keep_sparsity
589 REAL(kind=dp) :: alpha
590 TYPE(cp_fm_type), DIMENSION(:), POINTER :: mos_new, mos_occ
591 TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
592
593 CALL get_rtp(rtp=rtp, mos_new=mos_new)
594
595 my_keep_sparsity = .false.
596 IF (PRESENT(keep_sparsity)) my_keep_sparsity = keep_sparsity
597 CALL get_qs_env(qs_env, mos=mos)
598 ALLOCATE (mos_occ(SIZE(mos)*2))
599
600 DO i = 1, SIZE(mos_new)/2
601 re = 2*i - 1; im = 2*i
602 alpha = 3.0_dp - real(SIZE(matrix_p_im), dp)
603 CALL cp_fm_create(mos_occ(re), &
604 matrix_struct=mos(i)%mo_coeff%matrix_struct, &
605 name="mos_occ")
606 CALL cp_fm_create(mos_occ(im), &
607 matrix_struct=mos(i)%mo_coeff%matrix_struct, &
608 name="mos_occ")
609 CALL dbcsr_set(matrix_p_im(i)%matrix, 0.0_dp)
610 CALL cp_fm_get_info(mos_new(re), ncol_global=ncol)
611 CALL cp_fm_to_fm(mos_new(re), mos_occ(re))
612 CALL cp_fm_column_scale(mos_occ(re), mos(i)%occupation_numbers/alpha)
613 CALL cp_fm_to_fm(mos_new(im), mos_occ(im))
614 CALL cp_fm_column_scale(mos_occ(im), mos(i)%occupation_numbers/alpha)
615 CALL cp_dbcsr_plus_fm_fm_t(sparse_matrix=matrix_p_im(i)%matrix, &
616 matrix_v=mos_occ(im), &
617 matrix_g=mos_occ(re), &
618 ncol=ncol, &
619 keep_sparsity=my_keep_sparsity, &
620 alpha=2.0_dp*alpha, &
621 symmetry_mode=-1)
622 END DO
623 CALL cp_fm_release(mos_occ)
624
625 END SUBROUTINE calculate_p_imaginary
626
627! **************************************************************************************************
628!> \brief ...
629!> \param qs_env ...
630!> \param rtp ...
631! **************************************************************************************************
632 SUBROUTINE write_rtp_mos_to_output_unit(qs_env, rtp)
633 TYPE(qs_environment_type), POINTER :: qs_env
634 TYPE(rt_prop_type), POINTER :: rtp
635
636 CHARACTER(len=*), PARAMETER :: routinen = 'write_rtp_mos_to_output_unit'
637
638 CHARACTER(LEN=10) :: spin
639 CHARACTER(LEN=2*default_string_length) :: name
640 INTEGER :: handle, i, ispin, nao, nelectron, nmo, &
641 nspins
642 LOGICAL :: print_eigvecs, print_mo_info
643 REAL(kind=dp) :: flexible_electron_count, maxocc, n_el_f
644 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
645 TYPE(cp_fm_type), DIMENSION(:), POINTER :: mos_new
646 TYPE(cp_logger_type), POINTER :: logger
647 TYPE(mo_set_type) :: mo_set_rtp
648 TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
649 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
650 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
651 TYPE(section_vals_type), POINTER :: dft_section, input
652
653 NULLIFY (atomic_kind_set, particle_set, qs_kind_set, input, mos, dft_section)
654
655 CALL timeset(routinen, handle)
656
657 CALL get_qs_env(qs_env, &
658 atomic_kind_set=atomic_kind_set, &
659 qs_kind_set=qs_kind_set, &
660 particle_set=particle_set, &
661 input=input, &
662 mos=mos)
663 ! Quick return, if no printout of MO information is requested
664 dft_section => section_vals_get_subs_vals(input, "DFT")
665 CALL section_vals_val_get(dft_section, "PRINT%MO%EIGENVECTORS", l_val=print_eigvecs)
666
667 NULLIFY (logger)
668 logger => cp_get_default_logger()
669 print_mo_info = (cp_print_key_should_output(logger%iter_info, &
670 dft_section, "PRINT%MO") /= 0) .OR. &
671 (qs_env%sim_step == 1)
672
673 IF ((.NOT. print_mo_info) .OR. (.NOT. print_eigvecs)) THEN
674 CALL timestop(handle)
675 RETURN
676 END IF
677
678 CALL get_rtp(rtp=rtp, mos_new=mos_new)
679
680 nspins = SIZE(mos_new)/2
681
682 DO ispin = 1, nspins
683 ! initiate mo_set
684 CALL get_mo_set(mo_set=mos(ispin), nao=nao, nmo=nmo, nelectron=nelectron, &
685 n_el_f=n_el_f, maxocc=maxocc, flexible_electron_count=flexible_electron_count)
686
687 CALL allocate_mo_set(mo_set_rtp, &
688 nao=nao, &
689 nmo=nmo, &
690 nelectron=nelectron, &
691 n_el_f=n_el_f, &
692 maxocc=maxocc, &
693 flexible_electron_count=flexible_electron_count)
694
695 WRITE (name, fmt="(A,I6)") "RTP MO SET, SPIN ", ispin
696 CALL init_mo_set(mo_set_rtp, fm_ref=mos_new(2*ispin - 1), name=name)
697
698 IF (nspins > 1) THEN
699 IF (ispin == 1) THEN
700 spin = "ALPHA SPIN"
701 ELSE
702 spin = "BETA SPIN"
703 END IF
704 ELSE
705 spin = ""
706 END IF
707
708 mo_set_rtp%occupation_numbers = mos(ispin)%occupation_numbers
709
710 !loop for real (odd) and imaginary (even) parts
711 DO i = 1, 0, -1
712 CALL cp_fm_to_fm(mos_new(2*ispin - i), mo_set_rtp%mo_coeff)
713 CALL write_mo_set_to_output_unit(mo_set_rtp, qs_kind_set, particle_set, &
714 dft_section, 4, 0, rtp=.true., spin=trim(spin), &
715 cpart=mod(i, 2), sim_step=qs_env%sim_step)
716 END DO
717
718 CALL deallocate_mo_set(mo_set_rtp)
719 END DO
720
721 CALL timestop(handle)
722
723 END SUBROUTINE write_rtp_mos_to_output_unit
724
725! **************************************************************************************************
726!> \brief Write the time dependent amplitude of the MOs in real grid.
727!> Very close to qs_scf_post_gpw/qs_scf_post_occ_cubes subroutine.
728!> \param qs_env ...
729!> \param rtp ...
730!> \author Guillaume Le Breton (11.22)
731! **************************************************************************************************
732 SUBROUTINE write_rtp_mo_cubes(qs_env, rtp)
733 TYPE(qs_environment_type), POINTER :: qs_env
734 TYPE(rt_prop_type), POINTER :: rtp
735
736 CHARACTER(LEN=*), PARAMETER :: routinen = 'write_rtp_mo_cubes'
737
738 CHARACTER(LEN=default_path_length) :: filename, my_pos_cube, title
739 INTEGER :: handle, homo, i, ir, ispin, ivector, &
740 n_rep, nhomo, nlist, nspins, &
741 rt_time_step, unit_nr
742 INTEGER, DIMENSION(:), POINTER :: list, list_index
743 LOGICAL :: append_cube, do_kpoints, mpi_io
744 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
745 TYPE(cell_type), POINTER :: cell
746 TYPE(cp_blacs_env_type), POINTER :: blacs_env
747 TYPE(cp_fm_type), DIMENSION(:), POINTER :: mos_new
748 TYPE(cp_fm_type), POINTER :: mo_coeff
749 TYPE(cp_logger_type), POINTER :: logger
750 TYPE(dft_control_type), POINTER :: dft_control
751 TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
752 TYPE(mp_para_env_type), POINTER :: para_env
753 TYPE(particle_list_type), POINTER :: particles
754 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
755 TYPE(pw_c1d_gs_type) :: wf_g
756 TYPE(pw_env_type), POINTER :: pw_env
757 TYPE(pw_pool_p_type), DIMENSION(:), POINTER :: pw_pools
758 TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
759 TYPE(pw_r3d_rs_type) :: density_r, wf_r
760 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
761 TYPE(qs_subsys_type), POINTER :: subsys
762 TYPE(section_vals_type), POINTER :: dft_section, input
763
764 CALL timeset(routinen, handle)
765
766 NULLIFY (logger, auxbas_pw_pool, pw_pools, pw_env)
767
768 ! Get all the info from qs:
769 CALL get_qs_env(qs_env, do_kpoints=do_kpoints, &
770 input=input)
771
772 ! Kill the run in the case of K points
773 IF (do_kpoints) THEN
774 cpabort("K points not handled yet for printing MO_CUBE")
775 END IF
776
777 dft_section => section_vals_get_subs_vals(input, "DFT")
778 logger => cp_get_default_logger()
779
780 ! Quick return if no print required
781 IF (.NOT. btest(cp_print_key_should_output(logger%iter_info, dft_section, &
782 "PRINT%MO_CUBES"), cp_p_file)) THEN
783 CALL timestop(handle)
784 RETURN
785 END IF
786
787 CALL get_qs_env(qs_env, atomic_kind_set=atomic_kind_set, &
788 mos=mos, &
789 blacs_env=blacs_env, &
790 qs_kind_set=qs_kind_set, &
791 pw_env=pw_env, &
792 subsys=subsys, &
793 para_env=para_env, &
794 particle_set=particle_set, &
795 dft_control=dft_control)
796 CALL qs_subsys_get(subsys, particles=particles)
797
798 nspins = dft_control%nspins
799 rt_time_step = qs_env%sim_step
800
801 ! Setup the grids needed to compute a wavefunction given a vector
802 CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool, &
803 pw_pools=pw_pools)
804 CALL auxbas_pw_pool%create_pw(wf_r)
805 CALL auxbas_pw_pool%create_pw(wf_g)
806 CALL auxbas_pw_pool%create_pw(density_r)
807 CALL get_rtp(rtp=rtp, mos_new=mos_new)
808
809 DO ispin = 1, nspins
810 CALL get_mo_set(mo_set=mos(ispin), homo=homo)
811
812 nhomo = section_get_ival(dft_section, "PRINT%MO_CUBES%NHOMO")
813 append_cube = section_get_lval(dft_section, "PRINT%MO_CUBES%APPEND")
814 my_pos_cube = "REWIND"
815 IF (append_cube) THEN
816 my_pos_cube = "APPEND"
817 END IF
818 CALL section_vals_val_get(dft_section, "PRINT%MO_CUBES%HOMO_LIST", n_rep_val=n_rep)
819 IF (n_rep > 0) THEN ! write the cubes of the list
820 nlist = 0
821 DO ir = 1, n_rep
822 NULLIFY (list)
823 CALL section_vals_val_get(dft_section, "PRINT%MO_CUBES%HOMO_LIST", i_rep_val=ir, &
824 i_vals=list)
825 IF (ASSOCIATED(list)) THEN
826 CALL reallocate(list_index, 1, nlist + SIZE(list))
827 DO i = 1, SIZE(list)
828 list_index(i + nlist) = list(i)
829 END DO
830 nlist = nlist + SIZE(list)
831 END IF
832 END DO
833 ELSE
834
835 IF (nhomo == -1) nhomo = homo
836 nlist = homo - max(1, homo - nhomo + 1) + 1
837 ALLOCATE (list_index(nlist))
838 DO i = 1, nlist
839 list_index(i) = max(1, homo - nhomo + 1) + i - 1
840 END DO
841 END IF
842 DO i = 1, nlist
843 ivector = list_index(i)
844 CALL get_qs_env(qs_env=qs_env, &
845 atomic_kind_set=atomic_kind_set, &
846 qs_kind_set=qs_kind_set, &
847 cell=cell, &
848 particle_set=particle_set, &
849 pw_env=pw_env)
850
851 ! density_r contains the density of the MOs
852 CALL pw_zero(density_r)
853 mo_coeff => mos_new(2*ispin - 1)!Real coeff
854 CALL calculate_wavefunction(mo_coeff, ivector, wf_r, wf_g, atomic_kind_set, qs_kind_set, &
855 cell, dft_control, particle_set, pw_env)
856 ! Adding the real part
857 CALL pw_multiply(density_r, wf_r, wf_r, 1.0_dp)
858
859 mo_coeff => mos_new(2*ispin) !Im coeff
860 CALL calculate_wavefunction(mo_coeff, ivector, wf_r, wf_g, atomic_kind_set, qs_kind_set, &
861 cell, dft_control, particle_set, pw_env)
862 ! Adding the im part
863 CALL pw_multiply(density_r, wf_r, wf_r, 1.0_dp)
864
865 WRITE (filename, '(a4,I5.5,a1,I1.1)') "WFN_", ivector, "_", ispin
866 mpi_io = .true.
867 unit_nr = cp_print_key_unit_nr(logger, input, "DFT%PRINT%MO_CUBES", extension=".cube", &
868 middle_name=trim(filename), file_position=my_pos_cube, log_filename=.false., &
869 mpi_io=mpi_io)
870 WRITE (title, *) "DENSITY ", ivector, " spin ", ispin, " i.e. HOMO - ", ivector - homo
871 CALL cp_pw_to_cube(density_r, unit_nr, title, particles=particles, &
872 stride=section_get_ivals(dft_section, "PRINT%MO_CUBES%STRIDE"), mpi_io=mpi_io)
873 CALL cp_print_key_finished_output(unit_nr, logger, input, "DFT%PRINT%MO_CUBES", mpi_io=mpi_io)
874 END DO
875 IF (ASSOCIATED(list_index)) DEALLOCATE (list_index)
876 END DO
877
878 ! Deallocate grids needed to compute wavefunctions
879 CALL auxbas_pw_pool%give_back_pw(wf_r)
880 CALL auxbas_pw_pool%give_back_pw(wf_g)
881 CALL auxbas_pw_pool%give_back_pw(density_r)
882
883 CALL timestop(handle)
884
885 END SUBROUTINE write_rtp_mo_cubes
886
887! **************************************************************************************************
888!> \brief Warn about unused sections of the print section - only implemented for some of the methods
889!> \param section Parent section
890!> \param subsection_name Name of the subsection which is not used even when explicitly included
891!> \param error_message Message to display as warning
892!> \author Stepan Marek (01.25)
893! **************************************************************************************************
894 SUBROUTINE warn_section_unused(section, subsection_name, error_message)
895 TYPE(section_vals_type), POINTER :: section
896 CHARACTER(len=*) :: subsection_name, error_message
897
898 LOGICAL :: explicit
899 TYPE(section_vals_type), POINTER :: target_section
900
901 target_section => section_vals_get_subs_vals(section, subsection_name)
902 CALL section_vals_get(target_section, explicit=explicit)
903 IF (explicit) cpwarn(error_message)
904 END SUBROUTINE warn_section_unused
905! **************************************************************************************************
906!> \brief Attempt to read the moments from a previously written file
907!> \param moments_section Print key section for the moments trace
908!> \param orig_start Index from which the original calculation started (given by input file)
909!> \param current_start The current start index from which the calculation continues
910!> \param moments Moment trace, where the read entries are overwritten
911!> \param times Optional time trace, where again the read entries are overwritten
912!> \param mom_read Whether moments were actually read
913!> \author Stepan Marek (11.25)
914! **************************************************************************************************
915 SUBROUTINE read_moments(moments_section, orig_start, current_start, moments, times, mom_read)
916 TYPE(section_vals_type), POINTER :: moments_section
917 INTEGER :: orig_start, current_start
918 COMPLEX(kind=dp), DIMENSION(:, :, :), POINTER :: moments
919 REAL(kind=dp), DIMENSION(:), OPTIONAL, POINTER :: times
920 LOGICAL, OPTIONAL :: mom_read
921
922 CHARACTER(len=14), DIMENSION(4) :: file_extensions
923 CHARACTER(len=default_path_length) :: save_name_im, save_name_re
924 INTEGER :: i, k, max_n, n, n_spin, unit_im, unit_re
925 REAL(kind=dp) :: time
926 REAL(kind=dp), DIMENSION(3) :: moments_im, moments_re
927 TYPE(cp_logger_type), POINTER :: logger
928
929! Whether at least some moments where extracted from the files
930
931 file_extensions(1) = "_SPIN_A_RE.dat"
932 file_extensions(2) = "_SPIN_A_IM.dat"
933 file_extensions(3) = "_SPIN_B_RE.dat"
934 file_extensions(4) = "_SPIN_B_IM.dat"
935
936 IF (PRESENT(mom_read)) mom_read = .false.
937
938 n_spin = SIZE(moments, 1)
939
940 logger => cp_get_default_logger()
941 max_n = current_start - orig_start + 1
942
943 DO i = 1, n_spin
944 ! Open the relevant files
945 save_name_re = cp_print_key_generate_filename(logger, moments_section, &
946 extension=file_extensions(2*i - 1), my_local=.false.)
947 save_name_im = cp_print_key_generate_filename(logger, moments_section, &
948 extension=file_extensions(2*i), my_local=.false.)
949 IF (file_exists(save_name_re) .AND. file_exists(save_name_im)) THEN
950 CALL open_file(save_name_re, file_status="OLD", file_form="FORMATTED", file_action="READ", &
951 unit_number=unit_re)
952 CALL open_file(save_name_im, file_status="OLD", file_form="FORMATTED", file_action="READ", &
953 unit_number=unit_im)
954 DO n = 1, max_n
955 READ (unit_re, '(E20.8E3,E20.8E3,E20.8E3,E20.8E3,E20.8E3)') time, &
956 moments_re(1), moments_re(2), moments_re(3)
957 READ (unit_im, '(E20.8E3,E20.8E3,E20.8E3,E20.8E3,E20.8E3)') time, &
958 moments_im(1), moments_im(2), moments_im(3)
959 DO k = 1, 3
960 moments(i, k, n) = cmplx(moments_re(k), moments_im(k), kind=dp)
961 END DO
962 IF (PRESENT(times)) THEN
963 times(n) = time/femtoseconds
964 END IF
965 END DO
966 IF (PRESENT(mom_read)) mom_read = .true.
967 CALL close_file(unit_re)
968 CALL close_file(unit_im)
969 ELSE
970 moments(i, :, 1:max_n) = cmplx(0.0, 0.0, kind=dp)
971 cpwarn("Could not read at least one moments file - missing trace is set to zero.")
972 END IF
973 END DO
974 END SUBROUTINE read_moments
975
976! **************************************************************************************************
977!> \brief Recalculates the field for index range given by orig_start and i_start,
978!> with times taken from the times array
979!> \param fields Array, where relevant indices are recalculated
980!> \param times Array of input times
981!> \param orig_start original start (in case of restarts)
982!> \param i_start current start
983!> \param dft_control DFT parameters
984!> \date 11.2025
985!> \author Stepan Marek
986! **************************************************************************************************
987 SUBROUTINE recalculate_fields(fields, times, orig_start, i_start, dft_control)
988 COMPLEX(kind=dp), DIMENSION(:, :), POINTER :: fields
989 REAL(kind=dp), DIMENSION(:), POINTER :: times
990 INTEGER :: orig_start, i_start
991 TYPE(dft_control_type), POINTER :: dft_control
992
993 INTEGER :: i, max_n
994 REAL(kind=dp), DIMENSION(3) :: field
995
996 max_n = i_start - orig_start + 1
997
998 IF (dft_control%rtp_control%apply_delta_pulse .OR. &
999 dft_control%rtp_control%apply_delta_pulse_mag) THEN
1000 fields(:, 1:max_n) = 0.0_dp
1001 ELSE
1002 DO i = 1, max_n
1003 CALL make_field(dft_control, field, i + orig_start - 1, times(i))
1004 fields(:, i) = field(:)
1005 END DO
1006 END IF
1007 END SUBROUTINE recalculate_fields
1008
1009END MODULE rt_propagation_utils
Define the atomic kind types and their sub types.
Handles all functions related to the CELL.
Definition cell_types.F:15
methods related to the blacs parallel environment
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_scale(matrix, alpha_scalar)
...
subroutine, public dbcsr_deallocate_matrix(matrix)
...
subroutine, public dbcsr_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset)
...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_desymmetrize(matrix_a, matrix_b)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_get_info(matrix, nblkrows_total, nblkcols_total, nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, nfullrows_local, nfullcols_local, my_prow, my_pcol, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, matrix_type, group)
...
subroutine, public dbcsr_filter(matrix, eps)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_binary_read(filepath, distribution, matrix_new)
...
subroutine, public dbcsr_add(matrix_a, matrix_b, alpha_scalar, beta_scalar)
...
real(kind=dp) function, public dbcsr_checksum(matrix, pos)
Calculates the checksum of a DBCSR matrix.
DBCSR operations in CP2K.
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,...
Utility routines to open and close files. Tracking of preconnections.
Definition cp_files.F:16
subroutine, public open_file(file_name, file_status, file_form, file_action, file_position, file_pad, unit_number, debug, skip_get_unit_number, file_access)
Opens the requested file using a free unit number.
Definition cp_files.F:311
subroutine, public close_file(unit_number, file_status, keep_preconnection)
Close an open file given by its logical unit number. Optionally, keep the file and unit preconnected.
Definition cp_files.F:122
logical function, public file_exists(file_name)
Checks if file exists, considering also the file discovery mechanism.
Definition cp_files.F:504
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_column_scale(matrixa, scaling)
scales column i of matrix a with scaling(i)
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_set_all(matrix, alpha, beta)
set all elements of a matrix to the same value, and optionally the diagonal to a different one
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
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...
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)
...
character(len=default_path_length) function, public cp_print_key_generate_filename(logger, print_key, middle_name, extension, my_local)
Utility function that returns a unit number to write the print key. Might open a file with a unique f...
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
subroutine, public make_field(dft_control, field, sim_step, sim_time)
computes the amplitude of the efield within a given envelop
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public use_restart_wfn
integer, parameter, public use_rt_restart
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)
...
integer function, public section_get_ival(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
subroutine, public section_vals_get(section_vals, ref_count, n_repetition, n_subs_vals_rep, section, explicit)
returns various attributes about the section_vals
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
logical function, public section_get_lval(section_vals, keyword_name)
...
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
integer, parameter, public default_path_length
Definition kinds.F:58
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
Definition list.F:24
Definition of mathematical constants and functions.
real(kind=dp), parameter, public one
real(kind=dp), parameter, public zero
Utility routines for the memory handling.
Interface to the message passing library MPI.
Provides Cartesian and spherical orbital pointers and indices.
integer, dimension(:), allocatable, public ncoset
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
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
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Calculate the plane wave density by collocating the primitive Gaussian functions (pgf).
subroutine, public calculate_wavefunction(mo_vectors, ivector, rho, rho_gspace, atomic_kind_set, qs_kind_set, cell, dft_control, particle_set, pw_env, basis_type)
maps a given wavefunction on the grid
collects routines that calculate density matrices
Calculation of Overlap and Hamiltonian matrices in DFTB.
subroutine, public build_dftb_overlap(qs_env, nderivative, matrix_s)
...
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, 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, 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 qs_ks_did_change(ks_env, s_mstruct_changed, rho_changed, potential_changed, full_reset)
tells that some of the things relevant to the ks calculation did change. has to be called when change...
Definition and initialisation of the mo data type.
Definition qs_mo_io.F:21
subroutine, public read_mo_set_from_restart(mo_array, qs_kind_set, particle_set, para_env, id_nr, multiplicity, dft_section, natom_mismatch, cdft, out_unit)
...
Definition qs_mo_io.F:516
subroutine, public read_rt_mos_from_restart(mo_array, rt_mos, qs_kind_set, particle_set, para_env, id_nr, multiplicity, dft_section)
...
Definition qs_mo_io.F:609
subroutine, public write_mo_set_to_output_unit(mo_set, qs_kind_set, particle_set, dft_section, before, kpoint, final_mos, spin, solver_method, rtp, cpart, sim_step, umo_set)
Write MO information to output file (eigenvalues, occupation numbers, coefficients)
Definition qs_mo_io.F:1011
Definition and initialisation of the mo data type.
Definition qs_mo_types.F:22
subroutine, public allocate_mo_set(mo_set, nao, nmo, nelectron, n_el_f, maxocc, flexible_electron_count)
Allocates a mo set and partially initializes it (nao,nmo,nelectron, and flexible_electron_count are v...
subroutine, public deallocate_mo_set(mo_set)
Deallocate a wavefunction data structure.
subroutine, public get_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, mo_coeff, mo_coeff_b, uniform_occupation, kts, mu, flexible_electron_count)
Get the components of a MO set data structure.
subroutine, public init_mo_set(mo_set, fm_pool, fm_ref, fm_struct, name)
initializes an allocated mo_set. eigenvalues, mo_coeff, occupation_numbers are valid only after this ...
Define the neighbor list data types and the corresponding functionality.
Calculation of overlap matrix, its derivatives and forces.
Definition qs_overlap.F:19
subroutine, public build_overlap_matrix(ks_env, matrix_s, matrixkp_s, matrix_name, nderivative, basis_type_a, basis_type_b, sab_nl, calculate_forces, matrix_p, matrixkp_p)
Calculation of the overlap matrix over Cartesian Gaussian functions.
Definition qs_overlap.F:120
methods of the rho structure (defined in qs_rho_types)
subroutine, public qs_rho_update_rho(rho_struct, qs_env, rho_xc_external, local_rho_set, task_list_external, task_list_external_soft, pw_env_external, para_env_external)
updates rho_r and rho_g to the rhorho_ao. if use_kinetic_energy_density also computes tau_r and tau_g...
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public qs_rho_set(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)
...
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...
Does all kind of post scf calculations for GPW/GAPW.
subroutine, public wfn_mix(mos, particle_set, dft_section, qs_kind_set, para_env, output_unit, unoccupied_orbs, scf_env, matrix_s, marked_states, for_rtp)
writes a new 'mixed' set of mos to restart file, without touching the current MOs
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)
...
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 recalculate_fields(fields, times, orig_start, i_start, dft_control)
Recalculates the field for index range given by orig_start and i_start, with times taken from the tim...
subroutine, public calc_s_derivs(qs_env)
Calculates dS/dR respectily the velocity weighted derivatves only needed for ehrenfest MD.
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 read_moments(moments_section, orig_start, current_start, moments, times, mom_read)
Attempt to read the moments from a previously written file.
subroutine, public calculate_p_imaginary(qs_env, rtp, matrix_p_im, keep_sparsity)
...
subroutine, public calc_update_rho(qs_env)
calculates the density from the complex MOs and passes the density to qs_env.
subroutine, public warn_section_unused(section, subsection_name, error_message)
Warn about unused sections of the print section - only implemented for some of the methods.
subroutine, public get_restart_wfn(qs_env)
reads the restart file. At the moment only SCF (means only real)
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
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
contained for different pw related things
to create arrays of pools
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Provides all information about a quickstep kind.
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...
keeps the density in various representations, keeping track of which ones are valid.