62#include "./base/base_uses.f90"
68 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'population_analyses'
89 INTEGER,
INTENT(IN) :: output_unit, print_level
91 CHARACTER(LEN=*),
PARAMETER :: routinen =
'lowdin_population_analysis'
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
100 TYPE(
cp_fm_type) :: fm_s_half, fm_work1, fm_work2
101 TYPE(
cp_fm_type),
DIMENSION(:),
POINTER :: fmwork
103 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrixkp_p, matrixkp_s
110 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
114 CALL timeset(routinen, handle)
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, &
125 scf_control=scf_control, &
130 nspin =
SIZE(matrixkp_p, 1)
131 nimg =
SIZE(matrixkp_p, 2)
136 ALLOCATE (orbpop(nsgf, nspin))
137 orbpop(:, :) = 0.0_dp
140 IF (output_unit > 0)
THEN
141 WRITE (unit=output_unit, fmt=
"(/,/,T2,A)")
"LOWDIN POPULATION ANALYSIS"
146 CALL get_qs_env(qs_env=qs_env, kpoints=kpoints, sab_orb=sab_nl)
150 nrow_global=nsgf, ncol_global=nsgf)
153 CALL cp_fm_create(fmwork(i), matrix_struct=fmstruct, name=
"work")
166 sm_s => matrixkp_s(1, 1)%matrix
167 matrix_p => matrixkp_p(:, 1)
176 matrix_struct=fmstruct, &
177 name=
"S^(1/2) MATRIX")
179 matrix_struct=fmstruct, &
180 name=
"FULL WORK MATRIX 1")
181 headline =
"SYMMETRICALLY ORTHOGONALISED DENSITY MATRIX"
183 matrix_struct=fmstruct, &
189 CALL cp_fm_power(fm_s_half, fm_work1, 0.5_dp, scf_control%eps_eigval, ndep)
191 CALL cp_warn(__location__, &
192 "Overlap matrix exhibits linear dependencies. At least some "// &
193 "eigenvalues have been quenched.")
197 sm_p => matrix_p(ispin)%matrix
206 matrix_a=fm_s_half, &
210 IF (print_level > 2)
THEN
214 fm_work2%name = trim(headline)//
" FOR ALPHA SPIN"
216 fm_work2%name = trim(headline)//
" FOR BETA SPIN"
220 output_unit=output_unit)
233 IF (output_unit > 0)
THEN
234 print_gop = (print_level > 1)
235 CALL write_orbpop(orbpop, atomic_kind_set, qs_kind_set, particle_set, output_unit, print_gop)
240 CALL timestop(handle)
257 INTEGER,
INTENT(IN) :: output_unit, print_level
259 CHARACTER(LEN=*),
PARAMETER :: routinen =
'mulliken_population_analysis'
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
267 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: orbpop
268 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: p_block, ps_block, s_block
271 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_p, matrix_s
272 TYPE(
dbcsr_type),
POINTER :: sm_p, sm_ps, sm_s
275 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
278 CALL timeset(routinen, handle)
280 NULLIFY (sm_p, sm_ps, sm_s)
281 NULLIFY (p_block, s_block, ps_block)
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, &
292 nspin =
SIZE(matrix_p, 1)
297 ALLOCATE (first_sgf_atom(natom))
298 first_sgf_atom(:) = 0
303 ALLOCATE (orbpop(nsgf, nspin))
304 orbpop(:, :) = 0.0_dp
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"
314 IF (print_level > 2)
THEN
315 sm_s => matrix_s(1, 1)%matrix
317 headline =
"MULLIKEN NET ATOMIC ORBITAL AND OVERLAP POPULATION MATRIX"
320 headline = trim(headline)//
" For Alpha Spin"
322 headline = trim(headline)//
" For Beta Spin"
325 CALL dbcsr_copy(matrix_b=sm_ps, matrix_a=sm_s, name=trim(headline))
330 DO ic = 1,
SIZE(matrix_s, 2)
331 IF (print_level > 2)
THEN
334 sm_s => matrix_s(1, ic)%matrix
335 sm_p => matrix_p(ispin, ic)%matrix
341 IF (.NOT. (
ASSOCIATED(s_block))) cycle
347 IF (print_level > 2)
THEN
353 cpassert(
ASSOCIATED(ps_block))
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
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
379 IF (print_level > 2)
THEN
385 CALL para_env%sum(orbpop)
388 IF (output_unit > 0)
THEN
389 print_gop = (print_level > 1)
390 CALL write_orbpop(orbpop, atomic_kind_set, qs_kind_set, particle_set, output_unit, print_gop)
394 CALL save_mulliken_charges(orbpop, atomic_kind_set, qs_kind_set, particle_set, qs_env)
399 DEALLOCATE (first_sgf_atom)
401 IF (output_unit > 0)
THEN
402 WRITE (unit=output_unit, fmt=
"(T2,A)") &
403 '!-----------------------------------------------------------------------------!'
406 CALL timestop(handle)
425 SUBROUTINE save_mulliken_charges(orbpop, atomic_kind_set, qs_kind_set, particle_set, qs_env)
426 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: orbpop
428 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
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
443 nspin =
SIZE(orbpop, 2)
450 ALLOCATE (all_sumorbpop(natom))
451 ALLOCATE (charges_save(natom))
455 sumorbpop(:) = 0.0_dp
456 NULLIFY (orb_basis_set)
459 CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set, zeff=zeff)
460 IF (
ASSOCIATED(orb_basis_set))
THEN
467 DO ishell = 1, nshell(iset)
468 l = lshell(ishell, iset)
471 sumorbpop(1) = sumorbpop(1) + orbpop(iao, 1)
473 sumorbpop(1:2) = sumorbpop(1:2) + orbpop(iao, 1:2)
474 sumorbpop(3) = sumorbpop(3) + orbpop(iao, 1) - orbpop(iao, 2)
482 charges_save(iatom) = zeff - sumorbpop(1)
483 all_sumorbpop(iatom) = sumorbpop(1)
485 charges_save(iatom) = zeff - sumorbpop(1) - sumorbpop(2)
486 all_sumorbpop(iatom) = sumorbpop(1) + sumorbpop(2)
492 description =
"[MULLIKEN-ORBPOP]"
494 CALL put_results(results=results, description=description, &
498 description =
"[MULLIKEN-SUMORBPOP]"
500 CALL put_results(results=results, description=description, &
501 values=all_sumorbpop)
504 description =
"[MULLIKEN-CHARGES]"
506 CALL put_results(results=results, description=description, &
509 DEALLOCATE (all_sumorbpop)
510 DEALLOCATE (charges_save)
512 END SUBROUTINE save_mulliken_charges
527 SUBROUTINE write_orbpop(orbpop, atomic_kind_set, qs_kind_set, particle_set, output_unit, &
528 print_orbital_contributions)
530 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: orbpop
532 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
534 INTEGER,
INTENT(IN) :: output_unit
535 LOGICAL,
INTENT(IN) :: print_orbital_contributions
537 CHARACTER(LEN=*),
PARAMETER :: routinen =
'write_orbpop'
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, &
544 INTEGER,
DIMENSION(:),
POINTER :: nshell
545 INTEGER,
DIMENSION(:, :),
POINTER :: lshell
546 REAL(kind=
dp) :: zeff
547 REAL(kind=
dp),
DIMENSION(3) :: sumorbpop, totsumorbpop
550 CALL timeset(routinen, handle)
552 nspin =
SIZE(orbpop, 2)
559 IF (print_orbital_contributions)
THEN
560 WRITE (unit=output_unit, fmt=
"(/,T2,A)") &
561 "# Orbital AO symbol Orbital population Net charge"
563 WRITE (unit=output_unit, fmt=
"(/,T2,A)") &
564 "# Atom Element Kind Atomic population Net charge"
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"
571 WRITE (unit=output_unit, fmt=
"(/,T2,A)") &
572 "# Atom Element Kind Atomic population (alpha,beta) Net charge Spin moment"
576 totsumorbpop(:) = 0.0_dp
580 sumorbpop(:) = 0.0_dp
581 NULLIFY (orb_basis_set)
583 element_symbol=element_symbol, &
585 CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set, zeff=zeff)
586 IF (
ASSOCIATED(orb_basis_set))
THEN
591 sgf_symbol=sgf_symbol)
594 DO ishell = 1, nshell(iset)
595 l = lshell(ishell, iset)
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)
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)
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)
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)
639 IF (print_orbital_contributions)
WRITE (unit=output_unit, fmt=
"(A)")
""
641 WRITE (unit=output_unit, &
642 fmt=
"(T2,A,T42,F12.6,T68,F12.6,/)") &
643 "# Total charge", totsumorbpop(1), totsumorbpop(3)
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)
651 IF (output_unit > 0)
CALL m_flush(output_unit)
653 CALL timestop(handle)
655 END SUBROUTINE write_orbpop
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.
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...
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
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.
integer, parameter, public dp
integer, parameter, public default_string_length
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.
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
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
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.