71#include "./base/base_uses.f90"
77 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'ec_environment'
100 cpassert(.NOT.
ASSOCIATED(ec_env))
102 CALL init_ec_env(qs_env, ec_env, dft_section, ec_section)
116 SUBROUTINE init_ec_env(qs_env, ec_env, dft_section, ec_section)
122 CHARACTER(LEN=*),
PARAMETER :: routinen =
'init_ec_env'
124 INTEGER :: handle, ikind, maxlgto, nkind, unit_nr
125 LOGICAL :: explicit, gpw, gs_kpoints, paw_atom
126 REAL(kind=
dp) :: eps_pgf_orb, etemp, rc
136 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
140 pp_section, section1, section2, &
141 xc_fun_section, xc_section
143 CALL timeset(routinen, handle)
145 NULLIFY (atomic_kind_set, dispersion_env, ec_env%ls_env, para_env)
146 NULLIFY (ec_env%sab_orb, ec_env%sac_ae, ec_env%sac_ppl, ec_env%sap_ppnl)
147 NULLIFY (ec_env%matrix_ks, ec_env%matrix_h, ec_env%matrix_s)
148 NULLIFY (ec_env%matrix_t, ec_env%matrix_p, ec_env%matrix_w)
149 NULLIFY (ec_env%task_list)
150 NULLIFY (ec_env%mao_coef)
151 NULLIFY (ec_env%force)
152 NULLIFY (ec_env%dispersion_env)
153 NULLIFY (ec_env%xc_section)
154 NULLIFY (ec_env%matrix_z)
155 NULLIFY (ec_env%matrix_hz)
156 NULLIFY (ec_env%matrix_wz)
157 NULLIFY (ec_env%z_admm)
158 NULLIFY (ec_env%p_env)
159 NULLIFY (ec_env%vxc_rspace)
160 NULLIFY (ec_env%vtau_rspace)
161 NULLIFY (ec_env%vadmm_rspace)
162 NULLIFY (ec_env%rhoout_r, ec_env%rhoz_r)
163 NULLIFY (ec_env%x_data)
164 ec_env%should_update = .true.
166 ec_env%do_ec_admm = .false.
167 ec_env%do_kpoints = .false.
168 ec_env%do_ec_hfx = .false.
169 ec_env%reuse_hfx = .false.
171 IF (qs_env%energy_correction)
THEN
173 cpassert(
PRESENT(ec_section))
179 i_val=ec_env%ks_solver)
181 i_val=ec_env%energy_functional)
183 i_val=ec_env%factorization)
185 i_val=ec_env%ec_initial_guess)
187 r_val=ec_env%eps_default)
195 i_val=ec_env%mao_max_iter)
197 r_val=ec_env%mao_eps_grad)
199 r_val=ec_env%mao_eps1)
201 i_val=ec_env%mao_iolevel)
204 l_val=ec_env%skip_ec)
207 l_val=ec_env%debug_forces)
209 l_val=ec_env%debug_stress)
211 l_val=ec_env%debug_external)
216 c_val=ec_env%exresp_fn)
218 c_val=ec_env%exresperr_fn)
220 c_val=ec_env%exresult_fn)
222 l_val=ec_env%do_error)
224 c_val=ec_env%error_method)
227 r_val=ec_env%error_cutoff)
229 i_val=ec_env%error_subspace)
231 ec_env%do_skip = .false.
234 IF (etemp > 0.0_dp)
THEN
235 ec_env%smear%do_smear = .true.
237 ec_env%smear%electronic_temperature = etemp
238 ec_env%smear%eps_fermi_dirac = 1.0e-5_dp
239 ec_env%smear%fixed_mag_mom = -100.0_dp
241 ec_env%smear%do_smear = .false.
242 ec_env%smear%electronic_temperature = 0.0_dp
251 CALL get_qs_env(qs_env, do_kpoints=gs_kpoints)
254 ec_env%do_kpoints = gs_kpoints .OR. explicit
255 IF (ec_env%do_kpoints)
THEN
256 IF (.NOT. explicit)
THEN
259 CALL get_qs_env(qs_env, particle_set=particle_set, cell=cell)
263 CALL get_qs_env(qs_env, para_env=para_env, blacs_env=blacs_env)
266 NULLIFY (ec_env%kpoints)
270 CALL get_qs_env(qs_env, qs_kind_set=qs_kind_set, nkind=nkind)
272 SELECT CASE (ec_env%basis)
275 qs_kind => qs_kind_set(ikind)
276 CALL get_qs_kind(qs_kind=qs_kind, basis_set=basis_set, basis_type=
"ORB")
277 IF (
ASSOCIATED(basis_set))
THEN
278 NULLIFY (harris_basis)
279 CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type=
"HARRIS")
280 IF (
ASSOCIATED(harris_basis))
THEN
283 NULLIFY (harris_basis)
290 qs_kind => qs_kind_set(ikind)
291 CALL get_qs_kind(qs_kind=qs_kind, basis_set=basis_set, basis_type=
"ORB")
292 IF (
ASSOCIATED(basis_set))
THEN
293 NULLIFY (harris_basis)
294 CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type=
"HARRIS")
295 IF (
ASSOCIATED(harris_basis))
THEN
298 NULLIFY (harris_basis)
300 CALL get_qs_env(qs_env, dft_control=dft_control)
301 eps_pgf_orb = dft_control%qs_control%eps_pgf_orb
303 harris_basis%kind_radius = basis_set%kind_radius
309 qs_kind => qs_kind_set(ikind)
310 NULLIFY (harris_basis)
311 CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type=
"HARRIS")
312 IF (.NOT.
ASSOCIATED(harris_basis))
THEN
313 cpwarn(
"Harris Basis not defined for all types of atoms.")
317 cpabort(
"Unknown basis set for energy correction (Harris functional)")
320 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto, basis_type=
"HARRIS")
323 CALL get_qs_env(qs_env, dft_control=dft_control)
324 IF (dft_control%qs_control%gapw .OR. dft_control%qs_control%gapw_xc)
THEN
325 eps_pgf_orb = dft_control%qs_control%eps_pgf_orb
327 qs_kind => qs_kind_set(ikind)
328 NULLIFY (harris_basis)
329 CALL get_qs_kind(qs_kind, basis_set=harris_basis, basis_type=
"HARRIS")
330 CALL get_qs_kind(qs_kind, hard_radius=rc, gpw_type_forced=gpw)
331 NULLIFY (harris_soft_basis)
334 dft_control%qs_control%gapw_control%eps_fit, &
336 dft_control%qs_control%gapw_control%force_paw, gpw)
345 ec_env%basis_inconsistent = .false.
346 IF (ec_env%basis ==
"HARRIS")
THEN
348 qs_kind => qs_kind_set(ikind)
350 CALL get_qs_kind(qs_kind=qs_kind, basis_set=basis_set, basis_type=
"ORB")
352 CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type=
"HARRIS")
354 IF (basis_set%name /= harris_basis%name)
THEN
355 ec_env%basis_inconsistent = .true.
361 IF (ec_env%energy_functional ==
ec_functional_dc .AND. ec_env%basis_inconsistent)
THEN
362 CALL cp_abort(__location__, &
363 "DC-DFT: Correction and ground state need to use the same basis. "// &
364 "Checked by comparing basis set names only.")
366 IF (ec_env%energy_functional ==
ec_functional_ext .AND. ec_env%basis_inconsistent)
THEN
367 CALL cp_abort(__location__, &
368 "Exteranl Energy: Correction and ground state need to use the same basis. "// &
369 "Checked by comparing basis set names only.")
373 SELECT CASE (ec_env%energy_functional)
375 ec_env%ec_name =
"Harris"
377 ec_env%ec_name =
"DC-DFT"
379 ec_env%ec_name =
"External Energy"
381 cpabort(
"unknown energy correction")
391 ec_env%xc_section => section1
393 ec_env%xc_section => xc_section
396 CALL get_qs_env(qs_env, dft_control=dft_control, rho=rho)
398 dft_control%use_kinetic_energy_density = dft_control%use_kinetic_energy_density .OR. &
401 dft_control%drho_by_collocation = dft_control%drho_by_collocation .OR. &
405 ALLOCATE (dispersion_env)
407 xc_section => ec_env%xc_section
408 CALL get_qs_env(qs_env, atomic_kind_set=atomic_kind_set, para_env=para_env)
415 cpabort(
"nl-vdW functionals not available for EC calculations")
420 ec_env%dispersion_env => dispersion_env
427 ec_env%use_ls_solver = .false.
432 IF (ec_env%use_ls_solver)
THEN
433 CALL ec_ls_create(qs_env, ec_env)
439 cpabort(
"Harris functional with Fermi-Dirac smearing needs diagonalization solver.")
445 CALL timestop(handle)
447 END SUBROUTINE init_ec_env
458 SUBROUTINE ec_ls_create(qs_env, ec_env)
462 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ec_ls_create'
472 CALL timeset(routinen, handle)
474 ALLOCATE (ec_env%ls_env)
475 ls_env => ec_env%ls_env
477 NULLIFY (dft_control, input, ls_env%para_env)
480 dft_control=dft_control, &
482 molecule_set=molecule_set, &
483 particle_set=particle_set, &
484 para_env=ls_env%para_env, &
485 nelectron_spin=ls_env%nelectron_spin)
488 ls_env%nspins = dft_control%nspins
489 ls_env%natoms =
SIZE(particle_set, 1)
490 CALL ls_env%para_env%retain()
493 ALLOCATE (ls_env%ls_mstruct%atom_to_molecule(ls_env%natoms))
494 CALL molecule_of_atom(molecule_set, atom_to_mol=ls_env%ls_mstruct%atom_to_molecule)
496 ls_env%do_transport = .false.
497 ls_env%do_pao = .false.
498 ls_env%ls_mstruct%do_pao = ls_env%do_pao
499 ls_env%do_pexsi = .false.
500 ls_env%has_unit_metric = .false.
508 CALL section_vals_val_get(ec_section,
"MATRIX_CLUSTER_TYPE", i_val=ls_env%ls_mstruct%cluster_type)
511 CALL section_vals_val_get(ec_section,
"REPORT_ALL_SPARSITIES", l_val=ls_env%report_all_sparsities)
521 SELECT CASE (ec_env%ks_solver)
525 SELECT CASE (ls_env%s_inversion_type)
527 ls_env%needs_s_inv = .true.
528 ls_env%use_s_sqrt = .true.
530 ls_env%needs_s_inv = .true.
531 ls_env%use_s_sqrt = .false.
533 ls_env%needs_s_inv = .false.
534 ls_env%use_s_sqrt = .false.
539 ls_env%needs_s_inv = .false.
540 ls_env%use_s_sqrt = .true.
545 SELECT CASE (ls_env%s_preconditioner_type)
547 ls_env%has_s_preconditioner = .false.
549 ls_env%has_s_preconditioner = .true.
553 ls_env%extrapolation_order = 0
554 ls_env%scf_history%nstore = 0
555 ls_env%scf_history%istore = 0
556 ALLOCATE (ls_env%scf_history%matrix(ls_env%nspins, ls_env%scf_history%nstore))
558 NULLIFY (ls_env%mixing_store)
560 CALL timestop(handle)
562 END SUBROUTINE ec_ls_create
575 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ec_write_input'
577 INTEGER :: handle, unit_nr
580 CALL timeset(routinen, handle)
584 IF (unit_nr > 0)
THEN
586 WRITE (unit_nr,
'(T2,A)') &
587 "!"//repeat(
"-", 29)//
" Energy Correction "//repeat(
"-", 29)//
"!"
590 SELECT CASE (ec_env%energy_functional)
592 WRITE (unit_nr,
'(T2,A,T61,A20)')
"Energy Correction: ",
"HARRIS FUNCTIONAL"
594 WRITE (unit_nr,
'(T2,A,T61,A20)')
"Energy Correction: ",
"DC-DFT"
596 WRITE (unit_nr,
'(T2,A,T61,A20)')
"Energy Correction: ",
"External"
598 WRITE (unit_nr,
'()')
601 WRITE (unit_nr,
'(T2,A,T61,E20.3)')
"eps_default:", ec_env%eps_default
604 SELECT CASE (ec_env%basis)
606 WRITE (unit_nr,
'(T2,A,T61,A20)')
"EC basis: ",
"ORBITAL"
608 WRITE (unit_nr,
'(T2,A,T61,A20)')
"EC basis: ",
"PRIMITIVE"
610 WRITE (unit_nr,
'(T2,A,T61,A20)')
"EC Basis: ",
"HARRIS"
614 IF (ec_env%do_ec_hfx)
THEN
616 WRITE (unit_nr,
'(T2,A,T61,L20)')
"DC-DFT with HFX", ec_env%do_ec_hfx
617 WRITE (unit_nr,
'(T2,A,T61,L20)')
"Reuse HFX integrals", ec_env%reuse_hfx
618 WRITE (unit_nr,
'(T2,A,T61,L20)')
"DC-DFT HFX with ADMM", ec_env%do_ec_admm
622 IF (ec_env%do_kpoints)
THEN
630 SELECT CASE (ec_env%ks_solver)
632 WRITE (unit_nr,
'(T2,A,T61,A20)')
"Algorithm: ",
"DIAGONALIZATION"
634 WRITE (unit_nr,
'(T2,A,T61,A20)')
"Algorithm: ",
"OT DIAGONALIZATION"
636 WRITE (unit_nr,
'(T2,A,T61,A20)')
"Algorithm: ",
"MATRIX_SIGN"
638 WRITE (unit_nr,
'(T2,A,T61,A20)')
"Algorithm: ",
"TRS4"
641 WRITE (unit_nr,
'(T2,A,T61,A20)')
"Algorithm: ",
"TC2"
644 WRITE (unit_nr,
'()')
648 WRITE (unit_nr,
'(T2,A,T61,L20)')
"MAO:", ec_env%mao
649 WRITE (unit_nr,
'(T2,A,T61,L20)')
"MAO_IOLEVEL:", ec_env%mao_iolevel
650 WRITE (unit_nr,
'(T2,A,T61,I20)')
"MAO_MAX_ITER:", ec_env%mao_max_iter
651 WRITE (unit_nr,
'(T2,A,T61,E20.3)')
"MAO_EPS_GRAD:", ec_env%mao_eps_grad
652 WRITE (unit_nr,
'(T2,A,T61,E20.3)')
"MAO_EPS1:", ec_env%mao_eps1
653 WRITE (unit_nr,
'()')
657 IF (.NOT. ec_env%use_ls_solver)
THEN
659 WRITE (unit_nr,
'(T2,A)')
"MO Solver"
660 WRITE (unit_nr,
'()')
662 SELECT CASE (ec_env%ks_solver)
665 SELECT CASE (ec_env%factorization)
667 WRITE (unit_nr,
'(T2,A,T61,A20)')
"Factorization: ",
"CHOLESKY"
676 SELECT CASE (ec_env%ec_initial_guess)
678 WRITE (unit_nr,
'(T2,A,T61,A20)')
"OT Diag initial guess: ",
"ATOMIC"
680 WRITE (unit_nr,
'(T2,A,T61,A20)')
"OT Diag initial guess: ",
"GROUND STATE DM"
684 cpabort(
"Unknown Diagonalization algorithm for Harris functional")
689 WRITE (unit_nr,
'(T2,A)')
"AO Solver"
690 WRITE (unit_nr,
'()')
692 ls_env => ec_env%ls_env
693 WRITE (unit_nr,
'(T2,A,T61,E20.3)')
"eps_filter:", ls_env%eps_filter
694 WRITE (unit_nr,
'(T2,A,T61,L20)')
"fixed chemical potential (mu)", ls_env%fixed_mu
695 WRITE (unit_nr,
'(T2,A,T61,L20)')
"Computing inv(S):", ls_env%needs_s_inv
696 WRITE (unit_nr,
'(T2,A,T61,L20)')
"Computing sqrt(S):", ls_env%use_s_sqrt
697 WRITE (unit_nr,
'(T2,A,T61,L20)')
"Computing S preconditioner ", ls_env%has_s_preconditioner
699 IF (ls_env%use_s_sqrt)
THEN
700 SELECT CASE (ls_env%s_sqrt_method)
702 WRITE (unit_nr,
'(T2,A,T61,A20)')
"S sqrt method:",
"NEWTONSCHULZ"
704 WRITE (unit_nr,
'(T2,A,T61,A20)')
"S sqrt method:",
"PROOT"
706 cpabort(
"Unknown sqrt method.")
708 WRITE (unit_nr,
'(T2,A,T61,I20)')
"S sqrt order:", ls_env%s_sqrt_order
711 SELECT CASE (ls_env%s_preconditioner_type)
713 WRITE (unit_nr,
'(T2,A,T61,A20)')
"S preconditioner type ",
"NONE"
715 WRITE (unit_nr,
'(T2,A,T61,A20)')
"S preconditioner type ",
"ATOMIC"
717 WRITE (unit_nr,
'(T2,A,T61,A20)')
"S preconditioner type ",
"MOLECULAR"
720 SELECT CASE (ls_env%ls_mstruct%cluster_type)
722 WRITE (unit_nr,
'(T2,A,T61,A20)')
"Cluster type", adjustr(
"ATOMIC")
724 WRITE (unit_nr,
'(T2,A,T61,A20)')
"Cluster type", adjustr(
"MOLECULAR")
726 cpabort(
"Unknown cluster type")
733 WRITE (unit_nr,
'(T2,A)') repeat(
"-", 79)
734 WRITE (unit_nr,
'()')
738 CALL timestop(handle)
Define the atomic kind types and their sub types.
subroutine, public remove_basis_from_container(container, inum, basis_type)
...
subroutine, public add_basis_set_to_container(container, basis_set, basis_set_type)
...
subroutine, public allocate_gto_basis_set(gto_basis_set)
...
subroutine, public copy_gto_basis_set(basis_set_in, basis_set_out)
...
subroutine, public create_primitive_basis_set(basis_set, pbasis, lmax)
...
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public niklasson2014
integer, save, public niklasson2003
Handles all functions related to the CELL.
methods related to the blacs parallel environment
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
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
Types needed for a linear scaling quickstep SCF run based on the density matrix.
Types needed for a for a Energy Correction.
Energy correction environment setup and handling.
subroutine, public ec_write_input(ec_env)
Print out the energy correction input section.
subroutine, public ec_env_create(qs_env, ec_env, dft_section, ec_section)
Allocates and intitializes ec_env.
Defines the basic variable types.
integer, parameter, public dp
Routines needed for kpoint calculation.
subroutine, public kpoint_initialize(kpoint, particle_set, cell)
Generate the kpoints and initialize the kpoint environment.
subroutine, public kpoint_env_initialize(kpoint, para_env, blacs_env, with_aux_fit)
Initialize the kpoint environment.
Types and basic routines needed for a kpoint calculation.
subroutine, public read_kpoint_section(kpoint, kpoint_section, a_vec)
Read the kpoint input section.
subroutine, public write_kpoint_info(kpoint, iounit, dft_section)
Write information on the kpoints to output.
subroutine, public kpoint_create(kpoint)
Create a kpoint environment.
Interface to the message passing library MPI.
Define the data structure for the molecule information.
subroutine, public molecule_of_atom(molecule_set, atom_to_mol)
finds for each atom the molecule it belongs to
Provides Cartesian and spherical orbital pointers and indices.
subroutine, public init_orbital_pointers(maxl)
Initialize or update the orbital pointers.
Define the data structure for the particle information.
Calculation of non local dispersion functionals Some routines adapted from: Copyright (C) 2001-2009 Q...
subroutine, public qs_dispersion_nonloc_init(dispersion_env, para_env)
...
Calculation of dispersion using pair potentials.
subroutine, public qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, pp_section, para_env)
...
Definition of disperson types for DFT calculations.
Set disperson types for DFT calculations.
subroutine, public qs_dispersion_env_set(dispersion_env, xc_section)
...
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.
Calculate the interaction radii for the operator matrix calculation.
subroutine, public init_interaction_radii_orb_basis(orb_basis_set, eps_pgf_orb, eps_pgf_short)
...
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.
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public create_soft_basis(orb_basis, soft_basis, eps_fit, rc, paw_atom, paw_type_forced, gpw_r3d_rs_type_forced)
create the soft basis from a GTO basis
Utilities for string manipulations.
elemental subroutine, public uppercase(string)
Convert all lower case characters in a string to upper case.
Exchange and Correlation functional calculations.
logical function, public xc_uses_norm_drho(xc_fun_section, lsd)
...
logical function, public xc_uses_kinetic_energy_density(xc_fun_section, lsd)
...
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
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.