(git:419edc0)
Loading...
Searching...
No Matches
qs_ot_types.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 orbital transformations
10!> \par History
11!> Added Taylor expansion based computation of the matrix functions (01.2004)
12!> added additional rotation variables for non-equivalent occupied orbs (08.2004)
13!> \author Joost VandeVondele (06.2002)
14! **************************************************************************************************
17 weber2008,&
18 cite_reference
21 USE cp_dbcsr_api, ONLY: dbcsr_init_p,&
24 dbcsr_set,&
26 dbcsr_type_no_symmetry
33 USE input_constants, ONLY: &
42 USE kinds, ONLY: dp
46#include "./base/base_uses.f90"
47
48 IMPLICIT NONE
49
50 PRIVATE
51
52 PUBLIC :: qs_ot_type
53 PUBLIC :: qs_ot_settings_type
54 PUBLIC :: qs_ot_destroy
55 PUBLIC :: qs_ot_allocate
56 PUBLIC :: qs_ot_init
57 PUBLIC :: qs_ot_settings_init
58 PUBLIC :: ot_readwrite_input
59
60! **************************************************************************************************
61!> \brief notice, this variable needs to be copyable, needed for spins as e.g. in qs_ot_scf
62! **************************************************************************************************
64 LOGICAL :: do_rotation = .false., do_ener = .false.
65 LOGICAL :: ks = .false.
66 CHARACTER(LEN=4) :: ot_method = ""
67 CHARACTER(LEN=3) :: ot_algorithm = ""
68 CHARACTER(LEN=4) :: line_search_method = ""
69 CHARACTER(LEN=20) :: preconditioner_name = ""
70 INTEGER :: preconditioner_type = -1
71 INTEGER :: cholesky_type = -1
72 INTEGER :: ot_state = 0
73 CHARACTER(LEN=20) :: precond_solver_name = ""
74 INTEGER :: precond_solver_type = -1
75 LOGICAL :: safer_diis = .false.
76 REAL(kind=dp) :: ds_min = -1.0_dp
77 REAL(kind=dp) :: energy_gap = -1.0_dp
78 INTEGER :: diis_m = -1
79 INTEGER :: max_scf_diis = 0
80 REAL(kind=dp) :: gold_target = -1.0_dp
81 REAL(kind=dp) :: eps_taylor = -1.0_dp ! minimum accuracy of Taylor expansion
82 INTEGER :: max_taylor = -1 ! maximum order of Taylor expansion before switching to diagonalization
83 INTEGER :: irac_degree = -1 ! used to control the refinement polynomial degree
84 INTEGER :: max_irac = -1 ! maximum number of iteration for refinement
85 REAL(kind=dp) :: eps_irac = -1.0_dp ! target accuracy for refinement
86 REAL(kind=dp) :: eps_irac_quick_exit = -1.0_dp
87 REAL(kind=dp) :: eps_irac_filter_matrix = -1.0_dp
88 REAL(kind=dp) :: eps_irac_switch = -1.0_dp
89 LOGICAL :: on_the_fly_loc = .false.
90 CHARACTER(LEN=4) :: ortho_irac = ""
91 LOGICAL :: occupation_preconditioner = .false., add_nondiag_energy = .false.
92 REAL(kind=dp) :: nondiag_energy_strength = -1.0_dp
93 REAL(kind=dp) :: broyden_beta = -1.0_dp, broyden_gamma = -1.0_dp, broyden_sigma = -1.0_dp
94 REAL(kind=dp) :: broyden_eta = -1.0_dp, broyden_omega = -1.0_dp, broyden_sigma_decrease = -1.0_dp
95 REAL(kind=dp) :: broyden_sigma_min = -1.0_dp
96 LOGICAL :: broyden_forget_history = .false., broyden_adaptive_sigma = .false.
97 LOGICAL :: broyden_enable_flip = .false.
98 END TYPE qs_ot_settings_type
99
100! **************************************************************************************************
102 ! this sets the method to be used
104 LOGICAL :: restricted = .false.
105
106 ! first part of the variables, for occupied subspace invariant optimisation
107
108 ! add a preconditioner matrix. should be symmetric and positive definite
109 ! the type of this matrix might change in the future
110 TYPE(preconditioner_type), POINTER :: preconditioner => null()
111
112 ! these will/might change during iterations
113
114 ! OT / TOD
115 TYPE(dbcsr_type), POINTER :: matrix_p => null()
116 TYPE(dbcsr_type), POINTER :: matrix_r => null()
117 TYPE(dbcsr_type), POINTER :: matrix_sinp => null()
118 TYPE(dbcsr_type), POINTER :: matrix_cosp => null()
119 TYPE(dbcsr_type), POINTER :: matrix_sinp_b => null()
120 TYPE(dbcsr_type), POINTER :: matrix_cosp_b => null()
121 TYPE(dbcsr_type), POINTER :: matrix_buf1 => null()
122 TYPE(dbcsr_type), POINTER :: matrix_buf2 => null()
123 TYPE(dbcsr_type), POINTER :: matrix_buf3 => null()
124 TYPE(dbcsr_type), POINTER :: matrix_buf4 => null()
125 TYPE(dbcsr_type), POINTER :: matrix_os => null()
126 TYPE(dbcsr_type), POINTER :: matrix_buf1_ortho => null()
127 TYPE(dbcsr_type), POINTER :: matrix_buf2_ortho => null()
128
129 REAL(kind=dp), DIMENSION(:), POINTER :: evals => null()
130 REAL(kind=dp), DIMENSION(:), POINTER :: dum => null()
131
132 ! matrix os valid
133 LOGICAL :: os_valid = .false.
134
135 ! for efficient/parallel writing to the blacs_matrix
136 TYPE(mp_para_env_type), POINTER :: para_env => null()
137 TYPE(cp_blacs_env_type), POINTER :: blacs_env => null()
138
139 ! mo-like vectors
140 TYPE(dbcsr_type), POINTER :: matrix_c0 => null(), matrix_sc0 => null(), matrix_psc0 => null()
141
142 ! OT / IR
143 TYPE(dbcsr_type), POINTER :: buf1_k_k_nosym => null(), buf2_k_k_nosym => null(), &
144 buf3_k_k_nosym => null(), buf4_k_k_nosym => null(), &
145 buf1_k_k_sym => null(), buf2_k_k_sym => null(), &
146 buf3_k_k_sym => null(), buf4_k_k_sym => null(), &
147 p_k_k_sym => null(), buf1_n_k => null(), buf1_n_k_dp => null()
148
149 ! only here for the ease of programming. These will have to be supplied
150 ! explicitly at all times
151 TYPE(dbcsr_type), POINTER :: matrix_x => null(), matrix_sx => null(), matrix_gx => null()
152 TYPE(dbcsr_type), POINTER :: matrix_dx => null(), matrix_gx_old => null()
153
154 LOGICAL :: use_gx_old = .false., use_dx = .false.
155
156 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_h_e => null(), matrix_h_x => null()
157
158 REAL(kind=dp), DIMENSION(:, :), POINTER :: ls_diis => null()
159 REAL(kind=dp), DIMENSION(:, :), POINTER :: lss_diis => null()
160 REAL(kind=dp), DIMENSION(:), POINTER :: c_diis => null()
161 REAL(kind=dp), DIMENSION(:), POINTER :: c_broy => null()
162 REAL(kind=dp), DIMENSION(:), POINTER :: energy_h => null()
163 INTEGER, DIMENSION(:), POINTER :: ipivot => null()
164
165 REAL(kind=dp) :: ot_pos(53) = -1.0_dp, ot_energy(53) = -1.0_dp, ot_grad(53) = -1.0_dp ! HARD LIMIT FOR THE LS
166 INTEGER :: line_search_left = -1, line_search_right = -1, line_search_mid = -1
167 INTEGER :: line_search_count = -1
168 LOGICAL :: line_search_might_be_done = .false.
169 REAL(kind=dp) :: delta = -1.0_dp, gnorm = -1.0_dp, gnorm_old = -1.0_dp, etotal = -1.0_dp, gradient = -1.0_dp
170 LOGICAL :: energy_only = .false.
171 INTEGER :: diis_iter = -1
172 CHARACTER(LEN=8) :: ot_method_full = ""
173 INTEGER :: ot_count = -1
174 REAL(kind=dp) :: ds_min = -1.0_dp
175 REAL(kind=dp) :: broyden_adaptive_sigma = -1.0_dp
176
177 LOGICAL :: do_taylor = .false.
178 INTEGER :: taylor_order = -1
179 REAL(kind=dp) :: largest_eval_upper_bound = -1.0_dp
180
181 ! second part of the variables, if an explicit rotation is required as well
182 TYPE(dbcsr_type), POINTER :: rot_mat_u => null() ! rotation matrix
183 TYPE(dbcsr_type), POINTER :: rot_mat_x => null() ! antisymmetric matrix that parametrises rot_matrix_u
184 TYPE(dbcsr_type), POINTER :: rot_mat_dedu => null() ! derivative of the total energy wrt to u
185 TYPE(dbcsr_type), POINTER :: rot_mat_chc => null() ! for convencience, the matrix c^T H c
186
187 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rot_mat_h_e => null()
188 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: rot_mat_h_x => null()
189 TYPE(dbcsr_type), POINTER :: rot_mat_gx => null()
190 TYPE(dbcsr_type), POINTER :: rot_mat_gx_old => null()
191 TYPE(dbcsr_type), POINTER :: rot_mat_dx => null()
192
193 REAL(kind=dp), DIMENSION(:), POINTER :: rot_mat_evals => null()
194 TYPE(dbcsr_type), POINTER :: rot_mat_evec_re => null()
195 TYPE(dbcsr_type), POINTER :: rot_mat_evec_im => null()
196
197 ! third part of the variables, if we need to optimize orbital energies
198 REAL(kind=dp), POINTER, DIMENSION(:) :: ener_x => null()
199 REAL(kind=dp), POINTER, DIMENSION(:) :: ener_dx => null()
200 REAL(kind=dp), POINTER, DIMENSION(:) :: ener_gx => null()
201 REAL(kind=dp), POINTER, DIMENSION(:) :: ener_gx_old => null()
202 REAL(kind=dp), POINTER, DIMENSION(:, :) :: ener_h_e => null()
203 REAL(kind=dp), POINTER, DIMENSION(:, :) :: ener_h_x => null()
204 END TYPE qs_ot_type
205
206 CHARACTER(len=*), PARAMETER, PRIVATE :: modulen = 'qs_ot_types'
207
208CONTAINS
209
210! **************************************************************************************************
211!> \brief sets default values for the settings type
212!> \param settings ...
213!> \par History
214!> 10.2004 created [Joost VandeVondele]
215! **************************************************************************************************
216 SUBROUTINE qs_ot_settings_init(settings)
217 TYPE(qs_ot_settings_type) :: settings
218
219 settings%ot_method = "CG"
220 settings%ot_algorithm = "TOD"
221 settings%diis_m = 7
222 settings%preconditioner_name = "FULL_KINETIC"
223 settings%preconditioner_type = ot_precond_full_kinetic
224 settings%cholesky_type = cholesky_reduce
225 settings%precond_solver_name = "CHOLESKY_INVERSE"
226 settings%precond_solver_type = ot_precond_solver_inv_chol
227 settings%line_search_method = "2PNT"
228 settings%ds_min = 0.15_dp
229 settings%safer_diis = .true.
230 settings%energy_gap = 0.2_dp
231 settings%eps_taylor = 1.0e-16_dp
232 settings%max_taylor = 4
233 settings%gold_target = 0.01_dp
234 settings%do_rotation = .false.
235 settings%do_ener = .false.
236 settings%irac_degree = 4
237 settings%max_irac = 50
238 settings%max_scf_diis = 0
239 settings%eps_irac = 1.0e-10_dp
240 settings%eps_irac_quick_exit = 1.0e-5_dp
241 settings%eps_irac_switch = 1.0e-2
242 settings%eps_irac_filter_matrix = 0.0_dp
243 settings%on_the_fly_loc = .false.
244 settings%ortho_irac = "CHOL"
245 settings%ks = .true.
246 settings%occupation_preconditioner = .false.
247 settings%add_nondiag_energy = .false.
248 settings%nondiag_energy_strength = 0.0_dp
249
250 END SUBROUTINE qs_ot_settings_init
251
252! **************************************************************************************************
253!> \brief init matrices, needs c0 and sc0 so that c0*sc0=1
254!> \param qs_ot_env ...
255! **************************************************************************************************
256 SUBROUTINE qs_ot_init(qs_ot_env)
257 TYPE(qs_ot_type) :: qs_ot_env
258
259 qs_ot_env%OT_energy(:) = 0.0_dp
260 qs_ot_env%OT_pos(:) = 0.0_dp
261 qs_ot_env%OT_grad(:) = 0.0_dp
262 qs_ot_env%line_search_count = 0
263
264 qs_ot_env%energy_only = .false.
265 qs_ot_env%gnorm_old = 1.0_dp
266 qs_ot_env%diis_iter = 0
267 qs_ot_env%ds_min = qs_ot_env%settings%ds_min
268 qs_ot_env%os_valid = .false.
269
270 CALL dbcsr_set(qs_ot_env%matrix_gx, 0.0_dp)
271
272 IF (qs_ot_env%use_dx) &
273 CALL dbcsr_set(qs_ot_env%matrix_dx, 0.0_dp)
274
275 IF (qs_ot_env%use_gx_old) &
276 CALL dbcsr_set(qs_ot_env%matrix_gx_old, 0.0_dp)
277
278 IF (qs_ot_env%settings%do_rotation) THEN
279 CALL dbcsr_set(qs_ot_env%rot_mat_gx, 0.0_dp)
280 IF (qs_ot_env%use_dx) &
281 CALL dbcsr_set(qs_ot_env%rot_mat_dx, 0.0_dp)
282 IF (qs_ot_env%use_gx_old) &
283 CALL dbcsr_set(qs_ot_env%rot_mat_gx_old, 0.0_dp)
284 END IF
285 IF (qs_ot_env%settings%do_ener) THEN
286 qs_ot_env%ener_gx(:) = 0.0_dp
287 IF (qs_ot_env%use_dx) &
288 qs_ot_env%ener_dx(:) = 0.0_dp
289 IF (qs_ot_env%use_gx_old) &
290 qs_ot_env%ener_gx_old(:) = 0.0_dp
291 END IF
292
293 END SUBROUTINE qs_ot_init
294
295! **************************************************************************************************
296!> \brief allocates the data in qs_ot_env, for a calculation with fm_struct_ref
297!> ortho_k allows for specifying an additional orthogonal subspace (i.e. c will
298!> be kept orthogonal provided c0 was, used in qs_ot_eigensolver)
299!> \param qs_ot_env ...
300!> \param matrix_s ...
301!> \param fm_struct_ref ...
302!> \param ortho_k ...
303! **************************************************************************************************
304 SUBROUTINE qs_ot_allocate(qs_ot_env, matrix_s, fm_struct_ref, ortho_k)
305 TYPE(qs_ot_type) :: qs_ot_env
306 TYPE(dbcsr_type), POINTER :: matrix_s
307 TYPE(cp_fm_struct_type), POINTER :: fm_struct_ref
308 INTEGER, OPTIONAL :: ortho_k
309
310 INTEGER :: i, k, m_diis, my_ortho_k, n, ncoef
311 TYPE(cp_blacs_env_type), POINTER :: context
312 TYPE(mp_para_env_type), POINTER :: para_env
313
314 CALL cite_reference(vandevondele2003)
315
316 NULLIFY (qs_ot_env%preconditioner)
317 NULLIFY (qs_ot_env%matrix_psc0)
318 NULLIFY (qs_ot_env%para_env)
319 NULLIFY (qs_ot_env%blacs_env)
320
321 CALL cp_fm_struct_get(fm_struct_ref, nrow_global=n, ncol_global=k, &
322 para_env=para_env, context=context)
323
324 qs_ot_env%para_env => para_env
325 qs_ot_env%blacs_env => context
326 CALL para_env%retain()
327 CALL context%retain()
328
329 IF (PRESENT(ortho_k)) THEN
330 my_ortho_k = ortho_k
331 ELSE
332 my_ortho_k = k
333 END IF
334
335 m_diis = qs_ot_env%settings%diis_m
336
337 qs_ot_env%use_gx_old = .false.
338 qs_ot_env%use_dx = .false.
339
340 SELECT CASE (qs_ot_env%settings%ot_method)
341 CASE ("SD")
342 ! nothing
343 CASE ("CG")
344 qs_ot_env%use_gx_old = .true.
345 qs_ot_env%use_dx = .true.
346 CASE ("DIIS", "BROY")
347 IF (m_diis .LT. 1) cpabort("m_diis less than one")
348 CASE DEFAULT
349 cpabort("Unknown option")
350 END SELECT
351
352 IF (qs_ot_env%settings%ot_method .EQ. "DIIS" .OR. &
353 qs_ot_env%settings%ot_method .EQ. "BROY") THEN
354 ALLOCATE (qs_ot_env%ls_diis(m_diis + 1, m_diis + 1))
355 qs_ot_env%ls_diis = 0.0_dp
356 ALLOCATE (qs_ot_env%lss_diis(m_diis + 1, m_diis + 1))
357 ALLOCATE (qs_ot_env%c_diis(m_diis + 1))
358 ALLOCATE (qs_ot_env%c_broy(m_diis))
359 ALLOCATE (qs_ot_env%energy_h(m_diis))
360 ALLOCATE (qs_ot_env%ipivot(m_diis + 1))
361 END IF
362
363 ALLOCATE (qs_ot_env%evals(k))
364 ALLOCATE (qs_ot_env%dum(k))
365
366 NULLIFY (qs_ot_env%matrix_os)
367 NULLIFY (qs_ot_env%matrix_buf1_ortho)
368 NULLIFY (qs_ot_env%matrix_buf2_ortho)
369 NULLIFY (qs_ot_env%matrix_p)
370 NULLIFY (qs_ot_env%matrix_r)
371 NULLIFY (qs_ot_env%matrix_sinp)
372 NULLIFY (qs_ot_env%matrix_cosp)
373 NULLIFY (qs_ot_env%matrix_sinp_b)
374 NULLIFY (qs_ot_env%matrix_cosp_b)
375 NULLIFY (qs_ot_env%matrix_buf1)
376 NULLIFY (qs_ot_env%matrix_buf2)
377 NULLIFY (qs_ot_env%matrix_buf3)
378 NULLIFY (qs_ot_env%matrix_buf4)
379 NULLIFY (qs_ot_env%matrix_c0)
380 NULLIFY (qs_ot_env%matrix_sc0)
381 NULLIFY (qs_ot_env%matrix_x)
382 NULLIFY (qs_ot_env%matrix_sx)
383 NULLIFY (qs_ot_env%matrix_gx)
384 NULLIFY (qs_ot_env%matrix_gx_old)
385 NULLIFY (qs_ot_env%matrix_dx)
386 NULLIFY (qs_ot_env%buf1_k_k_nosym)
387 NULLIFY (qs_ot_env%buf2_k_k_nosym)
388 NULLIFY (qs_ot_env%buf3_k_k_nosym)
389 NULLIFY (qs_ot_env%buf4_k_k_nosym)
390 NULLIFY (qs_ot_env%buf1_k_k_sym)
391 NULLIFY (qs_ot_env%buf2_k_k_sym)
392 NULLIFY (qs_ot_env%buf3_k_k_sym)
393 NULLIFY (qs_ot_env%buf4_k_k_sym)
394 NULLIFY (qs_ot_env%buf1_n_k)
395 NULLIFY (qs_ot_env%buf1_n_k_dp)
396 NULLIFY (qs_ot_env%p_k_k_sym)
397
398 ! COMMON MATRICES
399 CALL dbcsr_init_p(qs_ot_env%matrix_c0)
400 CALL cp_dbcsr_m_by_n_from_row_template(qs_ot_env%matrix_c0, template=matrix_s, n=k, &
401 sym=dbcsr_type_no_symmetry)
402
403 CALL dbcsr_init_p(qs_ot_env%matrix_sc0)
404 CALL cp_dbcsr_m_by_n_from_row_template(qs_ot_env%matrix_sc0, template=matrix_s, n=my_ortho_k, &
405 sym=dbcsr_type_no_symmetry)
406
407 CALL dbcsr_init_p(qs_ot_env%matrix_x)
408 CALL cp_dbcsr_m_by_n_from_row_template(qs_ot_env%matrix_x, template=matrix_s, n=k, &
409 sym=dbcsr_type_no_symmetry)
410
411 CALL dbcsr_init_p(qs_ot_env%matrix_sx)
412 CALL cp_dbcsr_m_by_n_from_row_template(qs_ot_env%matrix_sx, template=matrix_s, n=k, &
413 sym=dbcsr_type_no_symmetry)
414
415 CALL dbcsr_init_p(qs_ot_env%matrix_gx)
416 CALL cp_dbcsr_m_by_n_from_row_template(qs_ot_env%matrix_gx, template=matrix_s, n=k, &
417 sym=dbcsr_type_no_symmetry)
418
419 IF (qs_ot_env%use_dx) THEN
420 CALL dbcsr_init_p(qs_ot_env%matrix_dx)
421 CALL cp_dbcsr_m_by_n_from_row_template(qs_ot_env%matrix_dx, template=matrix_s, n=k, &
422 sym=dbcsr_type_no_symmetry)
423 END IF
424
425 IF (qs_ot_env%use_gx_old) THEN
426 CALL dbcsr_init_p(qs_ot_env%matrix_gx_old)
427 CALL cp_dbcsr_m_by_n_from_row_template(qs_ot_env%matrix_gx_old, template=matrix_s, n=k, &
428 sym=dbcsr_type_no_symmetry)
429 END IF
430
431 SELECT CASE (qs_ot_env%settings%ot_algorithm)
432 CASE ("TOD")
433 CALL dbcsr_init_p(qs_ot_env%matrix_p)
434 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_p, template=matrix_s, m=k, n=k, &
435 sym=dbcsr_type_no_symmetry)
436
437 CALL dbcsr_init_p(qs_ot_env%matrix_r)
438 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_r, template=matrix_s, m=k, n=k, &
439 sym=dbcsr_type_no_symmetry)
440
441 CALL dbcsr_init_p(qs_ot_env%matrix_sinp)
442 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_sinp, template=matrix_s, m=k, n=k, &
443 sym=dbcsr_type_no_symmetry)
444
445 CALL dbcsr_init_p(qs_ot_env%matrix_cosp)
446 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_cosp, template=matrix_s, m=k, n=k, &
447 sym=dbcsr_type_no_symmetry)
448
449 CALL dbcsr_init_p(qs_ot_env%matrix_sinp_b)
450 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_sinp_b, template=matrix_s, m=k, n=k, &
451 sym=dbcsr_type_no_symmetry)
452
453 CALL dbcsr_init_p(qs_ot_env%matrix_cosp_b)
454 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_cosp_b, template=matrix_s, m=k, n=k, &
455 sym=dbcsr_type_no_symmetry)
456
457 CALL dbcsr_init_p(qs_ot_env%matrix_buf1)
458 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_buf1, template=matrix_s, m=k, n=k, &
459 sym=dbcsr_type_no_symmetry)
460
461 CALL dbcsr_init_p(qs_ot_env%matrix_buf2)
462 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_buf2, template=matrix_s, m=k, n=k, &
463 sym=dbcsr_type_no_symmetry)
464
465 CALL dbcsr_init_p(qs_ot_env%matrix_buf3)
466 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_buf3, template=matrix_s, m=k, n=k, &
467 sym=dbcsr_type_no_symmetry)
468
469 CALL dbcsr_init_p(qs_ot_env%matrix_buf4)
470 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_buf4, template=matrix_s, m=k, n=k, &
471 sym=dbcsr_type_no_symmetry)
472
473 CALL dbcsr_init_p(qs_ot_env%matrix_os)
474 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_os, template=matrix_s, m=my_ortho_k, n=my_ortho_k, &
475 sym=dbcsr_type_no_symmetry)
476
477 CALL dbcsr_init_p(qs_ot_env%matrix_buf1_ortho)
478 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_buf1_ortho, template=matrix_s, m=my_ortho_k, n=k, &
479 sym=dbcsr_type_no_symmetry)
480
481 CALL dbcsr_init_p(qs_ot_env%matrix_buf2_ortho)
482 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_buf2_ortho, template=matrix_s, m=my_ortho_k, n=k, &
483 sym=dbcsr_type_no_symmetry)
484
485 CASE ("REF")
486 CALL dbcsr_init_p(qs_ot_env%buf1_k_k_nosym)
487 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%buf1_k_k_nosym, template=matrix_s, m=k, n=k, &
488 sym=dbcsr_type_no_symmetry)
489
490 CALL dbcsr_init_p(qs_ot_env%buf2_k_k_nosym)
491 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%buf2_k_k_nosym, template=matrix_s, m=k, n=k, &
492 sym=dbcsr_type_no_symmetry)
493
494 CALL dbcsr_init_p(qs_ot_env%buf3_k_k_nosym)
495 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%buf3_k_k_nosym, template=matrix_s, m=k, n=k, &
496 sym=dbcsr_type_no_symmetry)
497
498 CALL dbcsr_init_p(qs_ot_env%buf4_k_k_nosym)
499 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%buf4_k_k_nosym, template=matrix_s, m=k, n=k, &
500 sym=dbcsr_type_no_symmetry)
501
502 ! It claims to be symmetric but to avoid dbcsr confusion nonsymmetric is kept
503 CALL dbcsr_init_p(qs_ot_env%buf1_k_k_sym)
504 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%buf1_k_k_sym, template=matrix_s, m=k, n=k, &
505 sym=dbcsr_type_no_symmetry)
506
507 CALL dbcsr_init_p(qs_ot_env%buf2_k_k_sym)
508 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%buf2_k_k_sym, template=matrix_s, m=k, n=k, &
509 sym=dbcsr_type_no_symmetry)
510
511 CALL dbcsr_init_p(qs_ot_env%buf3_k_k_sym)
512 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%buf3_k_k_sym, template=matrix_s, m=k, n=k, &
513 sym=dbcsr_type_no_symmetry)
514 !
515 CALL dbcsr_init_p(qs_ot_env%buf4_k_k_sym)
516 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%buf4_k_k_sym, template=matrix_s, m=k, n=k, &
517 sym=dbcsr_type_no_symmetry)
518 !
519 CALL dbcsr_init_p(qs_ot_env%p_k_k_sym)
520 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%p_k_k_sym, template=matrix_s, m=k, n=k, &
521 sym=dbcsr_type_no_symmetry)
522 !
523 CALL dbcsr_init_p(qs_ot_env%buf1_n_k)
524 CALL cp_dbcsr_m_by_n_from_row_template(qs_ot_env%buf1_n_k, template=matrix_s, n=k, &
525 sym=dbcsr_type_no_symmetry)
526 !
527 CALL dbcsr_init_p(qs_ot_env%matrix_buf1)
528 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%matrix_buf1, template=matrix_s, m=k, n=k, &
529 sym=dbcsr_type_no_symmetry)
530
531 END SELECT
532
533 IF (qs_ot_env%settings%ot_method .EQ. "DIIS" .OR. &
534 qs_ot_env%settings%ot_method .EQ. "BROY") THEN
535 NULLIFY (qs_ot_env%matrix_h_e)
536 NULLIFY (qs_ot_env%matrix_h_x)
537 CALL dbcsr_allocate_matrix_set(qs_ot_env%matrix_h_e, m_diis)
538 CALL dbcsr_allocate_matrix_set(qs_ot_env%matrix_h_x, m_diis)
539 DO i = 1, m_diis
540 CALL dbcsr_init_p(qs_ot_env%matrix_h_x(i)%matrix)
541 CALL cp_dbcsr_m_by_n_from_row_template(qs_ot_env%matrix_h_x(i)%matrix, template=matrix_s, n=k, &
542 sym=dbcsr_type_no_symmetry)
543
544 CALL dbcsr_init_p(qs_ot_env%matrix_h_e(i)%matrix)
545 CALL cp_dbcsr_m_by_n_from_row_template(qs_ot_env%matrix_h_e(i)%matrix, template=matrix_s, n=k, &
546 sym=dbcsr_type_no_symmetry)
547 END DO
548 END IF
549
550 NULLIFY (qs_ot_env%rot_mat_u, qs_ot_env%rot_mat_x, qs_ot_env%rot_mat_h_e, qs_ot_env%rot_mat_h_x, &
551 qs_ot_env%rot_mat_gx, qs_ot_env%rot_mat_gx_old, qs_ot_env%rot_mat_dx, &
552 qs_ot_env%rot_mat_evals, qs_ot_env%rot_mat_dedu, qs_ot_env%rot_mat_chc, &
553 qs_ot_env%rot_mat_evec_re, qs_ot_env%rot_mat_evec_im)
554
555 IF (qs_ot_env%settings%do_rotation) THEN
556 CALL dbcsr_init_p(qs_ot_env%rot_mat_u)
557 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%rot_mat_u, template=matrix_s, m=k, n=k, &
558 sym=dbcsr_type_no_symmetry)
559
560 CALL dbcsr_init_p(qs_ot_env%rot_mat_x)
561 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%rot_mat_x, template=matrix_s, m=k, n=k, &
562 sym=dbcsr_type_no_symmetry)
563
564 CALL dbcsr_init_p(qs_ot_env%rot_mat_dedu)
565 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%rot_mat_dedu, template=matrix_s, m=k, n=k, &
566 sym=dbcsr_type_no_symmetry)
567
568 CALL dbcsr_init_p(qs_ot_env%rot_mat_chc)
569 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%rot_mat_chc, template=matrix_s, m=k, n=k, &
570 sym=dbcsr_type_no_symmetry)
571
572 IF (qs_ot_env%settings%ot_method .EQ. "DIIS") THEN
573 CALL dbcsr_allocate_matrix_set(qs_ot_env%rot_mat_h_e, m_diis)
574 CALL dbcsr_allocate_matrix_set(qs_ot_env%rot_mat_h_x, m_diis)
575 DO i = 1, m_diis
576 CALL dbcsr_init_p(qs_ot_env%rot_mat_h_e(i)%matrix)
577 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%rot_mat_h_e(i)%matrix, template=matrix_s, m=k, n=k, &
578 sym=dbcsr_type_no_symmetry)
579
580 CALL dbcsr_init_p(qs_ot_env%rot_mat_h_x(i)%matrix)
581 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%rot_mat_h_x(i)%matrix, template=matrix_s, m=k, n=k, &
582 sym=dbcsr_type_no_symmetry)
583 END DO
584 END IF
585
586 ALLOCATE (qs_ot_env%rot_mat_evals(k))
587 CALL dbcsr_init_p(qs_ot_env%rot_mat_evec_re)
588 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%rot_mat_evec_re, template=matrix_s, m=k, n=k, &
589 sym=dbcsr_type_no_symmetry)
590 CALL dbcsr_init_p(qs_ot_env%rot_mat_evec_im)
591 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%rot_mat_evec_im, template=matrix_s, m=k, n=k, &
592 sym=dbcsr_type_no_symmetry)
593
594 CALL dbcsr_init_p(qs_ot_env%rot_mat_gx)
595 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%rot_mat_gx, template=matrix_s, m=k, n=k, &
596 sym=dbcsr_type_no_symmetry)
597
598 IF (qs_ot_env%use_gx_old) THEN
599 CALL dbcsr_init_p(qs_ot_env%rot_mat_gx_old)
600 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%rot_mat_gx_old, template=matrix_s, m=k, n=k, &
601 sym=dbcsr_type_no_symmetry)
602 END IF
603
604 IF (qs_ot_env%use_dx) THEN
605 CALL dbcsr_init_p(qs_ot_env%rot_mat_dx)
606 CALL cp_dbcsr_m_by_n_from_template(qs_ot_env%rot_mat_dx, template=matrix_s, m=k, n=k, &
607 sym=dbcsr_type_no_symmetry)
608 END IF
609
610 END IF
611
612 IF (qs_ot_env%settings%do_ener) THEN
613 ncoef = k
614 ALLOCATE (qs_ot_env%ener_x(ncoef))
615
616 IF (qs_ot_env%settings%ot_method .EQ. "DIIS") THEN
617 ALLOCATE (qs_ot_env%ener_h_e(m_diis, ncoef))
618 ALLOCATE (qs_ot_env%ener_h_x(m_diis, ncoef))
619 END IF
620
621 ALLOCATE (qs_ot_env%ener_gx(ncoef))
622
623 IF (qs_ot_env%use_gx_old) THEN
624 ALLOCATE (qs_ot_env%ener_gx_old(ncoef))
625 END IF
626
627 IF (qs_ot_env%use_dx) THEN
628 ALLOCATE (qs_ot_env%ener_dx(ncoef))
629 qs_ot_env%ener_dx = 0.0_dp
630 END IF
631 END IF
632
633 END SUBROUTINE qs_ot_allocate
634
635! **************************************************************************************************
636!> \brief deallocates data
637!> \param qs_ot_env ...
638! **************************************************************************************************
639 SUBROUTINE qs_ot_destroy(qs_ot_env)
640 TYPE(qs_ot_type) :: qs_ot_env
641
642 CALL mp_para_env_release(qs_ot_env%para_env)
643 CALL cp_blacs_env_release(qs_ot_env%blacs_env)
644
645 DEALLOCATE (qs_ot_env%evals)
646 DEALLOCATE (qs_ot_env%dum)
647
648 IF (ASSOCIATED(qs_ot_env%matrix_os)) CALL dbcsr_release_p(qs_ot_env%matrix_os)
649 IF (ASSOCIATED(qs_ot_env%matrix_p)) CALL dbcsr_release_p(qs_ot_env%matrix_p)
650 IF (ASSOCIATED(qs_ot_env%matrix_cosp)) CALL dbcsr_release_p(qs_ot_env%matrix_cosp)
651 IF (ASSOCIATED(qs_ot_env%matrix_sinp)) CALL dbcsr_release_p(qs_ot_env%matrix_sinp)
652 IF (ASSOCIATED(qs_ot_env%matrix_r)) CALL dbcsr_release_p(qs_ot_env%matrix_r)
653 IF (ASSOCIATED(qs_ot_env%matrix_cosp_b)) CALL dbcsr_release_p(qs_ot_env%matrix_cosp_b)
654 IF (ASSOCIATED(qs_ot_env%matrix_sinp_b)) CALL dbcsr_release_p(qs_ot_env%matrix_sinp_b)
655 IF (ASSOCIATED(qs_ot_env%matrix_buf1)) CALL dbcsr_release_p(qs_ot_env%matrix_buf1)
656 IF (ASSOCIATED(qs_ot_env%matrix_buf2)) CALL dbcsr_release_p(qs_ot_env%matrix_buf2)
657 IF (ASSOCIATED(qs_ot_env%matrix_buf3)) CALL dbcsr_release_p(qs_ot_env%matrix_buf3)
658 IF (ASSOCIATED(qs_ot_env%matrix_buf4)) CALL dbcsr_release_p(qs_ot_env%matrix_buf4)
659 IF (ASSOCIATED(qs_ot_env%matrix_buf1_ortho)) CALL dbcsr_release_p(qs_ot_env%matrix_buf1_ortho)
660 IF (ASSOCIATED(qs_ot_env%matrix_buf2_ortho)) CALL dbcsr_release_p(qs_ot_env%matrix_buf2_ortho)
661 IF (ASSOCIATED(qs_ot_env%matrix_c0)) CALL dbcsr_release_p(qs_ot_env%matrix_c0)
662 IF (ASSOCIATED(qs_ot_env%matrix_sc0)) CALL dbcsr_release_p(qs_ot_env%matrix_sc0)
663 IF (ASSOCIATED(qs_ot_env%matrix_psc0)) CALL dbcsr_release_p(qs_ot_env%matrix_psc0)
664 IF (ASSOCIATED(qs_ot_env%matrix_x)) CALL dbcsr_release_p(qs_ot_env%matrix_x)
665 IF (ASSOCIATED(qs_ot_env%matrix_sx)) CALL dbcsr_release_p(qs_ot_env%matrix_sx)
666 IF (ASSOCIATED(qs_ot_env%matrix_gx)) CALL dbcsr_release_p(qs_ot_env%matrix_gx)
667 IF (ASSOCIATED(qs_ot_env%matrix_dx)) CALL dbcsr_release_p(qs_ot_env%matrix_dx)
668 IF (ASSOCIATED(qs_ot_env%matrix_gx_old)) CALL dbcsr_release_p(qs_ot_env%matrix_gx_old)
669 IF (ASSOCIATED(qs_ot_env%buf1_k_k_nosym)) CALL dbcsr_release_p(qs_ot_env%buf1_k_k_nosym)
670 IF (ASSOCIATED(qs_ot_env%buf2_k_k_nosym)) CALL dbcsr_release_p(qs_ot_env%buf2_k_k_nosym)
671 IF (ASSOCIATED(qs_ot_env%buf3_k_k_nosym)) CALL dbcsr_release_p(qs_ot_env%buf3_k_k_nosym)
672 IF (ASSOCIATED(qs_ot_env%buf4_k_k_nosym)) CALL dbcsr_release_p(qs_ot_env%buf4_k_k_nosym)
673 IF (ASSOCIATED(qs_ot_env%p_k_k_sym)) CALL dbcsr_release_p(qs_ot_env%p_k_k_sym)
674 IF (ASSOCIATED(qs_ot_env%buf1_k_k_sym)) CALL dbcsr_release_p(qs_ot_env%buf1_k_k_sym)
675 IF (ASSOCIATED(qs_ot_env%buf2_k_k_sym)) CALL dbcsr_release_p(qs_ot_env%buf2_k_k_sym)
676 IF (ASSOCIATED(qs_ot_env%buf3_k_k_sym)) CALL dbcsr_release_p(qs_ot_env%buf3_k_k_sym)
677 IF (ASSOCIATED(qs_ot_env%buf4_k_k_sym)) CALL dbcsr_release_p(qs_ot_env%buf4_k_k_sym)
678 IF (ASSOCIATED(qs_ot_env%buf1_n_k)) CALL dbcsr_release_p(qs_ot_env%buf1_n_k)
679 IF (ASSOCIATED(qs_ot_env%buf1_n_k_dp)) CALL dbcsr_release_p(qs_ot_env%buf1_n_k_dp)
680
681 IF (qs_ot_env%settings%ot_method .EQ. "DIIS" .OR. &
682 qs_ot_env%settings%ot_method .EQ. "BROY") THEN
683 CALL dbcsr_deallocate_matrix_set(qs_ot_env%matrix_h_x)
684 CALL dbcsr_deallocate_matrix_set(qs_ot_env%matrix_h_e)
685 DEALLOCATE (qs_ot_env%ls_diis)
686 DEALLOCATE (qs_ot_env%lss_diis)
687 DEALLOCATE (qs_ot_env%c_diis)
688 DEALLOCATE (qs_ot_env%c_broy)
689 DEALLOCATE (qs_ot_env%energy_h)
690 DEALLOCATE (qs_ot_env%ipivot)
691 END IF
692
693 IF (qs_ot_env%settings%do_rotation) THEN
694
695 IF (ASSOCIATED(qs_ot_env%rot_mat_u)) CALL dbcsr_release_p(qs_ot_env%rot_mat_u)
696 IF (ASSOCIATED(qs_ot_env%rot_mat_x)) CALL dbcsr_release_p(qs_ot_env%rot_mat_x)
697 IF (ASSOCIATED(qs_ot_env%rot_mat_dedu)) CALL dbcsr_release_p(qs_ot_env%rot_mat_dedu)
698 IF (ASSOCIATED(qs_ot_env%rot_mat_chc)) CALL dbcsr_release_p(qs_ot_env%rot_mat_chc)
699
700 IF (qs_ot_env%settings%ot_method .EQ. "DIIS") THEN
701 CALL dbcsr_deallocate_matrix_set(qs_ot_env%rot_mat_h_x)
702 CALL dbcsr_deallocate_matrix_set(qs_ot_env%rot_mat_h_e)
703 END IF
704
705 DEALLOCATE (qs_ot_env%rot_mat_evals)
706 IF (ASSOCIATED(qs_ot_env%rot_mat_evec_re)) CALL dbcsr_release_p(qs_ot_env%rot_mat_evec_re)
707 IF (ASSOCIATED(qs_ot_env%rot_mat_evec_im)) CALL dbcsr_release_p(qs_ot_env%rot_mat_evec_im)
708 IF (ASSOCIATED(qs_ot_env%rot_mat_gx)) CALL dbcsr_release_p(qs_ot_env%rot_mat_gx)
709 IF (ASSOCIATED(qs_ot_env%rot_mat_gx_old)) CALL dbcsr_release_p(qs_ot_env%rot_mat_gx_old)
710 IF (ASSOCIATED(qs_ot_env%rot_mat_dx)) CALL dbcsr_release_p(qs_ot_env%rot_mat_dx)
711 END IF
712
713 IF (qs_ot_env%settings%do_ener) THEN
714 DEALLOCATE (qs_ot_env%ener_x)
715 DEALLOCATE (qs_ot_env%ener_gx)
716 IF (qs_ot_env%settings%ot_method .EQ. "DIIS") THEN
717 DEALLOCATE (qs_ot_env%ener_h_x)
718 DEALLOCATE (qs_ot_env%ener_h_e)
719 END IF
720 IF (qs_ot_env%use_dx) THEN
721 DEALLOCATE (qs_ot_env%ener_dx)
722 END IF
723 IF (qs_ot_env%use_gx_old) THEN
724 DEALLOCATE (qs_ot_env%ener_gx_old)
725 END IF
726 END IF
727
728 END SUBROUTINE qs_ot_destroy
729
730! **************************************************************************************************
731!> \brief ...
732!> \param settings ...
733!> \param ot_section ...
734!> \param output_unit ...
735! **************************************************************************************************
736 SUBROUTINE ot_readwrite_input(settings, ot_section, output_unit)
737 TYPE(qs_ot_settings_type) :: settings
738 TYPE(section_vals_type), POINTER :: ot_section
739 INTEGER, INTENT(IN) :: output_unit
740
741 CHARACTER(len=*), PARAMETER :: routinen = 'ot_readwrite_input'
742
743 INTEGER :: handle, ls_method, ot_algorithm, &
744 ot_method, ot_ortho_irac
745
746 CALL timeset(routinen, handle)
747
748 ! choose algorithm
749 CALL section_vals_val_get(ot_section, "ALGORITHM", i_val=ot_algorithm)
750 SELECT CASE (ot_algorithm)
752 settings%ot_algorithm = "TOD"
753 CASE (ot_algo_irac)
754 CALL cite_reference(weber2008)
755 settings%ot_algorithm = "REF"
756 CASE DEFAULT
757 cpabort("Value unknown")
758 END SELECT
759
760 ! irac input
761 CALL section_vals_val_get(ot_section, "IRAC_DEGREE", i_val=settings%irac_degree)
762 IF (settings%irac_degree < 2 .OR. settings%irac_degree > 4) THEN
763 cpabort("READ OT IRAC_DEGREE: Value unknown")
764 END IF
765 CALL section_vals_val_get(ot_section, "MAX_IRAC", i_val=settings%max_irac)
766 IF (settings%max_irac < 1) THEN
767 cpabort("READ OT MAX_IRAC: VALUE MUST BE GREATER THAN ZERO")
768 END IF
769 CALL section_vals_val_get(ot_section, "EPS_IRAC_FILTER_MATRIX", r_val=settings%eps_irac_filter_matrix)
770 CALL section_vals_val_get(ot_section, "EPS_IRAC", r_val=settings%eps_irac)
771 IF (settings%eps_irac < 0.0_dp) THEN
772 cpabort("READ OT EPS_IRAC: VALUE MUST BE GREATER THAN ZERO")
773 END IF
774 CALL section_vals_val_get(ot_section, "EPS_IRAC_QUICK_EXIT", r_val=settings%eps_irac_quick_exit)
775 IF (settings%eps_irac_quick_exit < 0.0_dp) THEN
776 cpabort("READ OT EPS_IRAC_QUICK_EXIT: VALUE MUST BE GREATER THAN ZERO")
777 END IF
778
779 CALL section_vals_val_get(ot_section, "EPS_IRAC_SWITCH", r_val=settings%eps_irac_switch)
780 IF (settings%eps_irac_switch < 0.0_dp) THEN
781 cpabort("READ OT EPS_IRAC_SWITCH: VALUE MUST BE GREATER THAN ZERO")
782 END IF
783
784 CALL section_vals_val_get(ot_section, "ORTHO_IRAC", i_val=ot_ortho_irac)
785 SELECT CASE (ot_ortho_irac)
786 CASE (ot_chol_irac)
787 settings%ortho_irac = "CHOL"
788 CASE (ot_poly_irac)
789 settings%ortho_irac = "POLY"
790 CASE (ot_lwdn_irac)
791 settings%ortho_irac = "LWDN"
792 CASE DEFAULT
793 cpabort("READ OT ORTHO_IRAC: Value unknown")
794 END SELECT
795
796 CALL section_vals_val_get(ot_section, "ON_THE_FLY_LOC", l_val=settings%on_the_fly_loc)
797
798 CALL section_vals_val_get(ot_section, "MINIMIZER", i_val=ot_method)
799 ! overwrite input if ot_state is set
800 IF (settings%ot_state == 1) THEN
801 ot_method = ot_mini_diis
802 END IF
803 ! compatibility
804 SELECT CASE (ot_method)
805 CASE (ot_mini_sd)
806 settings%ot_method = "SD"
807 CASE (ot_mini_cg)
808 settings%ot_method = "CG"
809 CASE (ot_mini_diis)
810 settings%ot_method = "DIIS"
811 CALL section_vals_val_get(ot_section, "N_HISTORY_VEC", i_val=settings%diis_m)
812 CASE (ot_mini_broyden)
813 CALL section_vals_val_get(ot_section, "N_HISTORY_VEC", i_val=settings%diis_m)
814 CALL section_vals_val_get(ot_section, "BROYDEN_BETA", r_val=settings%broyden_beta)
815 CALL section_vals_val_get(ot_section, "BROYDEN_GAMMA", r_val=settings%broyden_gamma)
816 CALL section_vals_val_get(ot_section, "BROYDEN_SIGMA", r_val=settings%broyden_sigma)
817 CALL section_vals_val_get(ot_section, "BROYDEN_ETA", r_val=settings%broyden_eta)
818 CALL section_vals_val_get(ot_section, "BROYDEN_OMEGA", r_val=settings%broyden_omega)
819 CALL section_vals_val_get(ot_section, "BROYDEN_SIGMA_DECREASE", r_val=settings%broyden_sigma_decrease)
820 CALL section_vals_val_get(ot_section, "BROYDEN_SIGMA_MIN", r_val=settings%broyden_sigma_min)
821 CALL section_vals_val_get(ot_section, "BROYDEN_FORGET_HISTORY", l_val=settings%broyden_forget_history)
822 CALL section_vals_val_get(ot_section, "BROYDEN_ADAPTIVE_SIGMA", l_val=settings%broyden_adaptive_sigma)
823 CALL section_vals_val_get(ot_section, "BROYDEN_ENABLE_FLIP", l_val=settings%broyden_enable_flip)
824 settings%ot_method = "BROY"
825 CASE DEFAULT
826 cpabort("READ OTSCF MINIMIZER: Value unknown")
827 END SELECT
828 CALL section_vals_val_get(ot_section, "SAFER_DIIS", l_val=settings%safer_diis)
829 CALL section_vals_val_get(ot_section, "LINESEARCH", i_val=ls_method)
830 SELECT CASE (ls_method)
831 CASE (ls_none)
832 settings%line_search_method = "NONE"
833 CASE (ls_2pnt)
834 settings%line_search_method = "2PNT"
835 CASE (ls_3pnt)
836 settings%line_search_method = "3PNT"
837 CASE (ls_adapt)
838 settings%line_search_method = "ADPT"
839 CASE (ls_gold)
840 settings%line_search_method = "GOLD"
841 CALL section_vals_val_get(ot_section, "GOLD_TARGET", r_val=settings%gold_target)
842 CASE DEFAULT
843 cpabort("READ OTSCF LS: Value unknown")
844 END SELECT
845
846 CALL section_vals_val_get(ot_section, "PRECOND_SOLVER", i_val=settings%precond_solver_type)
847 SELECT CASE (settings%precond_solver_type)
849 settings%precond_solver_name = "DEFAULT"
851 settings%precond_solver_name = "INVERSE_CHOLESKY"
853 settings%precond_solver_name = "DIRECT"
855 settings%precond_solver_name = "INVERSE_UPDATE"
856 CASE DEFAULT
857 cpabort("READ OTSCF SOLVER: Value unknown")
858 END SELECT
859
860 CALL section_vals_val_get(ot_section, "MAX_SCF_DIIS", i_val=settings%max_scf_diis)
861
862 !If these values are negative we will set them "optimal" for a given precondtioner below
863 CALL section_vals_val_get(ot_section, "STEPSIZE", r_val=settings%ds_min)
864 CALL section_vals_val_get(ot_section, "ENERGY_GAP", r_val=settings%energy_gap)
865
866 CALL section_vals_val_get(ot_section, "PRECONDITIONER", i_val=settings%preconditioner_type)
867 SELECT CASE (settings%preconditioner_type)
868 CASE (ot_precond_none)
869 settings%preconditioner_name = "NONE"
870 IF (settings%ds_min < 0.0_dp) settings%ds_min = 0.15_dp
871 IF (settings%energy_gap < 0.0_dp) settings%energy_gap = 0.2_dp
873 settings%preconditioner_name = "FULL_SINGLE"
874 IF (settings%ds_min < 0.0_dp) settings%ds_min = 0.15_dp
875 IF (settings%energy_gap < 0.0_dp) settings%energy_gap = 0.2_dp
877 settings%preconditioner_name = "FULL_SINGLE_INVERSE"
878 IF (settings%ds_min < 0.0_dp) settings%ds_min = 0.08_dp
879 IF (settings%energy_gap < 0.0_dp) settings%energy_gap = 0.08_dp
881 settings%preconditioner_name = "FULL_ALL"
882 IF (settings%ds_min < 0.0_dp) settings%ds_min = 0.15_dp
883 IF (settings%energy_gap < 0.0_dp) settings%energy_gap = 0.08_dp
885 settings%preconditioner_name = "FULL_KINETIC"
886 IF (settings%ds_min < 0.0_dp) settings%ds_min = 0.15_dp
887 IF (settings%energy_gap < 0.0_dp) settings%energy_gap = 0.2_dp
889 settings%preconditioner_name = "FULL_S_INVERSE"
890 IF (settings%ds_min < 0.0_dp) settings%ds_min = 0.15_dp
891 IF (settings%energy_gap < 0.0_dp) settings%energy_gap = 0.2_dp
892 CASE DEFAULT
893 cpabort("READ OTSCF PRECONDITIONER: Value unknown")
894 END SELECT
895 CALL section_vals_val_get(ot_section, "CHOLESKY", i_val=settings%cholesky_type)
896 CALL section_vals_val_get(ot_section, "EPS_TAYLOR", r_val=settings%eps_taylor)
897 CALL section_vals_val_get(ot_section, "MAX_TAYLOR", i_val=settings%max_taylor)
898 CALL section_vals_val_get(ot_section, "ROTATION", l_val=settings%do_rotation)
899 CALL section_vals_val_get(ot_section, "ENERGIES", l_val=settings%do_ener)
900 CALL section_vals_val_get(ot_section, "OCCUPATION_PRECONDITIONER", &
901 l_val=settings%occupation_preconditioner)
902 CALL section_vals_val_get(ot_section, "NONDIAG_ENERGY", l_val=settings%add_nondiag_energy)
903 CALL section_vals_val_get(ot_section, "NONDIAG_ENERGY_STRENGTH", &
904 r_val=settings%nondiag_energy_strength)
905 ! not yet fully implemented
906 cpassert(.NOT. settings%do_ener)
907
908 ! write OT output
909
910 IF (output_unit > 0) THEN
911 WRITE (output_unit, '(/,A)') " ----------------------------------- OT ---------------------------------------"
912 IF (settings%do_rotation) THEN
913 WRITE (output_unit, '(A)') " Allowing for rotations "
914 END IF
915 IF (settings%do_ener) THEN
916 WRITE (output_unit, '(A,L2)') " Optimizing orbital energies "
917 END IF
918 SELECT CASE (settings%OT_METHOD)
919 CASE ("SD")
920 WRITE (output_unit, '(A)') " Minimizer : SD : steepest descent"
921 CASE ("CG")
922 WRITE (output_unit, '(A)') " Minimizer : CG : conjugate gradient"
923 CASE ("DIIS")
924 WRITE (output_unit, '(A)') " Minimizer : DIIS : direct inversion"
925 WRITE (output_unit, '(A)') " in the iterative subspace"
926 WRITE (output_unit, '(A,I3,A)') " using ", settings%diis_m, " DIIS vectors"
927 IF (settings%safer_diis) THEN
928 WRITE (output_unit, '(A,I3,A)') " safer DIIS on"
929 ELSE
930 WRITE (output_unit, '(A,I3,A)') " safer DIIS off"
931 END IF
932 CASE ("BROY")
933 WRITE (output_unit, '(A)') " Minimizer : BROYDEN : Broyden "
934 WRITE (output_unit, '(A,F16.8)') " BETA : ", settings%broyden_beta
935 WRITE (output_unit, '(A,F16.8)') " GAMMA : ", settings%broyden_gamma
936 WRITE (output_unit, '(A,F16.8)') " SIGMA : ", settings%broyden_sigma
937 WRITE (output_unit, '(A,I3,A)') " using : - ", &
938 settings%diis_m, " BROYDEN vectors"
939 CASE DEFAULT
940 WRITE (output_unit, '(3A)') " Minimizer : ", settings%OT_METHOD, " : UNKNOWN"
941 END SELECT
942 SELECT CASE (settings%preconditioner_name)
943 CASE ("FULL_SINGLE")
944 WRITE (output_unit, '(A)') " Preconditioner : FULL_SINGLE : diagonalization based"
945 CASE ("FULL_SINGLE_INVERSE")
946 WRITE (output_unit, '(A,/,A)') " Preconditioner : FULL_SINGLE_INVERSE : inversion of ", &
947 " H + eS - 2*(Sc)(c^T*H*c+const)(Sc)^T"
948 CASE ("FULL_ALL")
949 WRITE (output_unit, '(A)') " Preconditioner : FULL_ALL : diagonalization, state selective"
950 CASE ("FULL_KINETIC")
951 WRITE (output_unit, '(A)') " Preconditioner : FULL_KINETIC : inversion of T + eS"
952 CASE ("FULL_S_INVERSE")
953 WRITE (output_unit, '(A)') " Preconditioner : FULL_S_INVERSE : cholesky inversion of S"
954 CASE ("SPARSE_DIAG")
955 WRITE (output_unit, '(A)') &
956 " Preconditioner : SPARSE_DIAG : diagonal atomic block diagonalization"
957 CASE ("SPARSE_KINETIC")
958 WRITE (output_unit, '(A)') " Preconditioner : SPARSE_KINETIC : sparse linear solver for T + eS"
959 CASE ("NONE")
960 WRITE (output_unit, '(A)') " Preconditioner : NONE"
961 CASE DEFAULT
962 WRITE (output_unit, '(3A)') " Preconditioner : ", settings%preconditioner_name, " : UNKNOWN"
963 END SELECT
964
965 WRITE (output_unit, '(A)') " Precond_solver : "//trim(settings%precond_solver_name)
966
967 IF (settings%OT_METHOD .EQ. "SD" .OR. settings%OT_METHOD .EQ. "CG") THEN
968 SELECT CASE (settings%line_search_method)
969 CASE ("2PNT")
970 WRITE (output_unit, '(A)') " Line search : 2PNT : 2 energies, one gradient"
971 CASE ("3PNT")
972 WRITE (output_unit, '(A)') " Line search : 3PNT : 3 energies"
973 CASE ("GOLD")
974 WRITE (output_unit, '(A)') " Line search : GOLD : bracketing and golden section search"
975 WRITE (output_unit, '(A,F14.8)') " target rel accuracy : ", settings%gold_target
976 CASE ("NONE")
977 WRITE (output_unit, '(A)') " Line search : NONE"
978 CASE DEFAULT
979 WRITE (output_unit, '(3A)') " Line search : ", settings%line_search_method, " : UNKNOWN"
980 END SELECT
981 END IF
982 WRITE (output_unit, '(A,F14.8,T49,A,F14.8)') " stepsize :", settings%ds_min, &
983 " energy_gap :", settings%energy_gap
984 IF (settings%ot_algorithm .EQ. 'TOD') THEN
985 WRITE (output_unit, '(A,E14.5,T49,A,I14)') " eps_taylor :", settings%eps_taylor, &
986 " max_taylor :", settings%max_taylor
987 END IF
988 IF (settings%ot_algorithm .EQ. 'REF') THEN
989 WRITE (output_unit, '(A,1X,A,T49,A,I14)') " ortho_irac :", settings%ortho_irac, &
990 " irac_degree :", settings%irac_degree
991 WRITE (output_unit, '(A,I14,T49,A,E14.5)') " max_irac :", settings%max_irac, &
992 " eps_irac :", settings%eps_irac
993 WRITE (output_unit, '(A,E14.5,T49,A,E10.3)') " eps_irac_switch:", settings%eps_irac_switch, &
994 " eps_irac_quick_exit:", settings%eps_irac_quick_exit
995 WRITE (output_unit, '(A,L2)') " on_the_fly_loc :", settings%on_the_fly_loc
996 END IF
997 WRITE (output_unit, '(A)') " ----------------------------------- OT ---------------------------------------"
998 WRITE (unit=output_unit, &
999 fmt="(/,T3,A,T12,A,T31,A,T39,A,T59,A,T75,A,/,T3,A)") &
1000 "Step", "Update method", "Time", "Convergence", "Total energy", "Change", &
1001 repeat("-", 78)
1002 END IF
1003
1004 CALL timestop(handle)
1005
1006 END SUBROUTINE ot_readwrite_input
1007
1008END MODULE qs_ot_types
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public vandevondele2003
integer, save, public weber2008
methods related to the blacs parallel environment
subroutine, public cp_blacs_env_release(blacs_env)
releases the given blacs_env
subroutine, public dbcsr_release_p(matrix)
...
subroutine, public dbcsr_init_p(matrix)
...
subroutine, public dbcsr_set(matrix, alpha)
...
DBCSR operations in CP2K.
subroutine, public cp_dbcsr_m_by_n_from_row_template(matrix, template, n, sym)
Utility function to create dbcsr matrix, m x n matrix (n arbitrary) with the same processor grid and ...
subroutine, public cp_dbcsr_m_by_n_from_template(matrix, template, m, n, sym)
Utility function to create an arbitrary shaped dbcsr matrix with the same processor grid as the templ...
represent the structure of a full matrix
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
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public ls_adapt
integer, parameter, public ot_chol_irac
integer, parameter, public ls_3pnt
integer, parameter, public ot_algo_irac
integer, parameter, public ot_algo_taylor_or_diag
integer, parameter, public ot_poly_irac
integer, parameter, public ot_mini_cg
integer, parameter, public ot_precond_full_kinetic
integer, parameter, public cholesky_reduce
integer, parameter, public ot_mini_diis
integer, parameter, public ot_precond_solver_default
integer, parameter, public ot_precond_full_single
integer, parameter, public ot_precond_solver_inv_chol
integer, parameter, public ot_precond_none
integer, parameter, public ls_2pnt
integer, parameter, public ls_none
integer, parameter, public ot_precond_full_single_inverse
integer, parameter, public ot_lwdn_irac
integer, parameter, public ls_gold
integer, parameter, public do_taylor
integer, parameter, public ot_precond_s_inverse
integer, parameter, public ot_mini_broyden
integer, parameter, public ot_precond_solver_update
integer, parameter, public ot_mini_sd
integer, parameter, public ot_precond_full_all
integer, parameter, public ot_precond_solver_direct
objects that represent the structure of input sections and the data contained in an input section
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
Interface to the message passing library MPI.
subroutine, public mp_para_env_release(para_env)
releases the para object (to be called when you don't want anymore the shared copy of this object)
types of preconditioners
computes preconditioners, and implements methods to apply them currently used in qs_ot
orbital transformations
Definition qs_ot_types.F:15
subroutine, public qs_ot_init(qs_ot_env)
init matrices, needs c0 and sc0 so that c0*sc0=1
subroutine, public qs_ot_allocate(qs_ot_env, matrix_s, fm_struct_ref, ortho_k)
allocates the data in qs_ot_env, for a calculation with fm_struct_ref ortho_k allows for specifying a...
subroutine, public qs_ot_settings_init(settings)
sets default values for the settings type
subroutine, public qs_ot_destroy(qs_ot_env)
deallocates data
subroutine, public ot_readwrite_input(settings, ot_section, output_unit)
...
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
keeps the information about the structure of a full matrix
stores all the informations relevant to an mpi environment
notice, this variable needs to be copyable, needed for spins as e.g. in qs_ot_scf
Definition qs_ot_types.F:63