(git:039bc03)
Loading...
Searching...
No Matches
ec_orth_solver.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief AO-based conjugate-gradient response solver routines
10!>
11!>
12!> \date 09.2019
13!> \author Fabian Belleflamme
14! **************************************************************************************************
16 USE admm_types, ONLY: admm_type,&
19 USE cp_dbcsr_api, ONLY: &
22 dbcsr_transposed, dbcsr_type, dbcsr_type_no_symmetry
29 USE ec_methods, ONLY: create_kernel
42 USE kinds, ONLY: dp
43 USE machine, ONLY: m_flush,&
45 USE mathlib, ONLY: abnormal_value
47 USE pw_env_types, ONLY: pw_env_get,&
49 USE pw_methods, ONLY: pw_axpy,&
50 pw_scale,&
55 USE pw_pool_types, ONLY: pw_pool_p_type,&
57 USE pw_types, ONLY: pw_c1d_gs_type,&
61 USE qs_integrate_potential, ONLY: integrate_v_rspace
63 USE qs_linres_kernel, ONLY: apply_hfx,&
70 USE qs_rho_types, ONLY: qs_rho_get,&
72 USE xc, ONLY: xc_prep_2nd_deriv
73#include "./base/base_uses.f90"
74
75 IMPLICIT NONE
76
77 PRIVATE
78
79 ! Global parameters
80
81 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'ec_orth_solver'
82
83 ! Public subroutines
84
85 PUBLIC :: ec_response_ao
86
87CONTAINS
88
89! **************************************************************************************************
90!> \brief Preconditioning of the AO-based CG linear response solver
91!> M * z_0 = r_0
92!> M(X) = [F,B], with B = [X,P]
93!> for M we need F and P in ortho basis
94!> Returns z_0, the preconditioned residual in orthonormal basis
95!>
96!> All matrices are in orthonormal Lowdin basis
97!>
98!> \param qs_env ...
99!> \param matrix_ks Ground-state Kohn-Sham matrix
100!> \param matrix_p Ground-state Density matrix
101!> \param matrix_rhs Unpreconditioned residual of linear response CG
102!> \param matrix_cg_z Preconditioned residual
103!> \param eps_filter ...
104!> \param iounit ...
105!>
106!> \param silent ...
107!> \date 01.2020
108!> \author Fabian Belleflamme
109! **************************************************************************************************
110 SUBROUTINE ec_preconditioner(qs_env, matrix_ks, matrix_p, matrix_rhs, &
111 matrix_cg_z, eps_filter, iounit, silent)
112
113 TYPE(qs_environment_type), POINTER :: qs_env
114 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN), &
115 POINTER :: matrix_ks, matrix_p, matrix_rhs
116 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(INOUT), &
117 POINTER :: matrix_cg_z
118 REAL(KIND=dp), INTENT(IN) :: eps_filter
119 INTEGER, INTENT(IN) :: iounit
120 LOGICAL, INTENT(IN), OPTIONAL :: silent
121
122 CHARACTER(len=*), PARAMETER :: routineN = 'ec_preconditioner'
123
124 INTEGER :: handle, i, ispin, max_iter, nao, nspins
125 LOGICAL :: converged, my_silent
126 REAL(KIND=dp) :: norm_res, t1, t2
127 REAL(KIND=dp), DIMENSION(:), POINTER :: alpha, beta, new_norm, norm_ca, norm_rr
128 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_Ax, matrix_b, matrix_cg, &
129 matrix_res
130 TYPE(dft_control_type), POINTER :: dft_control
131 TYPE(linres_control_type), POINTER :: linres_control
132
133 CALL timeset(routinen, handle)
134
135 my_silent = .false.
136 IF (PRESENT(silent)) my_silent = silent
137
138 cpassert(ASSOCIATED(qs_env))
139 cpassert(ASSOCIATED(matrix_ks))
140 cpassert(ASSOCIATED(matrix_p))
141 cpassert(ASSOCIATED(matrix_rhs))
142 cpassert(ASSOCIATED(matrix_cg_z))
143
144 NULLIFY (dft_control, linres_control)
145
146 t1 = m_walltime()
147
148 CALL get_qs_env(qs_env=qs_env, &
149 dft_control=dft_control, &
150 linres_control=linres_control)
151 nspins = dft_control%nspins
152 CALL dbcsr_get_info(matrix_ks(1)%matrix, nfullrows_total=nao)
153
154 ALLOCATE (alpha(nspins), beta(nspins), new_norm(nspins), norm_ca(nspins), norm_rr(nspins))
155
156 !----------------------------------------
157 ! Create non-symmetric matrices: Ax, B, cg, res
158 !----------------------------------------
159
160 NULLIFY (matrix_ax, matrix_b, matrix_cg, matrix_res)
161 CALL dbcsr_allocate_matrix_set(matrix_ax, nspins)
162 CALL dbcsr_allocate_matrix_set(matrix_b, nspins)
163 CALL dbcsr_allocate_matrix_set(matrix_cg, nspins)
164 CALL dbcsr_allocate_matrix_set(matrix_res, nspins)
165
166 DO ispin = 1, nspins
167 ALLOCATE (matrix_ax(ispin)%matrix)
168 ALLOCATE (matrix_b(ispin)%matrix)
169 ALLOCATE (matrix_cg(ispin)%matrix)
170 ALLOCATE (matrix_res(ispin)%matrix)
171 CALL dbcsr_create(matrix_ax(ispin)%matrix, name="linop MATRIX", &
172 template=matrix_ks(1)%matrix, &
173 matrix_type=dbcsr_type_no_symmetry)
174 CALL dbcsr_create(matrix_b(ispin)%matrix, name="MATRIX B", &
175 template=matrix_ks(1)%matrix, &
176 matrix_type=dbcsr_type_no_symmetry)
177 CALL dbcsr_create(matrix_cg(ispin)%matrix, name="TRIAL MATRIX", &
178 template=matrix_ks(1)%matrix, &
179 matrix_type=dbcsr_type_no_symmetry)
180 CALL dbcsr_create(matrix_res(ispin)%matrix, name="RESIDUE", &
181 template=matrix_ks(1)%matrix, &
182 matrix_type=dbcsr_type_no_symmetry)
183 END DO
184
185 !----------------------------------------
186 ! Get righ-hand-side operators
187 !----------------------------------------
188
189 ! Initial guess z_0
190 DO ispin = 1, nspins
191 CALL dbcsr_copy(matrix_cg_z(ispin)%matrix, matrix_rhs(ispin)%matrix)
192
193 ! r_0 = b
194 CALL dbcsr_copy(matrix_res(ispin)%matrix, matrix_rhs(ispin)%matrix)
195 END DO
196
197 ! Projector on trial matrix
198 ! Projector does not need to be applied here,
199 ! as matrix_rhs already had this done before entering preconditioner
200 !CALL projector(qs_env, matrix_p, matrix_cg_z, eps_filter)
201
202 ! Mz_0
203 CALL hessian_op1(matrix_ks, matrix_p, matrix_cg_z, matrix_b, matrix_ax, eps_filter)
204
205 ! r_0 = b - Ax_0
206 DO ispin = 1, nspins
207 CALL dbcsr_add(matrix_res(ispin)%matrix, matrix_ax(ispin)%matrix, 1.0_dp, -1.0_dp)
208 END DO
209
210 ! Matrix projector T
211 CALL projector(qs_env, matrix_p, matrix_res, eps_filter)
212
213 DO ispin = 1, nspins
214 ! cg = p_0 = z_0
215 CALL dbcsr_copy(matrix_cg(ispin)%matrix, matrix_res(ispin)%matrix)
216 END DO
217
218 ! header
219 IF (iounit > 0 .AND. .NOT. my_silent) THEN
220 WRITE (iounit, "(/,T10,A)") "Preconditioning of search direction"
221 WRITE (iounit, "(/,T10,A,T25,A,T42,A,T62,A,/,T10,A)") &
222 "Iteration", "Stepsize", "Convergence", "Time", &
223 repeat("-", 58)
224 END IF
225
226 alpha(:) = 0.0_dp
227 max_iter = 200
228 converged = .false.
229 norm_res = 0.0_dp
230
231 ! start iteration
232 iteration: DO i = 1, max_iter
233
234 ! Hessian Ax = [F,B] is updated preconditioner
235 CALL hessian_op1(matrix_ks, matrix_p, matrix_cg, matrix_b, matrix_ax, eps_filter)
236
237 ! Matrix projector
238 CALL projector(qs_env, matrix_p, matrix_ax, eps_filter)
239
240 DO ispin = 1, nspins
241
242 ! Tr(r_0 * r_0)
243 CALL dbcsr_dot(matrix_res(ispin)%matrix, matrix_res(ispin)%matrix, norm_rr(ispin))
244 IF (abnormal_value(norm_rr(ispin))) &
245 cpabort("Preconditioner: Tr[r_j*r_j] is an abnormal value (NaN/Inf)")
246
247 IF (norm_rr(ispin) < 0.0_dp) cpabort("norm_rr < 0")
248 norm_res = max(norm_res, abs(norm_rr(ispin)/real(nao, dp)))
249
250 ! norm_cA = tr(Ap_j * p_j)
251 CALL dbcsr_dot(matrix_cg(ispin)%matrix, matrix_ax(ispin)%matrix, norm_ca(ispin))
252
253 ! Determine step-size
254 IF (norm_ca(ispin) < linres_control%eps) THEN
255 alpha(ispin) = 1.0_dp
256 ELSE
257 alpha(ispin) = norm_rr(ispin)/norm_ca(ispin)
258 END IF
259
260 ! x_j+1 = x_j + alpha*p_j
261 ! save contribution of this iteration
262 CALL dbcsr_add(matrix_cg_z(ispin)%matrix, matrix_cg(ispin)%matrix, 1.0_dp, alpha(ispin))
263
264 ! r_j+1 = r_j - alpha * Ap_j
265 CALL dbcsr_add(matrix_res(ispin)%matrix, matrix_ax(ispin)%matrix, 1.0_dp, -alpha(ispin))
266
267 END DO
268
269 norm_res = 0.0_dp
270
271 DO ispin = 1, nspins
272 ! Tr[r_j+1*z_j+1]
273 CALL dbcsr_dot(matrix_res(ispin)%matrix, matrix_res(ispin)%matrix, new_norm(ispin))
274 IF (new_norm(ispin) < 0.0_dp) cpabort("tr(r_j+1*z_j+1) < 0")
275 IF (abnormal_value(new_norm(ispin))) &
276 cpabort("Preconditioner: Tr[r_j+1*z_j+1] is an abnormal value (NaN/Inf)")
277 norm_res = max(norm_res, new_norm(ispin)/real(nao, dp))
278
279 IF (norm_rr(ispin) < linres_control%eps*0.001_dp &
280 .OR. new_norm(ispin) < linres_control%eps*0.001_dp) THEN
281 beta(ispin) = 0.0_dp
282 converged = .true.
283 ELSE
284 beta(ispin) = new_norm(ispin)/norm_rr(ispin)
285 END IF
286
287 ! update new search vector (matrix cg)
288 ! cg_j+1 = z_j+1 + beta*cg_j
289 CALL dbcsr_add(matrix_cg(ispin)%matrix, matrix_res(ispin)%matrix, beta(ispin), 1.0_dp)
290 CALL dbcsr_filter(matrix_cg(ispin)%matrix, eps_filter)
291
292 norm_rr(ispin) = new_norm(ispin)
293 END DO
294
295 ! Convergence criteria
296 IF (norm_res < linres_control%eps) THEN
297 converged = .true.
298 END IF
299
300 t2 = m_walltime()
301 IF (i == 1 .OR. mod(i, 1) == 0 .OR. converged) THEN
302 IF (iounit > 0 .AND. .NOT. my_silent) THEN
303 WRITE (iounit, "(T10,I5,T25,1E8.2,T33,F25.14,T58,F8.2)") &
304 i, maxval(alpha), norm_res, t2 - t1
305 ! Convergence in scientific notation
306 !WRITE (iounit, "(T10,I5,T25,1E8.2,T42,1E14.8,T58,F8.2)") &
307 ! i, MAXVAL(alpha), norm_res, t2 - t1
308 CALL m_flush(iounit)
309 END IF
310 END IF
311 IF (converged) THEN
312 IF (iounit > 0 .AND. .NOT. my_silent) THEN
313 WRITE (iounit, "(/,T10,A,I4,A,/)") "The precon solver converged in ", i, " iterations."
314 CALL m_flush(iounit)
315 END IF
316 EXIT iteration
317 END IF
318
319 ! Max number of iteration reached
320 IF (i == max_iter) THEN
321 IF (iounit > 0) THEN
322 WRITE (iounit, "(/,T10,A/)") &
323 "The precon solver didnt converge! Maximum number of iterations reached."
324 CALL m_flush(iounit)
325 END IF
326 converged = .false.
327 END IF
328
329 END DO iteration
330
331 ! Matrix projector
332 CALL projector(qs_env, matrix_p, matrix_cg_z, eps_filter)
333
334 ! Release matrices
335 CALL dbcsr_deallocate_matrix_set(matrix_ax)
336 CALL dbcsr_deallocate_matrix_set(matrix_b)
337 CALL dbcsr_deallocate_matrix_set(matrix_res)
338 CALL dbcsr_deallocate_matrix_set(matrix_cg)
339
340 DEALLOCATE (alpha, beta, new_norm, norm_ca, norm_rr)
341
342 CALL timestop(handle)
343
344 END SUBROUTINE ec_preconditioner
345
346! **************************************************************************************************
347!> \brief AO-based conjugate gradient linear response solver.
348!> In goes the right hand side B of the equation AZ=B, and the linear transformation of the
349!> Hessian matrix A on trial matrices is iteratively solved. Result are
350!> the response density matrix_pz, and the energy-weighted response density matrix_wz.
351!>
352!> \param qs_env ...
353!> \param p_env ...
354!> \param matrix_hz Right hand-side of linear response equation
355!> \param matrix_pz Response density
356!> \param matrix_wz Energy-weighted response density matrix
357!> \param iounit ...
358!> \param should_stop ...
359!>
360!> \param silent ...
361!> \date 01.2020
362!> \author Fabian Belleflamme
363! **************************************************************************************************
364 SUBROUTINE ec_response_ao(qs_env, p_env, matrix_hz, matrix_pz, matrix_wz, iounit, &
365 should_stop, silent)
366
367 TYPE(qs_environment_type), POINTER :: qs_env
368 TYPE(qs_p_env_type), POINTER :: p_env
369 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN), &
370 POINTER :: matrix_hz
371 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(INOUT), &
372 POINTER :: matrix_pz, matrix_wz
373 INTEGER, INTENT(IN) :: iounit
374 LOGICAL, INTENT(OUT) :: should_stop
375 LOGICAL, INTENT(IN), OPTIONAL :: silent
376
377 CHARACTER(len=*), PARAMETER :: routinen = 'ec_response_ao'
378
379 INTEGER :: handle, i, ispin, max_iter_lanczos, nao, &
380 nspins, s_sqrt_method, s_sqrt_order
381 LOGICAL :: my_silent, restart
382 REAL(kind=dp) :: eps_filter, eps_lanczos, focc, &
383 min_shift, norm_res, old_conv, shift, &
384 t1, t2
385 REAL(kind=dp), DIMENSION(:), POINTER :: alpha, beta, new_norm, norm_ca, norm_rr, &
386 tr_rz00
387 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: ksmat, matrix_ax, matrix_cg, matrix_cg_z, &
388 matrix_ks, matrix_nsc, matrix_p, matrix_res, matrix_s, matrix_z, matrix_z0, rho_ao
389 TYPE(dbcsr_type) :: matrix_s_sqrt, matrix_s_sqrt_inv, &
390 matrix_tmp
391 TYPE(dft_control_type), POINTER :: dft_control
392 TYPE(linres_control_type), POINTER :: linres_control
393 TYPE(qs_rho_type), POINTER :: rho
394 TYPE(section_vals_type), POINTER :: solver_section
395
396 CALL timeset(routinen, handle)
397
398 my_silent = .false.
399 IF (PRESENT(silent)) my_silent = silent
400
401 cpassert(ASSOCIATED(qs_env))
402 cpassert(ASSOCIATED(matrix_hz))
403 cpassert(ASSOCIATED(matrix_pz))
404 cpassert(ASSOCIATED(matrix_wz))
405
406 NULLIFY (dft_control, ksmat, matrix_s, linres_control, rho)
407
408 t1 = m_walltime()
409
410 CALL get_qs_env(qs_env=qs_env, &
411 dft_control=dft_control, &
412 linres_control=linres_control, &
413 matrix_ks=ksmat, &
414 matrix_s=matrix_s, &
415 rho=rho)
416 nspins = dft_control%nspins
417
418 CALL dbcsr_get_info(matrix_s(1)%matrix, nfullrows_total=nao)
419
420 solver_section => section_vals_get_subs_vals(qs_env%input, "DFT%ENERGY_CORRECTION%RESPONSE_SOLVER")
421 CALL section_vals_val_get(solver_section, "S_SQRT_METHOD", i_val=s_sqrt_method)
422 CALL section_vals_val_get(solver_section, "S_SQRT_ORDER", i_val=s_sqrt_order)
423 CALL section_vals_val_get(solver_section, "EPS_LANCZOS", r_val=eps_lanczos)
424 CALL section_vals_val_get(solver_section, "MAX_ITER_LANCZOS", i_val=max_iter_lanczos)
425
426 eps_filter = linres_control%eps_filter
427
428 CALL qs_rho_get(rho, rho_ao=rho_ao)
429
430 ALLOCATE (alpha(nspins), beta(nspins), new_norm(nspins), norm_ca(nspins), norm_rr(nspins))
431 ALLOCATE (tr_rz00(nspins))
432
433 ! local matrix P, KS, and NSC
434 ! to bring into orthogonal basis
435 NULLIFY (matrix_p, matrix_ks, matrix_nsc)
436 CALL dbcsr_allocate_matrix_set(matrix_p, nspins)
437 CALL dbcsr_allocate_matrix_set(matrix_ks, nspins)
438 CALL dbcsr_allocate_matrix_set(matrix_nsc, nspins)
439 DO ispin = 1, nspins
440 ALLOCATE (matrix_p(ispin)%matrix)
441 ALLOCATE (matrix_ks(ispin)%matrix)
442 ALLOCATE (matrix_nsc(ispin)%matrix)
443 CALL dbcsr_create(matrix_p(ispin)%matrix, name="P_IN ORTHO", &
444 template=ksmat(1)%matrix, &
445 matrix_type=dbcsr_type_no_symmetry)
446 CALL dbcsr_create(matrix_ks(ispin)%matrix, name="KS_IN ORTHO", &
447 template=ksmat(1)%matrix, &
448 matrix_type=dbcsr_type_no_symmetry)
449 CALL dbcsr_create(matrix_nsc(ispin)%matrix, name="NSC IN ORTHO", &
450 template=ksmat(1)%matrix, &
451 matrix_type=dbcsr_type_no_symmetry)
452
453 CALL dbcsr_desymmetrize(rho_ao(ispin)%matrix, matrix_p(ispin)%matrix)
454 CALL dbcsr_desymmetrize(ksmat(ispin)%matrix, matrix_ks(ispin)%matrix)
455 CALL dbcsr_desymmetrize(matrix_hz(ispin)%matrix, matrix_nsc(ispin)%matrix)
456 END DO
457
458 ! Scale matrix_p by factor 1/2 in closed-shell
459 IF (nspins == 1) CALL dbcsr_scale(matrix_p(1)%matrix, 0.5_dp)
460
461 ! Transform P, KS, and Harris kernel matrix into Orthonormal basis
462 CALL dbcsr_create(matrix_s_sqrt, template=matrix_s(1)%matrix, &
463 matrix_type=dbcsr_type_no_symmetry)
464 CALL dbcsr_create(matrix_s_sqrt_inv, template=matrix_s(1)%matrix, &
465 matrix_type=dbcsr_type_no_symmetry)
466
467 SELECT CASE (s_sqrt_method)
468 CASE (ls_s_sqrt_proot)
469 CALL matrix_sqrt_proot(matrix_s_sqrt, matrix_s_sqrt_inv, &
470 matrix_s(1)%matrix, eps_filter, &
471 s_sqrt_order, eps_lanczos, max_iter_lanczos, symmetrize=.true.)
472 CASE (ls_s_sqrt_ns)
473 CALL matrix_sqrt_newton_schulz(matrix_s_sqrt, matrix_s_sqrt_inv, &
474 matrix_s(1)%matrix, eps_filter, &
475 s_sqrt_order, eps_lanczos, max_iter_lanczos)
476 CASE DEFAULT
477 cpabort("Unknown sqrt method.")
478 END SELECT
479
480 ! Transform into orthonormal Lowdin basis
481 DO ispin = 1, nspins
482 CALL transform_m_orth(matrix_p(ispin)%matrix, matrix_s_sqrt, eps_filter)
483 CALL transform_m_orth(matrix_ks(ispin)%matrix, matrix_s_sqrt_inv, eps_filter)
484 CALL transform_m_orth(matrix_nsc(ispin)%matrix, matrix_s_sqrt_inv, eps_filter)
485 END DO
486
487 !----------------------------------------
488 ! Create non-symmetric work matrices: Ax, cg, res
489 ! Content of Ax, cg, cg_z, res, z0 anti-symmetric
490 ! Content of z symmetric
491 !----------------------------------------
492
493 CALL dbcsr_create(matrix_tmp, template=matrix_s(1)%matrix, matrix_type=dbcsr_type_no_symmetry)
494
495 NULLIFY (matrix_ax, matrix_cg, matrix_cg_z, matrix_res, matrix_z, matrix_z0)
496 CALL dbcsr_allocate_matrix_set(matrix_ax, nspins)
497 CALL dbcsr_allocate_matrix_set(matrix_cg, nspins)
498 CALL dbcsr_allocate_matrix_set(matrix_cg_z, nspins)
499 CALL dbcsr_allocate_matrix_set(matrix_res, nspins)
500 CALL dbcsr_allocate_matrix_set(matrix_z, nspins)
501 CALL dbcsr_allocate_matrix_set(matrix_z0, nspins)
502
503 DO ispin = 1, nspins
504 ALLOCATE (matrix_ax(ispin)%matrix)
505 ALLOCATE (matrix_cg(ispin)%matrix)
506 ALLOCATE (matrix_cg_z(ispin)%matrix)
507 ALLOCATE (matrix_res(ispin)%matrix)
508 ALLOCATE (matrix_z(ispin)%matrix)
509 ALLOCATE (matrix_z0(ispin)%matrix)
510 CALL dbcsr_create(matrix_ax(ispin)%matrix, name="linop MATRIX", &
511 template=matrix_s(1)%matrix, &
512 matrix_type=dbcsr_type_no_symmetry)
513 CALL dbcsr_create(matrix_cg(ispin)%matrix, name="TRIAL MATRIX", &
514 template=matrix_s(1)%matrix, &
515 matrix_type=dbcsr_type_no_symmetry)
516 CALL dbcsr_create(matrix_cg_z(ispin)%matrix, name="MATRIX CG-Z", &
517 template=matrix_s(1)%matrix, &
518 matrix_type=dbcsr_type_no_symmetry)
519 CALL dbcsr_create(matrix_res(ispin)%matrix, name="RESIDUE", &
520 template=matrix_s(1)%matrix, &
521 matrix_type=dbcsr_type_no_symmetry)
522 CALL dbcsr_create(matrix_z(ispin)%matrix, name="Z-Matrix", &
523 template=matrix_s(1)%matrix, &
524 matrix_type=dbcsr_type_no_symmetry)
525 CALL dbcsr_create(matrix_z0(ispin)%matrix, name="p after precondi-Matrix", &
526 template=matrix_s(1)%matrix, &
527 matrix_type=dbcsr_type_no_symmetry)
528 END DO
529
530 !----------------------------------------
531 ! Get righ-hand-side operators
532 !----------------------------------------
533
534 ! Spin factor
535 focc = -2.0_dp
536 IF (nspins == 1) focc = -4.0_dp
537
538 ! E^[1]_Harris = -4*G[\delta P]*Pin - Pin*G[\delta P] = -4*[G[\delta P], Pin]
539 CALL commutator(matrix_nsc, matrix_p, matrix_res, eps_filter, .false., alpha=focc)
540
541 ! Initial guess cg_Z
542 DO ispin = 1, nspins
543 CALL dbcsr_copy(matrix_cg_z(ispin)%matrix, matrix_res(ispin)%matrix)
544 END DO
545
546 ! Projector on trial matrix
547 CALL projector(qs_env, matrix_p, matrix_cg_z, eps_filter)
548
549 ! Ax0
550 CALL build_hessian_op(qs_env=qs_env, &
551 p_env=p_env, &
552 matrix_ks=matrix_ks, &
553 matrix_p=matrix_p, & ! p
554 matrix_s_sqrt_inv=matrix_s_sqrt_inv, &
555 matrix_cg=matrix_cg_z, & ! cg
556 matrix_ax=matrix_ax, &
557 eps_filter=eps_filter)
558
559 ! r_0 = b - Ax0
560 DO ispin = 1, nspins
561 CALL dbcsr_add(matrix_res(ispin)%matrix, matrix_ax(ispin)%matrix, 1.0_dp, -1.0_dp)
562 END DO
563
564 ! Matrix projector T
565 CALL projector(qs_env, matrix_p, matrix_res, eps_filter)
566
567 ! Preconditioner
568 linres_control%flag = ""
569 IF (linres_control%preconditioner_type == precond_mlp) THEN
570 ! M * z_0 = r_0
571 ! Conjugate gradient returns z_0
572 CALL ec_preconditioner(qs_env=qs_env, &
573 matrix_ks=matrix_ks, &
574 matrix_p=matrix_p, &
575 matrix_rhs=matrix_res, &
576 matrix_cg_z=matrix_z0, &
577 eps_filter=eps_filter, &
578 iounit=iounit, silent=silent)
579 linres_control%flag = "PCG-AO"
580 ELSE
581 ! z_0 = r_0
582 DO ispin = 1, nspins
583 CALL dbcsr_copy(matrix_z0(ispin)%matrix, matrix_res(ispin)%matrix)
584 linres_control%flag = "CG-AO"
585 END DO
586 END IF
587
588 norm_res = 0.0_dp
589
590 DO ispin = 1, nspins
591 ! cg = p_0 = z_0
592 CALL dbcsr_copy(matrix_cg(ispin)%matrix, matrix_z0(ispin)%matrix)
593
594 ! Tr(r_0 * z_0)
595 CALL dbcsr_dot(matrix_res(ispin)%matrix, matrix_cg(ispin)%matrix, norm_rr(ispin))
596
597 IF (norm_rr(ispin) < 0.0_dp) cpabort("norm_rr < 0")
598 norm_res = max(norm_res, abs(norm_rr(ispin)/real(nao, dp)))
599 END DO
600
601 ! eigenvalue shifting
602 min_shift = 0.0_dp
603 old_conv = norm_rr(1)
604 shift = min(10.0_dp, max(min_shift, 0.05_dp*old_conv))
605 old_conv = 100.0_dp
606
607 ! header
608 IF (iounit > 0 .AND. .NOT. my_silent) THEN
609 WRITE (iounit, "(/,T3,A,T16,A,T25,A,T38,A,T52,A,/,T3,A)") &
610 "Iteration", "Method", "Stepsize", "Convergence", "Time", &
611 repeat("-", 80)
612 END IF
613
614 alpha(:) = 0.0_dp
615 restart = .false.
616 should_stop = .false.
617 linres_control%converged = .false.
618
619 ! start iteration
620 iteration: DO i = 1, linres_control%max_iter
621
622 ! Convergence criteria
623 ! default for eps 10E-6 in MO_linres
624 IF (norm_res < linres_control%eps) THEN
625 linres_control%converged = .true.
626 END IF
627
628 t2 = m_walltime()
629 IF (i == 1 .OR. mod(i, 1) == 0 .OR. linres_control%converged &
630 .OR. restart .OR. should_stop) THEN
631 IF (iounit > 0 .AND. .NOT. my_silent) THEN
632 WRITE (iounit, "(T5,I5,T18,A3,T28,L1,T38,1E8.2,T48,F16.10,T68,F8.2)") &
633 i, linres_control%flag, restart, maxval(alpha), norm_res, t2 - t1
634 CALL m_flush(iounit)
635 END IF
636 END IF
637 IF (linres_control%converged) THEN
638 IF (iounit > 0) THEN
639 WRITE (iounit, "(/,T2,A,I4,A,T73,F8.2,/)") "The linear solver converged in ", &
640 i, " iterations.", t2 - t1
641 CALL m_flush(iounit)
642 END IF
643 EXIT iteration
644 ELSE IF (should_stop) THEN
645 IF (iounit > 0) THEN
646 WRITE (iounit, "(/,T2,A,I4,A,/)") "The linear solver did NOT converge! External stop"
647 CALL m_flush(iounit)
648 END IF
649 EXIT iteration
650 END IF
651
652 ! Max number of iteration reached
653 IF (i == linres_control%max_iter) THEN
654 IF (iounit > 0) THEN
655 WRITE (iounit, "(/,T2,A/)") &
656 "The linear solver didnt converge! Maximum number of iterations reached."
657 CALL m_flush(iounit)
658 END IF
659 linres_control%converged = .false.
660 END IF
661
662 ! Hessian Ax = [F,B] + [G(B),P]
663 CALL build_hessian_op(qs_env=qs_env, &
664 p_env=p_env, &
665 matrix_ks=matrix_ks, &
666 matrix_p=matrix_p, & ! p
667 matrix_s_sqrt_inv=matrix_s_sqrt_inv, &
668 matrix_cg=matrix_cg, & ! cg
669 matrix_ax=matrix_ax, &
670 eps_filter=eps_filter)
671
672 ! Matrix projector T
673 CALL projector(qs_env, matrix_p, matrix_ax, eps_filter)
674
675 DO ispin = 1, nspins
676
677 CALL dbcsr_filter(matrix_ax(ispin)%matrix, eps_filter)
678 ! norm_cA = tr(Ap_j * p_j)
679 CALL dbcsr_dot(matrix_cg(ispin)%matrix, matrix_ax(ispin)%matrix, norm_ca(ispin))
680
681 IF (norm_ca(ispin) < 0.0_dp) THEN
682
683 ! Recalculate w/o preconditioner
684 IF (i > 1) THEN
685 ! p = -z + beta*p
686 CALL dbcsr_add(matrix_cg(ispin)%matrix, matrix_z0(ispin)%matrix, &
687 beta(ispin), -1.0_dp)
688 CALL dbcsr_dot(matrix_res(ispin)%matrix, matrix_res(ispin)%matrix, new_norm(ispin))
689 beta(ispin) = new_norm(ispin)/tr_rz00(ispin)
690 CALL dbcsr_add(matrix_cg(ispin)%matrix, matrix_res(ispin)%matrix, &
691 beta(ispin), 1.0_dp)
692 norm_rr(ispin) = new_norm(ispin)
693 ELSE
694 CALL dbcsr_copy(matrix_res(ispin)%matrix, matrix_cg(ispin)%matrix)
695 CALL dbcsr_dot(matrix_res(ispin)%matrix, matrix_res(ispin)%matrix, norm_rr(ispin))
696 END IF
697
698 CALL build_hessian_op(qs_env=qs_env, &
699 p_env=p_env, &
700 matrix_ks=matrix_ks, &
701 matrix_p=matrix_p, & ! p
702 matrix_s_sqrt_inv=matrix_s_sqrt_inv, &
703 matrix_cg=matrix_cg, & ! cg
704 matrix_ax=matrix_ax, &
705 eps_filter=eps_filter)
706
707 ! Matrix projector T
708 CALL projector(qs_env, matrix_p, matrix_ax, eps_filter)
709
710 CALL dbcsr_dot(matrix_cg(ispin)%matrix, matrix_ax(ispin)%matrix, norm_ca(ispin))
711
712 cpabort("tr(Ap_j*p_j) < 0")
713 IF (abnormal_value(norm_ca(ispin))) &
714 cpabort("Preconditioner: Tr[Ap_j*p_j] is an abnormal value (NaN/Inf)")
715
716 END IF
717
718 END DO
719
720 DO ispin = 1, nspins
721 ! Determine step-size
722 IF (norm_ca(ispin) < linres_control%eps) THEN
723 alpha(ispin) = 1.0_dp
724 ELSE
725 alpha(ispin) = norm_rr(ispin)/norm_ca(ispin)
726 END IF
727
728 ! x_j+1 = x_j + alpha*p_j
729 ! save response-denisty of this iteration
730 CALL dbcsr_add(matrix_cg_z(ispin)%matrix, matrix_cg(ispin)%matrix, 1.0_dp, alpha(ispin))
731 END DO
732
733 ! need to recompute the residue
734 restart = .false.
735 IF (mod(i, linres_control%restart_every) == 0) THEN
736 !
737 ! r_j+1 = b - A * x_j+1
738 CALL build_hessian_op(qs_env=qs_env, &
739 p_env=p_env, &
740 matrix_ks=matrix_ks, &
741 matrix_p=matrix_p, &
742 matrix_s_sqrt_inv=matrix_s_sqrt_inv, &
743 matrix_cg=matrix_cg_z, & ! cg
744 matrix_ax=matrix_ax, &
745 eps_filter=eps_filter)
746 ! b
747 CALL commutator(matrix_nsc, matrix_p, matrix_res, eps_filter, .false., alpha=focc)
748
749 DO ispin = 1, nspins
750 CALL dbcsr_add(matrix_res(ispin)%matrix, matrix_ax(ispin)%matrix, 1.0_dp, -1.0_dp)
751 END DO
752
753 CALL projector(qs_env, matrix_p, matrix_res, eps_filter)
754 !
755 restart = .true.
756 ELSE
757 ! proj Ap onto the virtual subspace
758 CALL projector(qs_env, matrix_p, matrix_ax, eps_filter)
759 !
760 ! r_j+1 = r_j - alpha * Ap_j
761 DO ispin = 1, nspins
762 CALL dbcsr_add(matrix_res(ispin)%matrix, matrix_ax(ispin)%matrix, 1.0_dp, -alpha(ispin))
763 END DO
764 restart = .false.
765 END IF
766
767 ! Preconditioner
768 linres_control%flag = ""
769 IF (linres_control%preconditioner_type == precond_mlp) THEN
770 ! M * z_j+1 = r_j+1
771 ! Conjugate gradient returns z_j+1
772 CALL ec_preconditioner(qs_env=qs_env, &
773 matrix_ks=matrix_ks, &
774 matrix_p=matrix_p, &
775 matrix_rhs=matrix_res, &
776 matrix_cg_z=matrix_z0, &
777 eps_filter=eps_filter, &
778 iounit=iounit, silent=silent)
779 linres_control%flag = "PCG-AO"
780 ELSE
781 DO ispin = 1, nspins
782 CALL dbcsr_copy(matrix_z0(ispin)%matrix, matrix_res(ispin)%matrix)
783 END DO
784 linres_control%flag = "CG-AO"
785 END IF
786
787 norm_res = 0.0_dp
788
789 DO ispin = 1, nspins
790 ! Tr[r_j+1*z_j+1]
791 CALL dbcsr_dot(matrix_res(ispin)%matrix, matrix_z0(ispin)%matrix, new_norm(ispin))
792 IF (new_norm(ispin) < 0.0_dp) cpabort("tr(r_j+1*z_j+1) < 0")
793 IF (abnormal_value(new_norm(ispin))) &
794 cpabort("Preconditioner: Tr[r_j+1*z_j+1] is an abnormal value (NaN/Inf)")
795 norm_res = max(norm_res, new_norm(ispin)/real(nao, dp))
796
797 IF (norm_rr(ispin) < linres_control%eps .OR. new_norm(ispin) < linres_control%eps) THEN
798 beta(ispin) = 0.0_dp
799 linres_control%converged = .true.
800 ELSE
801 beta(ispin) = new_norm(ispin)/norm_rr(ispin)
802 END IF
803
804 ! update new search vector (matrix cg)
805 ! Here: cg_j+1 = z_j+1 + beta*cg_j
806 CALL dbcsr_add(matrix_cg(ispin)%matrix, matrix_z0(ispin)%matrix, beta(ispin), 1.0_dp)
807 CALL dbcsr_filter(matrix_cg(ispin)%matrix, eps_filter)
808
809 tr_rz00(ispin) = norm_rr(ispin)
810 norm_rr(ispin) = new_norm(ispin)
811 END DO
812
813 ! Can we exit the loop?
814 CALL external_control(should_stop, "LS_SOLVER", target_time=qs_env%target_time, &
815 start_time=qs_env%start_time)
816
817 END DO iteration
818
819 ! Matrix projector
820 CALL projector(qs_env, matrix_p, matrix_cg_z, eps_filter)
821
822 ! Z = [cg_z,P]
823 CALL commutator(matrix_cg_z, matrix_p, matrix_z, eps_filter, .true., alpha=0.5_dp)
824
825 DO ispin = 1, nspins
826 ! Transform Z-matrix back into non-orthogonal basis
827 CALL transform_m_orth(matrix_z(ispin)%matrix, matrix_s_sqrt_inv, eps_filter)
828
829 ! Export Z-Matrix
830 CALL dbcsr_copy(matrix_pz(ispin)%matrix, matrix_z(ispin)%matrix, keep_sparsity=.true.)
831 END DO
832
833 ! Calculate energy-weighted response density matrix
834 ! AO: Wz = 0.5*(Z*KS*P + P*KS*Z)
835 CALL ec_wz_matrix(qs_env, matrix_pz, matrix_wz, eps_filter)
836
837 ! Release matrices
838 CALL dbcsr_release(matrix_tmp)
839
840 CALL dbcsr_release(matrix_s_sqrt)
841 CALL dbcsr_release(matrix_s_sqrt_inv)
842
843 CALL dbcsr_deallocate_matrix_set(matrix_p)
844 CALL dbcsr_deallocate_matrix_set(matrix_ks)
845 CALL dbcsr_deallocate_matrix_set(matrix_nsc)
846 CALL dbcsr_deallocate_matrix_set(matrix_z)
847 CALL dbcsr_deallocate_matrix_set(matrix_ax)
848 CALL dbcsr_deallocate_matrix_set(matrix_res)
849 CALL dbcsr_deallocate_matrix_set(matrix_cg)
850 CALL dbcsr_deallocate_matrix_set(matrix_cg_z)
851 CALL dbcsr_deallocate_matrix_set(matrix_z0)
852
853 DEALLOCATE (alpha, beta, new_norm, norm_ca, norm_rr)
854 DEALLOCATE (tr_rz00)
855
856 CALL timestop(handle)
857
858 END SUBROUTINE ec_response_ao
859
860! **************************************************************************************************
861!> \brief Compute matrix_wz as needed for the forces
862!> Wz = 0.5*(Z*KS*P + P*KS*Z) (closed-shell)
863!> \param qs_env ...
864!> \param matrix_z The response density we just calculated
865!> \param matrix_wz The energy weighted response-density matrix
866!> \param eps_filter ...
867!> \par History
868!> 2020.2 created [Fabian Belleflamme]
869!> \author Fabian Belleflamme
870! **************************************************************************************************
871 SUBROUTINE ec_wz_matrix(qs_env, matrix_z, matrix_wz, eps_filter)
872
873 TYPE(qs_environment_type), POINTER :: qs_env
874 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN), &
875 POINTER :: matrix_z
876 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(INOUT), &
877 POINTER :: matrix_wz
878 REAL(kind=dp), INTENT(IN) :: eps_filter
879
880 CHARACTER(len=*), PARAMETER :: routinen = 'ec_wz_matrix'
881
882 INTEGER :: handle, ispin, nspins
883 REAL(kind=dp) :: scaling
884 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_ks, matrix_p, matrix_s
885 TYPE(dbcsr_type) :: matrix_tmp, matrix_tmp2
886 TYPE(dft_control_type), POINTER :: dft_control
887 TYPE(qs_rho_type), POINTER :: rho
888
889 CALL timeset(routinen, handle)
890
891 cpassert(ASSOCIATED(qs_env))
892 cpassert(ASSOCIATED(matrix_z))
893 cpassert(ASSOCIATED(matrix_wz))
894
895 CALL get_qs_env(qs_env=qs_env, &
896 dft_control=dft_control, &
897 matrix_ks=matrix_ks, &
898 matrix_s=matrix_s, &
899 rho=rho)
900 nspins = dft_control%nspins
901
902 CALL qs_rho_get(rho, rho_ao=matrix_p)
903
904 ! Init temp matrices
905 CALL dbcsr_create(matrix_tmp, template=matrix_z(1)%matrix, &
906 matrix_type=dbcsr_type_no_symmetry)
907 CALL dbcsr_create(matrix_tmp2, template=matrix_z(1)%matrix, &
908 matrix_type=dbcsr_type_no_symmetry)
909
910 ! Scale matrix_p by factor 1/2 in closed-shell
911 scaling = 1.0_dp
912 IF (nspins == 1) scaling = 0.5_dp
913
914 ! Whz = ZFP + PFZ = Z(FP) + (Z(FP))^T
915 DO ispin = 1, nspins
916
917 ! tmp = FP
918 CALL dbcsr_multiply("N", "N", scaling, matrix_ks(ispin)%matrix, matrix_p(ispin)%matrix, &
919 0.0_dp, matrix_tmp, filter_eps=eps_filter, retain_sparsity=.false.)
920
921 ! tmp2 = ZFP
922 CALL dbcsr_multiply("N", "N", 1.0_dp, matrix_z(ispin)%matrix, matrix_tmp, &
923 0.0_dp, matrix_tmp2, filter_eps=eps_filter, retain_sparsity=.false.)
924
925 ! tmp = (ZFP)^T
926 CALL dbcsr_transposed(matrix_tmp, matrix_tmp2)
927
928 ! tmp = ZFP + (ZFP)^T
929 CALL dbcsr_add(matrix_tmp, matrix_tmp2, 1.0_dp, 1.0_dp)
930
931 CALL dbcsr_filter(matrix_tmp, eps_filter)
932
933 ! Whz = ZFP + PFZ
934 CALL dbcsr_copy(matrix_wz(ispin)%matrix, matrix_tmp, keep_sparsity=.true.)
935
936 END DO
937
938 ! Release matrices
939 CALL dbcsr_release(matrix_tmp)
940 CALL dbcsr_release(matrix_tmp2)
941
942 CALL timestop(handle)
943
944 END SUBROUTINE ec_wz_matrix
945
946! **************************************************************************************************
947!> \brief Calculate first term of electronic Hessian M = [F, B]
948!> acting as liner transformation on trial matrix (matrix_cg)
949!> with intermediate response density B = [cg,P] = cg*P - P*cg = cg*P + (cg*P)^T
950!>
951!> All matrices are in orthonormal basis
952!>
953!> \param matrix_ks Ground-state Kohn-Sham matrix
954!> \param matrix_p Ground-state Density matrix
955!> \param matrix_cg Trial matrix
956!> \param matrix_b Intermediate response density
957!> \param matrix_Ax First term of electronic Hessian applied on trial matrix (matrix_cg)
958!>
959!> \param eps_filter ...
960!> \date 12.2019
961!> \author Fabian Belleflamme
962! **************************************************************************************************
963 SUBROUTINE hessian_op1(matrix_ks, matrix_p, matrix_cg, matrix_b, matrix_Ax, eps_filter)
964
965 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN), &
966 POINTER :: matrix_ks, matrix_p, matrix_cg
967 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(INOUT), &
968 POINTER :: matrix_b, matrix_ax
969 REAL(kind=dp), INTENT(IN) :: eps_filter
970
971 CHARACTER(len=*), PARAMETER :: routinen = 'hessian_op1'
972
973 INTEGER :: handle
974
975 CALL timeset(routinen, handle)
976
977 cpassert(ASSOCIATED(matrix_ks))
978 cpassert(ASSOCIATED(matrix_p))
979 cpassert(ASSOCIATED(matrix_cg))
980 cpassert(ASSOCIATED(matrix_b))
981 cpassert(ASSOCIATED(matrix_ax))
982
983 ! Build intermediate density matrix
984 ! B = [cg, P] = cg*P - P*cg = cg*P + (cg*P)^T
985 CALL commutator(matrix_cg, matrix_p, matrix_b, eps_filter, .true.)
986
987 ! Build first part of operator
988 ! Ax = [F,[cg,P]] = [F,B]
989 CALL commutator(matrix_ks, matrix_b, matrix_ax, eps_filter, .false.)
990
991 CALL timestop(handle)
992
993 END SUBROUTINE hessian_op1
994
995! **************************************************************************************************
996!> \brief calculate linear transformation of Hessian matrix on a trial matrix matrix_cg
997!> which is stored in response density B = [cg,P] = cg*P - P*cg = cg*P + (cg*P)^T
998!> Ax = [F, B] + [G(B), Pin] in orthonormal basis
999!>
1000!> \param qs_env ...
1001!> \param p_env ...
1002!> \param matrix_ks Ground-state Kohn-Sham matrix
1003!> \param matrix_p Ground-state Density matrix
1004!> \param matrix_s_sqrt_inv S^(-1/2) needed for transformation to/from orthonormal basis
1005!> \param matrix_cg Trial matrix
1006!> \param matrix_Ax Electronic Hessian applied on trial matrix (matrix_cg)
1007!> \param eps_filter ...
1008!>
1009!> \date 12.2019
1010!> \author Fabian Belleflamme
1011! **************************************************************************************************
1012 SUBROUTINE build_hessian_op(qs_env, p_env, matrix_ks, matrix_p, matrix_s_sqrt_inv, &
1013 matrix_cg, matrix_Ax, eps_filter)
1014
1015 TYPE(qs_environment_type), POINTER :: qs_env
1016 TYPE(qs_p_env_type), POINTER :: p_env
1017 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN), &
1018 POINTER :: matrix_ks, matrix_p
1019 TYPE(dbcsr_type), INTENT(IN) :: matrix_s_sqrt_inv
1020 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN), &
1021 POINTER :: matrix_cg
1022 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(INOUT), &
1023 POINTER :: matrix_ax
1024 REAL(kind=dp), INTENT(IN) :: eps_filter
1025
1026 CHARACTER(len=*), PARAMETER :: routinen = 'build_hessian_op'
1027
1028 INTEGER :: handle, ispin, nspins
1029 REAL(kind=dp) :: chksum
1030 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_b, rho1_ao
1031 TYPE(dft_control_type), POINTER :: dft_control
1032 TYPE(mp_para_env_type), POINTER :: para_env
1033 TYPE(qs_rho_type), POINTER :: rho
1034
1035 CALL timeset(routinen, handle)
1036
1037 cpassert(ASSOCIATED(qs_env))
1038 cpassert(ASSOCIATED(matrix_ks))
1039 cpassert(ASSOCIATED(matrix_p))
1040 cpassert(ASSOCIATED(matrix_cg))
1041 cpassert(ASSOCIATED(matrix_ax))
1042
1043 CALL get_qs_env(qs_env=qs_env, &
1044 dft_control=dft_control, &
1045 para_env=para_env, &
1046 rho=rho)
1047 nspins = dft_control%nspins
1048
1049 NULLIFY (matrix_b)
1050 CALL dbcsr_allocate_matrix_set(matrix_b, nspins)
1051 DO ispin = 1, nspins
1052 ALLOCATE (matrix_b(ispin)%matrix)
1053 CALL dbcsr_create(matrix_b(ispin)%matrix, name="[X,P] RSP DNSTY", &
1054 template=matrix_p(1)%matrix, &
1055 matrix_type=dbcsr_type_no_symmetry)
1056 END DO
1057
1058 ! Build uncoupled term of Hessian linear transformation
1059 CALL hessian_op1(matrix_ks, matrix_p, matrix_cg, matrix_b, matrix_ax, eps_filter)
1060
1061 ! Avoid the buildup of noisy blocks
1062 DO ispin = 1, nspins
1063 CALL dbcsr_filter(matrix_b(ispin)%matrix, eps_filter)
1064 END DO
1065
1066 chksum = 0.0_dp
1067 DO ispin = 1, nspins
1068 chksum = chksum + dbcsr_checksum(matrix_b(ispin)%matrix)
1069 END DO
1070
1071 ! skip the kernel if the DM is very small
1072 IF (chksum > 1.0e-14_dp) THEN
1073
1074 ! Bring matrix B as density on grid
1075
1076 ! prepare perturbation environment
1077 CALL p_env_check_i_alloc(p_env, qs_env)
1078
1079 ! Get response density matrix
1080 CALL qs_rho_get(p_env%rho1, rho_ao=rho1_ao)
1081
1082 DO ispin = 1, nspins
1083 ! Transform B into NON-ortho basis for collocation
1084 CALL transform_m_orth(matrix_b(ispin)%matrix, matrix_s_sqrt_inv, eps_filter)
1085 ! Filter
1086 CALL dbcsr_filter(matrix_b(ispin)%matrix, eps_filter)
1087 ! Keep symmetry of density matrix
1088 CALL dbcsr_copy(rho1_ao(ispin)%matrix, matrix_b(ispin)%matrix, keep_sparsity=.true.)
1089 CALL dbcsr_copy(p_env%p1(ispin)%matrix, matrix_b(ispin)%matrix, keep_sparsity=.true.)
1090 END DO
1091
1092 ! Updates densities on grid wrt density matrix
1093 CALL p_env_update_rho(p_env, qs_env)
1094
1095 DO ispin = 1, nspins
1096 CALL dbcsr_set(p_env%kpp1(ispin)%matrix, 0.0_dp)
1097 IF (ASSOCIATED(p_env%kpp1_admm)) CALL dbcsr_set(p_env%kpp1_admm(ispin)%matrix, 0.0_dp)
1098 END DO
1099
1100 ! Calculate kernel
1101 ! Ax = F*B - B*F + G(B)*P - P*G(B)
1102 ! IN/OUT IN IN IN
1103 CALL hessian_op2(qs_env, p_env, matrix_ax, matrix_p, matrix_s_sqrt_inv, eps_filter)
1104
1105 END IF
1106
1107 CALL dbcsr_deallocate_matrix_set(matrix_b)
1108
1109 CALL timestop(handle)
1110
1111 END SUBROUTINE build_hessian_op
1112
1113! **************************************************************************************************
1114!> \brief Calculate lin transformation of Hessian matrix on a trial matrix matrix_cg
1115!> which is stored in response density B = [cg,P] = cg*P - P*cg = cg*P + (cg*P)^T
1116!> Ax = [F, B] + [G(B), Pin] in orthonormal basis
1117!>
1118!> \param qs_env ...
1119!> \param p_env p-environment with trial density environment
1120!> \param matrix_Ax contains first part of Hessian linear transformation, kernel contribution
1121!> is calculated and added in this routine
1122!> \param matrix_p Density matrix in orthogonal basis
1123!> \param matrix_s_sqrt_inv contains matrix S^(-1/2) for switching to orthonormal Lowdin basis
1124!> \param eps_filter ...
1125!>
1126!> \date 12.2019
1127!> \author Fabian Belleflamme
1128! **************************************************************************************************
1129 SUBROUTINE hessian_op2(qs_env, p_env, matrix_Ax, matrix_p, matrix_s_sqrt_inv, eps_filter)
1130
1131 TYPE(qs_environment_type), POINTER :: qs_env
1132 TYPE(qs_p_env_type), POINTER :: p_env
1133 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(INOUT), &
1134 POINTER :: matrix_ax
1135 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN), &
1136 POINTER :: matrix_p
1137 TYPE(dbcsr_type), INTENT(IN) :: matrix_s_sqrt_inv
1138 REAL(kind=dp), INTENT(IN) :: eps_filter
1139
1140 CHARACTER(len=*), PARAMETER :: routinen = 'hessian_op2'
1141
1142 INTEGER :: handle, ispin, nspins
1143 REAL(kind=dp) :: ekin_mol
1144 TYPE(admm_type), POINTER :: admm_env
1145 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_g, matrix_s, rho1_ao, rho_ao
1146 TYPE(dft_control_type), POINTER :: dft_control
1147 TYPE(mp_para_env_type), POINTER :: para_env
1148 TYPE(pw_c1d_gs_type) :: rho_tot_gspace, v_hartree_gspace
1149 TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER :: rho1_g
1150 TYPE(pw_env_type), POINTER :: pw_env
1151 TYPE(pw_poisson_type), POINTER :: poisson_env
1152 TYPE(pw_pool_p_type), DIMENSION(:), POINTER :: pw_pools
1153 TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
1154 TYPE(pw_r3d_rs_type) :: v_hartree_rspace
1155 TYPE(pw_r3d_rs_type), DIMENSION(:), POINTER :: rho1_r, rho_r, tau1_r, v_xc, v_xc_tau
1156 TYPE(pw_r3d_rs_type), POINTER :: weights
1157 TYPE(qs_kpp1_env_type), POINTER :: kpp1_env
1158 TYPE(qs_rho_type), POINTER :: rho, rho_aux
1159 TYPE(section_vals_type), POINTER :: input, xc_section, xc_section_aux
1160
1161 CALL timeset(routinen, handle)
1162
1163 NULLIFY (admm_env, dft_control, input, matrix_s, para_env, rho, rho_r, rho1_g, rho1_r)
1164
1165 CALL get_qs_env(qs_env=qs_env, &
1166 admm_env=admm_env, &
1167 dft_control=dft_control, &
1168 input=input, &
1169 matrix_s=matrix_s, &
1170 para_env=para_env, &
1171 rho=rho)
1172 nspins = dft_control%nspins
1173
1174 cpassert(ASSOCIATED(p_env%kpp1))
1175 cpassert(ASSOCIATED(p_env%kpp1_env))
1176 kpp1_env => p_env%kpp1_env
1177
1178 ! Get non-ortho input density matrix on grid
1179 CALL qs_rho_get(rho, rho_ao=rho_ao)
1180 ! Get non-ortho trial density stored in p_env
1181 CALL qs_rho_get(p_env%rho1, rho_g=rho1_g, rho_r=rho1_r, tau_r=tau1_r)
1182
1183 NULLIFY (pw_env)
1184 CALL get_qs_env(qs_env, pw_env=pw_env)
1185 cpassert(ASSOCIATED(pw_env))
1186
1187 NULLIFY (weights)
1188 CALL get_qs_env(qs_env, xcint_weights=weights)
1189
1190 NULLIFY (auxbas_pw_pool, poisson_env, pw_pools)
1191 ! gets the tmp grids
1192 CALL pw_env_get(pw_env=pw_env, &
1193 auxbas_pw_pool=auxbas_pw_pool, &
1194 pw_pools=pw_pools, &
1195 poisson_env=poisson_env)
1196
1197 ! Calculate the NSC Hartree potential
1198 CALL auxbas_pw_pool%create_pw(pw=v_hartree_gspace)
1199 CALL auxbas_pw_pool%create_pw(pw=rho_tot_gspace)
1200 CALL auxbas_pw_pool%create_pw(pw=v_hartree_rspace)
1201
1202 ! XC-Kernel
1203 NULLIFY (v_xc, v_xc_tau, xc_section)
1204
1205 IF (dft_control%do_admm) THEN
1206 xc_section => admm_env%xc_section_primary
1207 ELSE
1208 xc_section => section_vals_get_subs_vals(input, "DFT%XC")
1209 END IF
1210
1211 ! add xc-kernel
1212 CALL create_kernel(qs_env, &
1213 vxc=v_xc, &
1214 vxc_tau=v_xc_tau, &
1215 rho=rho, &
1216 rho1_r=rho1_r, &
1217 rho1_g=rho1_g, &
1218 tau1_r=tau1_r, &
1219 xc_section=xc_section)
1220
1221 DO ispin = 1, nspins
1222 CALL pw_scale(v_xc(ispin), v_xc(ispin)%pw_grid%dvol)
1223 IF (ASSOCIATED(v_xc_tau)) THEN
1224 CALL pw_scale(v_xc_tau(ispin), v_xc_tau(ispin)%pw_grid%dvol)
1225 END IF
1226 END DO
1227
1228 ! ADMM Correction
1229 IF (dft_control%do_admm) THEN
1230 IF (admm_env%aux_exch_func /= do_admm_aux_exch_func_none) THEN
1231 IF (.NOT. ASSOCIATED(kpp1_env%deriv_set_admm)) THEN
1232 xc_section_aux => admm_env%xc_section_aux
1233 CALL get_admm_env(qs_env%admm_env, rho_aux_fit=rho_aux)
1234 CALL qs_rho_get(rho_aux, rho_r=rho_r)
1235 ALLOCATE (kpp1_env%deriv_set_admm, kpp1_env%rho_set_admm)
1236 CALL xc_prep_2nd_deriv(kpp1_env%deriv_set_admm, kpp1_env%rho_set_admm, &
1237 rho_r, auxbas_pw_pool, weights, &
1238 xc_section=xc_section_aux)
1239 END IF
1240 END IF
1241 END IF
1242
1243 ! take trial density to build G^{H}[B]
1244 CALL pw_zero(rho_tot_gspace)
1245 DO ispin = 1, nspins
1246 CALL pw_axpy(rho1_g(ispin), rho_tot_gspace)
1247 END DO
1248
1249 ! get Hartree potential from rho_tot_gspace
1250 CALL pw_poisson_solve(poisson_env, rho_tot_gspace, &
1251 vhartree=v_hartree_gspace)
1252 CALL pw_transfer(v_hartree_gspace, v_hartree_rspace)
1253 CALL pw_scale(v_hartree_rspace, v_hartree_rspace%pw_grid%dvol)
1254
1255 ! Add v_xc + v_H
1256 DO ispin = 1, nspins
1257 CALL pw_axpy(v_hartree_rspace, v_xc(ispin))
1258 END DO
1259 IF (nspins == 1) THEN
1260 CALL pw_scale(v_xc(1), 2.0_dp)
1261 IF (ASSOCIATED(v_xc_tau)) CALL pw_scale(v_xc_tau(1), 2.0_dp)
1262 END IF
1263
1264 DO ispin = 1, nspins
1265 ! Integrate with ground-state density matrix, in non-orthogonal basis
1266 CALL integrate_v_rspace(v_rspace=v_xc(ispin), &
1267 pmat=rho_ao(ispin), &
1268 hmat=p_env%kpp1(ispin), &
1269 qs_env=qs_env, &
1270 calculate_forces=.false., &
1271 basis_type="ORB")
1272 IF (ASSOCIATED(v_xc_tau)) THEN
1273 CALL integrate_v_rspace(v_rspace=v_xc_tau(ispin), &
1274 pmat=rho_ao(ispin), &
1275 hmat=p_env%kpp1(ispin), &
1276 qs_env=qs_env, &
1277 compute_tau=.true., &
1278 calculate_forces=.false., &
1279 basis_type="ORB")
1280 END IF
1281 END DO
1282
1283 ! Hartree-Fock contribution
1284 CALL apply_hfx(qs_env, p_env)
1285 ! Calculate ADMM exchange correction to kernel
1286 CALL apply_xc_admm(qs_env, p_env)
1287 ! Add contribution from ADMM exchange correction to kernel
1288 CALL p_env_finish_kpp1(qs_env, p_env)
1289
1290 ! Calculate KG correction to kernel
1291 IF (dft_control%qs_control%do_kg) THEN
1292 IF (qs_env%kg_env%tnadd_method == kg_tnadd_embed .OR. &
1293 qs_env%kg_env%tnadd_method == kg_tnadd_embed_ri) THEN
1294
1295 cpassert(dft_control%nimages == 1)
1296 ekin_mol = 0.0_dp
1297 CALL qs_rho_get(p_env%rho1, rho_ao=rho1_ao)
1298 CALL kg_ekin_subset(qs_env=qs_env, &
1299 ks_matrix=p_env%kpp1, &
1300 ekin_mol=ekin_mol, &
1301 calc_force=.false., &
1302 do_kernel=.true., &
1303 pmat_ext=rho1_ao)
1304 END IF
1305 END IF
1306
1307 ! Init response kernel matrix
1308 ! matrix G(B)
1309 NULLIFY (matrix_g)
1310 CALL dbcsr_allocate_matrix_set(matrix_g, nspins)
1311 DO ispin = 1, nspins
1312 ALLOCATE (matrix_g(ispin)%matrix)
1313 CALL dbcsr_copy(matrix_g(ispin)%matrix, p_env%kpp1(ispin)%matrix, &
1314 name="MATRIX Kernel")
1315 END DO
1316
1317 ! Transforming G(B) into orthonormal basis
1318 ! Careful, this de-symmetrizes matrix_G
1319 DO ispin = 1, nspins
1320 CALL transform_m_orth(matrix_g(ispin)%matrix, matrix_s_sqrt_inv, eps_filter)
1321 CALL dbcsr_filter(matrix_g(ispin)%matrix, eps_filter)
1322 END DO
1323
1324 ! Hessian already contains Ax = [F,B] (ORTHO), now adding
1325 ! Ax = Ax + G(B)P - (G(B)P)^T
1326 CALL commutator(matrix_g, matrix_p, matrix_ax, eps_filter, .false., 1.0_dp, 1.0_dp)
1327
1328 ! release pw grids
1329 CALL auxbas_pw_pool%give_back_pw(v_hartree_gspace)
1330 CALL auxbas_pw_pool%give_back_pw(v_hartree_rspace)
1331 CALL auxbas_pw_pool%give_back_pw(rho_tot_gspace)
1332 DO ispin = 1, nspins
1333 CALL auxbas_pw_pool%give_back_pw(v_xc(ispin))
1334 END DO
1335 DEALLOCATE (v_xc)
1336 IF (ASSOCIATED(v_xc_tau)) THEN
1337 DO ispin = 1, nspins
1338 CALL auxbas_pw_pool%give_back_pw(v_xc_tau(ispin))
1339 END DO
1340 DEALLOCATE (v_xc_tau)
1341 END IF
1342
1343 CALL dbcsr_deallocate_matrix_set(matrix_g)
1344
1345 CALL timestop(handle)
1346
1347 END SUBROUTINE hessian_op2
1348
1349! **************************************************************************************************
1350!> \brief computes (anti-)commutator exploiting (anti-)symmetry:
1351!> A symmetric : RES = beta*RES + k*[A,B] = k*(AB-(AB)^T)
1352!> A anti-sym : RES = beta*RES + k*{A,B} = k*(AB+(AB)^T)
1353!>
1354!> \param a Matrix A
1355!> \param b Matrix B
1356!> \param res Commutator result
1357!> \param eps_filter filtering threshold for sparse matrices
1358!> \param anticomm Calculate anticommutator
1359!> \param alpha Scaling of anti-/commutator
1360!> \param beta Scaling of inital content of result matrix
1361!>
1362!> \par History
1363!> 2020.07 Fabian Belleflamme (based on commutator_symm)
1364! **************************************************************************************************
1365 SUBROUTINE commutator(a, b, res, eps_filter, anticomm, alpha, beta)
1366
1367 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN), &
1368 POINTER :: a, b, res
1369 REAL(kind=dp) :: eps_filter
1370 LOGICAL :: anticomm
1371 REAL(kind=dp), OPTIONAL :: alpha, beta
1372
1373 CHARACTER(LEN=*), PARAMETER :: routinen = 'commutator'
1374
1375 INTEGER :: handle, ispin
1376 REAL(kind=dp) :: facc, myalpha, mybeta
1377 TYPE(dbcsr_type) :: work, work2
1378
1379 CALL timeset(routinen, handle)
1380
1381 cpassert(ASSOCIATED(a))
1382 cpassert(ASSOCIATED(b))
1383 cpassert(ASSOCIATED(res))
1384
1385 CALL dbcsr_create(work, template=a(1)%matrix, matrix_type=dbcsr_type_no_symmetry)
1386 CALL dbcsr_create(work2, template=a(1)%matrix, matrix_type=dbcsr_type_no_symmetry)
1387
1388 ! Scaling of anti-/commutator
1389 myalpha = 1.0_dp
1390 IF (PRESENT(alpha)) myalpha = alpha
1391 ! Scaling of result matrix
1392 mybeta = 0.0_dp
1393 IF (PRESENT(beta)) mybeta = beta
1394 ! Add/subtract second term when calculating anti-/commutator
1395 facc = -1.0_dp
1396 IF (anticomm) facc = 1.0_dp
1397
1398 DO ispin = 1, SIZE(a)
1399
1400 CALL dbcsr_multiply("N", "N", myalpha, a(ispin)%matrix, b(ispin)%matrix, &
1401 0.0_dp, work, filter_eps=eps_filter)
1402 CALL dbcsr_transposed(work2, work)
1403
1404 ! RES= beta*RES + alpha*{A,B} = beta*RES + alpha*[AB+(AB)T]
1405 ! RES= beta*RES + alpha*[A,B] = beta*RES + alpha*[AB-(AB)T]
1406 CALL dbcsr_add(work, work2, 1.0_dp, facc)
1407
1408 CALL dbcsr_add(res(ispin)%matrix, work, mybeta, 1.0_dp)
1409
1410 END DO
1411
1412 CALL dbcsr_release(work)
1413 CALL dbcsr_release(work2)
1414
1415 CALL timestop(handle)
1416
1417 END SUBROUTINE commutator
1418
1419! **************************************************************************************************
1420!> \brief Projector P(M) = P*M*Q^T + Q*M*P^T
1421!> with P = D
1422!> with Q = (1-D)
1423!>
1424!> \param qs_env ...
1425!> \param matrix_p Ground-state density in orthonormal basis
1426!> \param matrix_io Matrix to which projector is applied.
1427!>
1428!> \param eps_filter ...
1429!> \date 06.2020
1430!> \author Fabian Belleflamme
1431! **************************************************************************************************
1432 SUBROUTINE projector(qs_env, matrix_p, matrix_io, eps_filter)
1433
1434 TYPE(qs_environment_type), POINTER :: qs_env
1435 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN), &
1436 POINTER :: matrix_p
1437 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(INOUT), &
1438 POINTER :: matrix_io
1439 REAL(kind=dp), INTENT(IN) :: eps_filter
1440
1441 CHARACTER(len=*), PARAMETER :: routinen = 'projector'
1442
1443 INTEGER :: handle, ispin, nspins
1444 TYPE(dbcsr_type) :: matrix_q, matrix_tmp
1445 TYPE(dft_control_type), POINTER :: dft_control
1446 TYPE(mp_para_env_type), POINTER :: para_env
1447
1448 CALL timeset(routinen, handle)
1449
1450 CALL get_qs_env(qs_env=qs_env, &
1451 dft_control=dft_control, &
1452 para_env=para_env)
1453 nspins = dft_control%nspins
1454
1455 CALL dbcsr_create(matrix_q, template=matrix_p(1)%matrix, &
1456 matrix_type=dbcsr_type_no_symmetry)
1457 CALL dbcsr_create(matrix_tmp, template=matrix_p(1)%matrix, &
1458 matrix_type=dbcsr_type_no_symmetry)
1459
1460 ! Q = (1 - P)
1461 CALL dbcsr_copy(matrix_q, matrix_p(1)%matrix)
1462 CALL dbcsr_scale(matrix_q, -1.0_dp)
1463 CALL dbcsr_add_on_diag(matrix_q, 1.0_dp)
1464 CALL dbcsr_finalize(matrix_q)
1465
1466 ! Proj(M) = P*M*Q + Q*M*P
1467 ! with P = D = CC^T
1468 ! and Q = (1 - P)
1469 DO ispin = 1, nspins
1470
1471 ! tmp1 = P*M
1472 CALL dbcsr_multiply("N", "N", 1.0_dp, matrix_p(ispin)%matrix, matrix_io(ispin)%matrix, &
1473 0.0_dp, matrix_tmp, filter_eps=eps_filter)
1474 ! m_io = P*M*Q
1475 CALL dbcsr_multiply("N", "N", 1.0_dp, matrix_tmp, matrix_q, &
1476 0.0_dp, matrix_io(ispin)%matrix, filter_eps=eps_filter)
1477
1478 ! tmp = (P^T*M^T*Q^T)^T = -(P*M*Q)^T
1479 CALL dbcsr_transposed(matrix_tmp, matrix_io(ispin)%matrix)
1480 CALL dbcsr_add(matrix_io(ispin)%matrix, matrix_tmp, 1.0_dp, -1.0_dp)
1481
1482 END DO
1483
1484 CALL dbcsr_release(matrix_tmp)
1485 CALL dbcsr_release(matrix_q)
1486
1487 CALL timestop(handle)
1488
1489 END SUBROUTINE projector
1490
1491! **************************************************************************************************
1492!> \brief performs a tranformation of a matrix back to/into orthonormal basis
1493!> in case of P a scaling of 0.5 has to be applied for closed shell case
1494!> \param matrix matrix to be transformed
1495!> \param matrix_trafo transformation matrix
1496!> \param eps_filter filtering threshold for sparse matrices
1497!> \par History
1498!> 2012.05 created [Florian Schiffmann]
1499!> \author Florian Schiffmann
1500!>
1501! **************************************************************************************************
1502
1503 SUBROUTINE transform_m_orth(matrix, matrix_trafo, eps_filter)
1504 TYPE(dbcsr_type) :: matrix, matrix_trafo
1505 REAL(kind=dp) :: eps_filter
1506
1507 CHARACTER(LEN=*), PARAMETER :: routinen = 'transform_m_orth'
1508
1509 INTEGER :: handle
1510 TYPE(dbcsr_type) :: matrix_tmp, matrix_work
1511
1512 CALL timeset(routinen, handle)
1513
1514 CALL dbcsr_create(matrix_work, template=matrix, matrix_type=dbcsr_type_no_symmetry)
1515 CALL dbcsr_create(matrix_tmp, template=matrix, matrix_type=dbcsr_type_no_symmetry)
1516
1517 CALL dbcsr_multiply("N", "N", 1.0_dp, matrix, matrix_trafo, &
1518 0.0_dp, matrix_work, filter_eps=eps_filter)
1519 CALL dbcsr_multiply("N", "N", 1.0_dp, matrix_trafo, matrix_work, &
1520 0.0_dp, matrix_tmp, filter_eps=eps_filter)
1521 ! symmetrize results (this is again needed to make sure everything is stable)
1522 CALL dbcsr_transposed(matrix_work, matrix_tmp)
1523 CALL dbcsr_add(matrix_tmp, matrix_work, 0.5_dp, 0.5_dp)
1524 CALL dbcsr_copy(matrix, matrix_tmp)
1525
1526 ! Avoid the buildup of noisy blocks
1527 CALL dbcsr_filter(matrix, eps_filter)
1528
1529 CALL dbcsr_release(matrix_tmp)
1530 CALL dbcsr_release(matrix_work)
1531 CALL timestop(handle)
1532
1533 END SUBROUTINE transform_m_orth
1534
1535END MODULE ec_orth_solver
Types and set/get functions for auxiliary density matrix methods.
Definition admm_types.F:15
subroutine, public get_admm_env(admm_env, mo_derivs_aux_fit, mos_aux_fit, sab_aux_fit, sab_aux_fit_asymm, sab_aux_fit_vs_orb, matrix_s_aux_fit, matrix_s_aux_fit_kp, matrix_s_aux_fit_vs_orb, matrix_s_aux_fit_vs_orb_kp, task_list_aux_fit, matrix_ks_aux_fit, matrix_ks_aux_fit_kp, matrix_ks_aux_fit_im, matrix_ks_aux_fit_dft, matrix_ks_aux_fit_hfx, matrix_ks_aux_fit_dft_kp, matrix_ks_aux_fit_hfx_kp, rho_aux_fit, rho_aux_fit_buffer, admm_dm)
Get routine for the ADMM env.
Definition admm_types.F:593
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_transposed(transposed, normal, shallow_data_copy, transpose_distribution, use_distribution)
...
subroutine, public dbcsr_scale(matrix, alpha_scalar)
...
subroutine, public dbcsr_desymmetrize(matrix_a, matrix_b)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_multiply(transa, transb, alpha, matrix_a, matrix_b, beta, matrix_c, first_row, last_row, first_column, last_column, first_k, last_k, retain_sparsity, filter_eps, flop)
...
subroutine, public dbcsr_get_info(matrix, nblkrows_total, nblkcols_total, nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, nfullrows_local, nfullcols_local, my_prow, my_pcol, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, matrix_type, group)
...
subroutine, public dbcsr_filter(matrix, eps)
...
subroutine, public dbcsr_finalize(matrix)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_release(matrix)
...
subroutine, public dbcsr_add(matrix_a, matrix_b, alpha_scalar, beta_scalar)
...
real(kind=dp) function, public dbcsr_checksum(matrix, pos)
Calculates the checksum of a DBCSR matrix.
subroutine, public dbcsr_add_on_diag(matrix, alpha)
Adds the given scalar to the diagonal of the matrix. Reserves any missing diagonal blocks.
subroutine, public dbcsr_dot(matrix_a, matrix_b, trace)
Computes the dot product of two matrices, also known as the trace of their matrix product.
DBCSR operations in CP2K.
Routines to handle the external control of CP2K.
subroutine, public external_control(should_stop, flag, globenv, target_time, start_time, force_check)
External manipulations during a run : when the <PROJECT_NAME>.EXIT_$runtype command is sent the progr...
Routines used for Harris functional Kohn-Sham calculation.
Definition ec_methods.F:15
subroutine, public create_kernel(qs_env, vxc, vxc_tau, rho, rho1_r, rho1_g, tau1_r, xc_section, compute_virial, virial_xc)
Creation of second derivative xc-potential.
Definition ec_methods.F:88
AO-based conjugate-gradient response solver routines.
subroutine, public ec_response_ao(qs_env, p_env, matrix_hz, matrix_pz, matrix_wz, iounit, should_stop, silent)
AO-based conjugate gradient linear response solver. In goes the right hand side B of the equation AZ=...
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public precond_mlp
integer, parameter, public kg_tnadd_embed_ri
integer, parameter, public kg_tnadd_embed
integer, parameter, public do_admm_aux_exch_func_none
integer, parameter, public ls_s_sqrt_proot
integer, parameter, public ls_s_sqrt_ns
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Routines useful for iterative matrix calculations.
subroutine, public matrix_sqrt_newton_schulz(matrix_sqrt, matrix_sqrt_inv, matrix, threshold, order, eps_lanczos, max_iter_lanczos, symmetrize, converged, iounit)
compute the sqrt of a matrix via the sign function and the corresponding Newton-Schulz iterations the...
subroutine, public matrix_sqrt_proot(matrix_sqrt, matrix_sqrt_inv, matrix, threshold, order, eps_lanczos, max_iter_lanczos, symmetrize, converged)
compute the sqrt of a matrix via the general algorithm for the p-th root of Richters et al....
Routines for a Kim-Gordon-like partitioning into molecular subunits.
subroutine, public kg_ekin_subset(qs_env, ks_matrix, ekin_mol, calc_force, do_kernel, pmat_ext)
Calculates the subsystem Hohenberg-Kohn kinetic energy and the forces.
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:136
real(kind=dp) function, public m_walltime()
returns time from a real-time clock, protected against rolling early/easily
Definition machine.F:153
Collection of simple mathematical functions and subroutines.
Definition mathlib.F:15
logical function, public abnormal_value(a)
determines if a value is not normal (e.g. for Inf and Nan) based on IO to work also under optimizatio...
Definition mathlib.F:159
Interface to the message passing library MPI.
container for various plainwaves related things
subroutine, public pw_env_get(pw_env, pw_pools, cube_info, gridlevel_info, auxbas_pw_pool, auxbas_grid, auxbas_rs_desc, auxbas_rs_grid, rs_descs, rs_grids, xc_pw_pool, vdw_pw_pool, poisson_env, interp_section)
returns the various attributes of the pw env
functions related to the poisson solver on regular grids
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, mimic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
Integrate single or product functions over a potential on a RS grid.
basis types for the calculation of the perturbation of density theory.
linres kernel functions
subroutine, public apply_xc_admm(qs_env, p_env)
...
subroutine, public apply_hfx(qs_env, p_env)
Update action of TDDFPT operator on trial vectors by adding exact-exchange term.
Type definitiona for linear response calculations.
Utility functions for the perturbation calculations.
subroutine, public p_env_finish_kpp1(qs_env, p_env)
...
subroutine, public p_env_update_rho(p_env, qs_env)
...
subroutine, public p_env_check_i_alloc(p_env, qs_env)
checks that the intenal storage is allocated, and allocs it if needed
basis types for the calculation of the perturbation of density theory.
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public qs_rho_get(rho_struct, rho_ao, rho_ao_im, rho_ao_kp, rho_ao_im_kp, rho_r, drho_r, rho_g, drho_g, tau_r, tau_g, rho_r_valid, drho_r_valid, rho_g_valid, drho_g_valid, tau_r_valid, tau_g_valid, tot_rho_r, tot_rho_g, rho_r_sccs, soft_valid, complex_rho_ao)
returns info about the density described by this object. If some representation is not available an e...
Exchange and Correlation functional calculations.
Definition xc.F:17
subroutine, public xc_prep_2nd_deriv(deriv_set, rho_set, rho_r, pw_pool, weights, xc_section, tau_r)
Prepare objects for the calculation of the 2nd derivatives of the density functional....
Definition xc.F:5533
stores some data used in wavefunction fitting
Definition admm_types.F:120
stores all the informations relevant to an mpi environment
contained for different pw related things
environment for the poisson solver
to create arrays of pools
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
environment that keeps the informations and temporary val to build the kpp1 kernel matrix
General settings for linear response calculations.
Represent a qs system that is perturbed. Can calculate the linear operator and the rhs of the system ...
keeps the density in various representations, keeping track of which ones are valid.