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