(git:e7e05ae)
qs_scf_initialization.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 Utility routines for qs_scf
10 ! **************************************************************************************************
12  USE cp_control_types, ONLY: dft_control_type
20  USE cp_fm_diag, ONLY: cp_fm_power
21  USE cp_fm_pool_types, ONLY: cp_fm_pool_p_type,&
26  cp_fm_struct_type
27  USE cp_fm_types, ONLY: cp_fm_create,&
30  cp_fm_to_fm,&
32  cp_fm_type
34  cp_logger_type,&
35  cp_to_string
36  USE cp_output_handling, ONLY: cp_p_file,&
40  USE dbcsr_api, ONLY: dbcsr_create,&
41  dbcsr_init_p,&
42  dbcsr_p_type,&
43  dbcsr_type,&
44  dbcsr_type_no_symmetry,&
45  dbcsr_type_real_default
46  USE input_constants, ONLY: &
52  section_vals_type,&
54  USE kinds, ONLY: dp
55  USE kpoint_types, ONLY: kpoint_type
56  USE message_passing, ONLY: mp_para_env_type
57  USE pw_types, ONLY: pw_c1d_gs_type
66  USE qs_environment_types, ONLY: get_qs_env,&
67  qs_environment_type,&
74  USE qs_fb_env_types, ONLY: fb_env_create,&
77  USE qs_kind_types, ONLY: get_qs_kind,&
78  qs_kind_type,&
80  USE qs_ks_types, ONLY: qs_ks_did_change
81  USE qs_matrix_pools, ONLY: mpools_get
85  USE qs_mo_occupation, ONLY: set_mo_occupation
86  USE qs_mo_types, ONLY: get_mo_set,&
87  init_mo_set,&
88  mo_set_type
92  USE qs_rho_atom_types, ONLY: rho_atom_type
95  USE qs_rho_types, ONLY: qs_rho_create,&
96  qs_rho_get,&
97  qs_rho_type
101  USE qs_scf_types, ONLY: &
108  wfi_update
109  USE scf_control_types, ONLY: scf_control_type
110  USE xas_env_types, ONLY: xas_environment_type
112 #include "./base/base_uses.f90"
113 
114  IMPLICIT NONE
115 
116  PRIVATE
117 
118  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_scf_initialization'
119 
121 
122 CONTAINS
123 
124 ! **************************************************************************************************
125 !> \brief initializes input parameters if needed or restores values from
126 !> previous runs to fill scf_env with the values required for scf
127 !> \param qs_env the qs_environment where to perform the scf procedure
128 !> \param scf_env ...
129 !> \param scf_control ...
130 !> \param scf_section ...
131 ! **************************************************************************************************
132  SUBROUTINE qs_scf_env_initialize(qs_env, scf_env, scf_control, scf_section)
133  TYPE(qs_environment_type), POINTER :: qs_env
134  TYPE(qs_scf_env_type), POINTER :: scf_env
135  TYPE(scf_control_type), OPTIONAL, POINTER :: scf_control
136  TYPE(section_vals_type), OPTIONAL, POINTER :: scf_section
137 
138  TYPE(dft_control_type), POINTER :: dft_control
139  TYPE(scf_control_type), POINTER :: my_scf_control
140  TYPE(section_vals_type), POINTER :: dft_section, input, my_scf_section
141 
142  CALL get_qs_env(qs_env, input=input, dft_control=dft_control)
143 
144  IF (PRESENT(scf_control)) THEN
145  my_scf_control => scf_control
146  ELSE
147  CALL get_qs_env(qs_env, scf_control=my_scf_control)
148  END IF
149 
150  dft_section => section_vals_get_subs_vals(input, "DFT")
151  IF (PRESENT(scf_section)) THEN
152  my_scf_section => scf_section
153  ELSE
154  my_scf_section => section_vals_get_subs_vals(dft_section, "SCF")
155  END IF
156 
157  CALL qs_scf_ensure_scf_env(qs_env, scf_env)
158 
159  CALL section_vals_val_get(my_scf_section, "CHOLESKY", i_val=scf_env%cholesky_method)
160 
161  CALL qs_scf_ensure_mos(qs_env)
162 
163  ! set flags for diagonalization
164  CALL qs_scf_ensure_diagonalization(scf_env, my_scf_section, qs_env, &
165  my_scf_control, qs_env%has_unit_metric)
166  ! set parameters for mixing/DIIS during scf
167  CALL qs_scf_ensure_mixing(my_scf_control, my_scf_section, scf_env, dft_control)
168 
169  CALL qs_scf_ensure_work_matrices(qs_env, scf_env)
170 
171  CALL qs_scf_ensure_mixing_store(qs_env, scf_env)
172 
173  ! Initialize outer loop variables: handle CDFT and regular outer loop separately
174  IF (dft_control%qs_control%cdft) THEN
175  CALL qs_scf_ensure_cdft_loop_vars(qs_env, scf_env, dft_control, &
176  scf_control=my_scf_control)
177  ELSE
178  CALL qs_scf_ensure_outer_loop_vars(scf_env, my_scf_control)
179  END IF
180 
181  CALL init_scf_run(scf_env, qs_env, my_scf_section, my_scf_control)
182 
183  END SUBROUTINE qs_scf_env_initialize
184 
185 ! **************************************************************************************************
186 !> \brief initializes input parameters if needed for non-scf calclulations using diagonalization
187 !> \param qs_env the qs_environment where to perform the scf procedure
188 !> \param scf_env ...
189 ! **************************************************************************************************
190  SUBROUTINE qs_scf_env_init_basic(qs_env, scf_env)
191  TYPE(qs_environment_type), POINTER :: qs_env
192  TYPE(qs_scf_env_type), POINTER :: scf_env
193 
194  TYPE(dft_control_type), POINTER :: dft_control
195  TYPE(scf_control_type), POINTER :: scf_control
196  TYPE(section_vals_type), POINTER :: dft_section, input, scf_section
197 
198  CALL get_qs_env(qs_env, input=input, dft_control=dft_control)
199 
200  CALL get_qs_env(qs_env, scf_control=scf_control)
201  dft_section => section_vals_get_subs_vals(input, "DFT")
202  scf_section => section_vals_get_subs_vals(dft_section, "SCF")
203 
204  CALL qs_scf_ensure_scf_env(qs_env, scf_env)
205 
206  CALL section_vals_val_get(scf_section, "CHOLESKY", i_val=scf_env%cholesky_method)
207  scf_control%use_diag = .true.
208  scf_control%diagonalization%method = diag_standard
209 
210  CALL qs_scf_ensure_mos(qs_env)
211 
212  ! set flags for diagonalization
213  CALL qs_scf_ensure_diagonalization(scf_env, scf_section, qs_env, &
214  scf_control, qs_env%has_unit_metric)
215  CALL qs_scf_ensure_work_matrices(qs_env, scf_env)
216 
217  CALL init_scf_run(scf_env, qs_env, scf_section, scf_control)
218 
219  END SUBROUTINE qs_scf_env_init_basic
220 
221 ! **************************************************************************************************
222 !> \brief makes sure scf_env is allocated (might already be from before)
223 !> in case it is present the g-space mixing storage is reset
224 !> \param qs_env ...
225 !> \param scf_env ...
226 ! **************************************************************************************************
227  SUBROUTINE qs_scf_ensure_scf_env(qs_env, scf_env)
228  TYPE(qs_environment_type), POINTER :: qs_env
229  TYPE(qs_scf_env_type), POINTER :: scf_env
230 
231  TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER :: rho_g
232  TYPE(qs_rho_type), POINTER :: rho
233 
234  NULLIFY (rho_g)
235 
236  IF (.NOT. ASSOCIATED(scf_env)) THEN ! i.e. for MD this is associated on the second step (it so seems)
237  ALLOCATE (scf_env)
238  CALL scf_env_create(scf_env)
239  ELSE
240  ! Reallocate mixing store, if the g space grid (cell) has changed
241  SELECT CASE (scf_env%mixing_method)
243  IF (ASSOCIATED(scf_env%mixing_store)) THEN
244  ! The current mixing_store data structure does not allow for an unique
245  ! grid comparison, but the probability that the 1d lengths of the old and
246  ! the new grid are accidentily equal is rather low
247  CALL get_qs_env(qs_env, rho=rho)
248  CALL qs_rho_get(rho, rho_g=rho_g)
249  IF (ASSOCIATED(scf_env%mixing_store%rhoin)) THEN
250  IF (SIZE(rho_g(1)%pw_grid%gsq) /= SIZE(scf_env%mixing_store%rhoin(1)%cc)) THEN
251  CALL mixing_storage_release(scf_env%mixing_store)
252  DEALLOCATE (scf_env%mixing_store)
253  END IF
254  END IF
255  END IF
256  END SELECT
257  END IF
258 
259  END SUBROUTINE qs_scf_ensure_scf_env
260 
261 ! **************************************************************************************************
262 !> \brief performs allocation of outer SCF variables
263 !> \param scf_env the SCF environment which contains the outer SCF variables
264 !> \param scf_control control settings for the outer SCF loop
265 !> \param nvar (optional) set number of outer SCF variables externally if CDFT SCF is active
266 ! **************************************************************************************************
267  SUBROUTINE qs_scf_ensure_outer_loop_vars(scf_env, scf_control, nvar)
268  TYPE(qs_scf_env_type), POINTER :: scf_env
269  TYPE(scf_control_type), POINTER :: scf_control
270  INTEGER, OPTIONAL :: nvar
271 
272  INTEGER :: nhistory, nvariables
273 
274  IF (scf_control%outer_scf%have_scf) THEN
275  nhistory = scf_control%outer_scf%max_scf + 1
276  IF (PRESENT(nvar)) THEN
277  IF (nvar > 0) THEN
278  nvariables = nvar
279  ELSE
280  nvariables = outer_loop_variables_count(scf_control)
281  END IF
282  ELSE
283  nvariables = outer_loop_variables_count(scf_control)
284  END IF
285  ALLOCATE (scf_env%outer_scf%variables(nvariables, nhistory))
286  ALLOCATE (scf_env%outer_scf%count(nhistory))
287  scf_env%outer_scf%count = 0
288  ALLOCATE (scf_env%outer_scf%gradient(nvariables, nhistory))
289  ALLOCATE (scf_env%outer_scf%energy(nhistory))
290  END IF
291 
292  END SUBROUTINE qs_scf_ensure_outer_loop_vars
293 
294 ! **************************************************************************************************
295 !> \brief performs allocation of CDFT SCF variables
296 !> \param qs_env the qs_env where to perform the allocation
297 !> \param scf_env the currently active scf_env
298 !> \param dft_control the dft_control that holds the cdft_control type
299 !> \param scf_control the currently active scf_control
300 ! **************************************************************************************************
301  SUBROUTINE qs_scf_ensure_cdft_loop_vars(qs_env, scf_env, dft_control, scf_control)
302  TYPE(qs_environment_type), POINTER :: qs_env
303  TYPE(qs_scf_env_type), POINTER :: scf_env
304  TYPE(dft_control_type), POINTER :: dft_control
305  TYPE(scf_control_type), POINTER :: scf_control
306 
307  INTEGER :: nhistory, nvariables
308  LOGICAL :: do_kpoints
309  REAL(kind=dp), DIMENSION(:, :), POINTER :: gradient_history, outer_scf_history, &
310  variable_history
311 
312  NULLIFY (outer_scf_history, gradient_history, variable_history)
313  CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints)
314  ! Test kpoints
315  IF (do_kpoints) &
316  cpabort("CDFT calculation not possible with kpoints")
317  ! Check that OUTER_SCF section in DFT&SCF is active
318  ! This section must always be active to facilitate
319  ! switching of the CDFT and SCF control parameters in outer_loop_switch
320  IF (.NOT. scf_control%outer_scf%have_scf) &
321  cpabort("Section SCF&OUTER_SCF must be active for CDFT calculations.")
322  ! Initialize CDFT and outer_loop variables (constraint settings active in scf_control)
323  IF (dft_control%qs_control%cdft_control%constraint_control%have_scf) THEN
324  nhistory = dft_control%qs_control%cdft_control%constraint_control%max_scf + 1
325  IF (scf_control%outer_scf%type /= outer_scf_none) THEN
326  nvariables = outer_loop_variables_count(scf_control, &
327  dft_control%qs_control%cdft_control)
328  ELSE
329  ! First iteration: scf_control has not yet been updated
330  nvariables = SIZE(dft_control%qs_control%cdft_control%target)
331  END IF
332  ALLOCATE (dft_control%qs_control%cdft_control%constraint%variables(nvariables, nhistory))
333  ALLOCATE (dft_control%qs_control%cdft_control%constraint%count(nhistory))
334  dft_control%qs_control%cdft_control%constraint%count = 0
335  ALLOCATE (dft_control%qs_control%cdft_control%constraint%gradient(nvariables, nhistory))
336  ALLOCATE (dft_control%qs_control%cdft_control%constraint%energy(nhistory))
337  CALL qs_scf_ensure_outer_loop_vars(scf_env, scf_control, nvariables)
338  END IF
339  ! Executed only on first call (OT settings active in scf_control)
340  ! Save OT settings and constraint initial values in CDFT control
341  ! Then switch to constraint outer_scf settings for proper initialization of history
342  IF (scf_control%outer_scf%have_scf) THEN
343  IF (scf_control%outer_scf%type == outer_scf_none) THEN
344  dft_control%qs_control%cdft_control%ot_control%have_scf = .true.
345  dft_control%qs_control%cdft_control%ot_control%max_scf = scf_control%outer_scf%max_scf
346  dft_control%qs_control%cdft_control%ot_control%eps_scf = scf_control%outer_scf%eps_scf
347  dft_control%qs_control%cdft_control%ot_control%step_size = scf_control%outer_scf%step_size
348  dft_control%qs_control%cdft_control%ot_control%type = scf_control%outer_scf%type
349  dft_control%qs_control%cdft_control%ot_control%optimizer = scf_control%outer_scf%optimizer
350  dft_control%qs_control%cdft_control%ot_control%diis_buffer_length = scf_control%outer_scf%diis_buffer_length
351  dft_control%qs_control%cdft_control%ot_control%bisect_trust_count = scf_control%outer_scf%bisect_trust_count
352  CALL cdft_opt_type_copy(dft_control%qs_control%cdft_control%ot_control%cdft_opt_control, &
353  scf_control%outer_scf%cdft_opt_control)
354  ! In case constraint and OT extrapolation orders are different, make sure to use former
355  nvariables = SIZE(dft_control%qs_control%cdft_control%target)
356  IF (scf_control%outer_scf%extrapolation_order /= &
357  dft_control%qs_control%cdft_control%constraint_control%extrapolation_order &
358  .OR. nvariables /= 1) THEN
359  DEALLOCATE (qs_env%outer_scf_history)
360  DEALLOCATE (qs_env%gradient_history)
361  DEALLOCATE (qs_env%variable_history)
362  nhistory = dft_control%qs_control%cdft_control%constraint_control%extrapolation_order
363  ALLOCATE (outer_scf_history(nvariables, nhistory))
364  ALLOCATE (gradient_history(nvariables, 2))
365  gradient_history = 0.0_dp
366  ALLOCATE (variable_history(nvariables, 2))
367  variable_history = 0.0_dp
368  CALL set_qs_env(qs_env, outer_scf_history=outer_scf_history, &
369  gradient_history=gradient_history, variable_history=variable_history)
370  END IF
371  CALL outer_loop_switch(scf_env, scf_control, dft_control%qs_control%cdft_control, ot2cdft)
372  END IF
373  END IF
374 
375  END SUBROUTINE qs_scf_ensure_cdft_loop_vars
376 
377 ! **************************************************************************************************
378 !> \brief performs allocation of the mixing storage
379 !> \param qs_env ...
380 !> \param scf_env ...
381 ! **************************************************************************************************
382  SUBROUTINE qs_scf_ensure_mixing_store(qs_env, scf_env)
383  TYPE(qs_environment_type), POINTER :: qs_env
384  TYPE(qs_scf_env_type), POINTER :: scf_env
385 
386  TYPE(dft_control_type), POINTER :: dft_control
387 
388  NULLIFY (dft_control)
389  CALL get_qs_env(qs_env=qs_env, dft_control=dft_control)
390 
391  IF (scf_env%mixing_method > 0) THEN
392  CALL mixing_allocate(qs_env, scf_env%mixing_method, scf_env%p_mix_new, &
393  scf_env%p_delta, dft_control%nspins, &
394  scf_env%mixing_store)
395  ELSE
396  NULLIFY (scf_env%p_mix_new)
397  END IF
398 
399  END SUBROUTINE qs_scf_ensure_mixing_store
400 
401 ! **************************************************************************************************
402 !> \brief Performs allocation of the SCF work matrices
403 !> In case of kpoints we probably don't need most of these matrices,
404 !> maybe we have to initialize some matrices in the fm_pool in kpoints
405 !> \param qs_env ...
406 !> \param scf_env ...
407 ! **************************************************************************************************
408  SUBROUTINE qs_scf_ensure_work_matrices(qs_env, scf_env)
409 
410  TYPE(qs_environment_type), POINTER :: qs_env
411  TYPE(qs_scf_env_type), POINTER :: scf_env
412 
413  CHARACTER(LEN=*), PARAMETER :: routinen = 'qs_scf_ensure_work_matrices'
414 
415  INTEGER :: handle, is, nao, nrow_block, nw
416  LOGICAL :: do_kpoints
417  TYPE(cp_fm_pool_p_type), DIMENSION(:), POINTER :: ao_mo_fm_pools
418  TYPE(cp_fm_struct_type), POINTER :: ao_ao_fmstruct, ao_mo_fmstruct
419  TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_s
420  TYPE(dbcsr_type), POINTER :: ref_matrix
421  TYPE(dft_control_type), POINTER :: dft_control
422  TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
423  TYPE(scf_control_type), POINTER :: scf_control
424 
425  CALL timeset(routinen, handle)
426 
427  NULLIFY (ao_mo_fm_pools, ao_mo_fmstruct, ao_ao_fmstruct, dft_control, matrix_s, mos)
428 
429  CALL get_qs_env(qs_env=qs_env, &
430  dft_control=dft_control, &
431  matrix_s_kp=matrix_s, &
432  mos=mos, &
433  scf_control=scf_control, &
434  do_kpoints=do_kpoints)
435  CALL mpools_get(qs_env%mpools, ao_mo_fm_pools=ao_mo_fm_pools)
436 
437  ! create an ao_ao parallel matrix structure
438  ao_mo_fmstruct => fm_pool_get_el_struct(ao_mo_fm_pools(1)%pool)
439  CALL cp_fm_struct_get(ao_mo_fmstruct, nrow_block=nrow_block)
440  CALL get_mo_set(mos(1), nao=nao)
441  CALL cp_fm_struct_create(fmstruct=ao_ao_fmstruct, &
442  nrow_block=nrow_block, &
443  ncol_block=nrow_block, &
444  nrow_global=nao, &
445  ncol_global=nao, &
446  template_fmstruct=ao_mo_fmstruct)
447 
448  IF ((scf_env%method /= ot_method_nr) .AND. &
449  (scf_env%method /= block_davidson_diag_method_nr)) THEN
450  IF (.NOT. ASSOCIATED(scf_env%scf_work1)) THEN
451  nw = dft_control%nspins
452  IF (do_kpoints) nw = 4
453  ALLOCATE (scf_env%scf_work1(nw))
454  DO is = 1, SIZE(scf_env%scf_work1)
455  CALL cp_fm_create(scf_env%scf_work1(is), &
456  matrix_struct=ao_ao_fmstruct, &
457  name="SCF-WORK_MATRIX-1-"//trim(adjustl(cp_to_string(is))))
458  END DO
459  END IF
460  IF ((.NOT. ASSOCIATED(scf_env%ortho)) .AND. &
461  (scf_env%method /= ot_diag_method_nr) .AND. &
462  (scf_env%method /= special_diag_method_nr)) THEN
463  ! Initialize fm matrix to store the Cholesky decomposition
464  ALLOCATE (scf_env%ortho)
465  CALL cp_fm_create(scf_env%ortho, &
466  matrix_struct=ao_ao_fmstruct, &
467  name="SCF-ORTHO_MATRIX")
468  ! Initialize dbcsr matrix to store the Cholesky decomposition
469  IF (scf_env%cholesky_method == cholesky_dbcsr) THEN
470  ref_matrix => matrix_s(1, 1)%matrix
471  CALL dbcsr_init_p(scf_env%ortho_dbcsr)
472  CALL dbcsr_create(scf_env%ortho_dbcsr, template=ref_matrix, &
473  matrix_type=dbcsr_type_no_symmetry)
474  CALL dbcsr_init_p(scf_env%buf1_dbcsr)
475  CALL dbcsr_create(scf_env%buf1_dbcsr, template=ref_matrix, &
476  matrix_type=dbcsr_type_no_symmetry)
477  CALL dbcsr_init_p(scf_env%buf2_dbcsr)
478  CALL dbcsr_create(scf_env%buf2_dbcsr, template=ref_matrix, &
479  matrix_type=dbcsr_type_no_symmetry)
480  ELSE IF (scf_env%cholesky_method == cholesky_inverse .OR. &
481  (scf_control%level_shift /= 0.0_dp .AND. &
482  scf_env%cholesky_method == cholesky_off)) THEN
483  ALLOCATE (scf_env%ortho_m1)
484  CALL cp_fm_create(scf_env%ortho_m1, &
485  matrix_struct=ao_ao_fmstruct, &
486  name="SCF-ORTHO_MATRIX-1")
487  END IF
488  END IF
489  IF (.NOT. ASSOCIATED(scf_env%scf_work2)) THEN
490  ALLOCATE (scf_env%scf_work2)
491  CALL cp_fm_create(scf_env%scf_work2, &
492  matrix_struct=ao_ao_fmstruct, &
493  name="SCF-WORK_MATRIX-2")
494  END IF
495  END IF
496 
497  IF (dft_control%dft_plus_u) THEN
498  IF (dft_control%plus_u_method_id == plus_u_lowdin) THEN
499  IF (.NOT. ASSOCIATED(scf_env%scf_work2)) THEN
500  ALLOCATE (scf_env%scf_work2)
501  CALL cp_fm_create(scf_env%scf_work2, &
502  matrix_struct=ao_ao_fmstruct, &
503  name="SCF-WORK_MATRIX-2")
504  END IF
505  IF (.NOT. ASSOCIATED(scf_env%s_half)) THEN
506  ALLOCATE (scf_env%s_half)
507  CALL cp_fm_create(scf_env%s_half, &
508  matrix_struct=ao_ao_fmstruct, &
509  name="S**(1/2) MATRIX")
510  END IF
511  END IF
512  END IF
513 
514  IF (do_kpoints) THEN
515  IF (.NOT. ASSOCIATED(scf_env%scf_work1)) THEN
516  nw = 4
517  ALLOCATE (scf_env%scf_work1(nw))
518  DO is = 1, SIZE(scf_env%scf_work1)
519  CALL cp_fm_create(scf_env%scf_work1(is), &
520  matrix_struct=ao_ao_fmstruct, &
521  name="SCF-WORK_MATRIX-1-"//trim(adjustl(cp_to_string(is))))
522  END DO
523  END IF
524  END IF
525 
526  CALL cp_fm_struct_release(ao_ao_fmstruct)
527 
528  CALL timestop(handle)
529 
530  END SUBROUTINE qs_scf_ensure_work_matrices
531 
532 ! **************************************************************************************************
533 !> \brief performs allocation of the MO matrices
534 !> \param qs_env ...
535 ! **************************************************************************************************
536  SUBROUTINE qs_scf_ensure_mos(qs_env)
537  TYPE(qs_environment_type), POINTER :: qs_env
538 
539  CHARACTER(len=*), PARAMETER :: routinen = 'qs_scf_ensure_mos'
540 
541  INTEGER :: handle, ic, ik, ikk, ispin, nmo, nmo_mat
542  TYPE(cp_fm_pool_p_type), DIMENSION(:), POINTER :: ao_mo_fm_pools
543  TYPE(cp_fm_type), POINTER :: mo_coeff, mo_coeff_last
544  TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: mo_derivs
545  TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_s
546  TYPE(dbcsr_type), POINTER :: mo_coeff_b
547  TYPE(dft_control_type), POINTER :: dft_control
548  TYPE(kpoint_type), POINTER :: kpoints
549  TYPE(mo_set_type), DIMENSION(:), POINTER :: mos, mos_last_converged
550  TYPE(mo_set_type), DIMENSION(:, :), POINTER :: mos_k
551  TYPE(xas_environment_type), POINTER :: xas_env
552 
553  CALL timeset(routinen, handle)
554 
555  NULLIFY (ao_mo_fm_pools, dft_control, mos, xas_env, matrix_s, mos_last_converged, mo_coeff_last)
556 
557  CALL get_qs_env(qs_env=qs_env, &
558  dft_control=dft_control, &
559  mos=mos, &
560  matrix_s_kp=matrix_s, &
561  xas_env=xas_env)
562  CALL mpools_get(qs_env%mpools, ao_mo_fm_pools=ao_mo_fm_pools)
563  IF (dft_control%switch_surf_dip) THEN
564  CALL get_qs_env(qs_env, mos_last_converged=mos_last_converged)
565  END IF
566 
567  nmo_mat = dft_control%nspins
568  IF (dft_control%restricted) nmo_mat = 1 ! right now, there might be more mos than needed derivs
569 
570 ! *** finish initialization of the MOs ***
571  cpassert(ASSOCIATED(mos))
572  DO ispin = 1, SIZE(mos)
573  CALL get_mo_set(mos(ispin), mo_coeff=mo_coeff, mo_coeff_b=mo_coeff_b)
574  IF (.NOT. ASSOCIATED(mo_coeff)) THEN
575  CALL init_mo_set(mos(ispin), &
576  fm_pool=ao_mo_fm_pools(ispin)%pool, &
577  name="qs_env%mo"//trim(adjustl(cp_to_string(ispin))))
578  END IF
579  IF (.NOT. ASSOCIATED(mo_coeff_b)) THEN
580  CALL cp_fm_get_info(mos(ispin)%mo_coeff, ncol_global=nmo)
581  CALL dbcsr_init_p(mos(ispin)%mo_coeff_b)
582  CALL cp_dbcsr_m_by_n_from_row_template(mos(ispin)%mo_coeff_b, template=matrix_s(1, 1)%matrix, n=nmo, &
583  sym=dbcsr_type_no_symmetry)
584  END IF
585  END DO
586 ! *** get the mo_derivs OK if needed ***
587  IF (qs_env%requires_mo_derivs) THEN
588  CALL get_qs_env(qs_env, mo_derivs=mo_derivs)
589  IF (.NOT. ASSOCIATED(mo_derivs)) THEN
590  ALLOCATE (mo_derivs(nmo_mat))
591  DO ispin = 1, nmo_mat
592  CALL get_mo_set(mos(ispin), mo_coeff_b=mo_coeff_b)
593  NULLIFY (mo_derivs(ispin)%matrix)
594  CALL dbcsr_init_p(mo_derivs(ispin)%matrix)
595  CALL dbcsr_create(mo_derivs(ispin)%matrix, template=mo_coeff_b, &
596  name="mo_derivs", matrix_type=dbcsr_type_no_symmetry, &
597  nze=0, data_type=dbcsr_type_real_default)
598  END DO
599  CALL set_qs_env(qs_env, mo_derivs=mo_derivs)
600  END IF
601 
602  ELSE
603  ! nothing should be done
604  END IF
605 
606 ! *** finish initialization of the MOs for ADMM and derivs if needed ***
607  IF (dft_control%do_admm) THEN
608  IF (dft_control%restricted) cpabort("ROKS with ADMM is not implemented")
609  END IF
610 
611 ! *** finish initialization of mos_last_converged *** [SGh]
612  IF (dft_control%switch_surf_dip) THEN
613  cpassert(ASSOCIATED(mos_last_converged))
614  DO ispin = 1, SIZE(mos_last_converged)
615  CALL get_mo_set(mos_last_converged(ispin), mo_coeff=mo_coeff_last)
616  IF (.NOT. ASSOCIATED(mo_coeff_last)) THEN
617  CALL init_mo_set(mos_last_converged(ispin), &
618  fm_ref=mos(ispin)%mo_coeff, &
619  name="qs_env%mos_last_converged"//trim(adjustl(cp_to_string(ispin))))
620  END IF
621  END DO
622  END IF
623  ! kpoints: we have to initialize all the k-point MOs
624  CALL get_qs_env(qs_env=qs_env, kpoints=kpoints)
625  IF (kpoints%nkp /= 0) THEN
626  ! check for some incompatible options
627  IF (qs_env%requires_mo_derivs) THEN
628  cpwarn("MO derivative methods flag has been switched off for kpoint calculation")
629  ! we switch it off to make band structure calculations
630  ! possible for OT gamma point calculations
631  qs_env%requires_mo_derivs = .false.
632  END IF
633  IF (dft_control%do_xas_calculation) &
634  cpabort("No XAS implemented with kpoints")
635  DO ik = 1, SIZE(kpoints%kp_env)
636  CALL mpools_get(kpoints%mpools, ao_mo_fm_pools=ao_mo_fm_pools)
637  mos_k => kpoints%kp_env(ik)%kpoint_env%mos
638  ikk = kpoints%kp_range(1) + ik - 1
639  cpassert(ASSOCIATED(mos_k))
640  DO ispin = 1, SIZE(mos_k, 2)
641  DO ic = 1, SIZE(mos_k, 1)
642  CALL get_mo_set(mos_k(ic, ispin), mo_coeff=mo_coeff, mo_coeff_b=mo_coeff_b)
643  IF (.NOT. ASSOCIATED(mo_coeff)) THEN
644  CALL init_mo_set(mos_k(ic, ispin), &
645  fm_pool=ao_mo_fm_pools(ispin)%pool, &
646  name="kpoints_"//trim(adjustl(cp_to_string(ikk)))// &
647  "%mo"//trim(adjustl(cp_to_string(ispin))))
648  END IF
649  ! no sparse matrix representation of kpoint MO vectors
650  cpassert(.NOT. ASSOCIATED(mo_coeff_b))
651  END DO
652  END DO
653  END DO
654  END IF
655 
656  CALL timestop(handle)
657 
658  END SUBROUTINE qs_scf_ensure_mos
659 
660 ! **************************************************************************************************
661 !> \brief sets flag for mixing/DIIS during scf
662 !> \param scf_control ...
663 !> \param scf_section ...
664 !> \param scf_env ...
665 !> \param dft_control ...
666 ! **************************************************************************************************
667  SUBROUTINE qs_scf_ensure_mixing(scf_control, scf_section, scf_env, dft_control)
668  TYPE(scf_control_type), POINTER :: scf_control
669  TYPE(section_vals_type), POINTER :: scf_section
670  TYPE(qs_scf_env_type), POINTER :: scf_env
671  TYPE(dft_control_type), POINTER :: dft_control
672 
673  TYPE(section_vals_type), POINTER :: mixing_section
674 
675  SELECT CASE (scf_control%mixing_method)
676  CASE (no_mix)
677  scf_env%mixing_method = no_mixing_nr
678  scf_env%p_mix_alpha = 1.0_dp
680  scf_env%mixing_method = scf_control%mixing_method
681  mixing_section => section_vals_get_subs_vals(scf_section, "MIXING")
682  IF (.NOT. ASSOCIATED(scf_env%mixing_store)) THEN
683  ALLOCATE (scf_env%mixing_store)
684  CALL mixing_storage_create(scf_env%mixing_store, mixing_section, scf_env%mixing_method, &
685  dft_control%qs_control%cutoff)
686  END IF
687  CASE DEFAULT
688  cpabort("Unknown mixing method")
689  END SELECT
690 
691  ! Disable DIIS for OT and g-space density mixing methods
692  IF (scf_env%method == ot_method_nr) THEN
693  ! No mixing is used with OT
694  scf_env%mixing_method = no_mixing_nr
695  scf_env%p_mix_alpha = 1.0_dp
696  scf_env%skip_diis = .true.
697  END IF
698 
699  IF (scf_control%use_diag .AND. scf_env%mixing_method == no_mixing_nr) THEN
700  cpabort("Diagonalization procedures without mixing are not recommendable")
701  END IF
702 
703  IF (scf_env%mixing_method > direct_mixing_nr) THEN
704  scf_env%skip_diis = .true.
705  scf_env%p_mix_alpha = scf_env%mixing_store%alpha
706  IF (scf_env%mixing_store%beta == 0.0_dp) THEN
707  cpabort("Mixing employing the Kerker damping factor needs BETA /= 0.0")
708  END IF
709  END IF
710 
711  IF (scf_env%mixing_method == direct_mixing_nr) THEN
712  scf_env%p_mix_alpha = scf_env%mixing_store%alpha
713  IF (scf_control%eps_diis < scf_control%eps_scf) THEN
714  scf_env%skip_diis = .true.
715  cpwarn("the DIIS scheme is disabled, since EPS_DIIS < EPS_SCF")
716  END IF
717  END IF
718 
719  END SUBROUTINE qs_scf_ensure_mixing
720 
721 ! **************************************************************************************************
722 !> \brief sets flags for diagonalization and ensure that everything is
723 !> allocated
724 !> \param scf_env ...
725 !> \param scf_section ...
726 !> \param qs_env ...
727 !> \param scf_control ...
728 !> \param has_unit_metric ...
729 ! **************************************************************************************************
730  SUBROUTINE qs_scf_ensure_diagonalization(scf_env, scf_section, qs_env, &
731  scf_control, has_unit_metric)
732  TYPE(qs_scf_env_type), POINTER :: scf_env
733  TYPE(section_vals_type), POINTER :: scf_section
734  TYPE(qs_environment_type), POINTER :: qs_env
735  TYPE(scf_control_type), POINTER :: scf_control
736  LOGICAL :: has_unit_metric
737 
738  INTEGER :: ispin, nao, nmo
739  LOGICAL :: do_kpoints, need_coeff_b, not_se_or_tb
740  TYPE(cp_fm_type), POINTER :: mo_coeff
741  TYPE(dft_control_type), POINTER :: dft_control
742  TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
743 
744  CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints, dft_control=dft_control, mos=mos)
745  not_se_or_tb = .NOT. (dft_control%qs_control%dftb .OR. dft_control%qs_control%xtb .OR. &
746  dft_control%qs_control%semi_empirical)
747  need_coeff_b = .false.
748  scf_env%needs_ortho = .false.
749 
750  IF (scf_control%use_diag) THEN
751  ! sanity check whether combinations are allowed
752  IF (dft_control%restricted) &
753  cpabort("OT only for restricted (ROKS)")
754  SELECT CASE (scf_control%diagonalization%method)
756  IF (.NOT. not_se_or_tb) &
757  cpabort("TB and SE not possible with OT diagonalization")
758  END SELECT
759  SELECT CASE (scf_control%diagonalization%method)
760  ! Diagonalization: additional check whether we are in an orthonormal basis
761  CASE (diag_standard)
762  scf_env%method = general_diag_method_nr
763  scf_env%needs_ortho = (.NOT. has_unit_metric) .AND. (.NOT. do_kpoints)
764  IF (has_unit_metric) THEN
765  scf_env%method = special_diag_method_nr
766  END IF
767  ! OT Diagonalization: not possible with ROKS
768  CASE (diag_ot)
769  IF (dft_control%roks) &
770  cpabort("ROKS with OT diagonalization not possible")
771  IF (do_kpoints) &
772  cpabort("OT diagonalization not possible with kpoint calculations")
773  scf_env%method = ot_diag_method_nr
774  need_coeff_b = .true.
775  ! Block Krylov diagonlization: not possible with ROKS,
776  ! allocation of additional matrices is needed
777  CASE (diag_block_krylov)
778  IF (dft_control%roks) &
779  cpabort("ROKS with block PF diagonalization not possible")
780  IF (do_kpoints) &
781  cpabort("Block Krylov diagonalization not possible with kpoint calculations")
782  scf_env%method = block_krylov_diag_method_nr
783  scf_env%needs_ortho = .true.
784  IF (.NOT. ASSOCIATED(scf_env%krylov_space)) &
785  CALL krylov_space_create(scf_env%krylov_space, scf_section)
786  CALL krylov_space_allocate(scf_env%krylov_space, scf_control, mos)
787  ! Block davidson diagonlization: allocation of additional matrices is needed
788  CASE (diag_block_davidson)
789  IF (do_kpoints) &
790  cpabort("Block Davidson diagonalization not possible with kpoint calculations")
791  scf_env%method = block_davidson_diag_method_nr
792  IF (.NOT. ASSOCIATED(scf_env%block_davidson_env)) &
793  CALL block_davidson_env_create(scf_env%block_davidson_env, dft_control%nspins, &
794  scf_section)
795  DO ispin = 1, dft_control%nspins
796  CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff, nao=nao, nmo=nmo)
797  CALL block_davidson_allocate(scf_env%block_davidson_env(ispin), mo_coeff, nao, nmo)
798  END DO
799  need_coeff_b = .true.
800  ! Filter matrix diagonalisation method
801  CASE (diag_filter_matrix)
802  scf_env%method = filter_matrix_diag_method_nr
803  IF (.NOT. fb_env_has_data(scf_env%filter_matrix_env)) THEN
804  CALL fb_env_create(scf_env%filter_matrix_env)
805  END IF
806  CALL fb_env_read_input(scf_env%filter_matrix_env, scf_section)
807  CALL fb_env_build_rcut_auto(scf_env%filter_matrix_env, qs_env)
808  CALL fb_env_write_info(scf_env%filter_matrix_env, qs_env, scf_section)
809  CALL fb_distribution_build(scf_env%filter_matrix_env, qs_env, scf_section)
810  CALL fb_env_build_atomic_halos(scf_env%filter_matrix_env, qs_env, scf_section)
811  CASE DEFAULT
812  cpabort("Unknown diagonalization method")
813  END SELECT
814  ! Check if subspace diagonlization is requested: allocation of additional matrices is needed
815  IF (scf_control%do_diag_sub) THEN
816  scf_env%needs_ortho = .true.
817  IF (.NOT. ASSOCIATED(scf_env%subspace_env)) &
818  CALL diag_subspace_env_create(scf_env%subspace_env, scf_section, &
819  dft_control%qs_control%cutoff)
820  CALL diag_subspace_allocate(scf_env%subspace_env, qs_env, mos)
821  IF (do_kpoints) &
822  cpabort("No subspace diagonlization with kpoint calculation")
823  END IF
824  ! OT: check if OT is used instead of diagonlization. Not possible with added MOS at the moment
825  ELSEIF (scf_control%use_ot) THEN
826  scf_env%method = ot_method_nr
827  need_coeff_b = .true.
828  IF (sum(abs(scf_control%added_mos)) > 0) &
829  cpabort("OT with ADDED_MOS/=0 not implemented")
830  IF (dft_control%restricted .AND. dft_control%nspins .NE. 2) &
831  cpabort("nspin must be 2 for restricted (ROKS)")
832  IF (do_kpoints) &
833  cpabort("OT not possible with kpoint calculations")
834  ELSE
835  cpabort("OT or DIAGONALIZATION have to be set")
836  END IF
837  DO ispin = 1, dft_control%nspins
838  mos(ispin)%use_mo_coeff_b = need_coeff_b
839  END DO
840 
841  END SUBROUTINE qs_scf_ensure_diagonalization
842 
843 ! **************************************************************************************************
844 !> \brief performs those initialisations that need to be done only once
845 !> (e.g. that only depend on the atomic positions)
846 !> this will be called in scf
847 !> \param scf_env ...
848 !> \param qs_env ...
849 !> \param scf_section ...
850 !> \param scf_control ...
851 !> \par History
852 !> 03.2006 created [Joost VandeVondele]
853 ! **************************************************************************************************
854  SUBROUTINE init_scf_run(scf_env, qs_env, scf_section, scf_control)
855 
856  TYPE(qs_scf_env_type), POINTER :: scf_env
857  TYPE(qs_environment_type), POINTER :: qs_env
858  TYPE(section_vals_type), POINTER :: scf_section
859  TYPE(scf_control_type), POINTER :: scf_control
860 
861  CHARACTER(LEN=*), PARAMETER :: routinen = 'init_scf_run'
862 
863  INTEGER :: after, handle, ikind, iw, nao, ndep, &
864  output_unit
865  LOGICAL :: dft_plus_u_atom, do_kpoints, &
866  init_u_ramping_each_scf, omit_headers, &
867  s_minus_half_available
868  REAL(kind=dp) :: u_ramping
869  TYPE(cp_logger_type), POINTER :: logger
870  TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s
871  TYPE(dft_control_type), POINTER :: dft_control
872  TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
873  TYPE(mp_para_env_type), POINTER :: para_env
874  TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
875  TYPE(qs_kind_type), POINTER :: qs_kind
876  TYPE(qs_rho_type), POINTER :: rho
877  TYPE(xas_environment_type), POINTER :: xas_env
878 
879  CALL timeset(routinen, handle)
880 
881  NULLIFY (qs_kind_set, matrix_s, dft_control, mos, qs_kind, rho, xas_env)
882 
883  logger => cp_get_default_logger()
884 
885  cpassert(ASSOCIATED(scf_env))
886  cpassert(ASSOCIATED(qs_env))
887  NULLIFY (para_env)
888 
889  s_minus_half_available = .false.
890  CALL get_qs_env(qs_env, &
891  dft_control=dft_control, &
892  qs_kind_set=qs_kind_set, &
893  mos=mos, &
894  rho=rho, &
895  nelectron_total=scf_env%nelectron, &
896  do_kpoints=do_kpoints, &
897  para_env=para_env, &
898  xas_env=xas_env)
899 
900  output_unit = cp_print_key_unit_nr(logger, scf_section, "PRINT%PROGRAM_RUN_INFO", &
901  extension=".scfLog")
902  CALL qs_scf_initial_info(output_unit, mos, dft_control)
903  CALL cp_print_key_finished_output(output_unit, logger, scf_section, &
904  "PRINT%PROGRAM_RUN_INFO")
905 
906  ! calc ortho matrix
907  ndep = 0
908  IF (scf_env%needs_ortho) THEN
909  CALL get_qs_env(qs_env, matrix_s=matrix_s)
910  CALL copy_dbcsr_to_fm(matrix_s(1)%matrix, scf_env%ortho)
911  IF (scf_env%cholesky_method > cholesky_off) THEN
912  CALL cp_fm_cholesky_decompose(scf_env%ortho)
913  IF (scf_env%cholesky_method == cholesky_dbcsr) THEN
914  CALL cp_fm_triangular_invert(scf_env%ortho)
915  CALL cp_fm_set_all(scf_env%scf_work2, 0.0_dp)
916  CALL cp_fm_to_fm_triangular(scf_env%ortho, scf_env%scf_work2, "U")
917  CALL copy_fm_to_dbcsr(scf_env%scf_work2, scf_env%ortho_dbcsr)
918  ELSE IF (scf_env%cholesky_method == cholesky_inverse) THEN
919  CALL cp_fm_to_fm(scf_env%ortho, scf_env%ortho_m1)
920  CALL cp_fm_triangular_invert(scf_env%ortho_m1)
921  END IF
922  ELSE
923  CALL cp_fm_power(scf_env%ortho, scf_env%scf_work2, -0.5_dp, &
924  scf_control%eps_eigval, ndep)
925  IF (scf_control%level_shift /= 0.0_dp) THEN
926  CALL copy_dbcsr_to_fm(matrix_s(1)%matrix, scf_env%ortho_m1)
927  CALL cp_fm_power(scf_env%ortho_m1, scf_env%scf_work2, 0.5_dp, &
928  scf_control%eps_eigval, ndep)
929  END IF
930  s_minus_half_available = .true.
931  END IF
932 
933  IF (btest(cp_print_key_should_output(logger%iter_info, &
934  qs_env%input, "DFT%PRINT%AO_MATRICES/ORTHO"), cp_p_file)) THEN
935  iw = cp_print_key_unit_nr(logger, qs_env%input, "DFT%PRINT%AO_MATRICES/ORTHO", &
936  extension=".Log")
937  CALL section_vals_val_get(qs_env%input, "DFT%PRINT%AO_MATRICES%NDIGITS", i_val=after)
938  CALL section_vals_val_get(qs_env%input, "DFT%PRINT%AO_MATRICES%OMIT_HEADERS", l_val=omit_headers)
939  after = min(max(after, 1), 16)
940  CALL write_fm_with_basis_info(scf_env%ortho, 4, after, qs_env, &
941  para_env, output_unit=iw, omit_headers=omit_headers)
942  CALL cp_print_key_finished_output(iw, logger, qs_env%input, &
943  "DFT%PRINT%AO_MATRICES/ORTHO")
944  END IF
945  END IF
946 
947  CALL get_mo_set(mo_set=mos(1), nao=nao)
948 
949  ! DFT+U methods based on Lowdin charges need S^(1/2)
950  IF (dft_control%dft_plus_u) THEN
951  CALL get_qs_env(qs_env, matrix_s=matrix_s)
952  IF (dft_control%plus_u_method_id == plus_u_lowdin) THEN
953  IF (s_minus_half_available) THEN
954  CALL cp_dbcsr_sm_fm_multiply(matrix_s(1)%matrix, scf_env%ortho, scf_env%s_half, &
955  nao)
956  ELSE
957  CALL copy_dbcsr_to_fm(matrix_s(1)%matrix, scf_env%s_half)
958  CALL cp_fm_power(scf_env%s_half, scf_env%scf_work2, 0.5_dp, &
959  scf_control%eps_eigval, ndep)
960  END IF
961  END IF
962  DO ikind = 1, SIZE(qs_kind_set)
963  qs_kind => qs_kind_set(ikind)
964  CALL get_qs_kind(qs_kind=qs_kind, &
965  dft_plus_u_atom=dft_plus_u_atom, &
966  u_ramping=u_ramping, &
967  init_u_ramping_each_scf=init_u_ramping_each_scf)
968  IF (dft_plus_u_atom .AND. (u_ramping /= 0.0_dp)) THEN
969  IF (init_u_ramping_each_scf) THEN
970  CALL set_qs_kind(qs_kind=qs_kind, u_minus_j=0.0_dp)
971  END IF
972  END IF
973  END DO
974  END IF
975 
976  output_unit = cp_print_key_unit_nr(logger, scf_section, "PRINT%PROGRAM_RUN_INFO", &
977  extension=".scfLog")
978  IF (output_unit > 0) THEN
979  WRITE (unit=output_unit, fmt="(T2,A,T71,I10)") &
980  "Number of independent orbital functions:", nao - ndep
981  END IF
982  CALL cp_print_key_finished_output(output_unit, logger, scf_section, &
983  "PRINT%PROGRAM_RUN_INFO")
984 
985  ! extrapolate outer loop variables
986  IF (scf_control%outer_scf%have_scf) THEN
987  CALL outer_loop_extrapolate(qs_env)
988  END IF
989 
990  ! initializes rho and the mos
991  IF (ASSOCIATED(qs_env%xas_env)) THEN
992  ! if just optimized wfn, e.g. ground state
993  ! changes come from a perturbation, e.g., the occupation numbers
994  ! it could be generalized for other cases, at the moment used only for core level spectroscopy
995  ! initialize the density with the localized mos
996  CALL xas_initialize_rho(qs_env, scf_env, scf_control)
997  ELSE
998  CALL scf_env_initial_rho_setup(scf_env, qs_env=qs_env, &
999  scf_section=scf_section, scf_control=scf_control)
1000  END IF
1001 
1002  ! Frozen density approximation
1003  IF (ASSOCIATED(qs_env%wf_history)) THEN
1004  IF (qs_env%wf_history%interpolation_method_nr == wfi_frozen_method_nr) THEN
1005  IF (.NOT. ASSOCIATED(qs_env%wf_history%past_states(1)%snapshot)) THEN
1006  CALL wfi_update(qs_env%wf_history, qs_env=qs_env, dt=1.0_dp)
1007  ALLOCATE (qs_env%wf_history%past_states(1)%snapshot%rho_frozen)
1008  CALL qs_rho_create(qs_env%wf_history%past_states(1)%snapshot%rho_frozen)
1009  CALL duplicate_rho_type(rho_input=rho, &
1010  rho_output=qs_env%wf_history%past_states(1)%snapshot%rho_frozen, &
1011  qs_env=qs_env)
1012  END IF
1013  END IF
1014  END IF
1015 
1016  !image charge method, calculate image_matrix if required
1017  IF (qs_env%qmmm) THEN
1018  IF (qs_env%qmmm .AND. qs_env%qmmm_env_qm%image_charge) THEN
1019  CALL conditional_calc_image_matrix(qs_env=qs_env, &
1020  qmmm_env=qs_env%qmmm_env_qm)
1021  END IF
1022  END IF
1023 
1024  CALL timestop(handle)
1025 
1026  END SUBROUTINE init_scf_run
1027 
1028 ! **************************************************************************************************
1029 !> \brief Initializes rho and the mos, so that an scf cycle can start
1030 !> \param scf_env the scf env in which to do the scf
1031 !> \param qs_env the qs env the scf_env lives in
1032 !> \param scf_section ...
1033 !> \param scf_control ...
1034 !> \par History
1035 !> 02.2003 created [fawzi]
1036 !> \author fawzi
1037 ! **************************************************************************************************
1038  SUBROUTINE scf_env_initial_rho_setup(scf_env, qs_env, scf_section, scf_control)
1039  TYPE(qs_scf_env_type), POINTER :: scf_env
1040  TYPE(qs_environment_type), POINTER :: qs_env
1041  TYPE(section_vals_type), POINTER :: scf_section
1042  TYPE(scf_control_type), POINTER :: scf_control
1043 
1044  CHARACTER(len=*), PARAMETER :: routinen = 'scf_env_initial_rho_setup'
1045 
1046  INTEGER :: extrapolation_method_nr, handle, ispin, &
1047  nmo, output_unit
1048  LOGICAL :: orthogonal_wf
1049  TYPE(cp_fm_type), POINTER :: mo_coeff
1050  TYPE(cp_logger_type), POINTER :: logger
1051  TYPE(dft_control_type), POINTER :: dft_control
1052  TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
1053  TYPE(mp_para_env_type), POINTER :: para_env
1054  TYPE(qs_rho_type), POINTER :: rho
1055  TYPE(rho_atom_type), DIMENSION(:), POINTER :: rho_atom
1056 
1057  CALL timeset(routinen, handle)
1058  NULLIFY (mo_coeff, rho, dft_control, para_env, mos)
1059  logger => cp_get_default_logger()
1060  cpassert(ASSOCIATED(scf_env))
1061  cpassert(ASSOCIATED(qs_env))
1062 
1063  CALL get_qs_env(qs_env, &
1064  rho=rho, &
1065  mos=mos, &
1066  dft_control=dft_control, &
1067  para_env=para_env)
1068 
1069  extrapolation_method_nr = wfi_use_guess_method_nr
1070  IF (ASSOCIATED(qs_env%wf_history)) THEN
1071  CALL wfi_extrapolate(qs_env%wf_history, &
1072  qs_env=qs_env, dt=1.0_dp, &
1073  extrapolation_method_nr=extrapolation_method_nr, &
1074  orthogonal_wf=orthogonal_wf)
1075  ! wfi_use_guess_method_nr the wavefunctions are not yet initialized
1076  IF ((.NOT. orthogonal_wf) .AND. &
1077  (scf_env%method == ot_method_nr) .AND. &
1078  (.NOT. (extrapolation_method_nr == wfi_use_guess_method_nr))) THEN
1079  DO ispin = 1, SIZE(mos)
1080  CALL get_mo_set(mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
1081  CALL reorthogonalize_vectors(qs_env, v_matrix=mo_coeff, n_col=nmo)
1082  CALL set_mo_occupation(mo_set=mos(ispin), &
1083  smear=scf_control%smear)
1084  END DO
1085  END IF
1086  END IF
1087  output_unit = cp_print_key_unit_nr(logger, scf_section, "PRINT%PROGRAM_RUN_INFO", &
1088  extension=".scfLog")
1089  IF (output_unit > 0) THEN
1090  WRITE (unit=output_unit, fmt="(/,T2,A,I0)") &
1091  "Extrapolation method: "// &
1092  trim(wfi_get_method_label(extrapolation_method_nr))
1093  IF (extrapolation_method_nr == wfi_ps_method_nr) THEN
1094  WRITE (unit=output_unit, fmt="(T2,A,I0,A)") &
1095  "Extrapolation order: ", &
1096  max((min(qs_env%wf_history%memory_depth, qs_env%wf_history%snapshot_count) - 1), 0)
1097  END IF
1098  END IF
1099 
1100  CALL cp_print_key_finished_output(output_unit, logger, scf_section, &
1101  "PRINT%PROGRAM_RUN_INFO")
1102  IF (extrapolation_method_nr == wfi_use_guess_method_nr) THEN
1103  CALL calculate_first_density_matrix(scf_env=scf_env, qs_env=qs_env)
1104  CALL qs_rho_update_rho(rho, qs_env=qs_env)
1105  CALL qs_ks_did_change(qs_env%ks_env, rho_changed=.true.)
1106  END IF
1107 
1108  ! Some preparation for the mixing
1109  IF (scf_env%mixing_method > 1) THEN
1110  IF (dft_control%qs_control%gapw) THEN
1111  CALL get_qs_env(qs_env=qs_env, rho_atom_set=rho_atom)
1112  CALL mixing_init(scf_env%mixing_method, rho, scf_env%mixing_store, &
1113  para_env, rho_atom=rho_atom)
1114  ELSEIF (dft_control%qs_control%dftb .OR. dft_control%qs_control%xtb) THEN
1115  CALL charge_mixing_init(scf_env%mixing_store)
1116  ELSEIF (dft_control%qs_control%semi_empirical) THEN
1117  cpabort('SE Code not possible')
1118  ELSE
1119  CALL mixing_init(scf_env%mixing_method, rho, scf_env%mixing_store, &
1120  para_env)
1121  END IF
1122  END IF
1123 
1124  DO ispin = 1, SIZE(mos) !fm->dbcsr
1125  IF (mos(ispin)%use_mo_coeff_b) THEN
1126  CALL copy_fm_to_dbcsr(mos(ispin)%mo_coeff, &
1127  mos(ispin)%mo_coeff_b) !fm->dbcsr
1128  END IF
1129  END DO !fm->dbcsr
1130 
1131  CALL timestop(handle)
1132 
1133  END SUBROUTINE scf_env_initial_rho_setup
1134 
1135 END MODULE qs_scf_initialization
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
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
subroutine, public cp_dbcsr_m_by_n_from_row_template(matrix, template, n, sym, data_type)
Utility function to create dbcsr matrix, m x n matrix (n arbitrary) with the same processor grid and ...
subroutine, public copy_fm_to_dbcsr(fm, matrix, keep_sparsity)
Copy a BLACS matrix to a dbcsr matrix.
DBCSR output in CP2K.
subroutine, public write_fm_with_basis_info(blacs_matrix, before, after, qs_env, para_env, first_row, last_row, first_col, last_col, output_unit, omit_headers)
Print a spherical matrix of blacs type.
basic linear algebra operations for full matrices
subroutine, public cp_fm_triangular_invert(matrix_a, uplo_tr)
inverts a triangular matrix
various cholesky decomposition related routines
subroutine, public cp_fm_cholesky_decompose(matrix, n, info_out)
used to replace a symmetric positive def. matrix M with its cholesky decomposition U: M = U^T * U,...
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
Definition: cp_fm_diag.F:17
subroutine, public cp_fm_power(matrix, work, exponent, threshold, n_dependent, verbose, eigvals)
...
Definition: cp_fm_diag.F:896
pool for for elements that are retained and released
type(cp_fm_struct_type) function, pointer, public fm_pool_get_el_struct(pool)
returns the structure of the elements in this pool
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_get(fmstruct, para_env, context, descriptor, ncol_block, nrow_block, nrow_global, ncol_global, first_p_pos, row_indices, col_indices, nrow_local, ncol_local, nrow_locals, ncol_locals, local_leading_dimension)
returns the values of various attributes of the matrix structure
Definition: cp_fm_struct.F:409
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_to_fm_triangular(msource, mtarget, uplo)
copy just a triangular matrix
Definition: cp_fm_types.F:1428
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 ...
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)
...
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...
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public wfi_frozen_method_nr
integer, parameter, public diag_block_krylov
integer, parameter, public broy_mix
integer, parameter, public plus_u_lowdin
integer, parameter, public pulay_mix
integer, parameter, public cholesky_dbcsr
integer, parameter, public direct_p_mix
integer, parameter, public wfi_use_guess_method_nr
integer, parameter, public cholesky_off
integer, parameter, public no_mix
integer, parameter, public kerker_mix
integer, parameter, public cholesky_inverse
integer, parameter, public diag_ot
integer, parameter, public diag_filter_matrix
integer, parameter, public multisec_mix
integer, parameter, public diag_block_davidson
integer, parameter, public wfi_ps_method_nr
integer, parameter, public outer_scf_none
integer, parameter, public diag_standard
integer, parameter, public ot2cdft
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
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
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
Types and basic routines needed for a kpoint calculation.
Definition: kpoint_types.F:15
Interface to the message passing library MPI.
Routines for image charge calculation within QM/MM.
subroutine, public conditional_calc_image_matrix(qs_env, qmmm_env)
calculate image matrix T depending on constraints on image atoms in case coefficients are estimated n...
module that contains the algorithms to perform an itrative diagonalization by the block-Davidson appr...
subroutine, public block_davidson_allocate(bdav_env, mo_coeff, nao, nmo)
...
subroutine, public block_davidson_env_create(bdav_env, nspins, scf_section)
...
Control parameters for optimizers that work with CDFT constraints.
subroutine, public cdft_opt_type_copy(new, old)
copies settings between two CDFT optimizer control objects retaining both
module that contains the definitions of the scf types
subroutine, public mixing_storage_release(mixing_store)
releases a mixing_storage
integer, parameter, public no_mixing_nr
integer, parameter, public direct_mixing_nr
subroutine, public mixing_storage_create(mixing_store, mixing_section, mixing_method, ecut)
creates a mixing_storage
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.
subroutine, public set_qs_env(qs_env, super_cell, mos, qmmm, qmmm_periodic, ewald_env, ewald_pw, mpools, rho_external, external_vxc, mask, scf_control, rel_control, qs_charges, ks_env, ks_qmmm_env, wf_history, scf_env, active_space, input, oce, rho_atom_set, rho0_atom_set, rho0_mpole, run_rtp, rtp, rhoz_set, rhoz_tot, ecoul_1c, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, efield, linres_control, xas_env, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, ls_scf_env, do_transport, transport_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, gcp_env, mp2_env, bs_env, kg_env, force, kpoints, WannierCentres, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, rhs)
Set the QUICKSTEP environment.
subroutine, public fb_distribution_build(fb_env, qs_env, scf_section)
Build local atoms associated to filter matrix algorithm for each MPI process, trying to balance the l...
subroutine, public fb_env_read_input(fb_env, scf_section)
Read input sections for filter matrix method.
subroutine, public fb_env_build_atomic_halos(fb_env, qs_env, scf_section)
Builds an fb_atomic_halo_list object using information from fb_env.
subroutine, public fb_env_write_info(fb_env, qs_env, scf_section)
Write out parameters used for the filter matrix method to output.
subroutine, public fb_env_build_rcut_auto(fb_env, qs_env)
Automatically generate the cutoff radii of atoms used for constructing the atomic halos,...
logical function, public fb_env_has_data(fb_env)
Checks if a fb_env object is associated with an actual data content or not.
subroutine, public fb_env_create(fb_env)
creates an empty fb_env object
Routines to somehow generate an initial guess.
subroutine, public calculate_first_density_matrix(scf_env, qs_env)
can use a variety of methods to come up with an initial density matrix and optionally an initial wave...
Define the quickstep kind type and their sub types.
Definition: qs_kind_types.F:23
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_r3d_rs_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, U_of_dft_plus_u, J_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, J0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, floating, name, element_symbol, pao_basis_size, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
subroutine, public set_qs_kind(qs_kind, paw_atom, ghost, floating, hard_radius, hard0_radius, covalent_radius, vdw_radius, lmax_rho0, zeff, no_optimize, dispersion, u_minus_j, reltmat, dftb_parameter, xtb_parameter, elec_conf, pao_basis_size)
Set the components of an atomic kind data set.
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: qs_ks_types.F:919
wrapper for the pools of matrixes
subroutine, public mpools_get(mpools, ao_mo_fm_pools, ao_ao_fm_pools, mo_mo_fm_pools, ao_mosub_fm_pools, mosub_mosub_fm_pools, maxao_maxmo_fm_pool, maxao_maxao_fm_pool, maxmo_maxmo_fm_pool)
returns various attributes of the mpools (notably the pools contained in it)
elemental subroutine, public charge_mixing_init(mixing_store)
initialiation needed when charge mixing is used
subroutine, public mixing_init(mixing_method, rho, mixing_store, para_env, rho_atom)
initialiation needed when gspace mixing is used
subroutine, public mixing_allocate(qs_env, mixing_method, p_mix_new, p_delta, nspins, mixing_store)
allocation needed when density mixing is used
Set occupation of molecular orbitals.
Definition and initialisation of the mo data type.
Definition: qs_mo_types.F:22
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.
Definition: qs_mo_types.F:397
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 ...
Definition: qs_mo_types.F:245
Routines for performing an outer scf loop.
Definition: qs_outer_scf.F:14
subroutine, public outer_loop_switch(scf_env, scf_control, cdft_control, dir)
switch between two outer_scf envs stored in cdft_control
Definition: qs_outer_scf.F:539
subroutine, public outer_loop_extrapolate(qs_env)
uses the outer_scf_history to extrapolate new values for the variables and updates their value in qs_...
Definition: qs_outer_scf.F:443
integer function, public outer_loop_variables_count(scf_control, cdft_control)
returns the number of variables that is employed in the outer loop. with a CDFT constraint this value...
Definition: qs_outer_scf.F:69
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...
subroutine, public duplicate_rho_type(rho_input, rho_output, qs_env)
Duplicates a pointer physically.
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
subroutine, public qs_rho_create(rho)
Allocates a new instance of rho.
Definition: qs_rho_types.F:99
Different diagonalization schemes that can be used for the iterative solution of the eigenvalue probl...
subroutine, public diag_subspace_allocate(subspace_env, qs_env, mos)
...
Utility routines for qs_scf.
subroutine, public qs_scf_env_init_basic(qs_env, scf_env)
initializes input parameters if needed for non-scf calclulations using diagonalization
subroutine, public qs_scf_env_initialize(qs_env, scf_env, scf_control, scf_section)
initializes input parameters if needed or restores values from previous runs to fill scf_env with the...
module that contains the algorithms to perform an itrative diagonalization by the block-Lanczos appro...
subroutine, public krylov_space_allocate(krylov_space, scf_control, mos)
allocates matrices and vectros used in the construction of the krylov space and for the lanczos refin...
subroutine, public qs_scf_initial_info(output_unit, mos, dft_control)
writes basic information at the beginning of an scf run
module that contains the definitions of the scf types
Definition: qs_scf_types.F:14
integer, parameter, public ot_diag_method_nr
Definition: qs_scf_types.F:51
subroutine, public krylov_space_create(krylov_space, scf_section)
creates krylov space
Definition: qs_scf_types.F:352
subroutine, public diag_subspace_env_create(subspace_env, scf_section, ecut)
creates subspace-rotation environment
Definition: qs_scf_types.F:447
integer, parameter, public filter_matrix_diag_method_nr
Definition: qs_scf_types.F:51
integer, parameter, public block_davidson_diag_method_nr
Definition: qs_scf_types.F:51
integer, parameter, public ot_method_nr
Definition: qs_scf_types.F:51
subroutine, public scf_env_create(scf_env)
allocates and initialize an scf_env
Definition: qs_scf_types.F:132
integer, parameter, public special_diag_method_nr
Definition: qs_scf_types.F:51
integer, parameter, public block_krylov_diag_method_nr
Definition: qs_scf_types.F:51
integer, parameter, public general_diag_method_nr
Definition: qs_scf_types.F:51
Storage of past states of the qs_env. Methods to interpolate (or actually normally extrapolate) the n...
character(len=30) function, public wfi_get_method_label(method_nr)
returns a string describing the interpolation method
subroutine, public reorthogonalize_vectors(qs_env, v_matrix, n_col)
reorthogonalizes the mos
subroutine, public wfi_update(wf_history, qs_env, dt)
updates the snapshot buffer, taking a new snapshot
subroutine, public wfi_extrapolate(wf_history, qs_env, dt, extrapolation_method_nr, orthogonal_wf)
calculates the new starting state for the scf for the next wf optimization
parameters that control an scf iteration
define create destroy get and put information in xas_env to calculate the x-ray absorption spectra
Definition: xas_env_types.F:15
Initialize the XAS orbitals for specific core excitations Either the GS orbitals are used as initial ...
Definition: xas_restart.F:20
subroutine, public xas_initialize_rho(qs_env, scf_env, scf_control)
Once the mos and the occupation numbers are initialized the electronic density of the excited state c...
Definition: xas_restart.F:387