30 ls_cluster_molecular,&
76#include "./base/base_uses.f90"
82 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'dm_ls_scf_qs'
104 CHARACTER(len=*),
PARAMETER :: routinen =
'matrix_ls_create'
106 CHARACTER(len=default_string_length) :: name
107 INTEGER :: handle, iatom, imol, jatom, natom, nmol
108 INTEGER,
ALLOCATABLE,
DIMENSION(:),
TARGET :: atom_to_cluster, atom_to_cluster_primus, &
109 clustered_blk_sizes, primus_of_mol
110 INTEGER,
DIMENSION(:),
POINTER :: clustered_col_dist, clustered_row_dist, &
111 ls_blk_sizes, ls_col_dist, ls_row_dist
114 CALL timeset(routinen, handle)
117 CALL dbcsr_get_info(matrix_qs, col_blk_size=ls_blk_sizes, distribution=ls_dist)
122 IF (ls_mstruct%do_pao)
THEN
123 CALL dbcsr_get_info(ls_mstruct%matrix_A, col_blk_size=ls_blk_sizes)
127 SELECT CASE (ls_mstruct%cluster_type)
133 nmol = maxval(ls_mstruct%atom_to_molecule)
134 ALLOCATE (atom_to_cluster_primus(natom))
135 ALLOCATE (atom_to_cluster(natom))
136 ALLOCATE (primus_of_mol(nmol))
138 atom_to_cluster(iatom) = ls_mstruct%atom_to_molecule(iatom)
142 DO jatom = iatom, 1, -1
143 IF (ls_mstruct%atom_to_molecule(jatom) == atom_to_cluster(iatom))
THEN
144 atom_to_cluster_primus(iatom) = jatom
149 primus_of_mol(atom_to_cluster(iatom)) = atom_to_cluster_primus(iatom)
153 ALLOCATE (clustered_row_dist(nmol))
155 clustered_row_dist(imol) = ls_row_dist(primus_of_mol(imol))
159 ALLOCATE (clustered_col_dist(nmol))
161 clustered_col_dist(imol) = ls_col_dist(primus_of_mol(imol))
164 ALLOCATE (clustered_blk_sizes(nmol))
165 clustered_blk_sizes = 0
167 clustered_blk_sizes(atom_to_cluster(iatom)) = clustered_blk_sizes(atom_to_cluster(iatom)) + &
170 ls_blk_sizes => clustered_blk_sizes
175 row_dist=clustered_row_dist, &
176 col_dist=clustered_col_dist, &
179 ls_dist = ls_dist_clustered
182 cpabort(
"Unknown LS cluster type")
191 row_blk_size=ls_blk_sizes, &
192 col_blk_size=ls_blk_sizes)
196 CALL timestop(handle)
216 LOGICAL,
INTENT(IN) :: covariant
218 CHARACTER(len=*),
PARAMETER :: routinen =
'matrix_qs_to_ls'
221 INTEGER,
DIMENSION(:),
POINTER :: pao_blk_sizes
225 CALL timeset(routinen, handle)
227 IF (.NOT. ls_mstruct%do_pao)
THEN
228 CALL matrix_cluster(matrix_ls, matrix_qs, ls_mstruct)
231 CALL dbcsr_get_info(ls_mstruct%matrix_A, col_blk_size=pao_blk_sizes)
234 template=matrix_qs, &
235 row_blk_size=pao_blk_sizes, &
236 col_blk_size=pao_blk_sizes)
238 matrix_trafo => ls_mstruct%matrix_A
239 IF (covariant) matrix_trafo => ls_mstruct%matrix_B
242 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_qs, matrix_trafo, 0.0_dp, matrix_tmp)
243 CALL dbcsr_multiply(
"T",
"N", 1.0_dp, matrix_trafo, matrix_tmp, 0.0_dp, matrix_pao)
246 CALL matrix_cluster(matrix_ls, matrix_pao, ls_mstruct)
250 CALL timestop(handle)
261 SUBROUTINE matrix_cluster(matrix_out, matrix_in, ls_mstruct)
265 CHARACTER(len=*),
PARAMETER :: routinen =
'matrix_cluster'
270 CALL timeset(routinen, handle)
272 SELECT CASE (ls_mstruct%cluster_type)
278 CALL dbcsr_create(matrix_in_nosym, template=matrix_in, matrix_type=
"N")
286 cpabort(
"Unknown LS cluster type")
289 CALL timestop(handle)
291 END SUBROUTINE matrix_cluster
307 SUBROUTINE matrix_ls_to_qs(matrix_qs, matrix_ls, ls_mstruct, covariant, keep_sparsity)
311 LOGICAL,
OPTIONAL :: keep_sparsity
313 CHARACTER(len=*),
PARAMETER :: routinen =
'matrix_ls_to_qs'
316 INTEGER,
DIMENSION(:),
POINTER :: pao_blk_sizes
317 LOGICAL :: my_keep_sparsity
318 TYPE(
dbcsr_type) :: matrix_declustered, matrix_tmp1, &
322 CALL timeset(routinen, handle)
324 my_keep_sparsity = .true.
325 IF (
PRESENT(keep_sparsity)) &
326 my_keep_sparsity = keep_sparsity
328 IF (.NOT. ls_mstruct%do_pao)
THEN
329 CALL dbcsr_create(matrix_declustered, template=matrix_qs)
331 CALL dbcsr_copy(matrix_qs, matrix_declustered, keep_sparsity=my_keep_sparsity)
335 CALL dbcsr_get_info(ls_mstruct%matrix_A, col_blk_size=pao_blk_sizes)
337 template=matrix_qs, &
338 row_blk_size=pao_blk_sizes, &
339 col_blk_size=pao_blk_sizes)
343 matrix_trafo => ls_mstruct%matrix_B
344 IF (covariant) matrix_trafo => ls_mstruct%matrix_A
347 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_trafo, matrix_declustered, 0.0_dp, matrix_tmp1)
348 CALL dbcsr_multiply(
"N",
"T", 1.0_dp, matrix_tmp1, matrix_trafo, 0.0_dp, matrix_tmp2)
349 CALL dbcsr_copy(matrix_qs, matrix_tmp2, keep_sparsity=my_keep_sparsity)
355 CALL timestop(handle)
370 CHARACTER(len=*),
PARAMETER :: routinen =
'matrix_decluster'
374 CALL timeset(routinen, handle)
376 SELECT CASE (ls_mstruct%cluster_type)
385 cpabort(
"Unknown LS cluster type")
388 CALL timestop(handle)
403 CHARACTER(len=*),
PARAMETER :: routinen =
'ls_scf_init_qs'
405 INTEGER :: handle, ispin, nspin, unit_nr
407 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_ks, matrix_s
414 CALL timeset(routinen, handle)
418 IF (logger%para_env%is_source())
THEN
425 CALL get_qs_env(qs_env, dft_control=dft_control, &
427 matrix_ks=matrix_ks, &
431 nspin = dft_control%nspins
434 IF (.NOT.
ASSOCIATED(matrix_ks))
THEN
437 ALLOCATE (matrix_ks(ispin)%matrix)
438 CALL dbcsr_create(matrix_ks(ispin)%matrix, template=matrix_s(1)%matrix)
440 CALL dbcsr_set(matrix_ks(ispin)%matrix, 0.0_dp)
445 CALL timestop(handle)
462 REAL(kind=
dp) :: energy
463 LOGICAL,
INTENT(IN),
OPTIONAL :: nonscf
465 CHARACTER(len=*),
PARAMETER :: routinen =
'ls_scf_qs_atomic_guess'
467 INTEGER :: handle, nspin, unit_nr
468 INTEGER,
DIMENSION(2) :: nelectron_spin
469 LOGICAL :: do_scf, has_unit_metric
472 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_ks, matrix_s, rho_ao
477 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
481 CALL timeset(routinen, handle)
482 NULLIFY (rho, rho_ao)
486 IF (logger%para_env%is_source())
THEN
493 CALL get_qs_env(qs_env, dft_control=dft_control, &
495 matrix_ks=matrix_ks, &
498 atomic_kind_set=atomic_kind_set, &
499 qs_kind_set=qs_kind_set, &
500 particle_set=particle_set, &
501 has_unit_metric=has_unit_metric, &
503 nelectron_spin=nelectron_spin, &
508 nspin = dft_control%nspins
511 IF (dft_control%qs_control%dftb .OR. dft_control%qs_control%semi_empirical .OR. &
512 dft_control%qs_control%xtb)
THEN
514 dft_control, particle_set, atomic_kind_set, qs_kind_set, &
515 nspin, nelectron_spin, para_env)
518 nspin, nelectron_spin, unit_nr, para_env)
522 IF (
PRESENT(nonscf)) do_scf = .NOT. nonscf
527 CALL ls_scf_tblite_energy(qs_env,
qs_energy)
533 CALL timestop(handle)
551 REAL(kind=
dp) :: energy_new
552 INTEGER,
INTENT(IN) :: iscf
554 CHARACTER(len=*),
PARAMETER :: routinen =
'ls_scf_dm_to_ks'
556 INTEGER :: handle, ispin, nspin, unit_nr
563 NULLIFY (energy, rho, rho_ao)
564 CALL timeset(routinen, handle)
567 IF (logger%para_env%is_source())
THEN
573 nspin = ls_scf_env%nspins
574 CALL get_qs_env(qs_env, para_env=para_env, energy=energy, rho=rho)
579 CALL matrix_ls_to_qs(rho_ao(ispin)%matrix, ls_scf_env%matrix_p(ispin), &
580 ls_scf_env%ls_mstruct, covariant=.false.)
585 IF (ls_scf_env%do_rho_mixing)
THEN
587 cpabort(
"Direct P mixing not implemented in linear scaling SCF. ")
589 IF (iscf > max(ls_scf_env%mixing_store%nskip_mixing, 1))
THEN
590 CALL gspace_mixing(qs_env, ls_scf_env%density_mixing_method, &
591 ls_scf_env%mixing_store, rho, para_env, &
593 IF (unit_nr > 0)
THEN
594 WRITE (unit_nr,
'(A57)') &
595 "*********************************************************"
596 WRITE (unit_nr,
'(A13,F5.3,A20,A6,A7,I3)') &
597 " Using ALPHA=", ls_scf_env%mixing_store%alpha, &
598 " to mix rho: method=", ls_scf_env%mixing_store%iter_method,
", iscf=", iscf
599 WRITE (unit_nr,
'(A8,F5.3,A6,F5.3,A8)') &
600 " rho_nw=", ls_scf_env%mixing_store%alpha,
"*rho + ", &
601 1.0_dp - ls_scf_env%mixing_store%alpha,
"*rho_old"
602 WRITE (unit_nr,
'(A57)') &
603 "*********************************************************"
611 just_energy=.false., print_active=.true.)
612 CALL ls_scf_tblite_energy(qs_env, energy)
613 energy_new = energy%total
615 CALL timestop(handle)
628 REAL(kind=
dp) :: energy_new
630 CHARACTER(len=*),
PARAMETER :: routinen =
'ls_nonscf_ks'
632 INTEGER :: handle, ispin, nspin
639 NULLIFY (energy, rho, rho_ao)
640 CALL timeset(routinen, handle)
642 nspin = ls_scf_env%nspins
643 CALL get_qs_env(qs_env, para_env=para_env, energy=energy, rho=rho)
648 CALL matrix_ls_to_qs(rho_ao(ispin)%matrix, ls_scf_env%matrix_p(ispin), &
649 ls_scf_env%ls_mstruct, covariant=.false.)
652 IF (qs_env%harris_method)
THEN
653 CALL get_qs_env(qs_env, harris_env=harris_env)
658 IF (ls_scf_env%do_rho_mixing)
THEN
659 cpabort(
"P mixing not implemented in linear scaling NONSCF. ")
664 just_energy=.false., print_active=.true.)
665 CALL ls_scf_tblite_energy(qs_env, energy)
666 energy_new = energy%total
668 CALL timestop(handle)
681 CHARACTER(len=*),
PARAMETER :: routinen =
'ls_nonscf_energy'
683 INTEGER :: handle, ispin, nspin
684 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_h, matrix_ks, rho_ao
689 NULLIFY (energy, rho, rho_ao)
690 CALL timeset(routinen, handle)
691 IF (qs_env%qmmm)
THEN
695 nspin = ls_scf_env%nspins
696 CALL get_qs_env(qs_env, para_env=para_env, energy=energy, rho=rho)
701 CALL matrix_ls_to_qs(rho_ao(ispin)%matrix, ls_scf_env%matrix_p(ispin), &
702 ls_scf_env%ls_mstruct, covariant=.false.)
711 energy%total = energy%total - energy%core
715 CALL timestop(handle)
724 SUBROUTINE ls_scf_tblite_energy(qs_env, energy)
731 NULLIFY (dft_control, tb)
732 CALL get_qs_env(qs_env, dft_control=dft_control, tb_tblite=tb)
734 IF (dft_control%qs_control%xtb .AND. dft_control%qs_control%xtb_control%do_tblite)
THEN
735 cpassert(
ASSOCIATED(tb))
739 END SUBROUTINE ls_scf_tblite_energy
754 INTEGER,
INTENT(IN) :: unit_nr
755 CHARACTER(LEN=*),
INTENT(IN) :: title
756 INTEGER,
DIMENSION(:),
POINTER :: stride
758 CHARACTER(len=*),
PARAMETER :: routinen =
'write_matrix_to_cube'
772 CALL timeset(routinen, handle)
774 NULLIFY (ks_env, pw_env, auxbas_pw_pool, pw_pools, particles, subsys, matrix_ks)
785 CALL dbcsr_copy(matrix_p_qs, matrix_ks(1)%matrix)
787 CALL matrix_ls_to_qs(matrix_p_qs, matrix_p_ls, ls_scf_env%ls_mstruct, covariant=.false.)
791 auxbas_pw_pool=auxbas_pw_pool, &
793 CALL auxbas_pw_pool%create_pw(pw=wf_r)
795 CALL auxbas_pw_pool%create_pw(pw=wf_g)
804 particles=particles, stride=stride)
807 CALL auxbas_pw_pool%give_back_pw(wf_r)
808 CALL auxbas_pw_pool%give_back_pw(wf_g)
811 CALL timestop(handle)
824 CHARACTER(len=*),
PARAMETER :: routinen =
'rho_mixing_ls_init'
831 CALL timeset(routinen, handle)
833 CALL get_qs_env(qs_env, dft_control=dft_control, rho=rho)
835 CALL mixing_allocate(qs_env, ls_scf_env%density_mixing_method, nspins=ls_scf_env%nspins, &
836 mixing_store=ls_scf_env%mixing_store)
838 IF (dft_control%qs_control%gapw)
THEN
839 CALL get_qs_env(qs_env, rho_atom_set=rho_atom)
840 CALL mixing_init(ls_scf_env%density_mixing_method, rho, ls_scf_env%mixing_store, &
841 ls_scf_env%para_env, rho_atom=rho_atom)
842 ELSEIF (dft_control%qs_control%dftb .OR. dft_control%qs_control%xtb)
THEN
844 ELSEIF (dft_control%qs_control%semi_empirical)
THEN
845 cpabort(
'SE Code not possible')
847 CALL mixing_init(ls_scf_env%density_mixing_method, rho, ls_scf_env%mixing_store, &
851 CALL timestop(handle)
Define the atomic kind types and their sub types.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_distribution_release(dist)
...
subroutine, public dbcsr_distribution_new(dist, template, group, pgrid, row_dist, col_dist, reuse_arrays)
...
subroutine, public dbcsr_distribution_hold(dist)
...
subroutine, public dbcsr_desymmetrize(matrix_a, matrix_b)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_multiply(transa, transb, alpha, matrix_a, matrix_b, beta, matrix_c, first_row, last_row, first_column, last_column, first_k, last_k, retain_sparsity, filter_eps, flop)
...
subroutine, public dbcsr_get_info(matrix, nblkrows_total, nblkcols_total, nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, nfullrows_local, nfullcols_local, my_prow, my_pcol, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, matrix_type, group)
...
subroutine, public dbcsr_finalize(matrix)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_release(matrix)
...
subroutine, public dbcsr_complete_redistribute(matrix, redist)
...
subroutine, public dbcsr_distribution_get(dist, row_dist, col_dist, nrows, ncols, has_threads, group, mynode, numnodes, nprows, npcols, myprow, mypcol, pgrid, subgroups_defined, prow_group, pcol_group)
...
Routines that link DBCSR and CP2K concepts together.
subroutine, public cp_dbcsr_alloc_block_from_nbl(matrix, sab_orb, desymmetrize)
allocate the blocks of a dbcsr based on the neighbor list
DBCSR operations in CP2K.
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
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
A wrapper around pw_to_cube() which accepts particle_list_type.
subroutine, public cp_pw_to_cube(pw, unit_nr, title, particles, zeff, stride, max_file_size_mb, zero_tails, silent, mpi_io)
...
Routines for a linear scaling quickstep SCF run based on the density matrix, with a focus on the inte...
subroutine, public ls_scf_qs_atomic_guess(qs_env, ls_scf_env, energy, nonscf)
get an atomic initial guess
subroutine, public ls_nonscf_ks(qs_env, ls_scf_env, energy_new)
use the external density in ls_scf_env to compute the new KS matrix
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_decluster(matrix_out, matrix_in, ls_mstruct)
Reverses molecular blocking and reduction to single precision if enabled.
subroutine, public ls_scf_dm_to_ks(qs_env, ls_scf_env, energy_new, iscf)
use the density matrix in ls_scf_env to compute the new energy and KS matrix
subroutine, public write_matrix_to_cube(qs_env, ls_scf_env, matrix_p_ls, unit_nr, title, stride)
...
subroutine, public rho_mixing_ls_init(qs_env, ls_scf_env)
Initialize g-space density mixing.
subroutine, public ls_nonscf_energy(qs_env, ls_scf_env)
use the new density matrix in ls_scf_env to compute the new energy
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 ...
subroutine, public ls_scf_init_qs(qs_env)
further required initialization of QS. Might be factored-out since this seems common code with the ot...
Types needed for a linear scaling quickstep SCF run based on the density matrix.
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
Interface to the message passing library MPI.
represent a simple array based list of the given type
Define the data structure for the particle information.
container for various plainwaves related things
subroutine, public pw_env_get(pw_env, pw_pools, cube_info, gridlevel_info, auxbas_pw_pool, auxbas_grid, auxbas_rs_desc, auxbas_rs_grid, rs_descs, rs_grids, xc_pw_pool, vdw_pw_pool, poisson_env, interp_section)
returns the various attributes of the pw env
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
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.
Calculate the plane wave density by collocating the primitive Gaussian functions (pgf).
subroutine, public calculate_rho_elec(matrix_p, matrix_p_kp, rho, rho_gspace, total_rho, ks_env, soft_valid, compute_tau, compute_grad, basis_type, der_type, idir, task_list_external, pw_env_external)
computes the density corresponding to a given density matrix on the grid
Calculation of the energies concerning the core charge distribution.
module that contains the definitions of the scf types
integer, parameter, public direct_mixing_nr
integer, parameter, public gspace_mixing_nr
Perform a QUICKSTEP wavefunction optimization (single point)
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.
subroutine, public gspace_mixing(qs_env, mixing_method, mixing_store, rho, para_env, iter_count)
Driver for the g-space mixing, calls the proper routine given the requested method.
Types needed for a for a Harris model calculation.
Harris method environment setup and handling.
subroutine, public harris_density_update(qs_env, harris_env)
...
Routines to somehow generate an initial guess.
subroutine, public calculate_mopac_dm(pmat, matrix_s, has_unit_metric, dft_control, particle_set, atomic_kind_set, qs_kind_set, nspin, nelectron_spin, para_env)
returns a block diagonal density matrix. Blocks correspond to the mopac initial guess.
Define the quickstep kind type and their sub types.
routines that build the Kohn-Sham matrix (i.e calculate the coulomb and xc parts
subroutine, public qs_ks_update_qs_env(qs_env, calculate_forces, just_energy, print_active)
updates the Kohn Sham matrix of the given qs_env (facility method)
subroutine, public qs_ks_did_change(ks_env, s_mstruct_changed, rho_changed, potential_changed, full_reset)
tells that some of the things relevant to the ks calculation did change. has to be called when change...
subroutine, public set_ks_env(ks_env, v_hartree_rspace, s_mstruct_changed, rho_changed, exc_accint, potential_changed, forces_up_to_date, complex_ks, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, kinetic, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_ks_im_kp, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, vee, neighbor_list_id, kpoints, sab_orb, sab_all, sac_ae, sac_ppl, sac_lri, sap_ppnl, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_vdw, sab_scp, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, task_list, task_list_soft, subsys, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env)
...
elemental subroutine, public charge_mixing_init(mixing_store)
initialiation needed when charge mixing is used
subroutine, public mixing_init(mixing_method, rho, mixing_store, para_env, rho_atom)
initialiation needed when gspace mixing is used
subroutine, public mixing_allocate(qs_env, mixing_method, p_mix_new, p_delta, nspins, mixing_store)
allocation needed when density mixing is used
Define the neighbor list data types and the corresponding functionality.
methods of the rho structure (defined in qs_rho_types)
subroutine, public qs_rho_update_rho(rho_struct, qs_env, rho_xc_external, local_rho_set, task_list_external, task_list_external_soft, pw_env_external, para_env_external)
updates rho_r and rho_g to the rhorho_ao. if use_kinetic_energy_density also computes tau_r and tau_g...
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...
types that represent a quickstep subsys
subroutine, public qs_subsys_get(subsys, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell, cell_ref, use_ref_cell, energy, force, qs_kind_set, cp_subsys, nelectron_total, nelectron_spin)
...
subroutine, public tb_get_energy(qs_env, tb, energy)
...
Provides all information about an atomic kind.
type of a logger, at the moment it contains just a print level starting at which level it should be l...
stores all the informations relevant to an mpi environment
represent a list of objects
contained for different pw related things
to create arrays of pools
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Contains information on the Harris method.
Provides all information about a quickstep kind.
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...
keeps the density in various representations, keeping track of which ones are valid.