85#include "./base/base_uses.f90"
93 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'ec_diag_solver'
117 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_ks, matrix_s, matrix_p, matrix_w
119 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ec_diag_solver_gamma'
121 INTEGER :: handle, ispin, nmo(2), nsize, nspins
122 REAL(kind=
dp) :: eps_filter, flexible_electron_count, &
124 REAL(kind=
dp),
DIMENSION(:),
POINTER :: eigvals
131 TYPE(
mo_set_type),
ALLOCATABLE,
DIMENSION(:) :: moset
134 CALL timeset(routinen, handle)
136 CALL get_qs_env(qs_env=qs_env, dft_control=dft_control)
137 eps_filter = dft_control%qs_control%eps_filter_matrix
138 nspins = dft_control%nspins
141 CALL get_qs_env(qs_env=qs_env, nelectron_spin=nmo)
143 IF (nspins == 1) focc = 2._dp
145 CALL dbcsr_get_info(matrix_ks(1, 1)%matrix, nfullrows_total=nsize)
147 NULLIFY (blacs_env, para_env)
148 CALL get_qs_env(qs_env=qs_env, blacs_env=blacs_env, para_env=para_env)
151 ncol_global=nsize, para_env=para_env)
158 flexible_electron_count = 0.0_dp
159 IF (ec_env%smear%do_smear) flexible_electron_count = 0.0001_dp
161 ALLOCATE (moset(nspins))
165 focc(ispin), flexible_electron_count)
166 CALL init_mo_set(moset(ispin), fm_ref=fm_w, name=
"MO")
170 ref_matrix => matrix_s(1, 1)%matrix
172 ref_matrix => matrix_ks(ispin, 1)%matrix
174 CALL get_mo_set(moset(ispin), eigenvalues=eigvals, mo_coeff=fm_mos)
175 CALL cp_fm_geeig(fm_k, fm_s, fm_mos, eigvals, fm_w)
181 ec_env%ekTS = ec_env%ekTS + moset(ispin)%kTS
194 CALL timestop(handle)
214 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_ks, matrix_s, matrix_p, matrix_w
216 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ec_diag_solver_kp'
218 INTEGER :: handle, i, ispin, nsize
221 TYPE(
cp_fm_type),
DIMENSION(:),
POINTER :: fmwork
226 CALL timeset(routinen, handle)
228 CALL dbcsr_get_info(matrix_ks(1, 1)%matrix, nfullrows_total=nsize)
230 NULLIFY (blacs_env, para_env)
231 CALL get_qs_env(qs_env=qs_env, blacs_env=blacs_env, para_env=para_env)
234 ncol_global=nsize, para_env=para_env)
242 CALL diag_kp_basic(matrix_ks, matrix_s, ec_env%kpoints, fmwork)
249 tempmat => matrix_s(1, 1)%matrix
251 ec_env%sab_orb, fmwork)
253 ec_env%sab_orb, fmwork)
256 mos => ec_env%kpoints%kp_env(1)%kpoint_env%mos
257 DO ispin = 1,
SIZE(mos, 2)
258 ec_env%ekTS = ec_env%ekTS + mos(1, ispin)%kTS
263 CALL timestop(handle)
287 POINTER :: matrix_ks, matrix_s
289 INTENT(INOUT),
POINTER :: matrix_p, matrix_w
291 CHARACTER(len=*),
PARAMETER :: routinen =
'ec_ot_diag_solver'
293 INTEGER :: handle, homo, ikind, iounit, ispin, &
294 max_iter, nao, nkind, nmo, nspins
295 INTEGER,
DIMENSION(2) :: nelectron_spin
296 REAL(kind=
dp),
DIMENSION(:),
POINTER :: eigenvalues
309 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
313 CALL timeset(routinen, handle)
317 cpassert(
ASSOCIATED(qs_env))
318 cpassert(
ASSOCIATED(ec_env))
319 cpassert(
ASSOCIATED(matrix_ks))
320 cpassert(
ASSOCIATED(matrix_s))
321 cpassert(
ASSOCIATED(matrix_p))
322 cpassert(
ASSOCIATED(matrix_w))
325 atomic_kind_set=atomic_kind_set, &
326 blacs_env=blacs_env, &
327 dft_control=dft_control, &
328 nelectron_spin=nelectron_spin, &
330 particle_set=particle_set, &
331 qs_kind_set=qs_kind_set)
332 nspins = dft_control%nspins
339 IF (dft_control%qs_control%do_ls_scf)
THEN
350 CALL get_qs_env(qs_env, qs_kind_set=qs_kind_set, nkind=nkind)
355 IF (ec_env%basis ==
"HARRIS")
THEN
357 qs_kind => qs_kind_set(ikind)
359 CALL get_qs_kind(qs_kind=qs_kind, basis_set=basis_set, basis_type=
"ORB")
361 CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type=
"HARRIS")
363 IF (basis_set%name /= harris_basis%name)
THEN
364 cpabort(
"OT-Diag initial guess: Harris and ground state need to use the same basis")
370 cpabort(
"OT-Diag initial guess: not implemented for MAOs")
374 SELECT CASE (ec_env%ec_initial_guess)
377 p_rmpv => matrix_p(:, 1)
380 nspins, nelectron_spin, iounit, para_env)
386 p_rmpv => rho_ao(:, 1)
389 cpabort(
"Unknown inital guess for OT-Diagonalization (Harris functional)")
415 NULLIFY (local_preconditioner)
416 ALLOCATE (local_preconditioner)
423 matrix_h=matrix_ks(ispin, 1)%matrix, &
424 matrix_s=matrix_s(ispin, 1)%matrix, &
425 convert_precond_to_dbcsr=.true., &
426 mo_set=mos(ispin), energy_gap=0.2_dp)
430 eigenvalues=eigenvalues, &
434 matrix_s=matrix_s(1, 1)%matrix, &
435 matrix_c_fm=mo_coeff, &
437 eps_gradient=ec_env%eps_default, &
441 evals_arg=eigenvalues, do_rotation=.true.)
445 DEALLOCATE (local_preconditioner)
449 mos(ispin)%mo_coeff_b)
460 IF (dft_control%qs_control%do_ls_scf)
THEN
464 IF (
ASSOCIATED(qs_env%mos))
THEN
465 DO ispin = 1,
SIZE(qs_env%mos)
468 DEALLOCATE (qs_env%mos)
472 CALL timestop(handle)
489 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_ks, matrix_s
491 CHARACTER(len=*),
PARAMETER :: routinen =
'ec_ls_init'
493 INTEGER :: handle, ispin, nspins
498 CALL timeset(routinen, handle)
501 dft_control=dft_control, &
503 nspins = dft_control%nspins
504 ls_env => ec_env%ls_env
507 CALL matrix_ls_create(matrix_ls=ls_env%matrix_s, matrix_qs=matrix_s(1, 1)%matrix, &
508 ls_mstruct=ls_env%ls_mstruct)
510 IF (
ALLOCATED(ls_env%matrix_p))
THEN
511 DO ispin = 1,
SIZE(ls_env%matrix_p)
515 ALLOCATE (ls_env%matrix_p(nspins))
519 CALL dbcsr_create(ls_env%matrix_p(ispin), template=ls_env%matrix_s, &
520 matrix_type=dbcsr_type_no_symmetry)
523 ALLOCATE (ls_env%matrix_ks(nspins))
525 CALL dbcsr_create(ls_env%matrix_ks(ispin), template=ls_env%matrix_s, &
526 matrix_type=dbcsr_type_no_symmetry)
536 matrix_qs=matrix_ks(ispin, 1)%matrix, &
537 ls_mstruct=ls_env%ls_mstruct, &
541 CALL timestop(handle)
561 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_p, matrix_w
562 INTEGER,
INTENT(IN) :: ec_ls_method
564 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ec_ls_solver'
566 INTEGER :: handle, ispin, nelectron_spin_real, &
568 INTEGER,
DIMENSION(2) :: nmo
570 TYPE(
dbcsr_type),
ALLOCATABLE,
DIMENSION(:) :: matrix_ks_deviation
576 CALL timeset(routinen, handle)
580 dft_control=dft_control, &
582 nspins = dft_control%nspins
583 ec_env => qs_env%ec_env
584 ls_env => ec_env%ls_env
587 CALL get_qs_env(qs_env=qs_env, nelectron_spin=nmo)
588 IF (nspins == 1) nmo(1) = nmo(1)/2
589 ls_env%homo_spin(:) = 0.0_dp
590 ls_env%lumo_spin(:) = 0.0_dp
592 ALLOCATE (matrix_ks_deviation(nspins))
594 CALL dbcsr_create(matrix_ks_deviation(ispin), template=ls_env%matrix_ks(ispin))
595 CALL dbcsr_set(matrix_ks_deviation(ispin), 0.0_dp)
599 IF (ls_env%has_s_preconditioner)
THEN
602 ls_env%matrix_bs_sqrt, ls_env%matrix_bs_sqrt_inv)
604 CALL dbcsr_filter(ls_env%matrix_ks(ispin), ls_env%eps_filter)
609 nelectron_spin_real = ls_env%nelectron_spin(ispin)
610 IF (ls_env%nspins == 1) nelectron_spin_real = nelectron_spin_real/2
612 SELECT CASE (ec_ls_method)
615 ls_env%mu_spin(ispin), &
617 ls_env%sign_method, &
619 ls_env%matrix_ks(ispin), &
621 ls_env%matrix_s_inv, &
622 nelectron_spin_real, &
627 ls_env%matrix_p(ispin), &
628 ls_env%matrix_ks(ispin), &
629 ls_env%matrix_s_sqrt_inv, &
630 nelectron_spin_real, &
631 ec_env%eps_default, &
632 ls_env%homo_spin(ispin), &
633 ls_env%lumo_spin(ispin), &
634 ls_env%mu_spin(ispin), &
635 matrix_ks_deviation=matrix_ks_deviation(ispin), &
636 dynamic_threshold=ls_env%dynamic_threshold, &
637 eps_lanczos=ls_env%eps_lanczos, &
638 max_iter_lanczos=ls_env%max_iter_lanczos)
642 ls_env%matrix_p(ispin), &
643 ls_env%matrix_ks(ispin), &
644 ls_env%matrix_s_sqrt_inv, &
645 nelectron_spin_real, &
646 ec_env%eps_default, &
647 ls_env%homo_spin(ispin), &
648 ls_env%lumo_spin(ispin), &
649 non_monotonic=ls_env%non_monotonic, &
650 eps_lanczos=ls_env%eps_lanczos, &
651 max_iter_lanczos=ls_env%max_iter_lanczos)
658 IF (ls_env%has_s_preconditioner)
THEN
662 ls_env%matrix_bs_sqrt, ls_env%matrix_bs_sqrt_inv)
664 CALL dbcsr_filter(ls_env%matrix_p(ispin), ls_env%eps_filter)
669 IF (nspins == 1)
CALL dbcsr_scale(ls_env%matrix_p(1), 2.0_dp)
677 matrix_ls=ls_env%matrix_p(ispin), &
678 ls_mstruct=ls_env%ls_mstruct, &
682 wmat => matrix_w(:, 1)
687 IF (ls_env%has_s_preconditioner)
THEN
691 IF (ls_env%needs_s_inv)
THEN
694 IF (ls_env%use_s_sqrt)
THEN
699 DO ispin = 1,
SIZE(ls_env%matrix_ks)
702 DEALLOCATE (ls_env%matrix_ks)
707 DEALLOCATE (matrix_ks_deviation)
709 CALL timestop(handle)
Define the atomic kind types and their sub types.
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_scale(matrix, alpha_scalar)
...
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_set(matrix, alpha)
...
subroutine, public dbcsr_release(matrix)
...
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 copy_fm_to_dbcsr(fm, matrix, keep_sparsity)
Copy a BLACS matrix to a dbcsr matrix.
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
subroutine, public cp_fm_geeig(amatrix, bmatrix, eigenvectors, eigenvalues, work)
General Eigenvalue Problem AX = BXE Single option version: Cholesky decomposition of B.
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_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
subroutine, public cp_fm_init_random(matrix, ncol, start_col)
fills a matrix with random numbers
various routines to log and control the output. The idea is that decisions about where to log should ...
recursive integer function, public cp_logger_get_default_unit_nr(logger, local, skip_not_ionode)
asks the default unit number of the given logger. try to use cp_logger_get_unit_nr
lower level routines for linear scaling SCF
subroutine, public density_matrix_trs4(matrix_p, matrix_ks, matrix_s_sqrt_inv, nelectron, threshold, e_homo, e_lumo, e_mu, dynamic_threshold, matrix_ks_deviation, max_iter_lanczos, eps_lanczos, converged, iounit)
compute the density matrix using a trace-resetting algorithm
subroutine, public ls_scf_init_matrix_s(matrix_s, ls_scf_env)
initialize S matrix related properties (sqrt, inverse...) Might be factored-out since this seems comm...
subroutine, public apply_matrix_preconditioner(matrix, direction, matrix_bs_sqrt, matrix_bs_sqrt_inv)
apply a preconditioner either forward (precondition) inv(sqrt(bs)) * A * inv(sqrt(bs)) backward (rest...
subroutine, public density_matrix_sign(matrix_p, mu, fixed_mu, sign_method, sign_order, matrix_ks, matrix_s, matrix_s_inv, nelectron, threshold, sign_symmetric, submatrix_sign_method, matrix_s_sqrt_inv)
compute the density matrix with a trace that is close to nelectron. take a mu as input,...
subroutine, public density_matrix_tc2(matrix_p, matrix_ks, matrix_s_sqrt_inv, nelectron, threshold, e_homo, e_lumo, non_monotonic, eps_lanczos, max_iter_lanczos, iounit)
compute the density matrix using a non monotonic trace conserving algorithm based on SIAM DOI....
Routines for a linear scaling quickstep SCF run based on the density matrix, with a focus on the inte...
subroutine, public matrix_ls_to_qs(matrix_qs, matrix_ls, ls_mstruct, covariant, keep_sparsity)
second link to QS, copy a LS matrix to QS matrix used to isolate QS style matrices from LS style will...
subroutine, public matrix_ls_create(matrix_ls, matrix_qs, ls_mstruct)
create a matrix for use (and as a template) in ls based on a qs template
subroutine, public matrix_qs_to_ls(matrix_ls, matrix_qs, ls_mstruct, covariant)
first link to QS, copy a QS matrix to LS matrix used to isolate QS style matrices from LS style will ...
Types needed for a linear scaling quickstep SCF run based on the density matrix.
Routines for a linear scaling quickstep SCF run based on the density matrix.
subroutine, public post_scf_sparsities(ls_scf_env)
Report on the sparsities of various interesting matrices.
subroutine, public calculate_w_matrix_ls(matrix_w, ls_scf_env)
Compute matrix_w as needed for the forces.
Routines for an energy correction on top of a Kohn-Sham calculation.
subroutine, public ec_diag_solver_gamma(qs_env, ec_env, matrix_ks, matrix_s, matrix_p, matrix_w)
Solve KS equation using diagonalization.
subroutine, public ec_ot_diag_solver(qs_env, ec_env, matrix_ks, matrix_s, matrix_p, matrix_w)
Use OT-diagonalziation to obtain density matrix from Harris Kohn-Sham matrix Initial guess of density...
subroutine, public ec_diag_solver_kp(qs_env, ec_env, matrix_ks, matrix_s, matrix_p, matrix_w)
Solve Kpoint-KS equation using diagonalization.
subroutine, public ec_ls_init(qs_env, matrix_ks, matrix_s)
Solve the Harris functional by linear scaling density purification scheme, instead of the diagonaliza...
Types needed for a for a Energy Correction.
Routines used for Harris functional Kohn-Sham calculation.
subroutine, public ec_mos_init(qs_env, matrix_s)
Allocate and initiate molecular orbitals environment.
Defines the basic variable types.
integer, parameter, public dp
Routines needed for kpoint calculation.
subroutine, public kpoint_density_transform(kpoint, denmat, wtype, tempmat, sab_nl, fmwork, for_aux_fit, pmat_ext)
generate real space density matrices in DBCSR format
subroutine, public kpoint_density_matrices(kpoint, energy_weighted, for_aux_fit)
Calculate kpoint density matrices (rho(k), owned by kpoint groups)
subroutine, public kpoint_set_mo_occupation(kpoint, smear, probe)
Given the eigenvalues of all kpoints, calculates the occupation numbers.
Interface to the message passing library MPI.
Define the data structure for the particle information.
subroutine, public init_preconditioner(preconditioner_env, para_env, blacs_env)
...
subroutine, public destroy_preconditioner(preconditioner_env)
...
computes preconditioners, and implements methods to apply them currently used in qs_ot
subroutine, public make_preconditioner(preconditioner_env, precon_type, solver_type, matrix_h, matrix_s, matrix_t, mo_set, energy_gap, convert_precond_to_dbcsr, chol_type)
...
Routine to return block diagonal density matrix. Blocks correspond to the atomic densities.
subroutine, public calculate_atomic_block_dm(pmatrix, matrix_s, atomic_kind_set, qs_kind_set, nspin, nelectron_spin, ounit, para_env)
returns a block diagonal density matrix. Blocks correspond to the atomic densities.
collects routines that calculate density matrices
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.
collects routines that perform operations directly related to MOs
subroutine, public make_basis_sm(vmatrix, ncol, matrix_s)
returns an S-orthonormal basis v (v^T S v ==1)
Set occupation of molecular orbitals.
Definition and initialisation of the mo data type.
subroutine, public allocate_mo_set(mo_set, nao, nmo, nelectron, n_el_f, maxocc, flexible_electron_count)
Allocates a mo set and partially initializes it (nao,nmo,nelectron, and flexible_electron_count are v...
subroutine, public deallocate_mo_set(mo_set)
Deallocate a wavefunction data structure.
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.
subroutine, public init_mo_set(mo_set, fm_pool, fm_ref, fm_struct, name)
initializes an allocated mo_set. eigenvalues, mo_coeff, occupation_numbers are valid only after this ...
an eigen-space solver for the generalised symmetric eigenvalue problem for sparse matrices,...
subroutine, public ot_eigensolver(matrix_h, matrix_s, matrix_orthogonal_space_fm, matrix_c_fm, preconditioner, eps_gradient, iter_max, size_ortho_space, silent, ot_settings)
...
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_basic(matrix_ks, matrix_s, kpoints, fmwork)
Kpoint diagonalization routine Transforms matrices to kpoint, distributes kpoint groups,...
Utilities for string manipulations.
elemental subroutine, public uppercase(string)
Convert all lower case characters in a string to upper case.
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 information on the energy correction functional for KG.
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.