(git:c2658d3)
Loading...
Searching...
No Matches
qs_mo_methods.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 collects routines that perform operations directly related to MOs
10!> \note
11!> first version : most routines imported
12!> \author Joost VandeVondele (2003-08)
13! **************************************************************************************************
15 USE admm_types, ONLY: admm_type
20 USE cp_dbcsr_api, ONLY: dbcsr_copy,&
27 dbcsr_type_no_symmetry
28 USE cp_dbcsr_diag, ONLY: cp_dbcsr_syevd,&
43 USE cp_fm_types, ONLY: cp_fm_create,&
50 USE kinds, ONLY: dp
53 USE physcon, ONLY: evolt
55 USE qs_mo_types, ONLY: get_mo_set,&
58#include "./base/base_uses.f90"
59
60 IMPLICIT NONE
61 PRIVATE
62 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_mo_methods'
63
67
69 MODULE PROCEDURE subspace_eigenvalues_ks_fm
70 MODULE PROCEDURE subspace_eigenvalues_ks_dbcsr
71 END INTERFACE
72
73 INTERFACE make_basis_sv
74 MODULE PROCEDURE make_basis_sv_fm
75 MODULE PROCEDURE make_basis_sv_dbcsr
76 END INTERFACE
77
78CONTAINS
79
80! **************************************************************************************************
81!> \brief returns an S-orthonormal basis v (v^T S v ==1)
82!> \param vmatrix ...
83!> \param ncol ...
84!> \param matrix_s ...
85!> \par History
86!> 03.2006 created [Joost VandeVondele]
87! **************************************************************************************************
88 SUBROUTINE make_basis_sm(vmatrix, ncol, matrix_s)
89 TYPE(cp_fm_type), INTENT(IN) :: vmatrix
90 INTEGER, INTENT(IN) :: ncol
91 TYPE(dbcsr_type), POINTER :: matrix_s
92
93 CHARACTER(LEN=*), PARAMETER :: routinen = 'make_basis_sm'
94 REAL(kind=dp), PARAMETER :: rone = 1.0_dp, rzero = 0.0_dp
95
96 INTEGER :: handle, n, ncol_global
97 TYPE(cp_fm_struct_type), POINTER :: fm_struct_tmp
98 TYPE(cp_fm_type) :: overlap_vv, svmatrix
99
100 IF (ncol == 0) RETURN
101
102 CALL timeset(routinen, handle)
103
104 CALL cp_fm_get_info(matrix=vmatrix, nrow_global=n, ncol_global=ncol_global)
105 IF (ncol > ncol_global) cpabort("Wrong ncol value")
106
107 CALL cp_fm_create(svmatrix, vmatrix%matrix_struct, "SV", set_zero=.true.)
108 CALL cp_dbcsr_sm_fm_multiply(matrix_s, vmatrix, svmatrix, ncol)
109
110 NULLIFY (fm_struct_tmp)
111 CALL cp_fm_struct_create(fm_struct_tmp, nrow_global=ncol, ncol_global=ncol, &
112 para_env=vmatrix%matrix_struct%para_env, &
113 context=vmatrix%matrix_struct%context)
114 CALL cp_fm_create(overlap_vv, fm_struct_tmp, "overlap_vv")
115 CALL cp_fm_struct_release(fm_struct_tmp)
116
117 CALL parallel_gemm('T', 'N', ncol, ncol, n, rone, vmatrix, svmatrix, rzero, overlap_vv)
118 CALL cp_fm_cholesky_decompose(overlap_vv)
119 CALL cp_fm_triangular_multiply(overlap_vv, vmatrix, n_cols=ncol, side='R', invert_tr=.true.)
120
121 CALL cp_fm_release(overlap_vv)
122 CALL cp_fm_release(svmatrix)
123
124 CALL timestop(handle)
125
126 END SUBROUTINE make_basis_sm
127
128! **************************************************************************************************
129!> \brief returns an S-orthonormal basis v and the corresponding matrix S*v as well
130!> \param vmatrix ...
131!> \param ncol ...
132!> \param svmatrix ...
133!> \par History
134!> 03.2006 created [Joost VandeVondele]
135! **************************************************************************************************
136 SUBROUTINE make_basis_sv_fm(vmatrix, ncol, svmatrix)
137
138 TYPE(cp_fm_type), INTENT(IN) :: vmatrix
139 INTEGER, INTENT(IN) :: ncol
140 TYPE(cp_fm_type), INTENT(IN) :: svmatrix
141
142 CHARACTER(LEN=*), PARAMETER :: routineN = 'make_basis_sv_fm'
143 REAL(KIND=dp), PARAMETER :: rone = 1.0_dp, rzero = 0.0_dp
144
145 INTEGER :: handle, n, ncol_global
146 TYPE(cp_fm_struct_type), POINTER :: fm_struct_tmp
147 TYPE(cp_fm_type) :: overlap_vv
148
149 IF (ncol == 0) RETURN
150
151 CALL timeset(routinen, handle)
152 NULLIFY (fm_struct_tmp)
153
154 CALL cp_fm_get_info(matrix=vmatrix, nrow_global=n, ncol_global=ncol_global)
155 IF (ncol > ncol_global) cpabort("Wrong ncol value")
156
157 CALL cp_fm_struct_create(fm_struct_tmp, nrow_global=ncol, ncol_global=ncol, &
158 para_env=vmatrix%matrix_struct%para_env, &
159 context=vmatrix%matrix_struct%context)
160 CALL cp_fm_create(overlap_vv, fm_struct_tmp, "overlap_vv")
161 CALL cp_fm_struct_release(fm_struct_tmp)
162
163 CALL parallel_gemm('T', 'N', ncol, ncol, n, rone, vmatrix, svmatrix, rzero, overlap_vv)
164 CALL cp_fm_cholesky_decompose(overlap_vv)
165 CALL cp_fm_triangular_multiply(overlap_vv, vmatrix, n_cols=ncol, side='R', invert_tr=.true.)
166 CALL cp_fm_triangular_multiply(overlap_vv, svmatrix, n_cols=ncol, side='R', invert_tr=.true.)
167
168 CALL cp_fm_release(overlap_vv)
169
170 CALL timestop(handle)
171
172 END SUBROUTINE make_basis_sv_fm
173
174! **************************************************************************************************
175!> \brief ...
176!> \param vmatrix ...
177!> \param ncol ...
178!> \param svmatrix ...
179!> \param para_env ...
180!> \param blacs_env ...
181! **************************************************************************************************
182 SUBROUTINE make_basis_sv_dbcsr(vmatrix, ncol, svmatrix, para_env, blacs_env)
183
184 TYPE(dbcsr_type) :: vmatrix
185 INTEGER, INTENT(IN) :: ncol
186 TYPE(dbcsr_type) :: svmatrix
187 TYPE(mp_para_env_type), POINTER :: para_env
188 TYPE(cp_blacs_env_type), POINTER :: blacs_env
189
190 CHARACTER(LEN=*), PARAMETER :: routineN = 'make_basis_sv_dbcsr'
191 REAL(KIND=dp), PARAMETER :: rone = 1.0_dp, rzero = 0.0_dp
192
193 INTEGER :: handle, n, ncol_global
194 TYPE(cp_fm_struct_type), POINTER :: fm_struct_tmp
195 TYPE(cp_fm_type) :: fm_svmatrix, fm_vmatrix, overlap_vv
196
197 IF (ncol == 0) RETURN
198
199 CALL timeset(routinen, handle)
200
201 !CALL cp_fm_get_info(matrix=vmatrix,nrow_global=n,ncol_global=ncol_global)
202 CALL dbcsr_get_info(vmatrix, nfullrows_total=n, nfullcols_total=ncol_global)
203 IF (ncol > ncol_global) cpabort("Wrong ncol value")
204
205 CALL cp_fm_struct_create(fm_struct_tmp, context=blacs_env, nrow_global=ncol, &
206 ncol_global=ncol, para_env=para_env)
207 CALL cp_fm_create(overlap_vv, fm_struct_tmp, name="fm_overlap_vv")
208 CALL cp_fm_struct_release(fm_struct_tmp)
209
210 CALL cp_fm_struct_create(fm_struct_tmp, context=blacs_env, nrow_global=n, &
211 ncol_global=ncol_global, para_env=para_env)
212 CALL cp_fm_create(fm_vmatrix, fm_struct_tmp, name="fm_vmatrix")
213 CALL cp_fm_create(fm_svmatrix, fm_struct_tmp, name="fm_svmatrix")
214 CALL cp_fm_struct_release(fm_struct_tmp)
215
216 CALL copy_dbcsr_to_fm(vmatrix, fm_vmatrix)
217 CALL copy_dbcsr_to_fm(svmatrix, fm_svmatrix)
218
219 CALL parallel_gemm('T', 'N', ncol, ncol, n, rone, fm_vmatrix, fm_svmatrix, rzero, overlap_vv)
220 CALL cp_fm_cholesky_decompose(overlap_vv)
221 CALL cp_fm_triangular_multiply(overlap_vv, fm_vmatrix, n_cols=ncol, side='R', invert_tr=.true.)
222 CALL cp_fm_triangular_multiply(overlap_vv, fm_svmatrix, n_cols=ncol, side='R', invert_tr=.true.)
223
224 CALL copy_fm_to_dbcsr(fm_vmatrix, vmatrix)
225 CALL copy_fm_to_dbcsr(fm_svmatrix, svmatrix)
226
227 CALL cp_fm_release(overlap_vv)
228 CALL cp_fm_release(fm_vmatrix)
229 CALL cp_fm_release(fm_svmatrix)
230
231 CALL timestop(handle)
232
233 END SUBROUTINE make_basis_sv_dbcsr
234
235! **************************************************************************************************
236!> \brief return a set of S orthonormal vectors (C^T S C == 1) where
237!> the cholesky decomposed form of S is passed as an argument
238!> \param vmatrix ...
239!> \param ncol ...
240!> \param ortho cholesky decomposed S matrix
241!> \par History
242!> 03.2006 created [Joost VandeVondele]
243!> \note
244!> if the cholesky decomposed S matrix is not available
245!> use make_basis_sm since this is much faster than computing the
246!> cholesky decomposition of S
247! **************************************************************************************************
248 SUBROUTINE make_basis_cholesky(vmatrix, ncol, ortho)
249
250 TYPE(cp_fm_type), INTENT(IN) :: vmatrix
251 INTEGER, INTENT(IN) :: ncol
252 TYPE(cp_fm_type), INTENT(IN) :: ortho
253
254 CHARACTER(LEN=*), PARAMETER :: routinen = 'make_basis_cholesky'
255 REAL(kind=dp), PARAMETER :: rone = 1.0_dp, rzero = 0.0_dp
256
257 INTEGER :: handle, n, ncol_global
258 TYPE(cp_fm_struct_type), POINTER :: fm_struct_tmp
259 TYPE(cp_fm_type) :: overlap_vv
260
261 IF (ncol == 0) RETURN
262
263 CALL timeset(routinen, handle)
264 NULLIFY (fm_struct_tmp)
265
266 CALL cp_fm_get_info(matrix=vmatrix, nrow_global=n, ncol_global=ncol_global)
267 IF (ncol > ncol_global) cpabort("Wrong ncol value")
268
269 CALL cp_fm_struct_create(fm_struct_tmp, nrow_global=ncol, ncol_global=ncol, &
270 para_env=vmatrix%matrix_struct%para_env, &
271 context=vmatrix%matrix_struct%context)
272 CALL cp_fm_create(overlap_vv, fm_struct_tmp, "overlap_vv")
273 CALL cp_fm_struct_release(fm_struct_tmp)
274
275 CALL cp_fm_triangular_multiply(ortho, vmatrix, n_cols=ncol)
276 CALL cp_fm_syrk('U', 'T', n, rone, vmatrix, 1, 1, rzero, overlap_vv)
277 CALL cp_fm_cholesky_decompose(overlap_vv)
278 CALL cp_fm_triangular_multiply(overlap_vv, vmatrix, n_cols=ncol, side='R', invert_tr=.true.)
279 CALL cp_fm_triangular_multiply(ortho, vmatrix, n_cols=ncol, invert_tr=.true.)
280
281 CALL cp_fm_release(overlap_vv)
282
283 CALL timestop(handle)
284
285 END SUBROUTINE make_basis_cholesky
286
287! **************************************************************************************************
288!> \brief return a set of S orthonormal vectors (C^T S C == 1) where
289!> a Loedwin transformation is applied to keep the rotated vectors as close
290!> as possible to the original ones
291!> \param vmatrix ...
292!> \param ncol ...
293!> \param matrix_s ...
294!> \param
295!> \par History
296!> 05.2009 created [MI]
297!> \note
298! **************************************************************************************************
299 SUBROUTINE make_basis_lowdin(vmatrix, ncol, matrix_s)
300
301 TYPE(cp_fm_type), INTENT(IN) :: vmatrix
302 INTEGER, INTENT(IN) :: ncol
303 TYPE(dbcsr_type) :: matrix_s
304
305 CHARACTER(LEN=*), PARAMETER :: routinen = 'make_basis_lowdin'
306 REAL(kind=dp), PARAMETER :: rone = 1.0_dp, rzero = 0.0_dp
307
308 INTEGER :: handle, n, ncol_global, ndep
309 REAL(dp) :: threshold
310 TYPE(cp_fm_struct_type), POINTER :: fm_struct_tmp
311 TYPE(cp_fm_type) :: csc, sc, work
312
313 IF (ncol == 0) RETURN
314
315 CALL timeset(routinen, handle)
316 NULLIFY (fm_struct_tmp)
317 threshold = 1.0e-7_dp
318 CALL cp_fm_get_info(matrix=vmatrix, nrow_global=n, ncol_global=ncol_global)
319 IF (ncol > ncol_global) cpabort("Wrong ncol value")
320
321 CALL cp_fm_create(sc, vmatrix%matrix_struct, "SC")
322 CALL cp_dbcsr_sm_fm_multiply(matrix_s, vmatrix, sc, ncol)
323
324 NULLIFY (fm_struct_tmp)
325 CALL cp_fm_struct_create(fm_struct_tmp, nrow_global=ncol, ncol_global=ncol, &
326 para_env=vmatrix%matrix_struct%para_env, &
327 context=vmatrix%matrix_struct%context)
328 CALL cp_fm_create(csc, fm_struct_tmp, "csc")
329 CALL cp_fm_set_all(csc, 0.0_dp)
330 CALL cp_fm_create(work, fm_struct_tmp, "work")
331 CALL cp_fm_set_all(work, 0.0_dp)
332 CALL cp_fm_struct_release(fm_struct_tmp)
333
334 CALL parallel_gemm('T', 'N', ncol, ncol, n, rone, vmatrix, sc, rzero, csc)
335 CALL cp_fm_power(csc, work, -0.5_dp, threshold, ndep)
336 CALL parallel_gemm('N', 'N', n, ncol, ncol, rone, vmatrix, csc, rzero, sc)
337 CALL cp_fm_to_fm(sc, vmatrix, ncol, 1, 1)
338
339 CALL cp_fm_release(csc)
340 CALL cp_fm_release(sc)
341 CALL cp_fm_release(work)
342
343 CALL timestop(handle)
344
345 END SUBROUTINE make_basis_lowdin
346
347! **************************************************************************************************
348!> \brief given a set of vectors, return an orthogonal (C^T C == 1) set
349!> spanning the same space (notice, only for cases where S==1)
350!> \param vmatrix ...
351!> \param ncol ...
352!> \par History
353!> 03.2006 created [Joost VandeVondele]
354! **************************************************************************************************
355 SUBROUTINE make_basis_simple(vmatrix, ncol)
356
357 TYPE(cp_fm_type), INTENT(IN) :: vmatrix
358 INTEGER, INTENT(IN) :: ncol
359
360 CHARACTER(LEN=*), PARAMETER :: routinen = 'make_basis_simple'
361 REAL(kind=dp), PARAMETER :: rone = 1.0_dp, rzero = 0.0_dp
362
363 INTEGER :: handle, n, ncol_global
364 TYPE(cp_fm_struct_type), POINTER :: fm_struct_tmp
365 TYPE(cp_fm_type) :: overlap_vv
366
367 IF (ncol == 0) RETURN
368
369 CALL timeset(routinen, handle)
370
371 NULLIFY (fm_struct_tmp)
372
373 CALL cp_fm_get_info(matrix=vmatrix, nrow_global=n, ncol_global=ncol_global)
374 IF (ncol > ncol_global) cpabort("Wrong ncol value")
375
376 CALL cp_fm_struct_create(fm_struct_tmp, nrow_global=ncol, ncol_global=ncol, &
377 para_env=vmatrix%matrix_struct%para_env, &
378 context=vmatrix%matrix_struct%context)
379 CALL cp_fm_create(overlap_vv, fm_struct_tmp, "overlap_vv")
380 CALL cp_fm_struct_release(fm_struct_tmp)
381
382 CALL parallel_gemm('T', 'N', ncol, ncol, n, rone, vmatrix, vmatrix, rzero, overlap_vv)
383 CALL cp_fm_cholesky_decompose(overlap_vv)
384 CALL cp_fm_triangular_multiply(overlap_vv, vmatrix, n_cols=ncol, side='R', invert_tr=.true.)
385
386 CALL cp_fm_release(overlap_vv)
387
388 CALL timestop(handle)
389
390 END SUBROUTINE make_basis_simple
391
392! **************************************************************************************************
393!> \brief computes ritz values of a set of orbitals given a ks_matrix
394!> rotates the orbitals into eigenstates depending on do_rotation
395!> writes the evals to the screen depending on ionode/scr
396!> \param orbitals S-orthonormal orbitals
397!> \param ks_matrix Kohn-Sham matrix
398!> \param evals_arg optional, filled with the evals
399!> \param ionode , scr if present write to unit scr where ionode
400!> \param scr ...
401!> \param do_rotation optional rotate orbitals (default=.TRUE.)
402!> note that rotating the orbitals is slower
403!> \param co_rotate an optional set of orbitals rotated by the same rotation matrix
404!> \param co_rotate_dbcsr ...
405!> \par History
406!> 08.2004 documented and added do_rotation [Joost VandeVondele]
407!> 09.2008 only compute eigenvalues if rotation is not needed
408! **************************************************************************************************
409 SUBROUTINE subspace_eigenvalues_ks_fm(orbitals, ks_matrix, evals_arg, ionode, scr, &
410 do_rotation, co_rotate, co_rotate_dbcsr)
411
412 TYPE(cp_fm_type), INTENT(IN) :: orbitals
413 TYPE(dbcsr_type), POINTER :: ks_matrix
414 REAL(KIND=dp), DIMENSION(:), OPTIONAL :: evals_arg
415 LOGICAL, INTENT(IN), OPTIONAL :: ionode
416 INTEGER, INTENT(IN), OPTIONAL :: scr
417 LOGICAL, INTENT(IN), OPTIONAL :: do_rotation
418 TYPE(cp_fm_type), INTENT(IN), OPTIONAL :: co_rotate
419 TYPE(dbcsr_type), OPTIONAL, POINTER :: co_rotate_dbcsr
420
421 CHARACTER(len=*), PARAMETER :: routineN = 'subspace_eigenvalues_ks_fm'
422
423 INTEGER :: handle, i, j, n, ncol_global, nrow_global
424 LOGICAL :: compute_evecs, do_rotation_local
425 REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: evals
426 TYPE(cp_fm_struct_type), POINTER :: fm_struct_tmp
427 TYPE(cp_fm_type) :: e_vectors, h_block, weighted_vectors, &
428 weighted_vectors2
429
430 CALL timeset(routinen, handle)
431
432 do_rotation_local = .true.
433 IF (PRESENT(do_rotation)) do_rotation_local = do_rotation
434
435 NULLIFY (fm_struct_tmp)
436 CALL cp_fm_get_info(matrix=orbitals, &
437 ncol_global=ncol_global, &
438 nrow_global=nrow_global)
439
440 IF (do_rotation_local) THEN
441 compute_evecs = .true.
442 ELSE
443 ! this would be the logical choice if syevx computing only evals were faster than syevd computing evecs and evals.
444 compute_evecs = .false.
445 ! this is not the case, so lets compute evecs always
446 compute_evecs = .true.
447 END IF
448
449 IF (ncol_global > 0) THEN
450
451 ALLOCATE (evals(ncol_global))
452
453 CALL cp_fm_create(weighted_vectors, orbitals%matrix_struct, "weighted_vectors")
454 CALL cp_fm_struct_create(fm_struct_tmp, nrow_global=ncol_global, ncol_global=ncol_global, &
455 para_env=orbitals%matrix_struct%para_env, &
456 context=orbitals%matrix_struct%context)
457 CALL cp_fm_create(h_block, fm_struct_tmp, name="h block")
458 IF (compute_evecs) THEN
459 CALL cp_fm_create(e_vectors, fm_struct_tmp, name="e vectors")
460 END IF
461 CALL cp_fm_struct_release(fm_struct_tmp)
462
463 ! h subblock and diag
464 CALL cp_dbcsr_sm_fm_multiply(ks_matrix, orbitals, weighted_vectors, ncol_global)
465
466 CALL parallel_gemm('T', 'N', ncol_global, ncol_global, nrow_global, 1.0_dp, &
467 orbitals, weighted_vectors, 0.0_dp, h_block)
468
469 ! if eigenvectors are required, go for syevd, otherwise only compute eigenvalues
470 IF (compute_evecs) THEN
471 CALL choose_eigv_solver(h_block, e_vectors, evals)
472 ELSE
473 CALL cp_fm_syevx(h_block, eigenvalues=evals)
474 END IF
475
476 ! rotate the orbitals
477 IF (do_rotation_local) THEN
478 CALL parallel_gemm('N', 'N', nrow_global, ncol_global, ncol_global, 1.0_dp, &
479 orbitals, e_vectors, 0.0_dp, weighted_vectors)
480 CALL cp_fm_to_fm(weighted_vectors, orbitals)
481 IF (PRESENT(co_rotate)) THEN
482 CALL parallel_gemm('N', 'N', nrow_global, ncol_global, ncol_global, 1.0_dp, &
483 co_rotate, e_vectors, 0.0_dp, weighted_vectors)
484 CALL cp_fm_to_fm(weighted_vectors, co_rotate)
485 END IF
486 IF (PRESENT(co_rotate_dbcsr)) THEN
487 IF (ASSOCIATED(co_rotate_dbcsr)) THEN
488 CALL cp_fm_create(weighted_vectors2, orbitals%matrix_struct, "weighted_vectors")
489 CALL copy_dbcsr_to_fm(co_rotate_dbcsr, weighted_vectors2)
490 CALL parallel_gemm('N', 'N', nrow_global, ncol_global, ncol_global, 1.0_dp, &
491 weighted_vectors2, e_vectors, 0.0_dp, weighted_vectors)
492 CALL copy_fm_to_dbcsr(weighted_vectors, co_rotate_dbcsr)
493 CALL cp_fm_release(weighted_vectors2)
494 END IF
495 END IF
496 END IF
497
498 ! give output
499 IF (PRESENT(evals_arg)) THEN
500 n = min(SIZE(evals_arg), SIZE(evals))
501 evals_arg(1:n) = evals(1:n)
502 END IF
503
504 IF (PRESENT(ionode) .OR. PRESENT(scr)) THEN
505 IF (.NOT. PRESENT(ionode)) cpabort("IONODE?")
506 IF (.NOT. PRESENT(scr)) cpabort("SCR?")
507 IF (ionode) THEN
508 DO i = 1, ncol_global, 4
509 j = min(3, ncol_global - i)
510 SELECT CASE (j)
511 CASE (3)
512 WRITE (scr, '(1X,4F16.8)') evals(i:i + j)
513 CASE (2)
514 WRITE (scr, '(1X,3F16.8)') evals(i:i + j)
515 CASE (1)
516 WRITE (scr, '(1X,2F16.8)') evals(i:i + j)
517 CASE (0)
518 WRITE (scr, '(1X,1F16.8)') evals(i:i + j)
519 END SELECT
520 END DO
521 END IF
522 END IF
523
524 CALL cp_fm_release(weighted_vectors)
525 CALL cp_fm_release(h_block)
526 IF (compute_evecs) THEN
527 CALL cp_fm_release(e_vectors)
528 END IF
529
530 DEALLOCATE (evals)
531
532 END IF
533
534 CALL timestop(handle)
535
536 END SUBROUTINE subspace_eigenvalues_ks_fm
537
538! **************************************************************************************************
539!> \brief ...
540!> \param orbitals ...
541!> \param ks_matrix ...
542!> \param evals_arg ...
543!> \param ionode ...
544!> \param scr ...
545!> \param do_rotation ...
546!> \param co_rotate ...
547!> \param para_env ...
548!> \param blacs_env ...
549! **************************************************************************************************
550 SUBROUTINE subspace_eigenvalues_ks_dbcsr(orbitals, ks_matrix, evals_arg, ionode, scr, &
551 do_rotation, co_rotate, para_env, blacs_env)
552
553 TYPE(dbcsr_type), POINTER :: orbitals, ks_matrix
554 REAL(KIND=dp), DIMENSION(:), OPTIONAL :: evals_arg
555 LOGICAL, INTENT(IN), OPTIONAL :: ionode
556 INTEGER, INTENT(IN), OPTIONAL :: scr
557 LOGICAL, INTENT(IN), OPTIONAL :: do_rotation
558 TYPE(dbcsr_type), OPTIONAL, POINTER :: co_rotate
559 TYPE(mp_para_env_type), POINTER :: para_env
560 TYPE(cp_blacs_env_type), POINTER :: blacs_env
561
562 CHARACTER(len=*), PARAMETER :: routineN = 'subspace_eigenvalues_ks_dbcsr'
563
564 INTEGER :: handle, i, j, ncol_global, nrow_global
565 LOGICAL :: compute_evecs, do_rotation_local
566 REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: evals
567 TYPE(dbcsr_type), POINTER :: e_vectors, h_block, weighted_vectors
568
569 CALL timeset(routinen, handle)
570
571 do_rotation_local = .true.
572 IF (PRESENT(do_rotation)) do_rotation_local = do_rotation
573
574 NULLIFY (e_vectors, h_block, weighted_vectors)
575
576 CALL dbcsr_get_info(matrix=orbitals, &
577 nfullcols_total=ncol_global, &
578 nfullrows_total=nrow_global)
579
580 IF (do_rotation_local) THEN
581 compute_evecs = .true.
582 ELSE
583 ! this would be the logical choice if syevx computing only evals were faster than syevd computing evecs and evals.
584 compute_evecs = .false.
585 ! this is not the case, so lets compute evecs always
586 compute_evecs = .true.
587 END IF
588
589 IF (ncol_global > 0) THEN
590
591 ALLOCATE (evals(ncol_global))
592
593 CALL dbcsr_init_p(weighted_vectors)
594 CALL dbcsr_copy(weighted_vectors, orbitals, name="weighted_vectors")
595
596 CALL dbcsr_init_p(h_block)
597 CALL cp_dbcsr_m_by_n_from_template(h_block, template=orbitals, m=ncol_global, n=ncol_global, &
598 sym=dbcsr_type_no_symmetry)
599
600!!!!!!!!!!!!!!XXXXXXXXXXXXXXXXXXX!!!!!!!!!!!!!
601 IF (compute_evecs) THEN
602 CALL dbcsr_init_p(e_vectors)
603 CALL cp_dbcsr_m_by_n_from_template(e_vectors, template=orbitals, m=ncol_global, n=ncol_global, &
604 sym=dbcsr_type_no_symmetry)
605 END IF
606
607 ! h subblock and diag
608 CALL dbcsr_multiply('N', 'N', 1.0_dp, ks_matrix, orbitals, &
609 0.0_dp, weighted_vectors)
610 !CALL cp_dbcsr_sm_fm_multiply(ks_matrix,orbitals,weighted_vectors, ncol_global)
611
612 CALL dbcsr_multiply('T', 'N', 1.0_dp, orbitals, weighted_vectors, 0.0_dp, h_block)
613 !CALL parallel_gemm('T','N',ncol_global,ncol_global,nrow_global,1.0_dp, &
614 ! orbitals,weighted_vectors,0.0_dp,h_block)
615
616 ! if eigenvectors are required, go for syevd, otherwise only compute eigenvalues
617 IF (compute_evecs) THEN
618 CALL cp_dbcsr_syevd(h_block, e_vectors, evals, para_env=para_env, blacs_env=blacs_env)
619 ELSE
620 CALL cp_dbcsr_syevx(h_block, eigenvalues=evals, para_env=para_env, blacs_env=blacs_env)
621 END IF
622
623 ! rotate the orbitals
624 IF (do_rotation_local) THEN
625 CALL dbcsr_multiply('N', 'N', 1.0_dp, orbitals, e_vectors, 0.0_dp, weighted_vectors)
626 !CALL parallel_gemm('N','N',nrow_global,ncol_global,ncol_global,1.0_dp, &
627 ! orbitals,e_vectors,0.0_dp,weighted_vectors)
628 CALL dbcsr_copy(orbitals, weighted_vectors)
629 !CALL cp_fm_to_fm(weighted_vectors,orbitals)
630 IF (PRESENT(co_rotate)) THEN
631 IF (ASSOCIATED(co_rotate)) THEN
632 CALL dbcsr_multiply('N', 'N', 1.0_dp, co_rotate, e_vectors, 0.0_dp, weighted_vectors)
633 !CALL parallel_gemm('N','N',nrow_global,ncol_global,ncol_global,1.0_dp, &
634 ! co_rotate,e_vectors,0.0_dp,weighted_vectors)
635 CALL dbcsr_copy(co_rotate, weighted_vectors)
636 !CALL cp_fm_to_fm(weighted_vectors,co_rotate)
637 END IF
638 END IF
639 END IF
640
641 ! give output
642 IF (PRESENT(evals_arg)) THEN
643 evals_arg(:) = evals(:)
644 END IF
645
646 IF (PRESENT(ionode) .OR. PRESENT(scr)) THEN
647 IF (.NOT. PRESENT(ionode)) cpabort("IONODE?")
648 IF (.NOT. PRESENT(scr)) cpabort("SCR?")
649 IF (ionode) THEN
650 DO i = 1, ncol_global, 4
651 j = min(3, ncol_global - i)
652 SELECT CASE (j)
653 CASE (3)
654 WRITE (scr, '(1X,4F16.8)') evals(i:i + j)
655 CASE (2)
656 WRITE (scr, '(1X,3F16.8)') evals(i:i + j)
657 CASE (1)
658 WRITE (scr, '(1X,2F16.8)') evals(i:i + j)
659 CASE (0)
660 WRITE (scr, '(1X,1F16.8)') evals(i:i + j)
661 END SELECT
662 END DO
663 END IF
664 END IF
665
666 CALL dbcsr_release_p(weighted_vectors)
667 CALL dbcsr_release_p(h_block)
668 IF (compute_evecs) THEN
669 CALL dbcsr_release_p(e_vectors)
670 END IF
671
672 DEALLOCATE (evals)
673
674 END IF
675
676 CALL timestop(handle)
677
678 END SUBROUTINE subspace_eigenvalues_ks_dbcsr
679
680! computes the effective orthonormality of a set of mos given an s-matrix
681! orthonormality is the max deviation from unity of the C^T S C
682! **************************************************************************************************
683!> \brief ...
684!> \param orthonormality ...
685!> \param mo_array ...
686!> \param matrix_s ...
687! **************************************************************************************************
688 SUBROUTINE calculate_orthonormality(orthonormality, mo_array, matrix_s)
689 REAL(kind=dp) :: orthonormality
690 TYPE(mo_set_type), DIMENSION(:), INTENT(IN) :: mo_array
691 TYPE(dbcsr_type), OPTIONAL, POINTER :: matrix_s
692
693 CHARACTER(len=*), PARAMETER :: routinen = 'calculate_orthonormality'
694
695 INTEGER :: handle, i, ispin, j, k, n, ncol_local, &
696 nrow_local, nspin
697 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
698 REAL(kind=dp) :: alpha, max_alpha
699 TYPE(cp_fm_struct_type), POINTER :: tmp_fm_struct
700 TYPE(cp_fm_type) :: overlap, svec
701
702 NULLIFY (tmp_fm_struct)
703
704 CALL timeset(routinen, handle)
705
706 nspin = SIZE(mo_array)
707 max_alpha = 0.0_dp
708
709 DO ispin = 1, nspin
710 IF (PRESENT(matrix_s)) THEN
711 ! get S*C
712 CALL cp_fm_create(svec, mo_array(ispin)%mo_coeff%matrix_struct)
713 CALL cp_fm_get_info(mo_array(ispin)%mo_coeff, &
714 nrow_global=n, ncol_global=k)
715 CALL cp_dbcsr_sm_fm_multiply(matrix_s, mo_array(ispin)%mo_coeff, &
716 svec, k)
717 ! get C^T (S*C)
718 CALL cp_fm_struct_create(tmp_fm_struct, nrow_global=k, ncol_global=k, &
719 para_env=mo_array(ispin)%mo_coeff%matrix_struct%para_env, &
720 context=mo_array(ispin)%mo_coeff%matrix_struct%context)
721 CALL cp_fm_create(overlap, tmp_fm_struct)
722 CALL cp_fm_struct_release(tmp_fm_struct)
723 CALL parallel_gemm('T', 'N', k, k, n, 1.0_dp, mo_array(ispin)%mo_coeff, &
724 svec, 0.0_dp, overlap)
725 CALL cp_fm_release(svec)
726 ELSE
727 ! orthogonal basis C^T C
728 CALL cp_fm_get_info(mo_array(ispin)%mo_coeff, &
729 nrow_global=n, ncol_global=k)
730 CALL cp_fm_struct_create(tmp_fm_struct, nrow_global=k, ncol_global=k, &
731 para_env=mo_array(ispin)%mo_coeff%matrix_struct%para_env, &
732 context=mo_array(ispin)%mo_coeff%matrix_struct%context)
733 CALL cp_fm_create(overlap, tmp_fm_struct)
734 CALL cp_fm_struct_release(tmp_fm_struct)
735 CALL parallel_gemm('T', 'N', k, k, n, 1.0_dp, mo_array(ispin)%mo_coeff, &
736 mo_array(ispin)%mo_coeff, 0.0_dp, overlap)
737 END IF
738 CALL cp_fm_get_info(overlap, nrow_local=nrow_local, ncol_local=ncol_local, &
739 row_indices=row_indices, col_indices=col_indices)
740 DO i = 1, nrow_local
741 DO j = 1, ncol_local
742 alpha = overlap%local_data(i, j)
743 IF (row_indices(i) == col_indices(j)) alpha = alpha - 1.0_dp
744 max_alpha = max(max_alpha, abs(alpha))
745 END DO
746 END DO
747 CALL cp_fm_release(overlap)
748 END DO
749 CALL mo_array(1)%mo_coeff%matrix_struct%para_env%max(max_alpha)
750 orthonormality = max_alpha
751 ! write(*,*) "max deviation from orthonormalization ",orthonormality
752
753 CALL timestop(handle)
754
755 END SUBROUTINE calculate_orthonormality
756
757! computes the minimum/maximum magnitudes of C^T C. This could be useful
758! to detect problems in the case of nearly singular overlap matrices.
759! in this case, we expect the ratio of min/max to be large
760! this routine is only similar to mo_orthonormality if S==1
761! **************************************************************************************************
762!> \brief ...
763!> \param mo_array ...
764!> \param mo_mag_min ...
765!> \param mo_mag_max ...
766! **************************************************************************************************
767 SUBROUTINE calculate_magnitude(mo_array, mo_mag_min, mo_mag_max)
768 TYPE(mo_set_type), DIMENSION(:), INTENT(IN) :: mo_array
769 REAL(kind=dp) :: mo_mag_min, mo_mag_max
770
771 CHARACTER(len=*), PARAMETER :: routinen = 'calculate_magnitude'
772
773 INTEGER :: handle, ispin, k, n, nspin
774 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: evals
775 TYPE(cp_fm_struct_type), POINTER :: tmp_fm_struct
776 TYPE(cp_fm_type) :: evecs, overlap
777
778 NULLIFY (tmp_fm_struct)
779
780 CALL timeset(routinen, handle)
781
782 nspin = SIZE(mo_array)
783 mo_mag_min = huge(0.0_dp)
784 mo_mag_max = -huge(0.0_dp)
785 DO ispin = 1, nspin
786 CALL cp_fm_get_info(mo_array(ispin)%mo_coeff, &
787 nrow_global=n, ncol_global=k)
788 ALLOCATE (evals(k))
789 CALL cp_fm_struct_create(tmp_fm_struct, nrow_global=k, ncol_global=k, &
790 para_env=mo_array(ispin)%mo_coeff%matrix_struct%para_env, &
791 context=mo_array(ispin)%mo_coeff%matrix_struct%context)
792 CALL cp_fm_create(overlap, tmp_fm_struct)
793 CALL cp_fm_create(evecs, tmp_fm_struct)
794 CALL cp_fm_struct_release(tmp_fm_struct)
795 CALL parallel_gemm('T', 'N', k, k, n, 1.0_dp, mo_array(ispin)%mo_coeff, &
796 mo_array(ispin)%mo_coeff, 0.0_dp, overlap)
797 CALL choose_eigv_solver(overlap, evecs, evals)
798 mo_mag_min = min(minval(evals), mo_mag_min)
799 mo_mag_max = max(maxval(evals), mo_mag_max)
800 CALL cp_fm_release(overlap)
801 CALL cp_fm_release(evecs)
802 DEALLOCATE (evals)
803 END DO
804 CALL timestop(handle)
805
806 END SUBROUTINE calculate_magnitude
807
808! **************************************************************************************************
809!> \brief Calculate KS eigenvalues starting from OF MOS
810!> \param mos ...
811!> \param nspins ...
812!> \param ks_rmpv ...
813!> \param scf_control ...
814!> \param mo_derivs ...
815!> \param admm_env ...
816!> \param hairy_probes ...
817!> \param probe ...
818!> \par History
819!> 02.2013 moved from qs_scf_post_gpw
820!>
821! **************************************************************************************************
822 SUBROUTINE make_mo_eig(mos, nspins, ks_rmpv, scf_control, mo_derivs, admm_env, &
823 hairy_probes, probe)
824
825 TYPE(mo_set_type), DIMENSION(:), INTENT(INOUT) :: mos
826 INTEGER, INTENT(IN) :: nspins
827 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: ks_rmpv
828 TYPE(scf_control_type), POINTER :: scf_control
829 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: mo_derivs
830 TYPE(admm_type), OPTIONAL, POINTER :: admm_env
831 LOGICAL, OPTIONAL :: hairy_probes
832 TYPE(hairy_probes_type), DIMENSION(:), OPTIONAL, &
833 POINTER :: probe
834
835 CHARACTER(len=*), PARAMETER :: routinen = 'make_mo_eig'
836
837 INTEGER :: handle, homo, ispin, nmo, output_unit
838 REAL(kind=dp), DIMENSION(:), POINTER :: mo_eigenvalues
839 TYPE(cp_fm_type), POINTER :: mo_coeff
840 TYPE(dbcsr_type), POINTER :: mo_coeff_deriv
841
842 CALL timeset(routinen, handle)
843
844 NULLIFY (mo_coeff_deriv, mo_coeff, mo_eigenvalues)
845 output_unit = cp_logger_get_default_io_unit()
846
847 DO ispin = 1, nspins
848 CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff, &
849 eigenvalues=mo_eigenvalues, homo=homo, nmo=nmo)
850 IF (output_unit > 0) WRITE (output_unit, *) " "
851 IF (output_unit > 0) WRITE (output_unit, *) " Eigenvalues of the occupied subspace spin ", ispin
852 ! IF (homo /= nmo) THEN
853 ! IF (output_unit>0) WRITE(output_unit,*)" and ",nmo-homo," added MO eigenvalues"
854 ! END IF
855 IF (output_unit > 0) WRITE (output_unit, *) "---------------------------------------------"
856
857 IF (scf_control%use_ot) THEN
858 ! Also rotate the OT derivs, since they are needed for force calculations
859 IF (ASSOCIATED(mo_derivs)) THEN
860 mo_coeff_deriv => mo_derivs(ispin)%matrix
861 ELSE
862 mo_coeff_deriv => null()
863 END IF
864
865 ! If we do ADMM, we add have to modify the Kohn-Sham matrix
866 IF (PRESENT(admm_env)) THEN
867 CALL admm_correct_for_eigenvalues(ispin, admm_env, ks_rmpv(ispin)%matrix)
868 END IF
869
870 CALL calculate_subspace_eigenvalues(mo_coeff, ks_rmpv(ispin)%matrix, mo_eigenvalues, &
871 scr=output_unit, ionode=output_unit > 0, do_rotation=.true., &
872 co_rotate_dbcsr=mo_coeff_deriv)
873
874 ! If we do ADMM, we restore the original Kohn-Sham matrix
875 IF (PRESENT(admm_env)) THEN
876 CALL admm_uncorrect_for_eigenvalues(ispin, admm_env, ks_rmpv(ispin)%matrix)
877 END IF
878 ELSE
879 IF (output_unit > 0) WRITE (output_unit, '(4(1X,1F16.8))') mo_eigenvalues(1:homo)
880 END IF
881 IF (.NOT. scf_control%diagonalization%mom) THEN
882 IF (PRESENT(hairy_probes) .EQV. .true.) THEN
883 scf_control%smear%do_smear = .false.
884 CALL set_mo_occupation(mo_set=mos(ispin), &
885 smear=scf_control%smear, &
886 probe=probe)
887 ELSE
888 CALL set_mo_occupation(mo_set=mos(ispin), smear=scf_control%smear)
889 END IF
890 END IF
891 IF (output_unit > 0) WRITE (output_unit, '(T2,A,F12.6)') &
892 "Fermi Energy [eV] :", mos(ispin)%mu*evolt
893 END DO
894
895 CALL timestop(handle)
896
897 END SUBROUTINE make_mo_eig
898
899END MODULE qs_mo_methods
Types and set/get functions for auxiliary density matrix methods.
Definition admm_types.F:15
Contains methods used in the context of density fitting.
Definition admm_utils.F:15
subroutine, public admm_uncorrect_for_eigenvalues(ispin, admm_env, ks_matrix)
...
Definition admm_utils.F:127
subroutine, public admm_correct_for_eigenvalues(ispin, admm_env, ks_matrix)
...
Definition admm_utils.F:53
methods related to the blacs parallel environment
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_release_p(matrix)
...
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_init_p(matrix)
...
Interface to (sca)lapack for the Cholesky based procedures.
subroutine, public cp_dbcsr_syevx(matrix, eigenvectors, eigenvalues, neig, work_syevx, para_env, blacs_env)
compute eigenvalues and optionally eigenvectors of a real symmetric matrix using scalapack....
subroutine, public cp_dbcsr_syevd(matrix, eigenvectors, eigenvalues, para_env, blacs_env)
...
DBCSR operations in CP2K.
subroutine, public cp_dbcsr_sm_fm_multiply(matrix, fm_in, fm_out, ncol, alpha, beta)
multiply a dbcsr with a fm matrix
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
subroutine, public cp_dbcsr_m_by_n_from_template(matrix, template, m, n, sym)
Utility function to create an arbitrary shaped dbcsr matrix with the same processor grid as the templ...
subroutine, public copy_fm_to_dbcsr(fm, matrix, keep_sparsity)
Copy a BLACS matrix to a dbcsr matrix.
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_syrk(uplo, trans, k, alpha, matrix_a, ia, ja, beta, matrix_c)
performs a rank-k update of a symmetric matrix_c matrix_c = beta * matrix_c + alpha * matrix_a * tran...
subroutine, public cp_fm_triangular_multiply(triangular_matrix, matrix_b, side, transpose_tr, invert_tr, uplo_tr, unit_diag_tr, n_rows, n_cols, alpha)
multiplies in place by a triangular matrix: matrix_b = alpha op(triangular_matrix) matrix_b or (if si...
various cholesky decomposition related routines
subroutine, public cp_fm_cholesky_decompose(matrix, n, info_out)
used to replace a symmetric positive def. matrix M with its cholesky decomposition U: M = U^T * U,...
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
Definition cp_fm_diag.F:17
subroutine, public cp_fm_power(matrix, work, exponent, threshold, n_dependent, verbose, eigvals)
...
subroutine, public choose_eigv_solver(matrix, eigenvectors, eigenvalues, info)
Choose the Eigensolver depending on which library is available ELPA seems to be unstable for small sy...
Definition cp_fm_diag.F:239
subroutine, public cp_fm_syevx(matrix, eigenvectors, eigenvalues, neig, work_syevx)
compute eigenvalues and optionally eigenvectors of a real symmetric matrix using scalapack....
Definition cp_fm_diag.F:672
represent the structure of a full matrix
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
subroutine, public cp_fm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
subroutine, public cp_fm_set_all(matrix, alpha, beta)
set all elements of a matrix to the same value, and optionally the diagonal to a different one
various routines to log and control the output. The idea is that decisions about where to log should ...
integer function, public cp_logger_get_default_io_unit(logger)
returns the unit nr for the ionode (-1 on all other processors) skips as well checks if the procs cal...
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Interface to the message passing library MPI.
basic linear algebra operations for full matrixes
Definition of physical constants:
Definition physcon.F:68
real(kind=dp), parameter, public evolt
Definition physcon.F:183
collects routines that perform operations directly related to MOs
subroutine, public make_basis_simple(vmatrix, ncol)
given a set of vectors, return an orthogonal (C^T C == 1) set spanning the same space (notice,...
subroutine, public make_basis_lowdin(vmatrix, ncol, matrix_s)
return a set of S orthonormal vectors (C^T S C == 1) where a Loedwin transformation is applied to kee...
subroutine, public make_basis_cholesky(vmatrix, ncol, ortho)
return a set of S orthonormal vectors (C^T S C == 1) where the cholesky decomposed form of S is passe...
subroutine, public calculate_magnitude(mo_array, mo_mag_min, mo_mag_max)
...
subroutine, public make_basis_sm(vmatrix, ncol, matrix_s)
returns an S-orthonormal basis v (v^T S v ==1)
subroutine, public make_mo_eig(mos, nspins, ks_rmpv, scf_control, mo_derivs, admm_env, hairy_probes, probe)
Calculate KS eigenvalues starting from OF MOS.
subroutine, public calculate_orthonormality(orthonormality, mo_array, matrix_s)
...
Set occupation of molecular orbitals.
Definition and initialisation of the mo data type.
Definition qs_mo_types.F:22
subroutine, public get_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, mo_coeff, mo_coeff_b, uniform_occupation, kts, mu, flexible_electron_count)
Get the components of a MO set data structure.
parameters that control an scf iteration
stores some data used in wavefunction fitting
Definition admm_types.F:120
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
keeps the information about the structure of a full matrix
represent a full matrix
stores all the informations relevant to an mpi environment