(git:495eafe)
Loading...
Searching...
No Matches
population_analyses.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!> \brief Provide various population analyses and print the requested output
9!> information
10!>
11!> \author Matthias Krack (MK)
12!> \date 09.07.2010
13!> \version 1.0
14! **************************************************************************************************
15
23 USE cp_dbcsr_api, ONLY: &
31 USE cp_fm_diag, ONLY: cp_fm_power
35 USE cp_fm_types, ONLY: cp_fm_create,&
42 USE kinds, ONLY: default_string_length,&
43 dp
45 USE kpoint_types, ONLY: kpoint_type
46 USE machine, ONLY: m_flush
48 USE orbital_pointers, ONLY: nso
54 USE qs_kind_types, ONLY: get_qs_kind,&
58 USE qs_rho_types, ONLY: qs_rho_get,&
62#include "./base/base_uses.f90"
63
64 IMPLICIT NONE
65
66 PRIVATE
67
68 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'population_analyses'
69
72
73CONTAINS
74
75! **************************************************************************************************
76!> \brief Perform a Lowdin population analysis based on a symmetric
77!> orthogonalisation of the density matrix using S^(1/2)
78!>
79!> \param qs_env ...
80!> \param output_unit ...
81!> \param print_level ...
82!> \date 06.07.2010
83!> \author Matthias Krack (MK)
84!> \version 1.0
85! **************************************************************************************************
86 SUBROUTINE lowdin_population_analysis(qs_env, output_unit, print_level)
87
88 TYPE(qs_environment_type), POINTER :: qs_env
89 INTEGER, INTENT(IN) :: output_unit, print_level
90
91 CHARACTER(LEN=*), PARAMETER :: routinen = 'lowdin_population_analysis'
92
93 CHARACTER(LEN=default_string_length) :: headline
94 INTEGER :: handle, i, ispin, ndep, nimg, nsgf, nspin
95 LOGICAL :: do_kpoints, print_gop
96 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: orbpop
97 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
98 TYPE(cp_blacs_env_type), POINTER :: blacs_env
99 TYPE(cp_fm_struct_type), POINTER :: fmstruct
100 TYPE(cp_fm_type) :: fm_s_half, fm_work1, fm_work2
101 TYPE(cp_fm_type), DIMENSION(:), POINTER :: fmwork
102 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_p
103 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrixkp_p, matrixkp_s
104 TYPE(dbcsr_type), POINTER :: sm_p, sm_s
105 TYPE(kpoint_type), POINTER :: kpoints
106 TYPE(mp_para_env_type), POINTER :: para_env
107 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
108 POINTER :: sab_nl
109 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
110 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
111 TYPE(qs_rho_type), POINTER :: rho
112 TYPE(scf_control_type), POINTER :: scf_control
113
114 CALL timeset(routinen, handle)
115
116 NULLIFY (sm_p, sm_s)
117
118 CALL get_qs_env(qs_env=qs_env, &
119 atomic_kind_set=atomic_kind_set, &
120 qs_kind_set=qs_kind_set, &
121 do_kpoints=do_kpoints, &
122 matrix_s_kp=matrixkp_s, &
123 particle_set=particle_set, &
124 rho=rho, &
125 scf_control=scf_control, &
126 para_env=para_env, &
127 blacs_env=blacs_env)
128
129 CALL qs_rho_get(rho, rho_ao_kp=matrixkp_p) ! Density matrices in sparse format
130 nspin = SIZE(matrixkp_p, 1)
131 nimg = SIZE(matrixkp_p, 2)
132
133 ! Get the total number of contracted spherical Gaussian basis functions
134 CALL get_qs_kind_set(qs_kind_set, nsgf=nsgf)
135 ! Provide an array to store the orbital populations for each spin
136 ALLOCATE (orbpop(nsgf, nspin))
137 orbpop(:, :) = 0.0_dp
138
139 ! Write headline
140 IF (output_unit > 0) THEN
141 WRITE (unit=output_unit, fmt="(/,/,T2,A)") "LOWDIN POPULATION ANALYSIS"
142 END IF
143
144 IF (do_kpoints) THEN
145
146 CALL get_qs_env(qs_env=qs_env, kpoints=kpoints, sab_orb=sab_nl)
147
148 ! work matrices
149 CALL cp_fm_struct_create(fmstruct, para_env=para_env, context=blacs_env, &
150 nrow_global=nsgf, ncol_global=nsgf)
151 ALLOCATE (fmwork(4))
152 DO i = 1, 4
153 CALL cp_fm_create(fmwork(i), matrix_struct=fmstruct, name="work")
154 END DO
155 CALL cp_fm_struct_release(fmstruct)
156
157 ! S^1/2
158 CALL diag_kp_smat(matrixkp_s, kpoints, fmwork)
159 ! Lowdin P Matrix Transform
160 CALL lowdin_kp_trans(kpoints, orbpop)
161
162 CALL cp_fm_release(fmwork)
163
164 ELSE
165
166 sm_s => matrixkp_s(1, 1)%matrix ! Overlap matrix in sparse format
167 matrix_p => matrixkp_p(:, 1)
168
169 ! Provide full size work matrices
170 CALL cp_fm_struct_create(fmstruct=fmstruct, &
171 para_env=para_env, &
172 context=blacs_env, &
173 nrow_global=nsgf, &
174 ncol_global=nsgf)
175 CALL cp_fm_create(matrix=fm_s_half, &
176 matrix_struct=fmstruct, &
177 name="S^(1/2) MATRIX")
178 CALL cp_fm_create(matrix=fm_work1, &
179 matrix_struct=fmstruct, &
180 name="FULL WORK MATRIX 1")
181 headline = "SYMMETRICALLY ORTHOGONALISED DENSITY MATRIX"
182 CALL cp_fm_create(matrix=fm_work2, &
183 matrix_struct=fmstruct, &
184 name=trim(headline))
185 CALL cp_fm_struct_release(fmstruct=fmstruct)
186
187 ! Build full S^(1/2) matrix (computationally expensive)
188 CALL copy_dbcsr_to_fm(sm_s, fm_s_half)
189 CALL cp_fm_power(fm_s_half, fm_work1, 0.5_dp, scf_control%eps_eigval, ndep)
190 IF (ndep /= 0) &
191 CALL cp_warn(__location__, &
192 "Overlap matrix exhibits linear dependencies. At least some "// &
193 "eigenvalues have been quenched.")
194
195 ! Build Lowdin population matrix for each spin
196 DO ispin = 1, nspin
197 sm_p => matrix_p(ispin)%matrix ! Density matrix for spin ispin in sparse format
198 ! Calculate S^(1/2)*P*S^(1/2) as a full matrix (Lowdin)
199 CALL cp_dbcsr_sm_fm_multiply(sm_p, fm_s_half, fm_work1, nsgf)
200 CALL parallel_gemm(transa="N", &
201 transb="N", &
202 m=nsgf, &
203 n=nsgf, &
204 k=nsgf, &
205 alpha=1.0_dp, &
206 matrix_a=fm_s_half, &
207 matrix_b=fm_work1, &
208 beta=0.0_dp, &
209 matrix_c=fm_work2)
210 IF (print_level > 2) THEN
211 ! Write the full Lowdin population matrix
212 IF (nspin > 1) THEN
213 IF (ispin == 1) THEN
214 fm_work2%name = trim(headline)//" FOR ALPHA SPIN"
215 ELSE
216 fm_work2%name = trim(headline)//" FOR BETA SPIN"
217 END IF
218 END IF
219 CALL write_fm_with_basis_info(fm_work2, 4, 6, qs_env, para_env, &
220 output_unit=output_unit)
221 END IF
222 CALL cp_fm_get_diag(fm_work2, orbpop(:, ispin))
223 END DO ! next spin ispin
224
225 ! Release local working storage
226 CALL cp_fm_release(matrix=fm_s_half)
227 CALL cp_fm_release(matrix=fm_work1)
228 CALL cp_fm_release(matrix=fm_work2)
229
230 END IF
231
232 ! Write atomic populations and charges
233 IF (output_unit > 0) THEN
234 print_gop = (print_level > 1) ! Print also orbital populations
235 CALL write_orbpop(orbpop, atomic_kind_set, qs_kind_set, particle_set, output_unit, print_gop)
236 END IF
237
238 DEALLOCATE (orbpop)
239
240 CALL timestop(handle)
241
242 END SUBROUTINE lowdin_population_analysis
243
244! **************************************************************************************************
245!> \brief Perform a Mulliken population analysis
246!>
247!> \param qs_env ...
248!> \param output_unit ...
249!> \param print_level ...
250!> \date 10.07.2010
251!> \author Matthias Krack (MK)
252!> \version 1.0
253! **************************************************************************************************
254 SUBROUTINE mulliken_population_analysis(qs_env, output_unit, print_level)
255
256 TYPE(qs_environment_type), POINTER :: qs_env
257 INTEGER, INTENT(IN) :: output_unit, print_level
258
259 CHARACTER(LEN=*), PARAMETER :: routinen = 'mulliken_population_analysis'
260
261 CHARACTER(LEN=default_string_length) :: headline
262 INTEGER :: handle, iatom, ic, isgf, ispin, jatom, &
263 jsgf, natom, nsgf, nspin, sgfa, sgfb
264 INTEGER, ALLOCATABLE, DIMENSION(:) :: first_sgf_atom
265 LOGICAL :: found, print_gop
266 REAL(kind=dp) :: ps
267 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: orbpop
268 REAL(kind=dp), DIMENSION(:, :), POINTER :: p_block, ps_block, s_block
269 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
270 TYPE(dbcsr_iterator_type) :: iter
271 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_p, matrix_s
272 TYPE(dbcsr_type), POINTER :: sm_p, sm_ps, sm_s
273 TYPE(mp_para_env_type), POINTER :: para_env
274 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
275 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
276 TYPE(qs_rho_type), POINTER :: rho
277
278 CALL timeset(routinen, handle)
279
280 NULLIFY (sm_p, sm_ps, sm_s)
281 NULLIFY (p_block, s_block, ps_block)
282
283 CALL get_qs_env(qs_env=qs_env, &
284 atomic_kind_set=atomic_kind_set, &
285 qs_kind_set=qs_kind_set, &
286 matrix_s_kp=matrix_s, &
287 particle_set=particle_set, &
288 rho=rho, &
289 para_env=para_env)
290
291 CALL qs_rho_get(rho, rho_ao_kp=matrix_p) ! Density matrices in sparse format
292 nspin = SIZE(matrix_p, 1)
293
294 ! Get the total number of contracted spherical Gaussian basis functions
295 CALL get_atomic_kind_set(atomic_kind_set, natom=natom)
296 CALL get_qs_kind_set(qs_kind_set, nsgf=nsgf)
297 ALLOCATE (first_sgf_atom(natom))
298 first_sgf_atom(:) = 0
299
300 CALL get_particle_set(particle_set, qs_kind_set, first_sgf=first_sgf_atom)
301
302 ! Provide an array to store the orbital populations for each spin
303 ALLOCATE (orbpop(nsgf, nspin))
304 orbpop(:, :) = 0.0_dp
305
306 ! Write headline
307 IF (output_unit > 0) THEN
308 WRITE (unit=output_unit, fmt="(/,T2,A)") &
309 '!-----------------------------------------------------------------------------!'
310 WRITE (unit=output_unit, fmt="(T22,A)") "Mulliken Population Analysis"
311 END IF
312
313 ! Create a DBCSR work matrix, if needed
314 IF (print_level > 2) THEN
315 sm_s => matrix_s(1, 1)%matrix ! Overlap matrix in sparse format
316 ALLOCATE (sm_ps)
317 headline = "MULLIKEN NET ATOMIC ORBITAL AND OVERLAP POPULATION MATRIX"
318 IF (nspin > 1) THEN
319 IF (ispin == 1) THEN
320 headline = trim(headline)//" For Alpha Spin"
321 ELSE
322 headline = trim(headline)//" For Beta Spin"
323 END IF
324 END IF
325 CALL dbcsr_copy(matrix_b=sm_ps, matrix_a=sm_s, name=trim(headline))
326 END IF
327
328 ! Build Mulliken population matrix for each spin
329 DO ispin = 1, nspin
330 DO ic = 1, SIZE(matrix_s, 2)
331 IF (print_level > 2) THEN
332 CALL dbcsr_set(sm_ps, 0.0_dp)
333 END IF
334 sm_s => matrix_s(1, ic)%matrix ! Overlap matrix in sparse format
335 sm_p => matrix_p(ispin, ic)%matrix ! Density matrix for spin ispin in sparse format
336 ! Calculate Hadamard product of P and S as sparse matrix (Mulliken)
337 ! CALL dbcsr_hadamard_product(sm_p,sm_s,sm_ps)
338 CALL dbcsr_iterator_start(iter, sm_s)
339 DO WHILE (dbcsr_iterator_blocks_left(iter))
340 CALL dbcsr_iterator_next_block(iter, iatom, jatom, s_block)
341 IF (.NOT. (ASSOCIATED(s_block))) cycle
342 CALL dbcsr_get_block_p(matrix=sm_p, &
343 row=iatom, &
344 col=jatom, &
345 block=p_block, &
346 found=found)
347 IF (print_level > 2) THEN
348 CALL dbcsr_get_block_p(matrix=sm_ps, &
349 row=iatom, &
350 col=jatom, &
351 block=ps_block, &
352 found=found)
353 cpassert(ASSOCIATED(ps_block))
354 END IF
355
356 sgfb = first_sgf_atom(jatom)
357 DO jsgf = 1, SIZE(s_block, 2)
358 DO isgf = 1, SIZE(s_block, 1)
359 ps = p_block(isgf, jsgf)*s_block(isgf, jsgf)
360 IF (ASSOCIATED(ps_block)) ps_block(isgf, jsgf) = ps_block(isgf, jsgf) + ps
361 orbpop(sgfb, ispin) = orbpop(sgfb, ispin) + ps
362 END DO
363 sgfb = sgfb + 1
364 END DO
365 IF (iatom /= jatom) THEN
366 sgfa = first_sgf_atom(iatom)
367 DO isgf = 1, SIZE(s_block, 1)
368 DO jsgf = 1, SIZE(s_block, 2)
369 ps = p_block(isgf, jsgf)*s_block(isgf, jsgf)
370 orbpop(sgfa, ispin) = orbpop(sgfa, ispin) + ps
371 END DO
372 sgfa = sgfa + 1
373 END DO
374 END IF
375 END DO
376 CALL dbcsr_iterator_stop(iter)
377 END DO
378
379 IF (print_level > 2) THEN
380 ! Write the full Mulliken net AO and overlap population matrix
381 CALL cp_dbcsr_write_sparse_matrix(sm_ps, 4, 6, qs_env, para_env, output_unit=output_unit)
382 END IF
383 END DO
384
385 CALL para_env%sum(orbpop)
386
387 ! Write atomic populations and charges
388 IF (output_unit > 0) THEN
389 print_gop = (print_level > 1) ! Print also orbital populations
390 CALL write_orbpop(orbpop, atomic_kind_set, qs_kind_set, particle_set, output_unit, print_gop)
391 END IF
392
393 ! Save the Mulliken charges to results
394 CALL save_mulliken_charges(orbpop, atomic_kind_set, qs_kind_set, particle_set, qs_env)
395
396 ! Release local working storage
397 IF (ASSOCIATED(sm_ps)) CALL dbcsr_deallocate_matrix(sm_ps)
398 DEALLOCATE (orbpop)
399 DEALLOCATE (first_sgf_atom)
400
401 IF (output_unit > 0) THEN
402 WRITE (unit=output_unit, fmt="(T2,A)") &
403 '!-----------------------------------------------------------------------------!'
404 END IF
405
406 CALL timestop(handle)
407
408 END SUBROUTINE mulliken_population_analysis
409
410! **************************************************************************************************
411!> \brief Save the Mulliken atomic orbital populations and charges in results
412!>
413!> \param orbpop ...
414!> \param atomic_kind_set ...
415!> \param qs_kind_set ...
416!> \param particle_set ...
417!> \param qs_env ...
418!> \par History
419!> 27.05.2022 BT
420!> 16.07.2025 RK
421!> \author Bo Thomsen (BT)
422!> Rangsiman Ketkaew (RK)
423!> \version 1.0
424! **************************************************************************************************
425 SUBROUTINE save_mulliken_charges(orbpop, atomic_kind_set, qs_kind_set, particle_set, qs_env)
426 REAL(kind=dp), DIMENSION(:, :), INTENT(IN) :: orbpop
427 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
428 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
429 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
430 TYPE(qs_environment_type), POINTER :: qs_env
431
432 CHARACTER(LEN=default_string_length) :: description
433 INTEGER :: iao, iatom, ikind, iset, isgf, ishell, &
434 iso, l, natom, nset, nsgf, nspin
435 INTEGER, DIMENSION(:), POINTER :: nshell
436 INTEGER, DIMENSION(:, :), POINTER :: lshell
437 REAL(kind=dp) :: zeff
438 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: all_sumorbpop, charges_save
439 REAL(kind=dp), DIMENSION(3) :: sumorbpop
440 TYPE(cp_result_type), POINTER :: results
441 TYPE(gto_basis_set_type), POINTER :: orb_basis_set
442
443 nspin = SIZE(orbpop, 2)
444
445 CALL get_atomic_kind_set(atomic_kind_set, natom=natom)
446 CALL get_qs_kind_set(qs_kind_set, nsgf=nsgf)
447 CALL get_atomic_kind_set(atomic_kind_set, natom=natom)
448 NULLIFY (results)
449 CALL get_qs_env(qs_env, results=results)
450 ALLOCATE (all_sumorbpop(natom))
451 ALLOCATE (charges_save(natom))
452
453 iao = 1
454 DO iatom = 1, natom
455 sumorbpop(:) = 0.0_dp
456 NULLIFY (orb_basis_set)
457 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
458 kind_number=ikind)
459 CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set, zeff=zeff)
460 IF (ASSOCIATED(orb_basis_set)) THEN
461 CALL get_gto_basis_set(gto_basis_set=orb_basis_set, &
462 nset=nset, &
463 nshell=nshell, &
464 l=lshell)
465 isgf = 1
466 DO iset = 1, nset
467 DO ishell = 1, nshell(iset)
468 l = lshell(ishell, iset)
469 DO iso = 1, nso(l)
470 IF (nspin == 1) THEN
471 sumorbpop(1) = sumorbpop(1) + orbpop(iao, 1)
472 ELSE
473 sumorbpop(1:2) = sumorbpop(1:2) + orbpop(iao, 1:2)
474 sumorbpop(3) = sumorbpop(3) + orbpop(iao, 1) - orbpop(iao, 2)
475 END IF
476 isgf = isgf + 1
477 iao = iao + 1
478 END DO
479 END DO
480 END DO
481 IF (nspin == 1) THEN
482 charges_save(iatom) = zeff - sumorbpop(1)
483 all_sumorbpop(iatom) = sumorbpop(1)
484 ELSE
485 charges_save(iatom) = zeff - sumorbpop(1) - sumorbpop(2)
486 all_sumorbpop(iatom) = sumorbpop(1) + sumorbpop(2)
487 END IF
488 END IF ! atom has an orbital basis
489 END DO ! next atom iatom
490
491 ! Store atomic orbital populations in results
492 description = "[MULLIKEN-ORBPOP]"
493 CALL cp_results_erase(results=results, description=description)
494 CALL put_results(results=results, description=description, &
495 values=orbpop)
496
497 ! Store sum orbital population in results
498 description = "[MULLIKEN-SUMORBPOP]"
499 CALL cp_results_erase(results=results, description=description)
500 CALL put_results(results=results, description=description, &
501 values=all_sumorbpop)
502
503 ! Store charges in results
504 description = "[MULLIKEN-CHARGES]"
505 CALL cp_results_erase(results=results, description=description)
506 CALL put_results(results=results, description=description, &
507 values=charges_save)
508
509 DEALLOCATE (all_sumorbpop)
510 DEALLOCATE (charges_save)
511
512 END SUBROUTINE save_mulliken_charges
513
514! **************************************************************************************************
515!> \brief Write atomic orbital populations and net atomic charges
516!>
517!> \param orbpop ...
518!> \param atomic_kind_set ...
519!> \param qs_kind_set ...
520!> \param particle_set ...
521!> \param output_unit ...
522!> \param print_orbital_contributions ...
523!> \date 07.07.2010
524!> \author Matthias Krack (MK)
525!> \version 1.0
526! **************************************************************************************************
527 SUBROUTINE write_orbpop(orbpop, atomic_kind_set, qs_kind_set, particle_set, output_unit, &
528 print_orbital_contributions)
529
530 REAL(kind=dp), DIMENSION(:, :), INTENT(IN) :: orbpop
531 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
532 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
533 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
534 INTEGER, INTENT(IN) :: output_unit
535 LOGICAL, INTENT(IN) :: print_orbital_contributions
536
537 CHARACTER(LEN=*), PARAMETER :: routinen = 'write_orbpop'
538
539 CHARACTER(LEN=2) :: element_symbol
540 CHARACTER(LEN=6), DIMENSION(:), POINTER :: sgf_symbol
541 INTEGER :: handle, iao, iatom, ikind, iset, isgf, &
542 ishell, iso, l, natom, nset, nsgf, &
543 nspin
544 INTEGER, DIMENSION(:), POINTER :: nshell
545 INTEGER, DIMENSION(:, :), POINTER :: lshell
546 REAL(kind=dp) :: zeff
547 REAL(kind=dp), DIMENSION(3) :: sumorbpop, totsumorbpop
548 TYPE(gto_basis_set_type), POINTER :: orb_basis_set
549
550 CALL timeset(routinen, handle)
551
552 nspin = SIZE(orbpop, 2)
553
554 CALL get_atomic_kind_set(atomic_kind_set, natom=natom)
555 CALL get_qs_kind_set(qs_kind_set, nsgf=nsgf)
556
557 ! Select and write headline
558 IF (nspin == 1) THEN
559 IF (print_orbital_contributions) THEN
560 WRITE (unit=output_unit, fmt="(/,T2,A)") &
561 "# Orbital AO symbol Orbital population Net charge"
562 ELSE
563 WRITE (unit=output_unit, fmt="(/,T2,A)") &
564 "# Atom Element Kind Atomic population Net charge"
565 END IF
566 ELSE
567 IF (print_orbital_contributions) THEN
568 WRITE (unit=output_unit, fmt="(/,T2,A)") &
569 "# Orbital AO symbol Orbital population (alpha,beta) Net charge Spin moment"
570 ELSE
571 WRITE (unit=output_unit, fmt="(/,T2,A)") &
572 "# Atom Element Kind Atomic population (alpha,beta) Net charge Spin moment"
573 END IF
574 END IF
575
576 totsumorbpop(:) = 0.0_dp
577
578 iao = 1
579 DO iatom = 1, natom
580 sumorbpop(:) = 0.0_dp
581 NULLIFY (orb_basis_set)
582 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
583 element_symbol=element_symbol, &
584 kind_number=ikind)
585 CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set, zeff=zeff)
586 IF (ASSOCIATED(orb_basis_set)) THEN
587 CALL get_gto_basis_set(gto_basis_set=orb_basis_set, &
588 nset=nset, &
589 nshell=nshell, &
590 l=lshell, &
591 sgf_symbol=sgf_symbol)
592 isgf = 1
593 DO iset = 1, nset
594 DO ishell = 1, nshell(iset)
595 l = lshell(ishell, iset)
596 DO iso = 1, nso(l)
597 IF (nspin == 1) THEN
598 sumorbpop(1) = sumorbpop(1) + orbpop(iao, 1)
599 IF (print_orbital_contributions) THEN
600 IF (isgf == 1) WRITE (unit=output_unit, fmt="(A)") ""
601 WRITE (unit=output_unit, &
602 fmt="(T2,I9,2X,A2,1X,A,T30,F12.6)") &
603 iao, element_symbol, sgf_symbol(isgf), orbpop(iao, 1)
604 END IF
605 ELSE
606 sumorbpop(1:2) = sumorbpop(1:2) + orbpop(iao, 1:2)
607 sumorbpop(3) = sumorbpop(3) + orbpop(iao, 1) - orbpop(iao, 2)
608 IF (print_orbital_contributions) THEN
609 IF (isgf == 1) WRITE (unit=output_unit, fmt="(A)") ""
610 WRITE (unit=output_unit, &
611 fmt="(T2,I9,2X,A2,1X,A,T29,2(1X,F12.6),T68,F12.6)") &
612 iao, element_symbol, sgf_symbol(isgf), orbpop(iao, 1:2), &
613 orbpop(iao, 1) - orbpop(iao, 2)
614 END IF
615 END IF
616 isgf = isgf + 1
617 iao = iao + 1
618 END DO
619 END DO
620 END DO
621 IF (nspin == 1) THEN
622 totsumorbpop(1) = totsumorbpop(1) + sumorbpop(1)
623 totsumorbpop(3) = totsumorbpop(3) + zeff - sumorbpop(1)
624 WRITE (unit=output_unit, &
625 fmt="(T2,I7,5X,A2,2X,I6,T30,F12.6,T68,F12.6)") &
626 iatom, element_symbol, ikind, sumorbpop(1), zeff - sumorbpop(1)
627 ELSE
628 totsumorbpop(1:2) = totsumorbpop(1:2) + sumorbpop(1:2)
629 totsumorbpop(3) = totsumorbpop(3) + zeff - sumorbpop(1) - sumorbpop(2)
630 WRITE (unit=output_unit, &
631 fmt="(T2,I7,5X,A2,2X,I6,T28,4(1X,F12.6))") &
632 iatom, element_symbol, ikind, sumorbpop(1:2), &
633 zeff - sumorbpop(1) - sumorbpop(2), sumorbpop(3)
634 END IF
635 END IF ! atom has an orbital basis
636 END DO ! next atom iatom
637
638 ! Write total sums
639 IF (print_orbital_contributions) WRITE (unit=output_unit, fmt="(A)") ""
640 IF (nspin == 1) THEN
641 WRITE (unit=output_unit, &
642 fmt="(T2,A,T42,F12.6,T68,F12.6,/)") &
643 "# Total charge", totsumorbpop(1), totsumorbpop(3)
644 ELSE
645 WRITE (unit=output_unit, &
646 fmt="(T2,A,T28,4(1X,F12.6),/)") &
647 "# Total charge and spin", totsumorbpop(1:2), totsumorbpop(3), &
648 totsumorbpop(1) - totsumorbpop(2)
649 END IF
650
651 IF (output_unit > 0) CALL m_flush(output_unit)
652
653 CALL timestop(handle)
654
655 END SUBROUTINE write_orbpop
656
657END MODULE population_analyses
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind_set(atomic_kind_set, atom_of_kind, kind_of, natom_of_kind, maxatom, natom, nshell, fist_potential_present, shell_present, shell_adiabatic, shell_check_distance, damping_present)
Get attributes of an atomic kind set.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
subroutine, public get_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, lmin, lx, ly, lz, m, ncgf_set, npgf, nsgf_set, nshell, cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, last_cgf, last_sgf, n, gcc, maxco, maxl, maxpgf, maxsgf_set, maxshell, maxso, nco_sum, npgf_sum, nshell_sum, maxder, short_kind_radius, npgf_seg_sum)
...
methods related to the blacs parallel environment
subroutine, public dbcsr_deallocate_matrix(matrix)
...
subroutine, public dbcsr_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset)
...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_get_block_p(matrix, row, col, block, found, row_size, col_size)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_set(matrix, alpha)
...
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.
DBCSR output in CP2K.
subroutine, public write_fm_with_basis_info(blacs_matrix, before, after, qs_env, para_env, first_row, last_row, first_col, last_col, output_unit, omit_headers)
Print a spherical matrix of blacs type.
subroutine, public cp_dbcsr_write_sparse_matrix(sparse_matrix, before, after, qs_env, para_env, first_row, last_row, first_col, last_col, scale, output_unit, omit_headers)
...
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)
...
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_diag(matrix, diag)
returns the diagonal elements of a fm
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
set of type/routines to handle the storage of results in force_envs
subroutine, public cp_results_erase(results, description, nval)
erase a part of result_list
set of type/routines to handle the storage of results in force_envs
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
Routines needed for kpoint calculation.
subroutine, public lowdin_kp_trans(kpoint, pmat_diag)
Calculate Lowdin transformation of density matrix S^1/2 P S^1/2 Integrate diagonal elements over k-po...
Types and basic routines needed for a kpoint calculation.
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
Interface to the message passing library MPI.
Provides Cartesian and spherical orbital pointers and indices.
integer, dimension(:), allocatable, public nso
basic linear algebra operations for full matrixes
Define methods related to particle_type.
subroutine, public get_particle_set(particle_set, qs_kind_set, first_sgf, last_sgf, nsgf, nmao, basis)
Get the components of a particle set.
Define the data structure for the particle information.
Provide various population analyses and print the requested output information.
subroutine, public lowdin_population_analysis(qs_env, output_unit, print_level)
Perform a Lowdin population analysis based on a symmetric orthogonalisation of the density matrix usi...
subroutine, public mulliken_population_analysis(qs_env, output_unit, print_level)
Perform a Mulliken population analysis.
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.
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, cneo_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zatom, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, monovalent, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
subroutine, public get_qs_kind_set(qs_kind_set, all_potential_present, tnadd_potential_present, gth_potential_present, sgp_potential_present, paw_atom_present, dft_plus_u_atom_present, maxcgf, maxsgf, maxco, maxco_proj, maxgtops, maxlgto, maxlprj, maxnset, maxsgf_set, ncgf, npgf, nset, nsgf, nshell, maxpol, maxlppl, maxlppnl, maxppnl, nelectron, maxder, max_ngrid_rad, max_sph_harm, maxg_iso_not0, lmax_rho0, basis_rcut, basis_type, total_zeff_corr, npgf_seg, cneo_potential_present, nkind_q, natom_q)
Get attributes of an atomic kind set.
Define the neighbor list data types and the corresponding functionality.
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...
Different diagonalization schemes that can be used for the iterative solution of the eigenvalue probl...
subroutine, public diag_kp_smat(matrix_s, kpoints, fmwork)
Kpoint diagonalization routine Transforms matrices to kpoint, distributes kpoint groups,...
parameters that control an scf iteration
Provides all information about an atomic kind.
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
contains arbitrary information which need to be stored
Contains information about kpoints.
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.
keeps the density in various representations, keeping track of which ones are valid.