221#include "./base/base_uses.f90"
228 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_environment'
252 SUBROUTINE qs_init(qs_env, para_env, root_section, globenv, cp_subsys, kpoint_env, &
253 qmmm, qmmm_env_qm, force_env_section, subsys_section, &
254 use_motion_section, silent)
262 LOGICAL,
INTENT(IN),
OPTIONAL :: qmmm
265 LOGICAL,
INTENT(IN) :: use_motion_section
266 LOGICAL,
INTENT(IN),
OPTIONAL :: silent
268 CHARACTER(LEN=default_string_length) :: basis_type
269 INTEGER :: ikind, method_id, nelectron_total, &
271 LOGICAL :: do_active_space, do_admm_rpa, do_bse, do_debug_fdiff, do_debug_forces, &
272 do_debug_stress_tensor, do_ec_hfx, do_et, do_exx, do_gw, do_hfx, do_kpoints, &
273 do_linear_response, do_mp2, do_ri_mp2, do_ri_rpa, do_ri_sos_mp2, do_tddfpt, &
274 do_wfc_low_scaling, do_wfc_low_scaling_kpoints, do_xtb_tblite, is_identical, is_semi, &
275 mp2_present, my_qmmm, qmmm_decoupl, same_except_frac, use_ref_cell
276 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: rtmat
278 TYPE(
cell_type),
POINTER :: my_cell, my_cell_ref
288 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
294 TYPE(
section_vals_type),
POINTER :: active_space_section, dft_section, ec_hfx_section, &
295 ec_section, et_coupling_section, gw_section, hfx_section, kpoint_section, mp2_section, &
296 rpa_hfx_section, tddfpt_section, transport_section
298 NULLIFY (my_cell, my_cell_ref, atomic_kind_set, particle_set, &
299 qs_kind_set, kpoint_section, dft_section, ec_section, &
300 subsys, ks_env, dft_control, blacs_env)
302 CALL set_qs_env(qs_env, input=force_env_section)
303 IF (.NOT.
ASSOCIATED(subsys_section))
THEN
309 IF (
PRESENT(qmmm)) my_qmmm = qmmm
310 qmmm_decoupl = .false.
311 IF (
PRESENT(qmmm_env_qm))
THEN
315 qmmm_decoupl = my_qmmm .AND. qmmm_env_qm%periodic .AND. qmmm_env_qm%multipole
317 qs_env%qmmm_env_qm => qmmm_env_qm
323 SELECT CASE (method_id)
336 force_env_section=force_env_section, &
337 subsys_section=subsys_section, &
338 use_motion_section=use_motion_section, &
339 root_section=root_section, &
340 cp_subsys=cp_subsys, &
341 elkind=is_semi, silent=silent)
350 cell_ref=my_cell_ref, &
351 use_ref_cell=use_ref_cell, &
352 atomic_kind_set=atomic_kind_set, &
353 qs_kind_set=qs_kind_set, &
354 particle_set=particle_set)
357 IF (
PRESENT(globenv))
THEN
359 globenv%blacs_repeatable)
368 force_env_section, subsys_section, para_env)
371 IF (
PRESENT(kpoint_env))
THEN
372 kpoints => kpoint_env
373 CALL set_qs_env(qs_env=qs_env, kpoints=kpoints)
378 CALL set_qs_env(qs_env=qs_env, kpoints=kpoints)
390 IF (.NOT. do_gw)
THEN
405 do_active_space = .false.
408 do_xtb_tblite = .false.
413 do_linear_response = .false.
415 do_debug_fdiff = .false.
416 IF (
PRESENT(globenv)) do_debug_fdiff = globenv%run_type_id ==
debug_run
417 IF (do_debug_fdiff .AND.
PRESENT(root_section))
THEN
419 l_val=do_debug_forces)
421 l_val=do_debug_stress_tensor)
422 do_debug_fdiff = do_debug_forces .OR. do_debug_stress_tensor
426 do_ri_sos_mp2 = .false.
428 do_wfc_low_scaling = .false.
429 do_wfc_low_scaling_kpoints = .false.
432 IF (mp2_present)
THEN
437 CALL section_vals_val_get(qs_env%input,
"DFT%XC%WF_CORRELATION%RI_SOS_MP2%_SECTION_PARAMETERS_", &
441 CALL section_vals_val_get(qs_env%input,
"DFT%XC%WF_CORRELATION%LOW_SCALING%_SECTION_PARAMETERS_", &
442 l_val=do_wfc_low_scaling)
444 l_val=do_wfc_low_scaling_kpoints)
447 "DFT%XC%WF_CORRELATION%RI_RPA%GW%BSE%_SECTION_PARAMETERS_", &
450 CALL restrict_unsupported_atomic_kpoint_symmetry(kpoints, method_id, do_hfx, do_exx, do_gw, &
451 do_tddfpt, do_active_space, do_linear_response, &
453 do_mp2 .OR. do_ri_mp2 .OR. do_ri_sos_mp2, &
454 do_ri_rpa .AND. .NOT. do_gw, do_bse, &
455 do_wfc_low_scaling, do_wfc_low_scaling_kpoints, &
462 CALL qs_init_subsys(qs_env, para_env, subsys, my_cell, my_cell_ref, use_ref_cell, &
463 subsys_section, silent=silent)
465 CALL get_qs_env(qs_env, dft_control=dft_control)
466 IF (method_id ==
do_method_lrigpw .OR. dft_control%qs_control%lri_optbas)
THEN
467 CALL get_qs_env(qs_env=qs_env, lri_env=lri_env)
470 CALL cp_warn(__location__,
"Experimental code: "// &
471 "RIGPW should only be used for testing.")
472 CALL get_qs_env(qs_env=qs_env, lri_env=lri_env)
476 IF (my_qmmm .AND.
PRESENT(qmmm_env_qm) .AND. .NOT. dft_control%qs_control%commensurate_mgrids)
THEN
478 CALL cp_abort(__location__,
"QM/MM with coupling GAUSS or S-WAVE requires "// &
479 "keyword FORCE_EVAL/DFT/MGRID/COMMENSURATE to be enabled.")
484 CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints, blacs_env=blacs_env)
488 CALL get_qs_env(qs_env=qs_env, wf_history=wf_history)
499 CALL get_qs_env(qs_env, dft_control=dft_control, scf_control=scf_control, nelectron_total=nelectron_total)
504 IF (dft_control%do_admm)
THEN
505 basis_type =
'AUX_FIT'
509 CALL hfx_create(qs_env%x_data, para_env, hfx_section, atomic_kind_set, &
510 qs_kind_set, particle_set, dft_control, my_cell, orb_basis=basis_type, &
511 nelectron_total=nelectron_total, nkp_grid=nkp_grid)
516 IF (mp2_present)
THEN
517 cpassert(
ASSOCIATED(qs_env%mp2_env))
530 qs_env%mp2_env%ri_rpa%reuse_hfx = .true.
531 IF (.NOT. do_hfx) qs_env%mp2_env%ri_rpa%reuse_hfx = .false.
533 IF (.NOT. (is_identical .OR. same_except_frac)) qs_env%mp2_env%ri_rpa%reuse_hfx = .false.
534 IF (dft_control%do_admm .AND. .NOT. do_admm_rpa) qs_env%mp2_env%ri_rpa%reuse_hfx = .false.
536 IF (.NOT. qs_env%mp2_env%ri_rpa%reuse_hfx)
THEN
537 IF (do_admm_rpa)
THEN
538 basis_type =
'AUX_FIT'
542 CALL hfx_create(qs_env%mp2_env%ri_rpa%x_data, para_env, rpa_hfx_section, atomic_kind_set, &
543 qs_kind_set, particle_set, dft_control, my_cell, orb_basis=basis_type, &
544 nelectron_total=nelectron_total)
546 qs_env%mp2_env%ri_rpa%x_data => qs_env%x_data
551 IF (dft_control%qs_control%do_kg)
THEN
553 CALL kg_env_create(qs_env, qs_env%kg_env, qs_kind_set, qs_env%input)
558 l_val=qs_env%excited_state)
559 NULLIFY (exstate_env)
560 CALL exstate_create(exstate_env, qs_env%excited_state, dft_section)
561 CALL set_qs_env(qs_env, exstate_env=exstate_env)
564 "PROPERTIES%ET_COUPLING")
570 IF (qs_env%do_transport)
THEN
574 CALL get_qs_env(qs_env, harris_env=harris_env)
575 IF (qs_env%harris_method)
THEN
577 CALL get_qs_env(qs_env, local_particles=local_particles)
578 CALL harris_rhoin_init(harris_env%rhoin,
"RHOIN", qs_kind_set, atomic_kind_set, &
579 local_particles, dft_control%nspins)
587 l_val=qs_env%energy_correction)
592 IF (qs_env%energy_correction)
THEN
597 IF (ec_env%do_ec_hfx)
THEN
600 IF (ec_env%do_kpoints)
THEN
601 CALL cp_abort(__location__, &
602 "Energy correction methods with hybrid functionals "// &
603 "and kpoints is not yet available.")
607 IF (ec_env%basis_inconsistent)
THEN
608 CALL cp_abort(__location__, &
609 "Energy correction methods with hybrid functionals: "// &
610 "correction and ground state need to use the same basis. "// &
611 "Checked by comparing basis set names only.")
615 IF (ec_env%do_ec_admm .AND. .NOT. dft_control%do_admm)
THEN
616 CALL cp_abort(__location__,
"Need an ADMM input section for ADMM EC to work")
619 ec_env%reuse_hfx = .true.
620 IF (.NOT. do_hfx) ec_env%reuse_hfx = .false.
622 IF (.NOT. (is_identical .OR. same_except_frac)) ec_env%reuse_hfx = .false.
623 IF (dft_control%do_admm .AND. .NOT. ec_env%do_ec_admm) ec_env%reuse_hfx = .false.
625 IF (.NOT. ec_env%reuse_hfx)
THEN
626 IF (ec_env%do_ec_admm)
THEN
627 basis_type =
'AUX_FIT'
631 CALL hfx_create(ec_env%x_data, para_env, ec_hfx_section, atomic_kind_set, &
632 qs_kind_set, particle_set, dft_control, my_cell, orb_basis=basis_type, &
633 nelectron_total=nelectron_total)
635 ec_env%x_data => qs_env%x_data
644 IF (dft_control%qs_control%do_almo_scf)
THEN
649 CALL get_qs_env(qs_env, rel_control=rel_control)
650 IF (rel_control%rel_method /=
rel_none)
THEN
652 nkind =
SIZE(atomic_kind_set)
656 IF (
ASSOCIATED(rtmat))
CALL set_qs_kind(qs_kind_set(ikind), reltmat=rtmat)
681 SUBROUTINE restrict_unsupported_atomic_kpoint_symmetry(kpoints, method_id, do_hfx, do_exx, do_gw, &
682 do_tddfpt, do_active_space, do_linear_response, &
684 do_mp2, do_rpa, do_bse, do_wfc_low_scaling, &
685 do_wfc_low_scaling_kpoints, do_xtb_tblite)
687 INTEGER,
INTENT(IN) :: method_id
688 LOGICAL,
INTENT(IN) :: do_hfx, do_exx, do_gw, do_tddfpt, do_active_space, &
689 do_linear_response, do_debug_fdiff, do_mp2, do_rpa, do_bse, do_wfc_low_scaling, &
690 do_wfc_low_scaling_kpoints, do_xtb_tblite
692 CHARACTER(LEN=default_string_length) :: kp_scheme, reason
693 LOGICAL :: full_grid, inversion_symmetry_only, &
696 reason = unsupported_kpoint_method_reason(method_id, do_gw, do_tddfpt, do_linear_response, &
697 do_mp2, do_bse, do_xtb_tblite)
698 IF (len_trim(reason) > 0)
THEN
700 IF (len_trim(kp_scheme) > 0 .AND. trim(kp_scheme) /=
"NONE")
THEN
701 IF (trim(reason) ==
"GW")
THEN
702 CALL cp_abort(__location__, &
703 "DFT%KPOINTS are not supported with GW; use "// &
704 "WF_CORRELATION%LOW_SCALING%KPOINTS and RI_RPA%GW%KPOINTS_SELF_ENERGY "// &
705 "for GW k-point sampling.")
707 CALL cp_abort(__location__, &
708 "DFT%KPOINTS are not supported with "//trim(reason)// &
709 "; remove DFT%KPOINTS for these calculations.")
713 IF (do_active_space)
THEN
715 IF (len_trim(kp_scheme) > 0 .AND. trim(kp_scheme) /=
"NONE" .AND. &
716 trim(kp_scheme) /=
"GAMMA")
THEN
717 CALL cp_abort(__location__, &
718 "Only Gamma-point DFT%KPOINTS are supported with ACTIVE_SPACE; "// &
719 "use SCHEME GAMMA, SCHEME NONE, or remove DFT%KPOINTS.")
723 CALL get_kpoint_info(kpoints, symmetry=kpoint_symmetry, full_grid=full_grid, &
724 inversion_symmetry_only=inversion_symmetry_only)
725 IF (.NOT. (kpoint_symmetry .AND. .NOT. full_grid .AND. .NOT. inversion_symmetry_only))
RETURN
727 reason = unsupported_atomic_kpoint_symmetry_reason(method_id, do_hfx, do_exx, do_gw, &
728 do_tddfpt, do_active_space, do_linear_response, &
730 do_mp2, do_rpa, do_bse, do_wfc_low_scaling, &
731 do_wfc_low_scaling_kpoints, do_xtb_tblite)
732 IF (len_trim(reason) == 0)
RETURN
734 CALL cp_warn(__location__, &
735 "Atomic k-point symmetry is currently not implemented for "//trim(reason)// &
736 "; restricting to inversion/time-reversal symmetry.")
739 END SUBROUTINE restrict_unsupported_atomic_kpoint_symmetry
752 FUNCTION unsupported_kpoint_method_reason(method_id, do_gw, do_tddfpt, do_linear_response, &
753 do_mp2, do_bse, do_xtb_tblite)
RESULT(reason)
754 INTEGER,
INTENT(IN) :: method_id
755 LOGICAL,
INTENT(IN) :: do_gw, do_tddfpt, do_linear_response, &
756 do_mp2, do_bse, do_xtb_tblite
757 CHARACTER(LEN=default_string_length) :: reason
762 mark_used(do_xtb_tblite)
769 reason =
"TDDFPT/TDDFT"
772 IF (do_linear_response)
THEN
773 reason =
"LINEAR_RESPONSE/DFPT"
776 SELECT CASE (method_id)
783 reason =
"semiempirical methods"
788 END FUNCTION unsupported_kpoint_method_reason
808 FUNCTION unsupported_atomic_kpoint_symmetry_reason(method_id, do_hfx, do_exx, do_gw, do_tddfpt, &
809 do_active_space, do_linear_response, do_debug_fdiff, &
810 do_mp2, do_rpa, do_bse, do_wfc_low_scaling, &
811 do_wfc_low_scaling_kpoints, do_xtb_tblite)
RESULT(reason)
812 INTEGER,
INTENT(IN) :: method_id
813 LOGICAL,
INTENT(IN) :: do_hfx, do_exx, do_gw, do_tddfpt, do_active_space, &
814 do_linear_response, do_debug_fdiff, do_mp2, do_rpa, do_bse, do_wfc_low_scaling, &
815 do_wfc_low_scaling_kpoints, do_xtb_tblite
816 CHARACTER(LEN=default_string_length) :: reason
819 mark_used(do_debug_fdiff)
820 mark_used(do_xtb_tblite)
822 SELECT CASE (method_id)
829 reason =
"semiempirical methods"
834 IF (len_trim(reason) > 0)
RETURN
839 ELSE IF (do_hfx .OR. do_exx)
THEN
841 ELSE IF (do_tddfpt)
THEN
842 reason =
"TDDFPT/TDDFT"
843 ELSE IF (do_active_space)
THEN
844 reason =
"ACTIVE_SPACE"
845 ELSE IF (do_linear_response)
THEN
846 reason =
"LINEAR_RESPONSE/DFPT"
847 ELSE IF (do_mp2)
THEN
849 ELSE IF (do_rpa .AND. do_wfc_low_scaling_kpoints)
THEN
850 reason =
"LOW_SCALING RPA"
851 ELSE IF (do_wfc_low_scaling)
THEN
852 reason =
"LOW_SCALING WF_CORRELATION"
853 ELSE IF (do_rpa)
THEN
857 END FUNCTION unsupported_atomic_kpoint_symmetry_reason
871 SUBROUTINE qs_init_subsys(qs_env, para_env, subsys, cell, cell_ref, use_ref_cell, subsys_section, &
877 TYPE(
cell_type),
POINTER :: cell, cell_ref
878 LOGICAL,
INTENT(in) :: use_ref_cell
880 LOGICAL,
INTENT(in),
OPTIONAL :: silent
882 CHARACTER(len=*),
PARAMETER :: routinen =
'qs_init_subsys'
884 CHARACTER(len=2) :: element_symbol
885 INTEGER :: gfn_type, handle, ikind, ispin, iw, lmax_sphere, maxl, maxlgto, maxlgto_lri, &
886 maxlgto_nuc, maxlppl, maxlppnl, method_id, multiplicity, my_ival, n_ao, n_mo_add, natom, &
887 nelectron, ngauss, nkind, output_unit, sort_basis, tnadd_method
888 INTEGER,
DIMENSION(2) :: n_mo, nelectron_spin
889 INTEGER,
DIMENSION(5) :: occ
890 INTEGER,
DIMENSION(:),
POINTER :: mo_index_range
891 LOGICAL :: all_potential_present, be_silent, cneo_potential_present, do_kpoints, do_ri_hfx, &
892 do_ri_mp2, do_ri_rpa, do_ri_sos_mp2, do_rpa_ri_exx, do_wfc_im_time, e1terms, &
893 has_unit_metric, lribas, mp2_present, orb_gradient, paw_atom
894 REAL(kind=
dp) :: alpha, ccore, ewald_rcut, fxx, maxocc, &
895 rc, rcut, total_zeff_corr, &
896 verlet_skin, zeff_correction
907 rhoin_basis, ri_aux_basis_set, &
908 ri_hfx_basis, ri_xas_basis, &
913 TYPE(
mo_set_type),
DIMENSION(:),
POINTER :: mos, mos_last_converged
922 POINTER :: dftb_potential
927 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
935 TYPE(
section_vals_type),
POINTER :: dft_section, et_coupling_section, et_ddapc_section, &
936 ewald_section, harris_section, lri_section, mp2_section, nl_section, poisson_section, &
937 pp_section, print_section, qs_section, rixs_section, se_section, tddfpt_section, &
943 CALL timeset(routinen, handle)
949 IF (
PRESENT(silent)) be_silent = silent
954 NULLIFY (mos, se_taper)
955 NULLIFY (dft_control)
958 NULLIFY (local_molecules)
959 NULLIFY (local_particles)
960 NULLIFY (scf_control)
961 NULLIFY (dft_section)
962 NULLIFY (et_coupling_section)
964 NULLIFY (mos_last_converged)
973 qs_kind_set=qs_kind_set, &
974 atomic_kind_set=atomic_kind_set, &
975 molecule_set=molecule_set, &
976 molecule_kind_set=molecule_kind_set)
982 dft_control%qs_control%periodicity = sum(cell%perd)
988 IF (.NOT. be_silent)
THEN
991 SELECT CASE (method_id)
1000 IF (dft_control%qs_control%xtb_control%do_tblite)
THEN
1001 CALL tblite_header(iw, dft_control%qs_control%xtb_control%tblite_method)
1003 gfn_type = dft_control%qs_control%xtb_control%gfn_type
1008 "PRINT%PROGRAM_BANNER")
1011 IF (dft_control%do_sccs .AND. dft_control%qs_control%gapw)
THEN
1012 cpabort(
"SCCS is not yet implemented with GAPW")
1014 CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints)
1015 IF (do_kpoints)
THEN
1017 SELECT CASE (dft_control%qs_control%wf_interpolation_method_nr)
1020 CALL cp_warn(__location__,
"Linear WFN-based extrapolation methods are not "// &
1021 "implemented for k-points. Switching to USE_PREV_WF.")
1028 dft_control%qs_control%et_coupling_calc = .false.
1031 dft_control%qs_control%et_coupling_calc = .true.
1032 dft_control%qs_control%ddapc_restraint = .true.
1033 CALL read_ddapc_section(dft_control%qs_control, ddapc_restraint_section=et_ddapc_section)
1043 IF (qs_env%do_rixs)
THEN
1044 CALL read_rixs_control(dft_control%rixs_control, rixs_section, dft_control%qs_control)
1050 ALLOCATE (rel_control)
1053 CALL set_qs_env(qs_env, rel_control=rel_control)
1058 NULLIFY (ewald_env, ewald_pw, dftb_potential)
1059 dftb_control => dft_control%qs_control%dftb_control
1060 CALL qs_dftb_param_init(atomic_kind_set, qs_kind_set, dftb_control, dftb_potential, &
1061 subsys_section=subsys_section, para_env=para_env)
1062 CALL set_qs_env(qs_env, dftb_potential=dftb_potential)
1064 IF (dftb_control%do_ewald)
THEN
1065 ALLOCATE (ewald_env)
1068 CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
1073 cell_periodic=cell%perd)
1075 CALL ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section=print_section)
1076 CALL set_qs_env(qs_env, ewald_env=ewald_env, ewald_pw=ewald_pw)
1078 ELSEIF (dft_control%qs_control%method_id ==
do_method_xtb)
THEN
1080 xtb_control => dft_control%qs_control%xtb_control
1082 IF (xtb_control%do_tblite)
THEN
1091 qs_kind => qs_kind_set(ikind)
1093 cpassert(.NOT.
ASSOCIATED(qs_kind%xtb_parameter))
1096 CALL get_qs_kind(qs_kind, element_symbol=element_symbol)
1098 NULLIFY (tmp_basis_set)
1099 CALL tb_get_basis(qs_env%tb_tblite, tmp_basis_set, element_symbol, qs_kind%xtb_parameter, occ)
1104 zeff_correction = 0.0_dp
1106 zeff=real(sum(occ),
dp), zeff_correction=zeff_correction)
1109 NULLIFY (ewald_env, ewald_pw)
1111 qs_kind => qs_kind_set(ikind)
1113 cpassert(.NOT.
ASSOCIATED(qs_kind%xtb_parameter))
1116 gfn_type = dft_control%qs_control%xtb_control%gfn_type
1117 CALL get_qs_kind(qs_kind, element_symbol=element_symbol)
1119 xtb_control%parameter_file_path, xtb_control%parameter_file_name, &
1124 NULLIFY (tmp_basis_set)
1125 IF (qs_kind%xtb_parameter%z == 1)
THEN
1127 ngauss = xtb_control%h_sto_ng
1129 ngauss = xtb_control%sto_ng
1131 IF (qs_kind%xtb_parameter%defined)
THEN
1132 CALL init_xtb_basis(qs_kind%xtb_parameter, tmp_basis_set, ngauss)
1136 IF (
ASSOCIATED(qs_kind%all_potential))
THEN
1137 DEALLOCATE (qs_kind%all_potential%elec_conf)
1138 DEALLOCATE (qs_kind%all_potential)
1142 IF (qs_kind%xtb_parameter%defined)
THEN
1143 zeff_correction = 0.0_dp
1145 zeff=qs_kind%xtb_parameter%zeff, zeff_correction=zeff_correction)
1146 CALL get_potential(qs_kind%all_potential, alpha_core_charge=alpha)
1147 ccore = qs_kind%xtb_parameter%zeff*sqrt((alpha/
pi)**3)
1148 CALL set_potential(qs_kind%all_potential, ccore_charge=ccore)
1149 qs_kind%xtb_parameter%zeff = qs_kind%xtb_parameter%zeff - zeff_correction
1155 ALLOCATE (xtb_control%rcpair(nkind, nkind))
1156 CALL xtb_pp_radius(qs_kind_set, xtb_control%rcpair, xtb_control%eps_pair, xtb_control%kf)
1158 IF (xtb_control%do_ewald)
THEN
1159 ALLOCATE (ewald_env)
1162 CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
1165 IF (gfn_type == 0)
THEN
1167 silent=silent, pset=
"EEQ", cell_periodic=cell%perd)
1170 silent=silent, cell_periodic=cell%perd)
1173 CALL ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section=print_section)
1174 CALL set_qs_env(qs_env, ewald_env=ewald_env, ewald_pw=ewald_pw)
1181 dft_control%qs_control%lri_optbas .OR. &
1191 CALL get_qs_kind_set(qs_kind_set, all_potential_present=all_potential_present)
1192 IF ((dft_control%qs_control%method_id ==
do_method_gpw) .OR. &
1195 IF (all_potential_present)
THEN
1196 cpabort(
"All-electron calculations with GPW, GAPW_XC, and OFGPW are not implemented")
1201 CALL get_qs_kind_set(qs_kind_set, cneo_potential_present=cneo_potential_present)
1202 IF (cneo_potential_present .AND. &
1204 cpabort(
"CNEO calculations require GAPW method")
1208 CALL get_qs_kind_set(qs_kind_set, dft_plus_u_atom_present=dft_control%dft_plus_u)
1210 IF (dft_control%do_admm)
THEN
1214 NULLIFY (aux_fit_basis)
1215 qs_kind => qs_kind_set(ikind)
1216 CALL get_qs_kind(qs_kind, basis_set=aux_fit_basis, basis_type=
"AUX_FIT")
1217 IF (.NOT. (
ASSOCIATED(aux_fit_basis)))
THEN
1219 cpabort(
"AUX_FIT basis set is not defined. ")
1229 e1terms = lri_env%exact_1c_terms
1231 IF (dft_control%qs_control%do_kg)
THEN
1239 NULLIFY (lri_aux_basis)
1240 qs_kind => qs_kind_set(ikind)
1241 CALL get_qs_kind(qs_kind, basis_set=lri_aux_basis, basis_type=
"LRI_AUX")
1242 IF (.NOT. (
ASSOCIATED(lri_aux_basis)))
THEN
1244 CALL cp_warn(__location__,
"Automatic Generation of LRI_AUX basis. "// &
1245 "This is experimental code.")
1254 CALL section_vals_val_get(qs_env%input,
"DFT%XC%WF_CORRELATION%RI_RPA%HF%RI%_SECTION_PARAMETERS_", &
1255 l_val=do_rpa_ri_exx)
1256 IF (do_ri_hfx .OR. do_rpa_ri_exx)
THEN
1260 NULLIFY (ri_hfx_basis)
1261 qs_kind => qs_kind_set(ikind)
1262 CALL get_qs_kind(qs_kind=qs_kind, basis_set=ri_hfx_basis, &
1263 basis_type=
"RI_HFX")
1264 IF (.NOT. (
ASSOCIATED(ri_hfx_basis)))
THEN
1266 IF (dft_control%do_admm)
THEN
1268 basis_type=
"AUX_FIT", basis_sort=sort_basis)
1271 basis_sort=sort_basis)
1282 NULLIFY (ri_hfx_basis)
1283 qs_kind => qs_kind_set(ikind)
1284 CALL get_qs_kind(qs_kind, basis_set=ri_hfx_basis, basis_type=
"RI_HXC")
1285 IF (.NOT. (
ASSOCIATED(ri_hfx_basis)))
THEN
1294 NULLIFY (harris_env)
1296 l_val=qs_env%harris_method)
1299 CALL set_qs_env(qs_env, harris_env=harris_env)
1301 IF (qs_env%harris_method)
THEN
1305 NULLIFY (tmp_basis_set)
1306 qs_kind => qs_kind_set(ikind)
1307 CALL get_qs_kind(qs_kind, basis_set=rhoin_basis, basis_type=
"RHOIN")
1308 IF (.NOT. (
ASSOCIATED(rhoin_basis)))
THEN
1311 IF (qs_env%harris_env%density_source ==
hden_atomic)
THEN
1315 rhoin_basis => tmp_basis_set
1324 IF (mp2_present)
THEN
1328 CALL section_vals_val_get(qs_env%input,
"DFT%XC%WF_CORRELATION%LOW_SCALING%_SECTION_PARAMETERS_", &
1329 l_val=do_wfc_im_time)
1332 CALL cp_warn(__location__, &
1333 "Low-scaling RPA requires SORT_BASIS EXP keyword (in DFT input section) for good performance")
1338 CALL get_qs_env(qs_env, mp2_env=mp2_env, nkind=nkind)
1339 CALL section_vals_val_get(qs_env%input,
"DFT%XC%WF_CORRELATION%RI_MP2%_SECTION_PARAMETERS_", l_val=do_ri_mp2)
1340 CALL section_vals_val_get(qs_env%input,
"DFT%XC%WF_CORRELATION%RI_SOS_MP2%_SECTION_PARAMETERS_", l_val=do_ri_sos_mp2)
1341 CALL section_vals_val_get(qs_env%input,
"DFT%XC%WF_CORRELATION%RI_RPA%_SECTION_PARAMETERS_", l_val=do_ri_rpa)
1342 IF (do_ri_mp2 .OR. do_ri_sos_mp2 .OR. do_ri_rpa)
THEN
1344 NULLIFY (ri_aux_basis_set)
1345 qs_kind => qs_kind_set(ikind)
1346 CALL get_qs_kind(qs_kind=qs_kind, basis_set=ri_aux_basis_set, &
1347 basis_type=
"RI_AUX")
1348 IF (.NOT. (
ASSOCIATED(ri_aux_basis_set)))
THEN
1355 qs_env%mp2_env%ri_aux_auto_generated = .true.
1362 IF (dft_control%do_xas_tdp_calculation .OR. qs_env%do_rixs)
THEN
1366 NULLIFY (ri_xas_basis)
1367 qs_kind => qs_kind_set(ikind)
1368 CALL get_qs_kind(qs_kind, basis_set=ri_xas_basis, basis_type=
"RI_XAS")
1369 IF (.NOT.
ASSOCIATED(ri_xas_basis))
THEN
1378 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto, maxlppl=maxlppl, maxlppnl=maxlppnl)
1381 IF (cneo_potential_present)
THEN
1382 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto_nuc, basis_type=
"NUC")
1383 maxlgto = max(maxlgto, maxlgto_nuc)
1385 lmax_sphere = dft_control%qs_control%gapw_control%lmax_sphere
1386 IF (lmax_sphere < 0)
THEN
1387 lmax_sphere = 2*maxlgto
1388 dft_control%qs_control%gapw_control%lmax_sphere = lmax_sphere
1390 IF (dft_control%qs_control%method_id ==
do_method_lrigpw .OR. dft_control%qs_control%lri_optbas)
THEN
1391 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto_lri, basis_type=
"LRI_AUX")
1393 maxlgto = max(maxlgto, maxlgto_lri)
1395 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto_lri, basis_type=
"RI_HXC")
1396 maxlgto = max(maxlgto, maxlgto_lri)
1398 IF (dft_control%do_xas_tdp_calculation .OR. qs_env%do_rixs)
THEN
1400 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto_lri, basis_type=
"RI_XAS")
1401 maxlgto = max(maxlgto, maxlgto_lri)
1403 maxl = max(2*maxlgto, maxlppl, maxlppnl, lmax_sphere) + 1
1415 qs_control => dft_control%qs_control
1420 IF (cneo_potential_present)
THEN
1427 maxl = max(3*maxlgto + 1, 0)
1434 IF (.NOT. dft_control%qs_control%xtb_control%do_tblite)
THEN
1438 qs_kind => qs_kind_set(ikind)
1439 IF (qs_kind%xtb_parameter%defined)
THEN
1440 CALL get_qs_kind(qs_kind, basis_set=tmp_basis_set)
1441 rcut = xtb_control%coulomb_sr_cut
1442 fxx = 2.0_dp*xtb_control%coulomb_sr_eps*qs_kind%xtb_parameter%eta**2
1443 fxx = 0.80_dp*(1.0_dp/fxx)**0.3333_dp
1444 rcut = min(rcut, xtb_control%coulomb_sr_cut)
1445 qs_kind%xtb_parameter%rcut = min(rcut, fxx)
1447 qs_kind%xtb_parameter%rcut = 0.0_dp
1453 IF (.NOT. be_silent)
THEN
1466 particle_set=particle_set, &
1467 local_particles=local_particles, &
1468 molecule_kind_set=molecule_kind_set, &
1469 molecule_set=molecule_set, &
1470 local_molecules=local_molecules, &
1471 force_env_section=qs_env%input)
1474 ALLOCATE (scf_control)
1476 IF (dft_control%qs_control%dftb)
THEN
1477 scf_control%non_selfconsistent = .NOT. dft_control%qs_control%dftb_control%self_consistent
1479 IF (dft_control%qs_control%xtb)
THEN
1480 IF (dft_control%qs_control%xtb_control%do_tblite)
THEN
1481 scf_control%non_selfconsistent = .false.
1483 scf_control%non_selfconsistent = (dft_control%qs_control%xtb_control%gfn_type == 0)
1486 IF (qs_env%harris_method)
THEN
1487 scf_control%non_selfconsistent = .true.
1496 has_unit_metric = .false.
1497 IF (dft_control%qs_control%semi_empirical)
THEN
1498 IF (dft_control%qs_control%se_control%orthogonal_basis) has_unit_metric = .true.
1500 IF (dft_control%qs_control%dftb)
THEN
1501 IF (dft_control%qs_control%dftb_control%orthogonal_basis) has_unit_metric = .true.
1503 CALL set_qs_env(qs_env, has_unit_metric=has_unit_metric)
1507 interpolation_method_nr= &
1508 dft_control%qs_control%wf_interpolation_method_nr, &
1509 extrapolation_order=dft_control%qs_control%wf_extrapolation_order, &
1510 has_unit_metric=has_unit_metric)
1514 scf_control=scf_control, &
1515 wf_history=wf_history)
1518 cell_ref=cell_ref, &
1519 use_ref_cell=use_ref_cell, &
1524 CALL set_ks_env(ks_env, dft_control=dft_control)
1526 CALL qs_subsys_set(subsys, local_molecules=local_molecules, &
1527 local_particles=local_particles, cell=cell)
1534 atomic_kind_set=atomic_kind_set, &
1535 dft_control=dft_control, &
1536 scf_control=scf_control)
1540 IF (dft_control%qs_control%do_ls_scf .OR. &
1541 dft_control%qs_control%do_almo_scf)
THEN
1542 CALL set_qs_env(qs_env=qs_env, requires_mo_derivs=.false.)
1544 IF (scf_control%use_ot)
THEN
1545 CALL set_qs_env(qs_env=qs_env, requires_mo_derivs=.true.)
1547 CALL set_qs_env(qs_env=qs_env, requires_mo_derivs=.false.)
1552 IF (dft_control%qs_control%xtb_control%do_tblite)
THEN
1553 IF (.NOT. scf_control%smear%do_smear)
THEN
1555 scf_control%smear%do_smear = .true.
1557 scf_control%smear%electronic_temperature = 300._dp/
kelvin
1558 scf_control%smear%eps_fermi_dirac = 1.e-6_dp
1561 dft_control%smear = scf_control%smear%do_smear
1564 IF (.NOT. (dft_control%qs_control%dftb .OR. dft_control%qs_control%xtb))
THEN
1565 IF (dft_control%apply_period_efield)
THEN
1566 CALL get_qs_env(qs_env=qs_env, requires_mo_derivs=orb_gradient)
1567 IF (.NOT. orb_gradient)
THEN
1568 CALL cp_abort(__location__,
"Periodic Efield needs orbital gradient and direct optimization."// &
1569 " Use the OT optimization method.")
1571 IF (dft_control%smear)
THEN
1572 CALL cp_abort(__location__,
"Periodic Efield needs equal occupation numbers."// &
1573 " Smearing option is not possible.")
1582 NULLIFY (rho_atom_set)
1583 gapw_control => dft_control%qs_control%gapw_control
1584 CALL init_rho_atom(rho_atom_set, atomic_kind_set, qs_kind_set, dft_control, para_env)
1585 CALL set_qs_env(qs_env=qs_env, rho_atom_set=rho_atom_set)
1587 CALL get_qs_env(qs_env=qs_env, local_rho_set=local_rho_set, natom=natom)
1589 CALL init_rho0(local_rho_set, qs_env, gapw_control)
1596 IF (gapw_control%accurate_xcint)
THEN
1597 cpassert(.NOT.
ASSOCIATED(gapw_control%aw))
1599 ALLOCATE (gapw_control%aw(nkind))
1600 alpha = gapw_control%aweights
1602 qs_kind => qs_kind_set(ikind)
1603 CALL get_qs_kind(qs_kind, hard_radius=rc, paw_atom=paw_atom)
1605 gapw_control%aw(ikind) = alpha*(1.2_dp/rc)**2
1607 gapw_control%aw(ikind) = 0.0_dp
1617 ELSE IF (dft_control%qs_control%semi_empirical)
THEN
1618 NULLIFY (se_store_int_env, se_nddo_mpole, se_nonbond_env)
1619 natom =
SIZE(particle_set)
1621 se_control => dft_control%qs_control%se_control
1626 SELECT CASE (dft_control%qs_control%method_id)
1632 CALL init_se_nlradius(se_control, atomic_kind_set, qs_kind_set, subsys_section)
1637 IF (se_control%do_ewald .OR. se_control%do_ewald_gks)
THEN
1638 ALLOCATE (ewald_env)
1641 CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
1644 "PRINT%GRID_INFORMATION")
1649 print_section=print_section)
1655 IF (se_control%do_ewald)
THEN
1656 CALL ewald_env_get(ewald_env, max_multipole=se_control%max_multipole)
1660 ALLOCATE (se_nonbond_env)
1662 do_electrostatics=.true., verlet_skin=verlet_skin, ewald_rcut=ewald_rcut, &
1663 ei_scale14=0.0_dp, vdw_scale14=0.0_dp, shift_cutoff=.false.)
1666 CALL set_qs_env(qs_env, ewald_env=ewald_env, ewald_pw=ewald_pw, &
1667 se_nonbond_env=se_nonbond_env, se_nddo_mpole=se_nddo_mpole)
1670 dft_control%qs_control%method_id)
1673 CALL se_taper_create(se_taper, se_control%integral_screening, se_control%do_ewald, &
1674 se_control%taper_cou, se_control%range_cou, &
1675 se_control%taper_exc, se_control%range_exc, &
1676 se_control%taper_scr, se_control%range_scr, &
1677 se_control%taper_lrc, se_control%range_lrc)
1681 CALL set_qs_env(qs_env, se_store_int_env=se_store_int_env)
1685 IF (dft_control%qs_control%method_id ==
do_method_gpw .OR. &
1691 ALLOCATE (dispersion_env)
1692 NULLIFY (xc_section)
1696 NULLIFY (pp_section)
1700 NULLIFY (nl_section)
1704 CALL set_qs_env(qs_env, dispersion_env=dispersion_env)
1705 ELSE IF (dft_control%qs_control%method_id ==
do_method_dftb)
THEN
1706 ALLOCATE (dispersion_env)
1708 dispersion_env%doabc = .false.
1709 dispersion_env%c9cnst = .false.
1710 dispersion_env%lrc = .false.
1711 dispersion_env%srb = .false.
1712 dispersion_env%verbose = .false.
1713 NULLIFY (dispersion_env%c6ab, dispersion_env%maxci, dispersion_env%r0ab, dispersion_env%rcov, &
1714 dispersion_env%r2r4, dispersion_env%cn, dispersion_env%cnkind, dispersion_env%cnlist, &
1715 dispersion_env%d3_exclude_pair)
1716 NULLIFY (dispersion_env%q_mesh, dispersion_env%kernel, dispersion_env%d2phi_dk2, &
1717 dispersion_env%d2y_dx2, dispersion_env%dftd_section)
1718 NULLIFY (dispersion_env%sab_vdw, dispersion_env%sab_cn)
1719 IF (dftb_control%dispersion .AND. dftb_control%dispersion_type ==
dispersion_d3)
THEN
1722 dispersion_env%eps_cn = dftb_control%epscn
1723 dispersion_env%s6 = dftb_control%sd3(1)
1724 dispersion_env%sr6 = dftb_control%sd3(2)
1725 dispersion_env%s8 = dftb_control%sd3(3)
1726 dispersion_env%domol = .false.
1727 dispersion_env%kgc8 = 0._dp
1728 dispersion_env%rc_disp = dftb_control%rcdisp
1729 dispersion_env%exp_pre = 0._dp
1730 dispersion_env%scaling = 0._dp
1731 dispersion_env%nd3_exclude_pair = 0
1732 dispersion_env%parameter_file_name = dftb_control%dispersion_parameter_file
1734 ELSEIF (dftb_control%dispersion .AND. dftb_control%dispersion_type ==
dispersion_d3bj)
THEN
1737 dispersion_env%eps_cn = dftb_control%epscn
1738 dispersion_env%s6 = dftb_control%sd3bj(1)
1739 dispersion_env%a1 = dftb_control%sd3bj(2)
1740 dispersion_env%s8 = dftb_control%sd3bj(3)
1741 dispersion_env%a2 = dftb_control%sd3bj(4)
1742 dispersion_env%domol = .false.
1743 dispersion_env%kgc8 = 0._dp
1744 dispersion_env%rc_disp = dftb_control%rcdisp
1745 dispersion_env%exp_pre = 0._dp
1746 dispersion_env%scaling = 0._dp
1747 dispersion_env%nd3_exclude_pair = 0
1748 dispersion_env%parameter_file_name = dftb_control%dispersion_parameter_file
1750 ELSEIF (dftb_control%dispersion .AND. dftb_control%dispersion_type ==
dispersion_d2)
THEN
1753 dispersion_env%exp_pre = dftb_control%exp_pre
1754 dispersion_env%scaling = dftb_control%scaling
1755 dispersion_env%parameter_file_name = dftb_control%dispersion_parameter_file
1756 dispersion_env%rc_disp = dftb_control%rcdisp
1761 CALL set_qs_env(qs_env, dispersion_env=dispersion_env)
1762 ELSE IF (dft_control%qs_control%method_id ==
do_method_xtb)
THEN
1763 IF (.NOT. (dft_control%qs_control%xtb_control%do_tblite))
THEN
1764 ALLOCATE (dispersion_env)
1766 dispersion_env%doabc = .false.
1767 dispersion_env%c9cnst = .false.
1768 dispersion_env%lrc = .false.
1769 dispersion_env%srb = .false.
1770 dispersion_env%verbose = .false.
1771 NULLIFY (dispersion_env%c6ab, dispersion_env%maxci, &
1772 dispersion_env%r0ab, dispersion_env%rcov, &
1773 dispersion_env%r2r4, dispersion_env%cn, &
1774 dispersion_env%cnkind, dispersion_env%cnlist, &
1775 dispersion_env%d3_exclude_pair)
1776 NULLIFY (dispersion_env%q_mesh, dispersion_env%kernel, dispersion_env%d2phi_dk2, &
1777 dispersion_env%d2y_dx2, dispersion_env%dftd_section)
1778 NULLIFY (dispersion_env%sab_vdw, dispersion_env%sab_cn)
1780 dispersion_env%eps_cn = xtb_control%epscn
1781 dispersion_env%s6 = xtb_control%s6
1782 dispersion_env%s8 = xtb_control%s8
1783 dispersion_env%a1 = xtb_control%a1
1784 dispersion_env%a2 = xtb_control%a2
1785 dispersion_env%domol = .false.
1786 dispersion_env%kgc8 = 0._dp
1787 dispersion_env%rc_disp = xtb_control%rcdisp
1788 dispersion_env%rc_d4 = xtb_control%rcdisp
1789 dispersion_env%exp_pre = 0._dp
1790 dispersion_env%scaling = 0._dp
1791 dispersion_env%nd3_exclude_pair = 0
1792 dispersion_env%parameter_file_name = xtb_control%dispersion_parameter_file
1794 SELECT CASE (xtb_control%vdw_type)
1801 dispersion_env%ref_functional =
"none"
1803 dispersion_env, para_env=para_env)
1804 dispersion_env%cnfun = 2
1808 CALL set_qs_env(qs_env, dispersion_env=dispersion_env)
1810 ELSE IF (dft_control%qs_control%semi_empirical)
THEN
1811 ALLOCATE (dispersion_env)
1813 dispersion_env%doabc = .false.
1814 dispersion_env%c9cnst = .false.
1815 dispersion_env%lrc = .false.
1816 dispersion_env%srb = .false.
1817 dispersion_env%verbose = .false.
1818 NULLIFY (dispersion_env%c6ab, dispersion_env%maxci, dispersion_env%r0ab, dispersion_env%rcov, &
1819 dispersion_env%r2r4, dispersion_env%cn, dispersion_env%cnkind, dispersion_env%cnlist, &
1820 dispersion_env%d3_exclude_pair)
1821 NULLIFY (dispersion_env%q_mesh, dispersion_env%kernel, dispersion_env%d2phi_dk2, &
1822 dispersion_env%d2y_dx2, dispersion_env%dftd_section)
1823 NULLIFY (dispersion_env%sab_vdw, dispersion_env%sab_cn)
1824 IF (se_control%dispersion)
THEN
1827 dispersion_env%eps_cn = se_control%epscn
1828 dispersion_env%s6 = se_control%sd3(1)
1829 dispersion_env%sr6 = se_control%sd3(2)
1830 dispersion_env%s8 = se_control%sd3(3)
1831 dispersion_env%domol = .false.
1832 dispersion_env%kgc8 = 0._dp
1833 dispersion_env%rc_disp = se_control%rcdisp
1834 dispersion_env%exp_pre = 0._dp
1835 dispersion_env%scaling = 0._dp
1836 dispersion_env%nd3_exclude_pair = 0
1837 dispersion_env%parameter_file_name = se_control%dispersion_parameter_file
1842 CALL set_qs_env(qs_env, dispersion_env=dispersion_env)
1846 IF (dft_control%qs_control%method_id ==
do_method_gpw .OR. &
1853 NULLIFY (xc_section)
1864 nelectron = nelectron - dft_control%charge
1866 IF (dft_control%multiplicity == 0)
THEN
1867 IF (
modulo(nelectron, 2) == 0)
THEN
1868 dft_control%multiplicity = 1
1870 dft_control%multiplicity = 2
1874 multiplicity = dft_control%multiplicity
1876 IF ((dft_control%nspins < 1) .OR. (dft_control%nspins > 2))
THEN
1877 cpabort(
"nspins should be 1 or 2 for the time being ...")
1880 IF ((
modulo(nelectron, 2) /= 0) .AND. (dft_control%nspins == 1))
THEN
1881 IF (.NOT. dft_control%qs_control%ofgpw .AND. .NOT. dft_control%smear)
THEN
1882 cpabort(
"Use the LSD option for an odd number of electrons")
1887 IF (dft_control%do_xas_calculation)
THEN
1888 IF (dft_control%nspins == 1)
THEN
1889 cpabort(
"Use the LSD option for XAS with transition potential")
1899 IF (dft_control%qs_control%ofgpw)
THEN
1901 IF (dft_control%nspins == 1)
THEN
1903 nelectron_spin(1) = nelectron
1904 nelectron_spin(2) = 0
1908 IF (
modulo(nelectron + multiplicity - 1, 2) /= 0)
THEN
1909 cpabort(
"LSD: try to use a different multiplicity")
1911 nelectron_spin(1) = (nelectron + multiplicity - 1)/2
1912 nelectron_spin(2) = (nelectron - multiplicity + 1)/2
1913 IF (nelectron_spin(1) < 0)
THEN
1914 cpabort(
"LSD: too few electrons for this multiplicity")
1916 maxocc = maxval(nelectron_spin)
1917 n_mo(1) = min(nelectron_spin(1), 1)
1918 n_mo(2) = min(nelectron_spin(2), 1)
1923 IF (dft_control%nspins == 1)
THEN
1925 nelectron_spin(1) = nelectron
1926 nelectron_spin(2) = 0
1927 IF (
modulo(nelectron, 2) == 0)
THEN
1928 n_mo(1) = nelectron/2
1930 n_mo(1) = int(nelectron/2._dp) + 1
1938 IF (
modulo(nelectron + multiplicity - 1, 2) /= 0)
THEN
1939 cpabort(
"LSD: try to use a different multiplicity")
1942 nelectron_spin(1) = (nelectron + multiplicity - 1)/2
1943 nelectron_spin(2) = (nelectron - multiplicity + 1)/2
1945 IF (nelectron_spin(2) < 0)
THEN
1946 cpabort(
"LSD: too few electrons for this multiplicity")
1949 n_mo(1) = nelectron_spin(1)
1950 n_mo(2) = nelectron_spin(2)
1959 qs_env%total_zeff_corr = total_zeff_corr
1963 nelectron_total=nelectron, &
1964 nelectron_spin=nelectron_spin)
1969 cpassert(
ASSOCIATED(mo_index_range))
1970 IF (all(mo_index_range > 0))
THEN
1971 IF (mo_index_range(1) > mo_index_range(2))
THEN
1972 CALL cp_abort(__location__, &
1973 "The upper orbital index ("// &
1975 ") of the MO_INDEX_RANGE should be equal or larger "// &
1976 "than the lower orbital index ("// &
1981 IF (.NOT. scf_control%use_ot)
THEN
1982 scf_control%added_mos(1) = min(max(scf_control%added_mos(1), &
1983 mo_index_range(2) - n_mo(1)), &
1985 IF (dft_control%nspins == 2)
THEN
1986 scf_control%added_mos(2) = min(max(scf_control%added_mos(2), &
1987 mo_index_range(2) - n_mo(2)), &
1991 ELSE IF (mo_index_range(2) < 0)
THEN
1992 IF (.NOT. scf_control%use_ot)
THEN
1994 scf_control%added_mos(1) = n_ao - n_mo(1)
1995 IF (dft_control%nspins == 2)
THEN
1997 scf_control%added_mos(2) = n_ao - n_mo(2)
2002 IF (dft_control%nspins == 2)
THEN
2004 IF (scf_control%added_mos(2) < 0)
THEN
2005 n_mo_add = n_ao - n_mo(2)
2006 ELSE IF (scf_control%added_mos(2) > 0)
THEN
2007 n_mo_add = scf_control%added_mos(2)
2009 n_mo_add = scf_control%added_mos(1)
2011 IF (n_mo_add > n_ao - n_mo(2))
THEN
2012 cpwarn(
"More ADDED_MOs requested for beta spin than available.")
2014 scf_control%added_mos(2) = min(n_mo_add, n_ao - n_mo(2))
2015 n_mo(2) = n_mo(2) + scf_control%added_mos(2)
2026 IF (dft_control%qs_control%xtb_control%do_tblite)
THEN
2027 scf_control%added_mos(1) = n_ao - n_mo(1)
2028 ELSE IF (scf_control%added_mos(1) < 0)
THEN
2029 scf_control%added_mos(1) = n_ao - n_mo(1)
2030 ELSE IF (scf_control%added_mos(1) > n_ao - n_mo(1))
THEN
2031 CALL cp_warn(__location__, &
2032 "More added MOs requested than available. "// &
2033 "The full set of unoccupied MOs will be used. "// &
2034 "Use 'ADDED_MOS -1' to always use all available MOs "// &
2035 "and to get rid of this warning.")
2037 scf_control%added_mos(1) = min(scf_control%added_mos(1), n_ao - n_mo(1))
2038 n_mo(1) = n_mo(1) + scf_control%added_mos(1)
2040 IF (dft_control%nspins == 2)
THEN
2041 IF (n_mo(2) > n_mo(1)) &
2042 CALL cp_warn(__location__, &
2043 "More beta than alpha MOs requested. "// &
2044 "The number of beta MOs will be reduced to the number alpha MOs.")
2045 n_mo(2) = min(n_mo(1), n_mo(2))
2046 cpassert(n_mo(1) >= nelectron_spin(1))
2047 cpassert(n_mo(2) >= nelectron_spin(2))
2051 CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints)
2052 IF (do_kpoints .AND. dft_control%nspins == 2)
THEN
2054 IF (n_mo(2) /= n_mo(1)) &
2055 CALL cp_warn(__location__, &
2056 "Kpoints: Different number of MOs requested. "// &
2057 "The number of beta MOs will be set to the number alpha MOs.")
2059 cpassert(n_mo(1) >= nelectron_spin(1))
2060 cpassert(n_mo(2) >= nelectron_spin(2))
2064 IF (scf_control%smear%do_smear)
THEN
2065 IF (scf_control%added_mos(1) == 0)
THEN
2066 cpabort(
"Extra MOs (ADDED_MOS) are required for smearing")
2072 "PRINT%MO/CARTESIAN"), &
2074 (scf_control%level_shift /= 0.0_dp) .OR. &
2075 (scf_control%diagonalization%eps_jacobi /= 0.0_dp) .OR. &
2076 (dft_control%roks .AND. (.NOT. scf_control%use_ot)))
THEN
2081 IF (dft_control%roks .AND. (.NOT. scf_control%use_ot))
THEN
2083 cpwarn(
"General ROKS scheme is not yet tested!")
2085 IF (scf_control%smear%do_smear)
THEN
2086 CALL cp_abort(__location__, &
2087 "The options ROKS and SMEAR are not compatible. "// &
2088 "Try UKS instead of ROKS")
2091 IF (dft_control%low_spin_roks)
THEN
2092 SELECT CASE (dft_control%qs_control%method_id)
2095 CALL cp_abort(__location__, &
2096 "xTB/DFTB methods are not compatible with low spin ROKS.")
2099 CALL cp_abort(__location__, &
2100 "SE methods are not compatible with low spin ROKS.")
2108 IF (dft_control%restricted .AND. (output_unit > 0))
THEN
2110 WRITE (output_unit, *)
""
2111 WRITE (output_unit, *)
" **************************************"
2112 WRITE (output_unit, *)
" restricted calculation cutting corners"
2113 WRITE (output_unit, *)
" experimental feature, check code "
2114 WRITE (output_unit, *)
" **************************************"
2118 IF (dft_control%qs_control%do_ls_scf)
THEN
2121 ALLOCATE (mos(dft_control%nspins))
2122 DO ispin = 1, dft_control%nspins
2126 nelectron=nelectron_spin(ispin), &
2127 n_el_f=real(nelectron_spin(ispin),
dp), &
2129 flexible_electron_count=dft_control%relax_multiplicity)
2136 IF (dft_control%switch_surf_dip)
THEN
2137 ALLOCATE (mos_last_converged(dft_control%nspins))
2138 DO ispin = 1, dft_control%nspins
2142 nelectron=nelectron_spin(ispin), &
2143 n_el_f=real(nelectron_spin(ispin),
dp), &
2145 flexible_electron_count=dft_control%relax_multiplicity)
2147 CALL set_qs_env(qs_env, mos_last_converged=mos_last_converged)
2150 IF (.NOT. be_silent)
THEN
2155 IF (dft_control%qs_control%method_id ==
do_method_gpw .OR. &
2162 (.NOT. dft_control%qs_control%xtb_control%do_tblite)) .OR. &
2164 CALL get_qs_env(qs_env, dispersion_env=dispersion_env)
2172 IF (dft_control%do_admm)
THEN
2177 IF (dft_control%do_xas_calculation)
THEN
2192 CALL write_total_numbers(qs_kind_set, particle_set, qs_env%input)
2207 IF ((.NOT. dft_control%qs_control%do_ls_scf) .AND. &
2208 (.NOT. dft_control%qs_control%do_almo_scf))
THEN
2218 CALL get_qs_env(qs_env=qs_env, pw_env=pw_env, rho0_mpole=rho0_mpole)
2222 IF (output_unit > 0)
CALL m_flush(output_unit)
2223 CALL timestop(handle)
2225 END SUBROUTINE qs_init_subsys
2235 SUBROUTINE write_total_numbers(qs_kind_set, particle_set, force_env_section)
2237 TYPE(qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
2238 TYPE(particle_type),
DIMENSION(:),
POINTER :: particle_set
2239 TYPE(section_vals_type),
POINTER :: force_env_section
2241 INTEGER :: maxlgto, maxlppl, maxlppnl, natom, &
2242 natom_q, ncgf, nkind, nkind_q, npgf, &
2243 nset, nsgf, nshell, output_unit
2244 TYPE(cp_logger_type),
POINTER :: logger
2247 logger => cp_get_default_logger()
2248 output_unit = cp_print_key_unit_nr(logger, force_env_section,
"PRINT%TOTAL_NUMBERS", &
2251 IF (output_unit > 0)
THEN
2252 natom =
SIZE(particle_set)
2253 nkind =
SIZE(qs_kind_set)
2255 CALL get_qs_kind_set(qs_kind_set, &
2265 WRITE (unit=output_unit, fmt=
"(/,/,T2,A)") &
2266 "TOTAL NUMBERS AND MAXIMUM NUMBERS"
2268 IF (nset + npgf + ncgf > 0)
THEN
2269 WRITE (unit=output_unit, fmt=
"(/,T3,A,(T30,A,T71,I10))") &
2270 "Total number of", &
2271 "- Atomic kinds: ", nkind, &
2272 "- Atoms: ", natom, &
2273 "- Shell sets: ", nset, &
2274 "- Shells: ", nshell, &
2275 "- Primitive Cartesian functions: ", npgf, &
2276 "- Cartesian basis functions: ", ncgf, &
2277 "- Spherical basis functions: ", nsgf
2278 ELSE IF (nshell + nsgf > 0)
THEN
2279 WRITE (unit=output_unit, fmt=
"(/,T3,A,(T30,A,T71,I10))") &
2280 "Total number of", &
2281 "- Atomic kinds: ", nkind, &
2282 "- Atoms: ", natom, &
2283 "- Shells: ", nshell, &
2284 "- Spherical basis functions: ", nsgf
2286 WRITE (unit=output_unit, fmt=
"(/,T3,A,(T30,A,T71,I10))") &
2287 "Total number of", &
2288 "- Atomic kinds: ", nkind, &
2292 IF ((maxlppl > -1) .AND. (maxlppnl > -1))
THEN
2293 WRITE (unit=output_unit, fmt=
"(/,T3,A,(T30,A,T75,I6))") &
2294 "Maximum angular momentum of the", &
2295 "- Orbital basis functions: ", maxlgto, &
2296 "- Local part of the GTH pseudopotential: ", maxlppl, &
2297 "- Non-local part of the GTH pseudopotential: ", maxlppnl
2298 ELSEIF (maxlppl > -1)
THEN
2299 WRITE (unit=output_unit, fmt=
"(/,T3,A,(T30,A,T75,I6))") &
2300 "Maximum angular momentum of the", &
2301 "- Orbital basis functions: ", maxlgto, &
2302 "- Local part of the GTH pseudopotential: ", maxlppl
2304 WRITE (unit=output_unit, fmt=
"(/,T3,A,T75,I6)") &
2305 "Maximum angular momentum of the orbital basis functions: ", maxlgto
2309 CALL get_qs_kind_set(qs_kind_set, &
2316 basis_type=
"LRI_AUX")
2317 IF (nset + npgf + ncgf > 0)
THEN
2318 WRITE (unit=output_unit, fmt=
"(/,T3,A,/,T3,A,(T30,A,T71,I10))") &
2319 "LRI_AUX Basis: ", &
2320 "Total number of", &
2321 "- Shell sets: ", nset, &
2322 "- Shells: ", nshell, &
2323 "- Primitive Cartesian functions: ", npgf, &
2324 "- Cartesian basis functions: ", ncgf, &
2325 "- Spherical basis functions: ", nsgf
2326 WRITE (unit=output_unit, fmt=
"(T30,A,T75,I6)") &
2327 " Maximum angular momentum ", maxlgto
2331 CALL get_qs_kind_set(qs_kind_set, &
2338 basis_type=
"RI_HXC")
2339 IF (nset + npgf + ncgf > 0)
THEN
2340 WRITE (unit=output_unit, fmt=
"(/,T3,A,/,T3,A,(T30,A,T71,I10))") &
2342 "Total number of", &
2343 "- Shell sets: ", nset, &
2344 "- Shells: ", nshell, &
2345 "- Primitive Cartesian functions: ", npgf, &
2346 "- Cartesian basis functions: ", ncgf, &
2347 "- Spherical basis functions: ", nsgf
2348 WRITE (unit=output_unit, fmt=
"(T30,A,T75,I6)") &
2349 " Maximum angular momentum ", maxlgto
2353 CALL get_qs_kind_set(qs_kind_set, &
2360 basis_type=
"AUX_FIT")
2361 IF (nset + npgf + ncgf > 0)
THEN
2362 WRITE (unit=output_unit, fmt=
"(/,T3,A,/,T3,A,(T30,A,T71,I10))") &
2363 "AUX_FIT ADMM-Basis: ", &
2364 "Total number of", &
2365 "- Shell sets: ", nset, &
2366 "- Shells: ", nshell, &
2367 "- Primitive Cartesian functions: ", npgf, &
2368 "- Cartesian basis functions: ", ncgf, &
2369 "- Spherical basis functions: ", nsgf
2370 WRITE (unit=output_unit, fmt=
"(T30,A,T75,I6)") &
2371 " Maximum angular momentum ", maxlgto
2375 CALL get_qs_kind_set(qs_kind_set, &
2385 IF (nset + npgf + ncgf > 0)
THEN
2386 WRITE (unit=output_unit, fmt=
"(/,T3,A,/,T3,A,(T30,A,T71,I10))") &
2387 "Nuclear Basis: ", &
2388 "Total number of", &
2389 "- Quantum atomic kinds: ", nkind_q, &
2390 "- Quantum atoms: ", natom_q, &
2391 "- Shell sets: ", nset, &
2392 "- Shells: ", nshell, &
2393 "- Primitive Cartesian functions: ", npgf, &
2394 "- Cartesian basis functions: ", ncgf, &
2395 "- Spherical basis functions: ", nsgf
2396 WRITE (unit=output_unit, fmt=
"(T30,A,T75,I6)") &
2397 " Maximum angular momentum ", maxlgto
2401 CALL cp_print_key_finished_output(output_unit, logger, force_env_section, &
2402 "PRINT%TOTAL_NUMBERS")
2404 END SUBROUTINE write_total_numbers
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
subroutine, public almo_scf_env_create(qs_env)
Creation and basic initialization of the almo environment.
calculate the orbitals for a given atomic kind type
subroutine, public calculate_atomic_relkin(atomic_kind, qs_kind, rel_control, rtmat)
...
Define the atomic kind types and their sub types.
Automatic generation of auxiliary basis sets of different kind.
subroutine, public create_lri_aux_basis_set(lri_aux_basis_set, qs_kind, basis_cntrl, exact_1c_terms, tda_kernel)
Create a LRI_AUX basis set using some heuristics.
subroutine, public create_ri_aux_basis_set(ri_aux_basis_set, qs_kind, basis_cntrl, basis_type, basis_sort)
Create a RI_AUX basis set using some heuristics.
subroutine, public add_basis_set_to_container(container, basis_set, basis_set_type)
...
integer, parameter, public basis_sort_zet
subroutine, public deallocate_gto_basis_set(gto_basis_set)
...
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 cp2kqs2020
integer, save, public iannuzzi2007
integer, save, public iannuzzi2006
Handles all functions related to the CELL.
methods related to the blacs parallel environment
subroutine, public cp_blacs_env_release(blacs_env)
releases the given blacs_env
subroutine, public cp_blacs_env_create(blacs_env, para_env, blacs_grid_layout, blacs_repeatable, row_major, grid_2d)
allocates and initializes a type that represent a blacs context
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
Utilities to set up the control types.
subroutine, public write_qs_control(qs_control, dft_section)
Purpose: Write the QS control parameters to the output unit.
subroutine, public read_rixs_control(rixs_control, rixs_section, qs_control)
Reads the input and stores in the rixs_control_type.
subroutine, public read_tddfpt2_control(t_control, t_section, qs_control)
Read TDDFPT-related input parameters.
subroutine, public read_dft_control(dft_control, dft_section)
...
subroutine, public read_qs_section(qs_control, qs_section)
...
subroutine, public write_admm_control(admm_control, dft_section)
Write the ADMM control parameters to the output unit.
subroutine, public write_dft_control(dft_control, dft_section)
Write the DFT control parameters to the output unit.
subroutine, public read_ddapc_section(qs_control, qs_section, ddapc_restraint_section)
reads the input parameters needed for ddapc.
subroutine, public read_mgrid_section(qs_control, dft_section)
...
contains information regarding the decoupling/recoupling method of Bloechl
subroutine, public cp_ddapc_ewald_create(cp_ddapc_ewald, qmmm_decoupl, qm_cell, force_env_section, subsys_section, para_env)
...
various routines to log and control the output. The idea is that decisions about where to log should ...
integer function, public cp_logger_get_default_io_unit(logger)
returns the unit nr for the ionode (-1 on all other processors) skips as well checks if the procs cal...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
integer, parameter, public cp_p_file
integer function, public cp_print_key_should_output(iteration_info, basis_section, print_key_path, used_print_key, first_time)
returns what should be done with the given property if btest(res,cp_p_store) then the property should...
types that represent a subsys, i.e. a part of the system
subroutine, public write_symmetry(particle_set, cell, input_section)
Write symmetry information to output.
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
subroutine, public distribution_1d_release(distribution_1d)
releases the given distribution_1d
Distribution methods for atoms, particles, or molecules.
subroutine, public distribute_molecules_1d(atomic_kind_set, particle_set, local_particles, molecule_kind_set, molecule_set, local_molecules, force_env_section, prev_molecule_kind_set, prev_local_molecules)
Distribute molecules and particles.
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.
Definition and initialisation of the et_coupling data type.
subroutine, public et_coupling_create(et_coupling)
...
subroutine, public ewald_env_set(ewald_env, ewald_type, alpha, epsilon, eps_pol, gmax, ns_max, precs, o_spline, para_env, poisson_section, interaction_cutoffs, cell_hmat)
Purpose: Set the EWALD environment.
subroutine, public ewald_env_create(ewald_env, para_env)
allocates and intitializes a ewald_env
subroutine, public read_ewald_section(ewald_env, ewald_section)
Purpose: read the EWALD section.
subroutine, public read_ewald_section_tb(ewald_env, ewald_section, hmat, silent, pset, cell_periodic)
Purpose: read the EWALD section for TB methods.
subroutine, public ewald_env_get(ewald_env, ewald_type, alpha, eps_pol, epsilon, gmax, ns_max, o_spline, group, para_env, poisson_section, precs, rcut, do_multipoles, max_multipole, do_ipol, max_ipol_iter, interaction_cutoffs, cell_hmat)
Purpose: Get the EWALD environment.
subroutine, public ewald_pw_grid_update(ewald_pw, ewald_env, cell_hmat)
Rescales pw_grids for given box, if necessary.
subroutine, public ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section)
creates the structure ewald_pw_type
Types for excited states potential energies.
subroutine, public exstate_create(ex_env, excited_state, dft_section)
Allocates and intitializes exstate_env.
Definition of the atomic potential types.
subroutine, public fist_nonbond_env_create(fist_nonbond_env, atomic_kind_set, potparm14, potparm, do_nonbonded, do_electrostatics, verlet_skin, ewald_rcut, ei_scale14, vdw_scale14, shift_cutoff)
allocates and intitializes a fist_nonbond_env
Calculation of the incomplete Gamma function F_n(t) for multi-center integrals over Cartesian Gaussia...
subroutine, public init_md_ftable(nmax)
Initialize a table of F_n(t) values in the range 0 <= t <= 12 with a stepsize of 0....
Define type storing the global information of a run. Keep the amount of stored data small....
subroutine, public init_coulomb_local(hartree_local, natom)
...
Types and set/get functions for HFX.
subroutine, public hfx_create(x_data, para_env, hfx_section, atomic_kind_set, qs_kind_set, particle_set, dft_control, cell, orb_basis, ri_basis, nelectron_total, nkp_grid)
This routine allocates and initializes all types in hfx_data
subroutine, public compare_hfx_sections(hfx_section1, hfx_section2, is_identical, same_except_frac)
Compares the non-technical parts of two HFX input section and check whether they are the same Ignore ...
Routines for a Kim-Gordon-like partitioning into molecular subunits.
subroutine, public kg_env_create(qs_env, kg_env, qs_kind_set, input)
Allocates and intitializes kg_env.
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
Routines needed for kpoint calculation.
subroutine, public kpoint_initialize_mos(kpoint, mos, added_mos, for_aux_fit)
Initialize a set of MOs and density matrix for each kpoint (kpoint group)
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 set_kpoint_info(kpoint, kp_scheme, nkp_grid, kp_shift, symmetry, verbose, full_grid, use_real_wfn, eps_geo, parallel_group_size, kp_range, nkp, xkp, wkp, para_env, blacs_env_all, para_env_kp, para_env_inter_kp, blacs_env, kp_env, kp_aux_env, mpools, iogrp, nkp_groups, kp_dist, cell_to_index, index_to_cell, sab_nl, sab_nl_nosym, inversion_symmetry_only, symmetry_backend, symmetry_reduction_method)
Set information in a kpoint environment.
subroutine, public kpoint_create(kpoint)
Create a kpoint environment.
subroutine, public get_kpoint_info(kpoint, kp_scheme, nkp_grid, kp_shift, symmetry, verbose, full_grid, use_real_wfn, eps_geo, parallel_group_size, kp_range, nkp, xkp, wkp, para_env, blacs_env_all, para_env_kp, para_env_inter_kp, blacs_env, kp_env, kp_aux_env, mpools, iogrp, nkp_groups, kp_dist, cell_to_index, index_to_cell, sab_nl, sab_nl_nosym, inversion_symmetry_only, symmetry_backend, symmetry_reduction_method)
Retrieve information from a kpoint environment.
initializes the environment for lri lri : local resolution of the identity
subroutine, public lri_env_init(lri_env, lri_section)
initializes the lri env
subroutine, public lri_env_basis(ri_type, qs_env, lri_env, qs_kind_set)
initializes the lri env
contains the types and subroutines for dealing with the lri_env lri : local resolution of the identit...
Machine interface based on Fortran 2003 and POSIX.
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
Interface to the message passing library MPI.
Define the molecule kind structure types and the corresponding functionality.
subroutine, public write_molecule_kind_set(molecule_kind_set, subsys_section)
Write a moleculeatomic kind set data set to the output unit.
Define the data structure for the molecule information.
Types needed for MP2 calculations.
subroutine, public read_mp2_section(input, mp2_env)
...
Types needed for MP2 calculations.
subroutine, public mp2_env_create(mp2_env)
...
Multipole structure: for multipole (fixed and induced) in FF based MD.
integer, parameter, public do_multipole_none
Provides Cartesian and spherical orbital pointers and indices.
subroutine, public init_orbital_pointers(maxl)
Initialize or update the orbital pointers.
Define methods related to particle_type.
subroutine, public write_qs_particle_coordinates(particle_set, qs_kind_set, subsys_section, label)
Write the atomic coordinates to the output unit.
subroutine, public write_structure_data(particle_set, cell, input_section)
Write structure data requested by a separate structure data input section to the output unit....
subroutine, public write_particle_distances(particle_set, cell, subsys_section)
Write the matrix of the particle distances to the output unit.
Define the data structure for the particle information.
Definition of physical constants:
real(kind=dp), parameter, public kelvin
container for various plainwaves related things
subroutine, public qs_basis_rotation(qs_env, kpoints, basis_type)
Construct basis set rotation matrices.
subroutine, public qs_dftb_param_init(atomic_kind_set, qs_kind_set, dftb_control, dftb_potential, subsys_section, para_env)
...
Definition of the DFTB parameter types.
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 qs_write_dispersion(qs_env, dispersion_env, ounit)
...
subroutine, public allocate_qs_energy(qs_energy)
Allocate and/or initialise a Quickstep energy data structure.
qs_environment methods that use many other modules
subroutine, public qs_env_setup(qs_env)
initializes various components of the qs_env, that need only atomic_kind_set, cell,...
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 set_qs_env(qs_env, super_cell, mos, qmmm, qmmm_periodic, mimic, ewald_env, ewald_pw, mpools, rho_external, external_vxc, mask, scf_control, rel_control, qs_charges, ks_env, ks_qmmm_env, wf_history, scf_env, active_space, input, oce, rho_atom_set, rho0_atom_set, rho0_mpole, run_rtp, rtp, rhoz_set, rhoz_tot, ecoul_1c, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, efield, rhoz_cneo_set, linres_control, xas_env, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, ls_scf_env, do_transport, transport_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, harris_env, gcp_env, mp2_env, bs_env, kg_env, force, kpoints, wanniercentres, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Set the QUICKSTEP environment.
subroutine, public qs_init(qs_env, para_env, root_section, globenv, cp_subsys, kpoint_env, qmmm, qmmm_env_qm, force_env_section, subsys_section, use_motion_section, silent)
Read the input and the database files for the setup of the QUICKSTEP environment.
Definition of gCP types for DFT calculations.
Set disperson types for DFT calculations.
subroutine, public qs_gcp_env_set(gcp_env, xc_section)
...
subroutine, public qs_gcp_init(qs_env, gcp_env)
...
Types needed for a for a Harris model calculation.
subroutine, public harris_rhoin_init(rhoin, basis_type, qs_kind_set, atomic_kind_set, local_particles, nspin)
...
Harris method environment setup and handling.
subroutine, public harris_write_input(harris_env)
Print out the Harris method input section.
subroutine, public harris_env_create(qs_env, harris_env, harris_section)
Allocates and intitializes harris_env.
Calculate the interaction radii for the operator matrix calculation.
subroutine, public write_pgf_orb_radii(basis, atomic_kind_set, qs_kind_set, subsys_section)
Write the orbital basis function radii to the output unit.
subroutine, public write_paw_radii(atomic_kind_set, qs_kind_set, subsys_section)
Write the radii of the one center projector.
subroutine, public write_ppnl_radii(atomic_kind_set, qs_kind_set, subsys_section)
Write the radii of the projector functions of the Goedecker pseudopotential (GTH, non-local part) to ...
subroutine, public init_se_nlradius(se_control, atomic_kind_set, qs_kind_set, subsys_section)
...
subroutine, public write_ppl_radii(atomic_kind_set, qs_kind_set, subsys_section)
Write the radii of the exponential functions of the Goedecker pseudopotential (GTH,...
subroutine, public init_interaction_radii(qs_control, qs_kind_set)
Initialize all the atomic kind radii for a given threshold value.
subroutine, public write_core_charge_radii(atomic_kind_set, qs_kind_set, subsys_section)
Write the radii of the core charge distributions to the output unit.
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.
subroutine, public init_gapw_nlcc(qs_kind_set)
...
subroutine, public init_cneo_basis_set(qs_kind_set, qs_control)
...
subroutine, public set_qs_kind(qs_kind, paw_atom, ghost, floating, hard_radius, hard0_radius, covalent_radius, vdw_radius, lmax_rho0, zeff, no_optimize, dispersion, u_minus_j, reltmat, dftb_parameter, xtb_parameter, elec_conf, pao_basis_size)
Set the components of an atomic kind data set.
subroutine, public write_qs_kind_set(qs_kind_set, subsys_section)
Write an atomic kind set data set to the output unit.
subroutine, public init_gapw_basis_set(qs_kind_set, qs_control, force_env_section, modify_qs_control)
...
subroutine, public init_qs_kind_set(qs_kind_set)
Initialise an atomic kind set data set.
subroutine, public check_qs_kind_set(qs_kind_set, dft_control, subsys_section)
...
subroutine, public write_gto_basis_sets(qs_kind_set, subsys_section)
Write all the GTO basis sets of an atomic kind set to the output unit (for the printing of the unnorm...
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)
...
subroutine, public qs_ks_env_create(ks_env)
Allocates a new instance of ks_env.
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 rho0_s_grid_create(pw_env, rho0_mpole)
...
subroutine, public init_rho0(local_rho_set, qs_env, gapw_control, zcore)
...
subroutine, public init_rho_atom(rho_atom_set, atomic_kind_set, qs_kind_set, dft_control, para_env)
...
Routines that work on qs_subsys_type.
subroutine, public qs_subsys_create(subsys, para_env, root_section, force_env_section, subsys_section, use_motion_section, cp_subsys, elkind, silent)
Creates a qs_subsys. Optionally an existsing cp_subsys is used.
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 qs_subsys_set(subsys, cp_subsys, local_particles, local_molecules, cell, cell_ref, use_ref_cell, energy, force, qs_kind_set, nelectron_total, nelectron_spin)
...
Storage of past states of the qs_env. Methods to interpolate (or actually normally extrapolate) the n...
subroutine, public wfi_create_for_kp(wf_history)
Adapts wf_history storage flags for k-point calculations. For ASPC, switches from Gamma WFN storage t...
subroutine, public wfi_create(wf_history, interpolation_method_nr, extrapolation_order, has_unit_metric)
...
interpolate the wavefunctions to speed up the convergence when doing MD
subroutine, public wfi_release(wf_history)
releases a wf_history of a wavefunction (see doc/ReferenceCounting.html)
parameters that control a relativistic calculation
subroutine, public rel_c_create(rel_control)
allocates and initializes an rel control object with the default values
subroutine, public rel_c_read_parameters(rel_control, dft_section)
reads the parameters of the relativistic section into the given rel_control
parameters that control an scf iteration
subroutine, public scf_c_read_parameters(scf_control, inp_section)
reads the parameters of the scf section into the given scf_control
subroutine, public scf_c_write_parameters(scf_control, dft_section)
writes out the scf parameters
subroutine, public scf_c_create(scf_control)
allocates and initializes an scf control object with the default values
Methods for handling the 1/R^3 residual integral part.
subroutine, public semi_empirical_expns3_setup(qs_kind_set, se_control, method_id)
Setup the quantity necessary to handle the slowly convergent residual integral term 1/R^3.
Arrays of parameters used in the semi-empirical calculations \References Everywhere in this module TC...
subroutine, public init_se_intd_array()
Initialize all arrays used for the evaluation of the integrals.
Setup and Methods for semi-empirical multipole types.
subroutine, public nddo_mpole_setup(nddo_mpole, natom)
Setup NDDO multipole type.
Definition of the semi empirical multipole integral expansions types.
Type to store integrals for semi-empirical calculations.
subroutine, public semi_empirical_si_create(store_int_env, se_section, compression)
Allocate semi-empirical store integrals type.
Definition of the semi empirical parameter types.
subroutine, public se_taper_create(se_taper, integral_screening, do_ewald, taper_cou, range_cou, taper_exc, range_exc, taper_scr, range_scr, taper_lrc, range_lrc)
Creates the taper type used in SE calculations.
Working with the semi empirical parameter types.
subroutine, public se_cutoff_compatible(se_control, se_section, cell, output_unit)
Reset cutoffs trying to be somehow a bit smarter.
subroutine, public tb_set_calculator(tb, typ)
...
subroutine, public tb_init_wf(tb)
initialize wavefunction ...
subroutine, public tb_init_geometry(qs_env, tb)
intialize geometry objects ...
subroutine, public tb_get_basis(tb, gto_basis_set, element_symbol, param, occ)
...
routines for DFT+NEGF calculations (coupling with the quantum transport code OMEN)
subroutine, public transport_env_create(qs_env)
creates the transport environment
subroutine, public xtb_parameters_set(param)
Read atom parameters for xTB Hamiltonian from input file.
subroutine, public xtb_parameters_init(param, gfn_type, element_symbol, parameter_file_path, parameter_file_name, para_env)
...
subroutine, public init_xtb_basis(param, gto_basis_set, ngauss)
...
xTB (repulsive) pair potentials Reference: Stefan Grimme, Christoph Bannwarth, Philip Shushkov JCTC 1...
subroutine, public xtb_pp_radius(qs_kind_set, ppradius, eps_pair, kfparam)
...
Definition of the xTB parameter types.
subroutine, public allocate_xtb_atom_param(xtb_parameter)
...
subroutine, public set_xtb_atom_param(xtb_parameter, aname, typ, defined, z, zeff, natorb, lmax, nao, lao, rcut, rcov, kx, eta, xgamma, alpha, zneff, nshell, nval, lval, kpoly, kappa, hen, zeta, xi, kappa0, alpg, electronegativity, occupation, chmax, en, kqat2, kcn, kq)
...
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...
type of a logger, at the moment it contains just a print level starting at which level it should be l...
represents a system: atoms, molecules, their pos,vel,...
structure to store local (to a processor) ordered lists of integers.
Contains information on the energy correction functional for KG.
to build arrays of pointers
Contains information on the excited states energy.
contains the initially parsed file and the initial parallel environment
Contains information about kpoints.
stores all the informations relevant to an mpi environment
contained for different pw related things
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 track of the previous wavefunctions and can extrapolate them for the next step of md
contains the parameters needed by a relativistic calculation
Global Multipolar NDDO information type.
Semi-empirical store integrals type.
Taper type use in semi-empirical calculations.