Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \par History
10 : !> - Merged with the Quickstep MODULE method_specification (17.01.2002,MK)
11 : !> - USE statements cleaned, added
12 : !> (25.09.2002,MK)
13 : !> - Added more LSD structure (01.2003,Joost VandeVondele)
14 : !> - New molecule data types introduced (Sep. 2003,MK)
15 : !> - Cleaning; getting rid of pnode (02.10.2003,MK)
16 : !> - Sub-system setup added (08.10.2003,MK)
17 : !> \author MK (18.05.2000)
18 : ! **************************************************************************************************
19 : MODULE qs_environment
20 : USE almo_scf_env_methods, ONLY: almo_scf_env_create
21 : USE atom_kind_orbitals, ONLY: calculate_atomic_relkin
22 : USE atomic_kind_types, ONLY: atomic_kind_type
23 : USE auto_basis, ONLY: create_lri_aux_basis_set,&
24 : create_ri_aux_basis_set
25 : USE basis_set_container_types, ONLY: add_basis_set_to_container
26 : USE basis_set_types, ONLY: basis_sort_zet,&
27 : create_primitive_basis_set,&
28 : deallocate_gto_basis_set,&
29 : gto_basis_set_type
30 : USE bibliography, ONLY: Iannuzzi2006,&
31 : Iannuzzi2007,&
32 : cite_reference,&
33 : cp2kqs2020
34 : USE cell_types, ONLY: cell_type
35 : USE cp_blacs_env, ONLY: cp_blacs_env_create,&
36 : cp_blacs_env_release,&
37 : cp_blacs_env_type
38 : USE cp_control_types, ONLY: dft_control_type,&
39 : dftb_control_type,&
40 : gapw_control_type,&
41 : qs_control_type,&
42 : semi_empirical_control_type,&
43 : xtb_control_type
44 : USE cp_control_utils, ONLY: &
45 : read_ddapc_section, read_dft_control, read_mgrid_section, read_qs_section, &
46 : read_rixs_control, read_tddfpt2_control, write_admm_control, write_dft_control, &
47 : write_qs_control
48 : USE cp_ddapc_types, ONLY: cp_ddapc_ewald_create
49 : USE cp_log_handling, ONLY: cp_get_default_logger,&
50 : cp_logger_get_default_io_unit,&
51 : cp_logger_type
52 : USE cp_output_handling, ONLY: cp_p_file,&
53 : cp_print_key_finished_output,&
54 : cp_print_key_should_output,&
55 : cp_print_key_unit_nr
56 : USE cp_subsys_types, ONLY: cp_subsys_type
57 : USE cp_symmetry, ONLY: write_symmetry
58 : USE distribution_1d_types, ONLY: distribution_1d_release,&
59 : distribution_1d_type
60 : USE distribution_methods, ONLY: distribute_molecules_1d
61 : USE ec_env_types, ONLY: energy_correction_type
62 : USE ec_environment, ONLY: ec_env_create,&
63 : ec_write_input
64 : USE et_coupling_types, ONLY: et_coupling_create
65 : USE ewald_environment_types, ONLY: ewald_env_create,&
66 : ewald_env_get,&
67 : ewald_env_set,&
68 : ewald_environment_type,&
69 : read_ewald_section,&
70 : read_ewald_section_tb
71 : USE ewald_pw_methods, ONLY: ewald_pw_grid_update
72 : USE ewald_pw_types, ONLY: ewald_pw_create,&
73 : ewald_pw_type
74 : USE exstates_types, ONLY: excited_energy_type,&
75 : exstate_create
76 : USE external_potential_types, ONLY: get_potential,&
77 : init_potential,&
78 : set_potential
79 : USE fist_nonbond_env_types, ONLY: fist_nonbond_env_create,&
80 : fist_nonbond_env_type
81 : USE gamma, ONLY: init_md_ftable
82 : USE global_types, ONLY: global_environment_type
83 : USE hartree_local_methods, ONLY: init_coulomb_local
84 : USE header, ONLY: dftb_header,&
85 : qs_header,&
86 : se_header,&
87 : tblite_header,&
88 : xtb_header
89 : USE hfx_types, ONLY: compare_hfx_sections,&
90 : hfx_create
91 : USE input_constants, ONLY: &
92 : dispersion_d2, dispersion_d3, dispersion_d3bj, do_et_ddapc, do_method_am1, do_method_dftb, &
93 : do_method_gapw, do_method_gapw_xc, do_method_gpw, do_method_lrigpw, do_method_mndo, &
94 : do_method_mndod, do_method_ofgpw, do_method_pdg, do_method_pm3, do_method_pm6, &
95 : do_method_pm6fm, do_method_pnnl, do_method_rigpw, do_method_rm1, do_method_xtb, &
96 : do_qmmm_gauss, do_qmmm_swave, general_roks, hden_atomic, kg_tnadd_embed_ri, rel_none, &
97 : rel_trans_atom, vdw_pairpot_dftd2, vdw_pairpot_dftd3, vdw_pairpot_dftd3bj, &
98 : vdw_pairpot_dftd4, wfi_aspc_nr, wfi_linear_ps_method_nr, wfi_linear_wf_method_nr, &
99 : wfi_ps_method_nr, wfi_use_guess_method_nr, xc_vdw_fun_none, xc_vdw_fun_nonloc, &
100 : xc_vdw_fun_pairpot, xtb_vdw_type_d3, xtb_vdw_type_d4, xtb_vdw_type_none
101 : USE input_section_types, ONLY: section_vals_get,&
102 : section_vals_get_subs_vals,&
103 : section_vals_type,&
104 : section_vals_val_get
105 : USE kg_environment, ONLY: kg_env_create
106 : USE kinds, ONLY: default_string_length,&
107 : dp
108 : USE kpoint_methods, ONLY: kpoint_env_initialize,&
109 : kpoint_initialize,&
110 : kpoint_initialize_mos
111 : USE kpoint_types, ONLY: get_kpoint_info,&
112 : kpoint_create,&
113 : kpoint_type,&
114 : read_kpoint_section,&
115 : write_kpoint_info
116 : USE lri_environment_init, ONLY: lri_env_basis,&
117 : lri_env_init
118 : USE lri_environment_types, ONLY: lri_environment_type
119 : USE machine, ONLY: m_flush
120 : USE mathconstants, ONLY: pi
121 : USE message_passing, ONLY: mp_para_env_type
122 : USE molecule_kind_types, ONLY: molecule_kind_type,&
123 : write_molecule_kind_set
124 : USE molecule_types, ONLY: molecule_type
125 : USE mp2_setup, ONLY: read_mp2_section
126 : USE mp2_types, ONLY: mp2_env_create,&
127 : mp2_type
128 : USE multipole_types, ONLY: do_multipole_none
129 : USE orbital_pointers, ONLY: init_orbital_pointers
130 : USE orbital_transformation_matrices, ONLY: init_spherical_harmonics
131 : USE particle_methods, ONLY: write_particle_distances,&
132 : write_qs_particle_coordinates,&
133 : write_structure_data
134 : USE particle_types, ONLY: particle_type
135 : USE pw_env_types, ONLY: pw_env_type
136 : USE qmmm_types_low, ONLY: qmmm_env_qm_type
137 : USE qs_basis_rotation_methods, ONLY: qs_basis_rotation
138 : USE qs_dftb_parameters, ONLY: qs_dftb_param_init
139 : USE qs_dftb_types, ONLY: qs_dftb_pairpot_type
140 : USE qs_dispersion_nonloc, ONLY: qs_dispersion_nonloc_init
141 : USE qs_dispersion_pairpot, ONLY: qs_dispersion_pairpot_init
142 : USE qs_dispersion_types, ONLY: qs_dispersion_type
143 : USE qs_dispersion_utils, ONLY: qs_dispersion_env_set,&
144 : qs_write_dispersion
145 : USE qs_energy_types, ONLY: allocate_qs_energy,&
146 : qs_energy_type
147 : USE qs_environment_methods, ONLY: qs_env_setup
148 : USE qs_environment_types, ONLY: get_qs_env,&
149 : qs_environment_type,&
150 : set_qs_env
151 : USE qs_force_types, ONLY: qs_force_type
152 : USE qs_gcp_types, ONLY: qs_gcp_type
153 : USE qs_gcp_utils, ONLY: qs_gcp_env_set,&
154 : qs_gcp_init
155 : USE qs_harris_types, ONLY: harris_rhoin_init,&
156 : harris_type
157 : USE qs_harris_utils, ONLY: harris_env_create,&
158 : harris_write_input
159 : USE qs_interactions, ONLY: init_interaction_radii,&
160 : init_se_nlradius,&
161 : write_core_charge_radii,&
162 : write_paw_radii,&
163 : write_pgf_orb_radii,&
164 : write_ppl_radii,&
165 : write_ppnl_radii
166 : USE qs_kind_types, ONLY: &
167 : check_qs_kind_set, get_qs_kind, get_qs_kind_set, init_gapw_basis_set, init_gapw_nlcc, &
168 : init_qs_kind_set, qs_kind_type, set_qs_kind, write_gto_basis_sets, write_qs_kind_set
169 : USE qs_ks_types, ONLY: qs_ks_env_create,&
170 : qs_ks_env_type,&
171 : set_ks_env
172 : USE qs_local_rho_types, ONLY: local_rho_type
173 : USE qs_mo_types, ONLY: allocate_mo_set,&
174 : mo_set_type
175 : USE qs_rho0_ggrid, ONLY: rho0_s_grid_create
176 : USE qs_rho0_methods, ONLY: init_rho0
177 : USE qs_rho0_types, ONLY: rho0_mpole_type
178 : USE qs_rho_atom_methods, ONLY: init_rho_atom
179 : USE qs_rho_atom_types, ONLY: rho_atom_type
180 : USE qs_subsys_methods, ONLY: qs_subsys_create
181 : USE qs_subsys_types, ONLY: qs_subsys_get,&
182 : qs_subsys_set,&
183 : qs_subsys_type
184 : USE qs_wf_history_methods, ONLY: wfi_create,&
185 : wfi_create_for_kp
186 : USE qs_wf_history_types, ONLY: qs_wf_history_type,&
187 : wfi_release
188 : USE rel_control_types, ONLY: rel_c_create,&
189 : rel_c_read_parameters,&
190 : rel_control_type
191 : USE scf_control_types, ONLY: scf_c_create,&
192 : scf_c_read_parameters,&
193 : scf_c_write_parameters,&
194 : scf_control_type
195 : USE semi_empirical_expns3_methods, ONLY: semi_empirical_expns3_setup
196 : USE semi_empirical_int_arrays, ONLY: init_se_intd_array
197 : USE semi_empirical_mpole_methods, ONLY: nddo_mpole_setup
198 : USE semi_empirical_mpole_types, ONLY: nddo_mpole_type
199 : USE semi_empirical_store_int_types, ONLY: semi_empirical_si_create,&
200 : semi_empirical_si_type
201 : USE semi_empirical_types, ONLY: se_taper_create,&
202 : se_taper_type
203 : USE semi_empirical_utils, ONLY: se_cutoff_compatible
204 : USE tblite_interface, ONLY: tb_get_basis,&
205 : tb_init_geometry,&
206 : tb_init_wf,&
207 : tb_set_calculator
208 : USE transport, ONLY: transport_env_create
209 : USE xtb_parameters, ONLY: init_xtb_basis,&
210 : xtb_parameters_init,&
211 : xtb_parameters_set
212 : USE xtb_potentials, ONLY: xtb_pp_radius
213 : USE xtb_types, ONLY: allocate_xtb_atom_param,&
214 : set_xtb_atom_param
215 : #include "./base/base_uses.f90"
216 :
217 : IMPLICIT NONE
218 :
219 : PRIVATE
220 :
221 : ! *** Global parameters ***
222 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_environment'
223 :
224 : ! *** Public subroutines ***
225 : PUBLIC :: qs_init
226 :
227 : CONTAINS
228 :
229 : ! **************************************************************************************************
230 : !> \brief Read the input and the database files for the setup of the
231 : !> QUICKSTEP environment.
232 : !> \param qs_env ...
233 : !> \param para_env ...
234 : !> \param root_section ...
235 : !> \param globenv ...
236 : !> \param cp_subsys ...
237 : !> \param kpoint_env ...
238 : !> \param cell ...
239 : !> \param cell_ref ...
240 : !> \param qmmm ...
241 : !> \param qmmm_env_qm ...
242 : !> \param force_env_section ...
243 : !> \param subsys_section ...
244 : !> \param use_motion_section ...
245 : !> \param silent ...
246 : !> \author Creation (22.05.2000,MK)
247 : ! **************************************************************************************************
248 51828 : SUBROUTINE qs_init(qs_env, para_env, root_section, globenv, cp_subsys, kpoint_env, cell, cell_ref, &
249 : qmmm, qmmm_env_qm, force_env_section, subsys_section, &
250 : use_motion_section, silent)
251 :
252 : TYPE(qs_environment_type), POINTER :: qs_env
253 : TYPE(mp_para_env_type), POINTER :: para_env
254 : TYPE(section_vals_type), OPTIONAL, POINTER :: root_section
255 : TYPE(global_environment_type), OPTIONAL, POINTER :: globenv
256 : TYPE(cp_subsys_type), OPTIONAL, POINTER :: cp_subsys
257 : TYPE(kpoint_type), OPTIONAL, POINTER :: kpoint_env
258 : TYPE(cell_type), OPTIONAL, POINTER :: cell, cell_ref
259 : LOGICAL, INTENT(IN), OPTIONAL :: qmmm
260 : TYPE(qmmm_env_qm_type), OPTIONAL, POINTER :: qmmm_env_qm
261 : TYPE(section_vals_type), POINTER :: force_env_section, subsys_section
262 : LOGICAL, INTENT(IN) :: use_motion_section
263 : LOGICAL, INTENT(IN), OPTIONAL :: silent
264 :
265 : CHARACTER(LEN=default_string_length) :: basis_type
266 : INTEGER :: ikind, method_id, nelectron_total, &
267 : nkind, nkp_grid(3)
268 : LOGICAL :: do_admm_rpa, do_ec_hfx, do_et, do_exx, do_hfx, do_kpoints, is_identical, is_semi, &
269 : mp2_present, my_qmmm, qmmm_decoupl, same_except_frac, use_ref_cell
270 7404 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: rtmat
271 7404 : TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
272 : TYPE(cell_type), POINTER :: my_cell, my_cell_ref
273 : TYPE(cp_blacs_env_type), POINTER :: blacs_env
274 : TYPE(dft_control_type), POINTER :: dft_control
275 : TYPE(distribution_1d_type), POINTER :: local_particles
276 : TYPE(energy_correction_type), POINTER :: ec_env
277 : TYPE(excited_energy_type), POINTER :: exstate_env
278 : TYPE(harris_type), POINTER :: harris_env
279 : TYPE(kpoint_type), POINTER :: kpoints
280 : TYPE(lri_environment_type), POINTER :: lri_env
281 7404 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
282 7404 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
283 : TYPE(qs_ks_env_type), POINTER :: ks_env
284 : TYPE(qs_subsys_type), POINTER :: subsys
285 : TYPE(qs_wf_history_type), POINTER :: wf_history
286 : TYPE(rel_control_type), POINTER :: rel_control
287 : TYPE(scf_control_type), POINTER :: scf_control
288 : TYPE(section_vals_type), POINTER :: dft_section, ec_hfx_section, ec_section, &
289 : et_coupling_section, hfx_section, kpoint_section, mp2_section, rpa_hfx_section, &
290 : transport_section
291 :
292 7404 : NULLIFY (my_cell, my_cell_ref, atomic_kind_set, particle_set, &
293 7404 : qs_kind_set, kpoint_section, dft_section, ec_section, &
294 7404 : subsys, ks_env, dft_control, blacs_env)
295 :
296 7404 : CALL set_qs_env(qs_env, input=force_env_section)
297 7404 : IF (.NOT. ASSOCIATED(subsys_section)) THEN
298 108 : subsys_section => section_vals_get_subs_vals(force_env_section, "SUBSYS")
299 : END IF
300 :
301 : ! QMMM
302 7404 : my_qmmm = .FALSE.
303 7404 : IF (PRESENT(qmmm)) my_qmmm = qmmm
304 7404 : qmmm_decoupl = .FALSE.
305 7404 : IF (PRESENT(qmmm_env_qm)) THEN
306 394 : IF (qmmm_env_qm%qmmm_coupl_type == do_qmmm_gauss .OR. &
307 : qmmm_env_qm%qmmm_coupl_type == do_qmmm_swave) THEN
308 : ! For GAUSS/SWAVE methods there could be a DDAPC decoupling requested
309 458 : qmmm_decoupl = my_qmmm .AND. qmmm_env_qm%periodic .AND. qmmm_env_qm%multipole
310 : END IF
311 394 : qs_env%qmmm_env_qm => qmmm_env_qm
312 : END IF
313 7404 : CALL set_qs_env(qs_env=qs_env, qmmm=my_qmmm)
314 :
315 : ! Possibly initialize arrays for SE
316 7404 : CALL section_vals_val_get(force_env_section, "DFT%QS%METHOD", i_val=method_id)
317 998 : SELECT CASE (method_id)
318 : CASE (do_method_rm1, do_method_am1, do_method_mndo, do_method_pdg, &
319 : do_method_pm3, do_method_pm6, do_method_pm6fm, do_method_mndod, do_method_pnnl)
320 998 : CALL init_se_intd_array()
321 998 : is_semi = .TRUE.
322 : CASE (do_method_xtb, do_method_dftb)
323 1206 : is_semi = .TRUE.
324 : CASE DEFAULT
325 7404 : is_semi = .FALSE.
326 : END SELECT
327 :
328 29616 : ALLOCATE (subsys)
329 : CALL qs_subsys_create(subsys, para_env, &
330 : force_env_section=force_env_section, &
331 : subsys_section=subsys_section, &
332 : use_motion_section=use_motion_section, &
333 : root_section=root_section, &
334 : cp_subsys=cp_subsys, cell=cell, cell_ref=cell_ref, &
335 7404 : elkind=is_semi, silent=silent)
336 :
337 7404 : ALLOCATE (ks_env)
338 7404 : CALL qs_ks_env_create(ks_env)
339 7404 : CALL set_ks_env(ks_env, subsys=subsys)
340 7404 : CALL set_qs_env(qs_env, ks_env=ks_env)
341 :
342 : CALL qs_subsys_get(subsys, &
343 : cell=my_cell, &
344 : cell_ref=my_cell_ref, &
345 : use_ref_cell=use_ref_cell, &
346 : atomic_kind_set=atomic_kind_set, &
347 : qs_kind_set=qs_kind_set, &
348 7404 : particle_set=particle_set)
349 :
350 7404 : CALL set_ks_env(ks_env, para_env=para_env)
351 7404 : IF (PRESENT(globenv)) THEN
352 : CALL cp_blacs_env_create(blacs_env, para_env, globenv%blacs_grid_layout, &
353 7398 : globenv%blacs_repeatable)
354 : ELSE
355 6 : CALL cp_blacs_env_create(blacs_env, para_env)
356 : END IF
357 7404 : CALL set_ks_env(ks_env, blacs_env=blacs_env)
358 7404 : CALL cp_blacs_env_release(blacs_env)
359 :
360 : ! *** Setup the grids for the G-space Interpolation if any
361 : CALL cp_ddapc_ewald_create(qs_env%cp_ddapc_ewald, qmmm_decoupl, my_cell, &
362 7404 : force_env_section, subsys_section, para_env)
363 :
364 : ! kpoints
365 7404 : IF (PRESENT(kpoint_env)) THEN
366 2 : kpoints => kpoint_env
367 2 : CALL set_qs_env(qs_env=qs_env, kpoints=kpoints)
368 2 : CALL kpoint_initialize(kpoints, particle_set, my_cell)
369 : ELSE
370 7402 : NULLIFY (kpoints)
371 7402 : CALL kpoint_create(kpoints)
372 7402 : CALL set_qs_env(qs_env=qs_env, kpoints=kpoints)
373 7402 : kpoint_section => section_vals_get_subs_vals(qs_env%input, "DFT%KPOINTS")
374 7402 : CALL read_kpoint_section(kpoints, kpoint_section, my_cell%hmat)
375 7402 : CALL kpoint_initialize(kpoints, particle_set, my_cell)
376 7402 : dft_section => section_vals_get_subs_vals(qs_env%input, "DFT")
377 7402 : CALL write_kpoint_info(kpoints, dft_section)
378 : END IF
379 :
380 : CALL qs_init_subsys(qs_env, para_env, subsys, my_cell, my_cell_ref, use_ref_cell, &
381 7404 : subsys_section, silent=silent)
382 :
383 7404 : CALL get_qs_env(qs_env, dft_control=dft_control)
384 7404 : IF (method_id == do_method_lrigpw .OR. dft_control%qs_control%lri_optbas) THEN
385 46 : CALL get_qs_env(qs_env=qs_env, lri_env=lri_env)
386 46 : CALL lri_env_basis("LRI", qs_env, lri_env, qs_kind_set)
387 7358 : ELSE IF (method_id == do_method_rigpw) THEN
388 : CALL cp_warn(__LOCATION__, "Experimental code: "// &
389 0 : "RIGPW should only be used for testing.")
390 0 : CALL get_qs_env(qs_env=qs_env, lri_env=lri_env)
391 0 : CALL lri_env_basis("RI", qs_env, lri_env, qs_kind_set)
392 : END IF
393 :
394 7404 : IF (my_qmmm .AND. PRESENT(qmmm_env_qm) .AND. .NOT. dft_control%qs_control%commensurate_mgrids) THEN
395 132 : IF (qmmm_env_qm%qmmm_coupl_type == do_qmmm_gauss .OR. qmmm_env_qm%qmmm_coupl_type == do_qmmm_swave) THEN
396 : CALL cp_abort(__LOCATION__, "QM/MM with coupling GAUSS or S-WAVE requires "// &
397 0 : "keyword FORCE_EVAL/DFT/MGRID/COMMENSURATE to be enabled.")
398 : END IF
399 : END IF
400 :
401 : ! more kpoint stuff
402 7404 : CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints, blacs_env=blacs_env)
403 7404 : IF (do_kpoints) THEN
404 170 : CALL kpoint_env_initialize(kpoints, para_env, blacs_env, with_aux_fit=dft_control%do_admm)
405 170 : CALL kpoint_initialize_mos(kpoints, qs_env%mos)
406 170 : CALL get_qs_env(qs_env=qs_env, wf_history=wf_history)
407 170 : CALL wfi_create_for_kp(wf_history)
408 : END IF
409 : ! basis set symmetry rotations
410 7404 : IF (do_kpoints) THEN
411 170 : CALL qs_basis_rotation(qs_env, kpoints)
412 : END IF
413 :
414 : do_hfx = .FALSE.
415 7404 : hfx_section => section_vals_get_subs_vals(qs_env%input, "DFT%XC%HF")
416 7404 : CALL section_vals_get(hfx_section, explicit=do_hfx)
417 7404 : CALL get_qs_env(qs_env, dft_control=dft_control, scf_control=scf_control, nelectron_total=nelectron_total)
418 7404 : IF (do_hfx) THEN
419 : ! Retrieve particle_set and atomic_kind_set (needed for both kinds of initialization)
420 4664 : nkp_grid = 1
421 1166 : IF (do_kpoints) CALL get_kpoint_info(kpoints, nkp_grid=nkp_grid)
422 1166 : IF (dft_control%do_admm) THEN
423 442 : basis_type = 'AUX_FIT'
424 : ELSE
425 724 : basis_type = 'ORB'
426 : END IF
427 : CALL hfx_create(qs_env%x_data, para_env, hfx_section, atomic_kind_set, &
428 : qs_kind_set, particle_set, dft_control, my_cell, orb_basis=basis_type, &
429 1166 : nelectron_total=nelectron_total, nkp_grid=nkp_grid)
430 : END IF
431 :
432 7404 : mp2_section => section_vals_get_subs_vals(qs_env%input, "DFT%XC%WF_CORRELATION")
433 7404 : CALL section_vals_get(mp2_section, explicit=mp2_present)
434 7404 : IF (mp2_present) THEN
435 470 : CPASSERT(ASSOCIATED(qs_env%mp2_env))
436 470 : CALL read_mp2_section(qs_env%input, qs_env%mp2_env)
437 : ! create the EXX section if necessary
438 : do_exx = .FALSE.
439 470 : rpa_hfx_section => section_vals_get_subs_vals(qs_env%input, "DFT%XC%WF_CORRELATION%RI_RPA%HF")
440 470 : CALL section_vals_get(rpa_hfx_section, explicit=do_exx)
441 470 : IF (do_exx) THEN
442 :
443 : ! do_exx in call of hfx_create decides whether to go without ADMM (do_exx=.TRUE.) or with
444 : ! ADMM (do_exx=.FALSE.)
445 142 : CALL section_vals_val_get(mp2_section, "RI_RPA%ADMM", l_val=do_admm_rpa)
446 :
447 : ! Reuse the HFX integrals from the qs_env if applicable
448 142 : qs_env%mp2_env%ri_rpa%reuse_hfx = .TRUE.
449 142 : IF (.NOT. do_hfx) qs_env%mp2_env%ri_rpa%reuse_hfx = .FALSE.
450 142 : CALL compare_hfx_sections(hfx_section, rpa_hfx_section, is_identical, same_except_frac)
451 142 : IF (.NOT. (is_identical .OR. same_except_frac)) qs_env%mp2_env%ri_rpa%reuse_hfx = .FALSE.
452 142 : IF (dft_control%do_admm .AND. .NOT. do_admm_rpa) qs_env%mp2_env%ri_rpa%reuse_hfx = .FALSE.
453 :
454 142 : IF (.NOT. qs_env%mp2_env%ri_rpa%reuse_hfx) THEN
455 124 : IF (do_admm_rpa) THEN
456 10 : basis_type = 'AUX_FIT'
457 : ELSE
458 114 : basis_type = 'ORB'
459 : END IF
460 : CALL hfx_create(qs_env%mp2_env%ri_rpa%x_data, para_env, rpa_hfx_section, atomic_kind_set, &
461 : qs_kind_set, particle_set, dft_control, my_cell, orb_basis=basis_type, &
462 124 : nelectron_total=nelectron_total)
463 : ELSE
464 18 : qs_env%mp2_env%ri_rpa%x_data => qs_env%x_data
465 : END IF
466 : END IF
467 : END IF
468 :
469 7404 : IF (dft_control%qs_control%do_kg) THEN
470 66 : CALL cite_reference(Iannuzzi2006)
471 66 : CALL kg_env_create(qs_env, qs_env%kg_env, qs_kind_set, qs_env%input)
472 : END IF
473 :
474 7404 : dft_section => section_vals_get_subs_vals(qs_env%input, "DFT")
475 : CALL section_vals_val_get(dft_section, "EXCITED_STATES%_SECTION_PARAMETERS_", &
476 7404 : l_val=qs_env%excited_state)
477 7404 : NULLIFY (exstate_env)
478 7404 : CALL exstate_create(exstate_env, qs_env%excited_state, dft_section)
479 7404 : CALL set_qs_env(qs_env, exstate_env=exstate_env)
480 :
481 : et_coupling_section => section_vals_get_subs_vals(qs_env%input, &
482 7404 : "PROPERTIES%ET_COUPLING")
483 7404 : CALL section_vals_get(et_coupling_section, explicit=do_et)
484 7404 : IF (do_et) CALL et_coupling_create(qs_env%et_coupling)
485 :
486 7404 : transport_section => section_vals_get_subs_vals(qs_env%input, "DFT%TRANSPORT")
487 7404 : CALL section_vals_get(transport_section, explicit=qs_env%do_transport)
488 7404 : IF (qs_env%do_transport) THEN
489 0 : CALL transport_env_create(qs_env)
490 : END IF
491 :
492 7404 : CALL get_qs_env(qs_env, harris_env=harris_env)
493 7404 : IF (qs_env%harris_method) THEN
494 : ! initialize the Harris input density and potential integrals
495 6 : CALL get_qs_env(qs_env, local_particles=local_particles)
496 : CALL harris_rhoin_init(harris_env%rhoin, "RHOIN", qs_kind_set, atomic_kind_set, &
497 6 : local_particles, dft_control%nspins)
498 : ! Print information of the HARRIS section
499 6 : CALL harris_write_input(harris_env)
500 : END IF
501 :
502 7404 : NULLIFY (ec_env)
503 7404 : dft_section => section_vals_get_subs_vals(qs_env%input, "DFT")
504 : CALL section_vals_val_get(dft_section, "ENERGY_CORRECTION%_SECTION_PARAMETERS_", &
505 7404 : l_val=qs_env%energy_correction)
506 7404 : ec_section => section_vals_get_subs_vals(qs_env%input, "DFT%ENERGY_CORRECTION")
507 7404 : CALL ec_env_create(qs_env, ec_env, dft_section, ec_section)
508 7404 : CALL set_qs_env(qs_env, ec_env=ec_env)
509 :
510 7404 : IF (qs_env%energy_correction) THEN
511 : ! Energy correction with Hartree-Fock exchange
512 236 : ec_hfx_section => section_vals_get_subs_vals(ec_section, "XC%HF")
513 236 : CALL section_vals_get(ec_hfx_section, explicit=do_ec_hfx)
514 :
515 236 : IF (ec_env%do_ec_hfx) THEN
516 :
517 : ! Hybrid functionals require same basis
518 16 : IF (ec_env%basis_inconsistent) THEN
519 : CALL cp_abort(__LOCATION__, &
520 : "Energy correction methods with hybrid functionals: "// &
521 : "correction and ground state need to use the same basis. "// &
522 0 : "Checked by comparing basis set names only.")
523 : END IF
524 :
525 : ! Similar to RPA_HFX we can check if HFX integrals from the qs_env can be reused
526 16 : IF (ec_env%do_ec_admm .AND. .NOT. dft_control%do_admm) THEN
527 0 : CALL cp_abort(__LOCATION__, "Need an ADMM input section for ADMM EC to work")
528 : END IF
529 :
530 16 : ec_env%reuse_hfx = .TRUE.
531 16 : IF (.NOT. do_hfx) ec_env%reuse_hfx = .FALSE.
532 16 : CALL compare_hfx_sections(hfx_section, ec_hfx_section, is_identical, same_except_frac)
533 16 : IF (.NOT. (is_identical .OR. same_except_frac)) ec_env%reuse_hfx = .FALSE.
534 16 : IF (dft_control%do_admm .AND. .NOT. ec_env%do_ec_admm) ec_env%reuse_hfx = .FALSE.
535 :
536 16 : IF (.NOT. ec_env%reuse_hfx) THEN
537 6 : IF (ec_env%do_ec_admm) THEN
538 2 : basis_type = 'AUX_FIT'
539 : ELSE
540 4 : basis_type = 'ORB'
541 : END IF
542 : CALL hfx_create(ec_env%x_data, para_env, ec_hfx_section, atomic_kind_set, &
543 : qs_kind_set, particle_set, dft_control, my_cell, orb_basis=basis_type, &
544 6 : nelectron_total=nelectron_total)
545 : ELSE
546 10 : ec_env%x_data => qs_env%x_data
547 : END IF
548 : END IF
549 :
550 : ! Print information of the EC section
551 236 : CALL ec_write_input(ec_env)
552 :
553 : END IF
554 :
555 7404 : IF (dft_control%qs_control%do_almo_scf) THEN
556 66 : CALL almo_scf_env_create(qs_env)
557 : END IF
558 :
559 : ! see if we have atomic relativistic corrections
560 7404 : CALL get_qs_env(qs_env, rel_control=rel_control)
561 7404 : IF (rel_control%rel_method /= rel_none) THEN
562 16 : IF (rel_control%rel_transformation == rel_trans_atom) THEN
563 16 : nkind = SIZE(atomic_kind_set)
564 42 : DO ikind = 1, nkind
565 26 : NULLIFY (rtmat)
566 26 : CALL calculate_atomic_relkin(atomic_kind_set(ikind), qs_kind_set(ikind), rel_control, rtmat)
567 42 : IF (ASSOCIATED(rtmat)) CALL set_qs_kind(qs_kind_set(ikind), reltmat=rtmat)
568 : END DO
569 : END IF
570 : END IF
571 :
572 7404 : END SUBROUTINE qs_init
573 :
574 : ! **************************************************************************************************
575 : !> \brief Initialize the qs environment (subsys)
576 : !> \param qs_env ...
577 : !> \param para_env ...
578 : !> \param subsys ...
579 : !> \param cell ...
580 : !> \param cell_ref ...
581 : !> \param use_ref_cell ...
582 : !> \param subsys_section ...
583 : !> \param silent ...
584 : !> \author Creation (22.05.2000,MK)
585 : ! **************************************************************************************************
586 7404 : SUBROUTINE qs_init_subsys(qs_env, para_env, subsys, cell, cell_ref, use_ref_cell, subsys_section, &
587 : silent)
588 :
589 : TYPE(qs_environment_type), POINTER :: qs_env
590 : TYPE(mp_para_env_type), POINTER :: para_env
591 : TYPE(qs_subsys_type), POINTER :: subsys
592 : TYPE(cell_type), POINTER :: cell, cell_ref
593 : LOGICAL, INTENT(in) :: use_ref_cell
594 : TYPE(section_vals_type), POINTER :: subsys_section
595 : LOGICAL, INTENT(in), OPTIONAL :: silent
596 :
597 : CHARACTER(len=*), PARAMETER :: routineN = 'qs_init_subsys'
598 :
599 : CHARACTER(len=2) :: element_symbol
600 : INTEGER :: gfn_type, handle, ikind, ispin, iw, lmax_sphere, maxl, maxlgto, maxlgto_lri, &
601 : maxlppl, maxlppnl, method_id, multiplicity, my_ival, n_ao, n_mo_add, natom, nelectron, &
602 : ngauss, nkind, output_unit, sort_basis, tnadd_method
603 : INTEGER, DIMENSION(2) :: n_mo, nelectron_spin
604 : INTEGER, DIMENSION(5) :: occ
605 : LOGICAL :: all_potential_present, be_silent, do_kpoints, do_ri_hfx, do_ri_mp2, do_ri_rpa, &
606 : do_ri_sos_mp2, do_rpa_ri_exx, do_wfc_im_time, e1terms, has_unit_metric, lribas, &
607 : mp2_present, orb_gradient
608 : REAL(KIND=dp) :: alpha, ccore, ewald_rcut, fxx, maxocc, &
609 : rcut, total_zeff_corr, verlet_skin, &
610 : zeff_correction
611 7404 : TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
612 : TYPE(cp_logger_type), POINTER :: logger
613 : TYPE(dft_control_type), POINTER :: dft_control
614 : TYPE(dftb_control_type), POINTER :: dftb_control
615 : TYPE(distribution_1d_type), POINTER :: local_molecules, local_particles
616 : TYPE(ewald_environment_type), POINTER :: ewald_env
617 : TYPE(ewald_pw_type), POINTER :: ewald_pw
618 : TYPE(fist_nonbond_env_type), POINTER :: se_nonbond_env
619 : TYPE(gapw_control_type), POINTER :: gapw_control
620 : TYPE(gto_basis_set_type), POINTER :: aux_fit_basis, lri_aux_basis, &
621 : rhoin_basis, ri_aux_basis_set, &
622 : ri_hfx_basis, ri_xas_basis, &
623 : tmp_basis_set
624 : TYPE(harris_type), POINTER :: harris_env
625 : TYPE(local_rho_type), POINTER :: local_rho_set
626 : TYPE(lri_environment_type), POINTER :: lri_env
627 7404 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos, mos_last_converged
628 7404 : TYPE(molecule_kind_type), DIMENSION(:), POINTER :: molecule_kind_set
629 7404 : TYPE(molecule_type), DIMENSION(:), POINTER :: molecule_set
630 : TYPE(mp2_type), POINTER :: mp2_env
631 : TYPE(nddo_mpole_type), POINTER :: se_nddo_mpole
632 7404 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
633 : TYPE(pw_env_type), POINTER :: pw_env
634 : TYPE(qs_control_type), POINTER :: qs_control
635 : TYPE(qs_dftb_pairpot_type), DIMENSION(:, :), &
636 7404 : POINTER :: dftb_potential
637 : TYPE(qs_dispersion_type), POINTER :: dispersion_env
638 : TYPE(qs_energy_type), POINTER :: energy
639 7404 : TYPE(qs_force_type), DIMENSION(:), POINTER :: force
640 : TYPE(qs_gcp_type), POINTER :: gcp_env
641 7404 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
642 : TYPE(qs_kind_type), POINTER :: qs_kind
643 : TYPE(qs_ks_env_type), POINTER :: ks_env
644 : TYPE(qs_wf_history_type), POINTER :: wf_history
645 : TYPE(rho0_mpole_type), POINTER :: rho0_mpole
646 7404 : TYPE(rho_atom_type), DIMENSION(:), POINTER :: rho_atom_set
647 : TYPE(scf_control_type), POINTER :: scf_control
648 : TYPE(se_taper_type), POINTER :: se_taper
649 : TYPE(section_vals_type), POINTER :: dft_section, et_coupling_section, et_ddapc_section, &
650 : ewald_section, harris_section, lri_section, mp2_section, nl_section, poisson_section, &
651 : pp_section, print_section, qs_section, rixs_section, se_section, tddfpt_section, &
652 : xc_section
653 : TYPE(semi_empirical_control_type), POINTER :: se_control
654 : TYPE(semi_empirical_si_type), POINTER :: se_store_int_env
655 : TYPE(xtb_control_type), POINTER :: xtb_control
656 :
657 7404 : CALL timeset(routineN, handle)
658 7404 : NULLIFY (logger)
659 7404 : logger => cp_get_default_logger()
660 7404 : output_unit = cp_logger_get_default_io_unit(logger)
661 :
662 7404 : be_silent = .FALSE.
663 7404 : IF (PRESENT(silent)) be_silent = silent
664 :
665 7404 : CALL cite_reference(cp2kqs2020)
666 :
667 : ! Initialise the Quickstep environment
668 7404 : NULLIFY (mos, se_taper)
669 7404 : NULLIFY (dft_control)
670 7404 : NULLIFY (energy)
671 7404 : NULLIFY (force)
672 7404 : NULLIFY (local_molecules)
673 7404 : NULLIFY (local_particles)
674 7404 : NULLIFY (scf_control)
675 7404 : NULLIFY (dft_section)
676 7404 : NULLIFY (et_coupling_section)
677 7404 : NULLIFY (ks_env)
678 7404 : NULLIFY (mos_last_converged)
679 7404 : dft_section => section_vals_get_subs_vals(qs_env%input, "DFT")
680 7404 : qs_section => section_vals_get_subs_vals(dft_section, "QS")
681 7404 : et_coupling_section => section_vals_get_subs_vals(qs_env%input, "PROPERTIES%ET_COUPLING")
682 : ! reimplemented TDDFPT
683 7404 : tddfpt_section => section_vals_get_subs_vals(qs_env%input, "PROPERTIES%TDDFPT")
684 7404 : rixs_section => section_vals_get_subs_vals(qs_env%input, "PROPERTIES%RIXS")
685 :
686 : CALL qs_subsys_get(subsys, particle_set=particle_set, &
687 : qs_kind_set=qs_kind_set, &
688 : atomic_kind_set=atomic_kind_set, &
689 : molecule_set=molecule_set, &
690 7404 : molecule_kind_set=molecule_kind_set)
691 :
692 : ! *** Read the input section with the DFT control parameters ***
693 7404 : CALL read_dft_control(dft_control, dft_section)
694 :
695 : ! set periodicity flag
696 29616 : dft_control%qs_control%periodicity = SUM(cell%perd)
697 :
698 : ! Read the input section with the Quickstep control parameters
699 7404 : CALL read_qs_section(dft_control%qs_control, qs_section)
700 :
701 : ! *** Print the Quickstep program banner (copyright and version number) ***
702 7404 : IF (.NOT. be_silent) THEN
703 7398 : iw = cp_print_key_unit_nr(logger, dft_section, "PRINT%PROGRAM_BANNER", extension=".Log")
704 7398 : CALL section_vals_val_get(qs_section, "METHOD", i_val=method_id)
705 5198 : SELECT CASE (method_id)
706 : CASE DEFAULT
707 5198 : CALL qs_header(iw)
708 : CASE (do_method_rm1, do_method_am1, do_method_mndo, do_method_pdg, &
709 : do_method_pm3, do_method_pm6, do_method_pm6fm, do_method_mndod, do_method_pnnl)
710 998 : CALL se_header(iw)
711 : CASE (do_method_dftb)
712 222 : CALL dftb_header(iw)
713 : CASE (do_method_xtb)
714 7398 : IF (dft_control%qs_control%xtb_control%do_tblite) THEN
715 40 : CALL tblite_header(iw, dft_control%qs_control%xtb_control%tblite_method)
716 : ELSE
717 940 : gfn_type = dft_control%qs_control%xtb_control%gfn_type
718 940 : CALL xtb_header(iw, gfn_type)
719 : END IF
720 : END SELECT
721 : CALL cp_print_key_finished_output(iw, logger, dft_section, &
722 7398 : "PRINT%PROGRAM_BANNER")
723 : END IF
724 :
725 7404 : IF (dft_control%do_sccs .AND. dft_control%qs_control%gapw) THEN
726 0 : CPABORT("SCCS is not yet implemented with GAPW")
727 : END IF
728 7404 : CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints)
729 7404 : IF (do_kpoints) THEN
730 : ! reset some of the settings for wfn extrapolation for kpoints
731 274 : SELECT CASE (dft_control%qs_control%wf_interpolation_method_nr)
732 : CASE (wfi_linear_wf_method_nr, wfi_linear_ps_method_nr, wfi_ps_method_nr, wfi_aspc_nr)
733 170 : dft_control%qs_control%wf_interpolation_method_nr = wfi_use_guess_method_nr
734 : END SELECT
735 : END IF
736 :
737 : ! ******* check if any kind of electron transfer calculation has to be performed
738 7404 : CALL section_vals_val_get(et_coupling_section, "TYPE_OF_CONSTRAINT", i_val=my_ival)
739 7404 : dft_control%qs_control%et_coupling_calc = .FALSE.
740 7404 : IF (my_ival == do_et_ddapc) THEN
741 0 : et_ddapc_section => section_vals_get_subs_vals(et_coupling_section, "DDAPC_RESTRAINT_A")
742 0 : dft_control%qs_control%et_coupling_calc = .TRUE.
743 0 : dft_control%qs_control%ddapc_restraint = .TRUE.
744 0 : CALL read_ddapc_section(dft_control%qs_control, ddapc_restraint_section=et_ddapc_section)
745 : END IF
746 :
747 7404 : CALL read_mgrid_section(dft_control%qs_control, dft_section)
748 :
749 : ! reimplemented TDDFPT
750 7404 : CALL read_tddfpt2_control(dft_control%tddfpt2_control, tddfpt_section, dft_control%qs_control)
751 :
752 : ! rixs
753 7404 : CALL section_vals_get(rixs_section, explicit=qs_env%do_rixs)
754 7404 : IF (qs_env%do_rixs) THEN
755 10 : CALL read_rixs_control(dft_control%rixs_control, rixs_section, dft_control%qs_control)
756 : END IF
757 :
758 : ! Create relativistic control section
759 : BLOCK
760 : TYPE(rel_control_type), POINTER :: rel_control
761 7404 : ALLOCATE (rel_control)
762 7404 : CALL rel_c_create(rel_control)
763 7404 : CALL rel_c_read_parameters(rel_control, dft_section)
764 7404 : CALL set_qs_env(qs_env, rel_control=rel_control)
765 : END BLOCK
766 :
767 : ! *** Read DFTB parameter files ***
768 7404 : IF (dft_control%qs_control%method_id == do_method_dftb) THEN
769 222 : NULLIFY (ewald_env, ewald_pw, dftb_potential)
770 222 : dftb_control => dft_control%qs_control%dftb_control
771 : CALL qs_dftb_param_init(atomic_kind_set, qs_kind_set, dftb_control, dftb_potential, &
772 222 : subsys_section=subsys_section, para_env=para_env)
773 222 : CALL set_qs_env(qs_env, dftb_potential=dftb_potential)
774 : ! check for Ewald
775 222 : IF (dftb_control%do_ewald) THEN
776 1888 : ALLOCATE (ewald_env)
777 118 : CALL ewald_env_create(ewald_env, para_env)
778 118 : poisson_section => section_vals_get_subs_vals(dft_section, "POISSON")
779 118 : CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
780 118 : ewald_section => section_vals_get_subs_vals(poisson_section, "EWALD")
781 118 : print_section => section_vals_get_subs_vals(qs_env%input, "PRINT%GRID_INFORMATION")
782 118 : CALL get_qs_kind_set(qs_kind_set, basis_rcut=ewald_rcut)
783 118 : CALL read_ewald_section_tb(ewald_env, ewald_section, cell_ref%hmat)
784 118 : ALLOCATE (ewald_pw)
785 118 : CALL ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section=print_section)
786 118 : CALL set_qs_env(qs_env, ewald_env=ewald_env, ewald_pw=ewald_pw)
787 : END IF
788 7182 : ELSEIF (dft_control%qs_control%method_id == do_method_xtb) THEN
789 : ! *** Read xTB parameter file ***
790 984 : xtb_control => dft_control%qs_control%xtb_control
791 984 : CALL get_qs_env(qs_env, nkind=nkind)
792 984 : IF (xtb_control%do_tblite) THEN
793 : ! put geometry to tblite
794 40 : CALL tb_init_geometry(qs_env, qs_env%tb_tblite)
795 : ! select tblite method
796 40 : CALL tb_set_calculator(qs_env%tb_tblite, xtb_control%tblite_method)
797 : !set up wave function
798 40 : CALL tb_init_wf(qs_env%tb_tblite)
799 : !get basis set
800 152 : DO ikind = 1, nkind
801 112 : qs_kind => qs_kind_set(ikind)
802 : ! Setup proper xTB parameters
803 112 : CPASSERT(.NOT. ASSOCIATED(qs_kind%xtb_parameter))
804 112 : CALL allocate_xtb_atom_param(qs_kind%xtb_parameter)
805 : ! Set default parameters
806 112 : CALL get_qs_kind(qs_kind, element_symbol=element_symbol)
807 :
808 112 : NULLIFY (tmp_basis_set)
809 112 : CALL tb_get_basis(qs_env%tb_tblite, tmp_basis_set, element_symbol, qs_kind%xtb_parameter, occ)
810 112 : CALL add_basis_set_to_container(qs_kind%basis_sets, tmp_basis_set, "ORB")
811 112 : CALL set_xtb_atom_param(qs_kind%xtb_parameter, occupation=occ)
812 :
813 : !setting the potential for the computation
814 112 : zeff_correction = 0.0_dp
815 : CALL init_potential(qs_kind%all_potential, itype="BARE", &
816 712 : zeff=REAL(SUM(occ), dp), zeff_correction=zeff_correction)
817 : END DO
818 : ELSE
819 944 : NULLIFY (ewald_env, ewald_pw)
820 3040 : DO ikind = 1, nkind
821 2096 : qs_kind => qs_kind_set(ikind)
822 : ! Setup proper xTB parameters
823 2096 : CPASSERT(.NOT. ASSOCIATED(qs_kind%xtb_parameter))
824 2096 : CALL allocate_xtb_atom_param(qs_kind%xtb_parameter)
825 : ! Set default parameters
826 2096 : gfn_type = dft_control%qs_control%xtb_control%gfn_type
827 2096 : CALL get_qs_kind(qs_kind, element_symbol=element_symbol)
828 : CALL xtb_parameters_init(qs_kind%xtb_parameter, gfn_type, element_symbol, &
829 : xtb_control%parameter_file_path, xtb_control%parameter_file_name, &
830 2096 : para_env)
831 : ! set dependent parameters
832 2096 : CALL xtb_parameters_set(qs_kind%xtb_parameter)
833 : ! Generate basis set
834 2096 : NULLIFY (tmp_basis_set)
835 2096 : IF (qs_kind%xtb_parameter%z == 1) THEN
836 : ! special case hydrogen
837 456 : ngauss = xtb_control%h_sto_ng
838 : ELSE
839 1640 : ngauss = xtb_control%sto_ng
840 : END IF
841 2096 : IF (qs_kind%xtb_parameter%defined) THEN
842 2094 : CALL init_xtb_basis(qs_kind%xtb_parameter, tmp_basis_set, ngauss)
843 2094 : CALL add_basis_set_to_container(qs_kind%basis_sets, tmp_basis_set, "ORB")
844 : ELSE
845 2 : CALL set_qs_kind(qs_kind, ghost=.TRUE.)
846 2 : IF (ASSOCIATED(qs_kind%all_potential)) THEN
847 2 : DEALLOCATE (qs_kind%all_potential%elec_conf)
848 2 : DEALLOCATE (qs_kind%all_potential)
849 : END IF
850 : END IF
851 : ! potential
852 3040 : IF (qs_kind%xtb_parameter%defined) THEN
853 2094 : zeff_correction = 0.0_dp
854 : CALL init_potential(qs_kind%all_potential, itype="BARE", &
855 2094 : zeff=qs_kind%xtb_parameter%zeff, zeff_correction=zeff_correction)
856 2094 : CALL get_potential(qs_kind%all_potential, alpha_core_charge=alpha)
857 2094 : ccore = qs_kind%xtb_parameter%zeff*SQRT((alpha/pi)**3)
858 2094 : CALL set_potential(qs_kind%all_potential, ccore_charge=ccore)
859 2094 : qs_kind%xtb_parameter%zeff = qs_kind%xtb_parameter%zeff - zeff_correction
860 : END IF
861 : END DO
862 : !
863 : ! set repulsive potential range
864 : !
865 3776 : ALLOCATE (xtb_control%rcpair(nkind, nkind))
866 944 : CALL xtb_pp_radius(qs_kind_set, xtb_control%rcpair, xtb_control%eps_pair, xtb_control%kf)
867 : ! check for Ewald
868 944 : IF (xtb_control%do_ewald) THEN
869 2848 : ALLOCATE (ewald_env)
870 178 : CALL ewald_env_create(ewald_env, para_env)
871 178 : poisson_section => section_vals_get_subs_vals(dft_section, "POISSON")
872 178 : CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
873 178 : ewald_section => section_vals_get_subs_vals(poisson_section, "EWALD")
874 178 : print_section => section_vals_get_subs_vals(qs_env%input, "PRINT%GRID_INFORMATION")
875 178 : IF (gfn_type == 0) THEN
876 : CALL read_ewald_section_tb(ewald_env, ewald_section, cell_ref%hmat, &
877 34 : silent=silent, pset="EEQ")
878 : ELSE
879 : CALL read_ewald_section_tb(ewald_env, ewald_section, cell_ref%hmat, &
880 144 : silent=silent)
881 : END IF
882 178 : ALLOCATE (ewald_pw)
883 178 : CALL ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section=print_section)
884 178 : CALL set_qs_env(qs_env, ewald_env=ewald_env, ewald_pw=ewald_pw)
885 : END IF
886 : END IF
887 : END IF
888 : ! lri or ri env initialization
889 7404 : lri_section => section_vals_get_subs_vals(qs_section, "LRIGPW")
890 : IF (dft_control%qs_control%method_id == do_method_lrigpw .OR. &
891 7404 : dft_control%qs_control%lri_optbas .OR. &
892 : dft_control%qs_control%method_id == do_method_rigpw) THEN
893 46 : CALL lri_env_init(lri_env, lri_section)
894 46 : CALL set_qs_env(qs_env, lri_env=lri_env)
895 : END IF
896 :
897 : ! *** Check basis and fill in missing parts ***
898 7404 : CALL check_qs_kind_set(qs_kind_set, dft_control, subsys_section=subsys_section)
899 :
900 : ! *** Check that no all-electron potential is present if GPW or GAPW_XC
901 7404 : CALL get_qs_kind_set(qs_kind_set, all_potential_present=all_potential_present)
902 : IF ((dft_control%qs_control%method_id == do_method_gpw) .OR. &
903 7404 : (dft_control%qs_control%method_id == do_method_gapw_xc) .OR. &
904 : (dft_control%qs_control%method_id == do_method_ofgpw)) THEN
905 4332 : IF (all_potential_present) THEN
906 0 : CPABORT("All-electron calculations with GPW, GAPW_XC, and OFGPW are not implemented")
907 : END IF
908 : END IF
909 :
910 : ! DFT+U
911 7404 : CALL get_qs_kind_set(qs_kind_set, dft_plus_u_atom_present=dft_control%dft_plus_u)
912 :
913 7404 : IF (dft_control%do_admm) THEN
914 : ! Check if ADMM basis is available
915 450 : CALL get_qs_env(qs_env, nkind=nkind)
916 1276 : DO ikind = 1, nkind
917 826 : NULLIFY (aux_fit_basis)
918 826 : qs_kind => qs_kind_set(ikind)
919 826 : CALL get_qs_kind(qs_kind, basis_set=aux_fit_basis, basis_type="AUX_FIT")
920 1276 : IF (.NOT. (ASSOCIATED(aux_fit_basis))) THEN
921 : ! AUX_FIT basis set is not available
922 0 : CPABORT("AUX_FIT basis set is not defined. ")
923 : END IF
924 : END DO
925 : END IF
926 :
927 7404 : lribas = .FALSE.
928 7404 : e1terms = .FALSE.
929 7404 : IF (dft_control%qs_control%method_id == do_method_lrigpw) THEN
930 40 : lribas = .TRUE.
931 40 : CALL get_qs_env(qs_env, lri_env=lri_env)
932 40 : e1terms = lri_env%exact_1c_terms
933 : END IF
934 7404 : IF (dft_control%qs_control%do_kg) THEN
935 66 : CALL section_vals_val_get(dft_section, "KG_METHOD%TNADD_METHOD", i_val=tnadd_method)
936 66 : IF (tnadd_method == kg_tnadd_embed_ri) lribas = .TRUE.
937 : END IF
938 7402 : IF (lribas) THEN
939 : ! Check if LRI_AUX basis is available, auto-generate if needed
940 42 : CALL get_qs_env(qs_env, nkind=nkind)
941 122 : DO ikind = 1, nkind
942 80 : NULLIFY (lri_aux_basis)
943 80 : qs_kind => qs_kind_set(ikind)
944 80 : CALL get_qs_kind(qs_kind, basis_set=lri_aux_basis, basis_type="LRI_AUX")
945 122 : IF (.NOT. (ASSOCIATED(lri_aux_basis))) THEN
946 : ! LRI_AUX basis set is not yet loaded
947 : CALL cp_warn(__LOCATION__, "Automatic Generation of LRI_AUX basis. "// &
948 18 : "This is experimental code.")
949 : ! Generate a default basis
950 18 : CALL create_lri_aux_basis_set(lri_aux_basis, qs_kind, dft_control%auto_basis_lri_aux, e1terms)
951 18 : CALL add_basis_set_to_container(qs_kind%basis_sets, lri_aux_basis, "LRI_AUX")
952 : END IF
953 : END DO
954 : END IF
955 :
956 7404 : CALL section_vals_val_get(qs_env%input, "DFT%XC%HF%RI%_SECTION_PARAMETERS_", l_val=do_ri_hfx)
957 : CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%RI_RPA%HF%RI%_SECTION_PARAMETERS_", &
958 7404 : l_val=do_rpa_ri_exx)
959 7404 : IF (do_ri_hfx .OR. do_rpa_ri_exx) THEN
960 108 : CALL get_qs_env(qs_env, nkind=nkind)
961 108 : CALL section_vals_val_get(qs_env%input, "DFT%SORT_BASIS", i_val=sort_basis)
962 290 : DO ikind = 1, nkind
963 182 : NULLIFY (ri_hfx_basis)
964 182 : qs_kind => qs_kind_set(ikind)
965 : CALL get_qs_kind(qs_kind=qs_kind, basis_set=ri_hfx_basis, &
966 182 : basis_type="RI_HFX")
967 7586 : IF (.NOT. (ASSOCIATED(ri_hfx_basis))) THEN
968 178 : CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto)
969 178 : IF (dft_control%do_admm) THEN
970 : CALL create_ri_aux_basis_set(ri_hfx_basis, qs_kind, dft_control%auto_basis_ri_hfx, &
971 58 : basis_type="AUX_FIT", basis_sort=sort_basis)
972 : ELSE
973 : CALL create_ri_aux_basis_set(ri_hfx_basis, qs_kind, dft_control%auto_basis_ri_hfx, &
974 120 : basis_sort=sort_basis)
975 : END IF
976 178 : CALL add_basis_set_to_container(qs_kind%basis_sets, ri_hfx_basis, "RI_HFX")
977 : END IF
978 : END DO
979 : END IF
980 :
981 7404 : IF (dft_control%qs_control%method_id == do_method_rigpw) THEN
982 : ! Check if RI_HXC basis is available, auto-generate if needed
983 0 : CALL get_qs_env(qs_env, nkind=nkind)
984 0 : DO ikind = 1, nkind
985 0 : NULLIFY (ri_hfx_basis)
986 0 : qs_kind => qs_kind_set(ikind)
987 0 : CALL get_qs_kind(qs_kind, basis_set=ri_hfx_basis, basis_type="RI_HXC")
988 0 : IF (.NOT. (ASSOCIATED(ri_hfx_basis))) THEN
989 : ! Generate a default basis
990 0 : CALL create_ri_aux_basis_set(ri_hfx_basis, qs_kind, dft_control%auto_basis_ri_hxc)
991 0 : CALL add_basis_set_to_container(qs_kind%basis_sets, ri_hfx_basis, "RI_HXC")
992 : END IF
993 : END DO
994 : END IF
995 :
996 : ! Harris method
997 7404 : NULLIFY (harris_env)
998 : CALL section_vals_val_get(dft_section, "HARRIS_METHOD%_SECTION_PARAMETERS_", &
999 7404 : l_val=qs_env%harris_method)
1000 7404 : harris_section => section_vals_get_subs_vals(dft_section, "HARRIS_METHOD")
1001 7404 : CALL harris_env_create(qs_env, harris_env, harris_section)
1002 7404 : CALL set_qs_env(qs_env, harris_env=harris_env)
1003 : !
1004 7404 : IF (qs_env%harris_method) THEN
1005 6 : CALL get_qs_env(qs_env, nkind=nkind)
1006 : ! Check if RI_HXC basis is available, auto-generate if needed
1007 22 : DO ikind = 1, nkind
1008 16 : NULLIFY (tmp_basis_set)
1009 16 : qs_kind => qs_kind_set(ikind)
1010 16 : CALL get_qs_kind(qs_kind, basis_set=rhoin_basis, basis_type="RHOIN")
1011 22 : IF (.NOT. (ASSOCIATED(rhoin_basis))) THEN
1012 : ! Generate a default basis
1013 16 : CALL create_ri_aux_basis_set(tmp_basis_set, qs_kind, dft_control%auto_basis_ri_hxc)
1014 16 : IF (qs_env%harris_env%density_source == hden_atomic) THEN
1015 16 : CALL create_primitive_basis_set(tmp_basis_set, rhoin_basis, lmax=0)
1016 16 : CALL deallocate_gto_basis_set(tmp_basis_set)
1017 : ELSE
1018 0 : rhoin_basis => tmp_basis_set
1019 : END IF
1020 16 : CALL add_basis_set_to_container(qs_kind%basis_sets, rhoin_basis, "RHOIN")
1021 : END IF
1022 : END DO
1023 : END IF
1024 :
1025 7404 : mp2_section => section_vals_get_subs_vals(qs_env%input, "DFT%XC%WF_CORRELATION")
1026 7404 : CALL section_vals_get(mp2_section, explicit=mp2_present)
1027 7404 : IF (mp2_present) THEN
1028 :
1029 : ! basis should be sorted for imaginary time RPA/GW
1030 470 : CALL section_vals_val_get(qs_env%input, "DFT%SORT_BASIS", i_val=sort_basis)
1031 : CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%LOW_SCALING%_SECTION_PARAMETERS_", &
1032 470 : l_val=do_wfc_im_time)
1033 :
1034 470 : IF (do_wfc_im_time .AND. sort_basis /= basis_sort_zet) THEN
1035 : CALL cp_warn(__LOCATION__, &
1036 10 : "Low-scaling RPA requires SORT_BASIS EXP keyword (in DFT input section) for good performance")
1037 : END IF
1038 :
1039 : ! Check if RI_AUX basis (for MP2/RPA) is given, auto-generate if not
1040 470 : CALL mp2_env_create(qs_env%mp2_env)
1041 470 : CALL get_qs_env(qs_env, mp2_env=mp2_env, nkind=nkind)
1042 470 : CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%RI_MP2%_SECTION_PARAMETERS_", l_val=do_ri_mp2)
1043 470 : CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%RI_SOS_MP2%_SECTION_PARAMETERS_", l_val=do_ri_sos_mp2)
1044 470 : CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%RI_RPA%_SECTION_PARAMETERS_", l_val=do_ri_rpa)
1045 470 : IF (do_ri_mp2 .OR. do_ri_sos_mp2 .OR. do_ri_rpa) THEN
1046 1264 : DO ikind = 1, nkind
1047 832 : NULLIFY (ri_aux_basis_set)
1048 832 : qs_kind => qs_kind_set(ikind)
1049 : CALL get_qs_kind(qs_kind=qs_kind, basis_set=ri_aux_basis_set, &
1050 832 : basis_type="RI_AUX")
1051 1302 : IF (.NOT. (ASSOCIATED(ri_aux_basis_set))) THEN
1052 : ! RI_AUX basis set is not yet loaded
1053 : ! Generate a default basis
1054 8 : CALL create_ri_aux_basis_set(ri_aux_basis_set, qs_kind, dft_control%auto_basis_ri_aux, basis_sort=sort_basis)
1055 8 : CALL add_basis_set_to_container(qs_kind%basis_sets, ri_aux_basis_set, "RI_AUX")
1056 : ! Add a flag, which allows to check if the basis was generated
1057 : ! when applying ERI_METHOD OS to mp2, ri-rpa, gw etc
1058 8 : qs_env%mp2_env%ri_aux_auto_generated = .TRUE.
1059 : END IF
1060 : END DO
1061 : END IF
1062 :
1063 : END IF
1064 :
1065 7404 : IF (dft_control%do_xas_tdp_calculation .OR. qs_env%do_rixs) THEN
1066 : ! Check if RI_XAS basis is given, auto-generate if not
1067 60 : CALL get_qs_env(qs_env, nkind=nkind)
1068 154 : DO ikind = 1, nkind
1069 94 : NULLIFY (ri_xas_basis)
1070 94 : qs_kind => qs_kind_set(ikind)
1071 94 : CALL get_qs_kind(qs_kind, basis_Set=ri_xas_basis, basis_type="RI_XAS")
1072 7498 : IF (.NOT. ASSOCIATED(ri_xas_basis)) THEN
1073 : ! Generate a default basis
1074 90 : CALL create_ri_aux_basis_set(ri_xas_basis, qs_kind, dft_control%auto_basis_ri_xas)
1075 90 : CALL add_basis_set_to_container(qs_kind%basis_sets, ri_xas_basis, "RI_XAS")
1076 : END IF
1077 : END DO
1078 : END IF
1079 :
1080 : ! *** Initialize the spherical harmonics and ***
1081 : ! *** the orbital transformation matrices ***
1082 7404 : CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto, maxlppl=maxlppl, maxlppnl=maxlppnl)
1083 :
1084 7404 : lmax_sphere = dft_control%qs_control%gapw_control%lmax_sphere
1085 7404 : IF (lmax_sphere .LT. 0) THEN
1086 7286 : lmax_sphere = 2*maxlgto
1087 7286 : dft_control%qs_control%gapw_control%lmax_sphere = lmax_sphere
1088 : END IF
1089 7404 : IF (dft_control%qs_control%method_id == do_method_lrigpw .OR. dft_control%qs_control%lri_optbas) THEN
1090 46 : CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto_lri, basis_type="LRI_AUX")
1091 : !take maxlgto from lri basis if larger (usually)
1092 46 : maxlgto = MAX(maxlgto, maxlgto_lri)
1093 7358 : ELSE IF (dft_control%qs_control%method_id == do_method_rigpw) THEN
1094 0 : CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto_lri, basis_type="RI_HXC")
1095 0 : maxlgto = MAX(maxlgto, maxlgto_lri)
1096 : END IF
1097 7404 : IF (dft_control%do_xas_tdp_calculation .OR. qs_env%do_rixs) THEN
1098 : !done as a precaution
1099 60 : CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto_lri, basis_type="RI_XAS")
1100 60 : maxlgto = MAX(maxlgto, maxlgto_lri)
1101 : END IF
1102 7404 : maxl = MAX(2*maxlgto, maxlppl, maxlppnl, lmax_sphere) + 1
1103 :
1104 7404 : CALL init_orbital_pointers(maxl)
1105 :
1106 7404 : CALL init_spherical_harmonics(maxl, 0)
1107 :
1108 : ! *** Initialise the qs_kind_set ***
1109 7404 : CALL init_qs_kind_set(qs_kind_set)
1110 :
1111 : ! *** Initialise GAPW soft basis and projectors
1112 7404 : IF (dft_control%qs_control%method_id == do_method_gapw .OR. &
1113 : dft_control%qs_control%method_id == do_method_gapw_xc) THEN
1114 938 : qs_control => dft_control%qs_control
1115 938 : CALL init_gapw_basis_set(qs_kind_set, qs_control, qs_env%input)
1116 : END IF
1117 :
1118 : ! *** Initialize the pretabulation for the calculation of the ***
1119 : ! *** incomplete Gamma function F_n(t) after McMurchie-Davidson ***
1120 7404 : CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto)
1121 7404 : maxl = MAX(3*maxlgto + 1, 0)
1122 7404 : CALL init_md_ftable(maxl)
1123 :
1124 : ! *** Initialize the atomic interaction radii ***
1125 7404 : CALL init_interaction_radii(dft_control%qs_control, qs_kind_set)
1126 : !
1127 7404 : IF (dft_control%qs_control%method_id == do_method_xtb) THEN
1128 984 : IF (.NOT. dft_control%qs_control%xtb_control%do_tblite) THEN
1129 : ! cutoff radius
1130 3040 : DO ikind = 1, nkind
1131 2096 : qs_kind => qs_kind_set(ikind)
1132 3040 : IF (qs_kind%xtb_parameter%defined) THEN
1133 2094 : CALL get_qs_kind(qs_kind, basis_set=tmp_basis_set)
1134 2094 : rcut = xtb_control%coulomb_sr_cut
1135 2094 : fxx = 2.0_dp*xtb_control%coulomb_sr_eps*qs_kind%xtb_parameter%eta**2
1136 2094 : fxx = 0.80_dp*(1.0_dp/fxx)**0.3333_dp
1137 2094 : rcut = MIN(rcut, xtb_control%coulomb_sr_cut)
1138 2094 : qs_kind%xtb_parameter%rcut = MIN(rcut, fxx)
1139 : ELSE
1140 2 : qs_kind%xtb_parameter%rcut = 0.0_dp
1141 : END IF
1142 : END DO
1143 : END IF
1144 : END IF
1145 :
1146 7404 : IF (.NOT. be_silent) THEN
1147 7398 : CALL write_pgf_orb_radii("orb", atomic_kind_set, qs_kind_set, subsys_section)
1148 7398 : CALL write_pgf_orb_radii("aux", atomic_kind_set, qs_kind_set, subsys_section)
1149 7398 : CALL write_pgf_orb_radii("lri", atomic_kind_set, qs_kind_set, subsys_section)
1150 7398 : CALL write_core_charge_radii(atomic_kind_set, qs_kind_set, subsys_section)
1151 7398 : CALL write_ppl_radii(atomic_kind_set, qs_kind_set, subsys_section)
1152 7398 : CALL write_ppnl_radii(atomic_kind_set, qs_kind_set, subsys_section)
1153 7398 : CALL write_paw_radii(atomic_kind_set, qs_kind_set, subsys_section)
1154 : END IF
1155 :
1156 : ! *** Distribute molecules and atoms using the new data structures ***
1157 : CALL distribute_molecules_1d(atomic_kind_set=atomic_kind_set, &
1158 : particle_set=particle_set, &
1159 : local_particles=local_particles, &
1160 : molecule_kind_set=molecule_kind_set, &
1161 : molecule_set=molecule_set, &
1162 : local_molecules=local_molecules, &
1163 7404 : force_env_section=qs_env%input)
1164 :
1165 : ! *** SCF parameters ***
1166 214716 : ALLOCATE (scf_control)
1167 : ! set (non)-self consistency
1168 7404 : IF (dft_control%qs_control%dftb) THEN
1169 222 : scf_control%non_selfconsistent = .NOT. dft_control%qs_control%dftb_control%self_consistent
1170 : END IF
1171 7404 : IF (dft_control%qs_control%xtb) THEN
1172 984 : IF (dft_control%qs_control%xtb_control%do_tblite) THEN
1173 40 : scf_control%non_selfconsistent = .FALSE.
1174 : ELSE
1175 944 : scf_control%non_selfconsistent = (dft_control%qs_control%xtb_control%gfn_type == 0)
1176 : END IF
1177 : END IF
1178 7404 : IF (qs_env%harris_method) THEN
1179 6 : scf_control%non_selfconsistent = .TRUE.
1180 : END IF
1181 7404 : CALL scf_c_create(scf_control)
1182 7404 : CALL scf_c_read_parameters(scf_control, dft_section)
1183 : ! *** Allocate the data structure for Quickstep energies ***
1184 7404 : CALL allocate_qs_energy(energy)
1185 :
1186 : ! check for orthogonal basis
1187 7404 : has_unit_metric = .FALSE.
1188 7404 : IF (dft_control%qs_control%semi_empirical) THEN
1189 998 : IF (dft_control%qs_control%se_control%orthogonal_basis) has_unit_metric = .TRUE.
1190 : END IF
1191 7404 : IF (dft_control%qs_control%dftb) THEN
1192 222 : IF (dft_control%qs_control%dftb_control%orthogonal_basis) has_unit_metric = .TRUE.
1193 : END IF
1194 7404 : CALL set_qs_env(qs_env, has_unit_metric=has_unit_metric)
1195 :
1196 : ! *** Activate the interpolation ***
1197 : CALL wfi_create(wf_history, &
1198 : interpolation_method_nr= &
1199 : dft_control%qs_control%wf_interpolation_method_nr, &
1200 : extrapolation_order=dft_control%qs_control%wf_extrapolation_order, &
1201 7404 : has_unit_metric=has_unit_metric)
1202 :
1203 : ! *** Set the current Quickstep environment ***
1204 : CALL set_qs_env(qs_env=qs_env, &
1205 : scf_control=scf_control, &
1206 7404 : wf_history=wf_history)
1207 :
1208 : CALL qs_subsys_set(subsys, &
1209 : cell_ref=cell_ref, &
1210 : use_ref_cell=use_ref_cell, &
1211 : energy=energy, &
1212 7404 : force=force)
1213 :
1214 7404 : CALL get_qs_env(qs_env, ks_env=ks_env)
1215 7404 : CALL set_ks_env(ks_env, dft_control=dft_control)
1216 :
1217 : CALL qs_subsys_set(subsys, local_molecules=local_molecules, &
1218 7404 : local_particles=local_particles, cell=cell)
1219 :
1220 7404 : CALL distribution_1d_release(local_particles)
1221 7404 : CALL distribution_1d_release(local_molecules)
1222 7404 : CALL wfi_release(wf_history)
1223 :
1224 : CALL get_qs_env(qs_env=qs_env, &
1225 : atomic_kind_set=atomic_kind_set, &
1226 : dft_control=dft_control, &
1227 7404 : scf_control=scf_control)
1228 :
1229 : ! decide what conditions need mo_derivs
1230 : ! right now, this only appears to be OT
1231 7404 : IF (dft_control%qs_control%do_ls_scf .OR. &
1232 : dft_control%qs_control%do_almo_scf) THEN
1233 404 : CALL set_qs_env(qs_env=qs_env, requires_mo_derivs=.FALSE.)
1234 : ELSE
1235 7000 : IF (scf_control%use_ot) THEN
1236 2018 : CALL set_qs_env(qs_env=qs_env, requires_mo_derivs=.TRUE.)
1237 : ELSE
1238 4982 : CALL set_qs_env(qs_env=qs_env, requires_mo_derivs=.FALSE.)
1239 : END IF
1240 : END IF
1241 :
1242 : ! XXXXXXX this is backwards XXXXXXXX
1243 7404 : dft_control%smear = scf_control%smear%do_smear
1244 :
1245 : ! Periodic efield needs equal occupation and orbital gradients
1246 7404 : IF (.NOT. (dft_control%qs_control%dftb .OR. dft_control%qs_control%xtb)) THEN
1247 6198 : IF (dft_control%apply_period_efield) THEN
1248 30 : CALL get_qs_env(qs_env=qs_env, requires_mo_derivs=orb_gradient)
1249 30 : IF (.NOT. orb_gradient) THEN
1250 : CALL cp_abort(__LOCATION__, "Periodic Efield needs orbital gradient and direct optimization."// &
1251 0 : " Use the OT optimization method.")
1252 : END IF
1253 30 : IF (dft_control%smear) THEN
1254 : CALL cp_abort(__LOCATION__, "Periodic Efield needs equal occupation numbers."// &
1255 0 : " Smearing option is not possible.")
1256 : END IF
1257 : END IF
1258 : END IF
1259 :
1260 : ! Initialize the GAPW local densities and potentials
1261 7404 : IF (dft_control%qs_control%method_id == do_method_gapw .OR. &
1262 : dft_control%qs_control%method_id == do_method_gapw_xc) THEN
1263 : ! *** Allocate and initialize the set of atomic densities ***
1264 938 : NULLIFY (rho_atom_set)
1265 938 : gapw_control => dft_control%qs_control%gapw_control
1266 938 : CALL init_rho_atom(rho_atom_set, atomic_kind_set, qs_kind_set, dft_control, para_env)
1267 938 : CALL set_qs_env(qs_env=qs_env, rho_atom_set=rho_atom_set)
1268 938 : IF (dft_control%qs_control%method_id /= do_method_gapw_xc) THEN
1269 828 : CALL get_qs_env(qs_env=qs_env, local_rho_set=local_rho_set, natom=natom)
1270 : ! *** Allocate and initialize the compensation density rho0 ***
1271 828 : CALL init_rho0(local_rho_set, qs_env, gapw_control)
1272 : ! *** Allocate and Initialize the local coulomb term ***
1273 828 : CALL init_coulomb_local(qs_env%hartree_local, natom)
1274 : END IF
1275 : ! NLCC
1276 938 : CALL init_gapw_nlcc(qs_kind_set)
1277 6466 : ELSE IF (dft_control%qs_control%method_id == do_method_lrigpw) THEN
1278 : ! allocate local ri environment
1279 : ! nothing to do here?
1280 6426 : ELSE IF (dft_control%qs_control%method_id == do_method_rigpw) THEN
1281 : ! allocate ri environment
1282 : ! nothing to do here?
1283 6426 : ELSE IF (dft_control%qs_control%semi_empirical) THEN
1284 998 : NULLIFY (se_store_int_env, se_nddo_mpole, se_nonbond_env)
1285 998 : natom = SIZE(particle_set)
1286 998 : se_section => section_vals_get_subs_vals(qs_section, "SE")
1287 998 : se_control => dft_control%qs_control%se_control
1288 :
1289 : ! Make the cutoff radii choice a bit smarter
1290 998 : CALL se_cutoff_compatible(se_control, se_section, cell, output_unit)
1291 :
1292 1994 : SELECT CASE (dft_control%qs_control%method_id)
1293 : CASE DEFAULT
1294 : CASE (do_method_rm1, do_method_am1, do_method_mndo, do_method_pm3, &
1295 : do_method_pm6, do_method_pm6fm, do_method_mndod, do_method_pnnl)
1296 : ! Neighbor lists have to be MAX(interaction range, orbital range)
1297 : ! set new kind radius
1298 998 : CALL init_se_nlradius(se_control, atomic_kind_set, qs_kind_set, subsys_section)
1299 : END SELECT
1300 : ! Initialize to zero the max multipole to treat in the EWALD scheme..
1301 998 : se_control%max_multipole = do_multipole_none
1302 : ! check for Ewald
1303 998 : IF (se_control%do_ewald .OR. se_control%do_ewald_gks) THEN
1304 512 : ALLOCATE (ewald_env)
1305 32 : CALL ewald_env_create(ewald_env, para_env)
1306 32 : poisson_section => section_vals_get_subs_vals(dft_section, "POISSON")
1307 32 : CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
1308 32 : ewald_section => section_vals_get_subs_vals(poisson_section, "EWALD")
1309 : print_section => section_vals_get_subs_vals(qs_env%input, &
1310 32 : "PRINT%GRID_INFORMATION")
1311 32 : CALL read_ewald_section(ewald_env, ewald_section)
1312 : ! Create ewald grids
1313 32 : ALLOCATE (ewald_pw)
1314 : CALL ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, &
1315 32 : print_section=print_section)
1316 : ! Initialize ewald grids
1317 32 : CALL ewald_pw_grid_update(ewald_pw, ewald_env, cell%hmat)
1318 : ! Setup the nonbond environment (real space part of Ewald)
1319 32 : CALL ewald_env_get(ewald_env, rcut=ewald_rcut)
1320 : ! Setup the maximum level of multipoles to be treated in the periodic SE scheme
1321 32 : IF (se_control%do_ewald) THEN
1322 30 : CALL ewald_env_get(ewald_env, max_multipole=se_control%max_multipole)
1323 : END IF
1324 : CALL section_vals_val_get(se_section, "NEIGHBOR_LISTS%VERLET_SKIN", &
1325 32 : r_val=verlet_skin)
1326 32 : ALLOCATE (se_nonbond_env)
1327 : CALL fist_nonbond_env_create(se_nonbond_env, atomic_kind_set, do_nonbonded=.TRUE., &
1328 : do_electrostatics=.TRUE., verlet_skin=verlet_skin, ewald_rcut=ewald_rcut, &
1329 32 : ei_scale14=0.0_dp, vdw_scale14=0.0_dp, shift_cutoff=.FALSE.)
1330 : ! Create and Setup NDDO multipole environment
1331 32 : CALL nddo_mpole_setup(se_nddo_mpole, natom)
1332 : CALL set_qs_env(qs_env, ewald_env=ewald_env, ewald_pw=ewald_pw, &
1333 32 : se_nonbond_env=se_nonbond_env, se_nddo_mpole=se_nddo_mpole)
1334 : ! Handle the residual integral part 1/R^3
1335 : CALL semi_empirical_expns3_setup(qs_kind_set, se_control, &
1336 32 : dft_control%qs_control%method_id)
1337 : END IF
1338 : ! Taper function
1339 : CALL se_taper_create(se_taper, se_control%integral_screening, se_control%do_ewald, &
1340 : se_control%taper_cou, se_control%range_cou, &
1341 : se_control%taper_exc, se_control%range_exc, &
1342 : se_control%taper_scr, se_control%range_scr, &
1343 998 : se_control%taper_lrc, se_control%range_lrc)
1344 998 : CALL set_qs_env(qs_env, se_taper=se_taper)
1345 : ! Store integral environment
1346 998 : CALL semi_empirical_si_create(se_store_int_env, se_section)
1347 998 : CALL set_qs_env(qs_env, se_store_int_env=se_store_int_env)
1348 : END IF
1349 :
1350 : ! Initialize possible dispersion parameters
1351 : IF (dft_control%qs_control%method_id == do_method_gpw .OR. &
1352 : dft_control%qs_control%method_id == do_method_gapw .OR. &
1353 : dft_control%qs_control%method_id == do_method_gapw_xc .OR. &
1354 : dft_control%qs_control%method_id == do_method_lrigpw .OR. &
1355 7404 : dft_control%qs_control%method_id == do_method_rigpw .OR. &
1356 : dft_control%qs_control%method_id == do_method_ofgpw) THEN
1357 26000 : ALLOCATE (dispersion_env)
1358 5200 : NULLIFY (xc_section)
1359 5200 : xc_section => section_vals_get_subs_vals(dft_section, "XC")
1360 5200 : CALL qs_dispersion_env_set(dispersion_env, xc_section)
1361 5200 : IF (dispersion_env%type == xc_vdw_fun_pairpot) THEN
1362 112 : NULLIFY (pp_section)
1363 112 : pp_section => section_vals_get_subs_vals(xc_section, "VDW_POTENTIAL%PAIR_POTENTIAL")
1364 112 : CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, pp_section, para_env)
1365 5088 : ELSE IF (dispersion_env%type == xc_vdw_fun_nonloc) THEN
1366 46 : NULLIFY (nl_section)
1367 46 : nl_section => section_vals_get_subs_vals(xc_section, "VDW_POTENTIAL%NON_LOCAL")
1368 46 : CALL qs_dispersion_nonloc_init(dispersion_env, para_env)
1369 : END IF
1370 5200 : CALL set_qs_env(qs_env, dispersion_env=dispersion_env)
1371 2204 : ELSE IF (dft_control%qs_control%method_id == do_method_dftb) THEN
1372 1110 : ALLOCATE (dispersion_env)
1373 : ! set general defaults
1374 : dispersion_env%doabc = .FALSE.
1375 : dispersion_env%c9cnst = .FALSE.
1376 : dispersion_env%lrc = .FALSE.
1377 : dispersion_env%srb = .FALSE.
1378 : dispersion_env%verbose = .FALSE.
1379 : NULLIFY (dispersion_env%c6ab, dispersion_env%maxci, dispersion_env%r0ab, dispersion_env%rcov, &
1380 : dispersion_env%r2r4, dispersion_env%cn, dispersion_env%cnkind, dispersion_env%cnlist, &
1381 : dispersion_env%d3_exclude_pair)
1382 : NULLIFY (dispersion_env%q_mesh, dispersion_env%kernel, dispersion_env%d2phi_dk2, &
1383 : dispersion_env%d2y_dx2, dispersion_env%dftd_section)
1384 : NULLIFY (dispersion_env%sab_vdw, dispersion_env%sab_cn)
1385 222 : IF (dftb_control%dispersion .AND. dftb_control%dispersion_type == dispersion_d3) THEN
1386 14 : dispersion_env%type = xc_vdw_fun_pairpot
1387 14 : dispersion_env%pp_type = vdw_pairpot_dftd3
1388 14 : dispersion_env%eps_cn = dftb_control%epscn
1389 14 : dispersion_env%s6 = dftb_control%sd3(1)
1390 14 : dispersion_env%sr6 = dftb_control%sd3(2)
1391 14 : dispersion_env%s8 = dftb_control%sd3(3)
1392 : dispersion_env%domol = .FALSE.
1393 14 : dispersion_env%kgc8 = 0._dp
1394 14 : dispersion_env%rc_disp = dftb_control%rcdisp
1395 14 : dispersion_env%exp_pre = 0._dp
1396 14 : dispersion_env%scaling = 0._dp
1397 14 : dispersion_env%nd3_exclude_pair = 0
1398 14 : dispersion_env%parameter_file_name = dftb_control%dispersion_parameter_file
1399 14 : CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, para_env=para_env)
1400 208 : ELSEIF (dftb_control%dispersion .AND. dftb_control%dispersion_type == dispersion_d3bj) THEN
1401 2 : dispersion_env%type = xc_vdw_fun_pairpot
1402 2 : dispersion_env%pp_type = vdw_pairpot_dftd3bj
1403 2 : dispersion_env%eps_cn = dftb_control%epscn
1404 2 : dispersion_env%s6 = dftb_control%sd3bj(1)
1405 2 : dispersion_env%a1 = dftb_control%sd3bj(2)
1406 2 : dispersion_env%s8 = dftb_control%sd3bj(3)
1407 2 : dispersion_env%a2 = dftb_control%sd3bj(4)
1408 : dispersion_env%domol = .FALSE.
1409 2 : dispersion_env%kgc8 = 0._dp
1410 2 : dispersion_env%rc_disp = dftb_control%rcdisp
1411 2 : dispersion_env%exp_pre = 0._dp
1412 2 : dispersion_env%scaling = 0._dp
1413 2 : dispersion_env%nd3_exclude_pair = 0
1414 2 : dispersion_env%parameter_file_name = dftb_control%dispersion_parameter_file
1415 2 : CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, para_env=para_env)
1416 206 : ELSEIF (dftb_control%dispersion .AND. dftb_control%dispersion_type == dispersion_d2) THEN
1417 2 : dispersion_env%type = xc_vdw_fun_pairpot
1418 2 : dispersion_env%pp_type = vdw_pairpot_dftd2
1419 2 : dispersion_env%exp_pre = dftb_control%exp_pre
1420 2 : dispersion_env%scaling = dftb_control%scaling
1421 2 : dispersion_env%parameter_file_name = dftb_control%dispersion_parameter_file
1422 2 : dispersion_env%rc_disp = dftb_control%rcdisp
1423 2 : CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, para_env=para_env)
1424 : ELSE
1425 204 : dispersion_env%type = xc_vdw_fun_none
1426 : END IF
1427 222 : CALL set_qs_env(qs_env, dispersion_env=dispersion_env)
1428 1982 : ELSE IF (dft_control%qs_control%method_id == do_method_xtb) THEN
1429 984 : IF (.NOT. (dft_control%qs_control%xtb_control%do_tblite)) THEN
1430 4720 : ALLOCATE (dispersion_env)
1431 : ! set general defaults
1432 : dispersion_env%doabc = .FALSE.
1433 : dispersion_env%c9cnst = .FALSE.
1434 : dispersion_env%lrc = .FALSE.
1435 : dispersion_env%srb = .FALSE.
1436 : dispersion_env%verbose = .FALSE.
1437 : NULLIFY (dispersion_env%c6ab, dispersion_env%maxci, &
1438 : dispersion_env%r0ab, dispersion_env%rcov, &
1439 : dispersion_env%r2r4, dispersion_env%cn, &
1440 : dispersion_env%cnkind, dispersion_env%cnlist, &
1441 : dispersion_env%d3_exclude_pair)
1442 : NULLIFY (dispersion_env%q_mesh, dispersion_env%kernel, dispersion_env%d2phi_dk2, &
1443 : dispersion_env%d2y_dx2, dispersion_env%dftd_section)
1444 : NULLIFY (dispersion_env%sab_vdw, dispersion_env%sab_cn)
1445 944 : dispersion_env%type = xc_vdw_fun_pairpot
1446 944 : dispersion_env%eps_cn = xtb_control%epscn
1447 944 : dispersion_env%s6 = xtb_control%s6
1448 944 : dispersion_env%s8 = xtb_control%s8
1449 944 : dispersion_env%a1 = xtb_control%a1
1450 944 : dispersion_env%a2 = xtb_control%a2
1451 : dispersion_env%domol = .FALSE.
1452 944 : dispersion_env%kgc8 = 0._dp
1453 944 : dispersion_env%rc_disp = xtb_control%rcdisp
1454 944 : dispersion_env%rc_d4 = xtb_control%rcdisp
1455 944 : dispersion_env%exp_pre = 0._dp
1456 944 : dispersion_env%scaling = 0._dp
1457 944 : dispersion_env%nd3_exclude_pair = 0
1458 944 : dispersion_env%parameter_file_name = xtb_control%dispersion_parameter_file
1459 : !
1460 1248 : SELECT CASE (xtb_control%vdw_type)
1461 : CASE (xtb_vdw_type_none, xtb_vdw_type_d3)
1462 304 : dispersion_env%pp_type = vdw_pairpot_dftd3bj
1463 304 : CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, para_env=para_env)
1464 304 : IF (xtb_control%vdw_type == xtb_vdw_type_none) dispersion_env%type = xc_vdw_fun_none
1465 : CASE (xtb_vdw_type_d4)
1466 640 : dispersion_env%pp_type = vdw_pairpot_dftd4
1467 640 : dispersion_env%ref_functional = "none"
1468 : CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, &
1469 640 : dispersion_env, para_env=para_env)
1470 640 : dispersion_env%cnfun = 2
1471 : CASE DEFAULT
1472 944 : CPABORT("vdw type")
1473 : END SELECT
1474 944 : CALL set_qs_env(qs_env, dispersion_env=dispersion_env)
1475 : END IF
1476 998 : ELSE IF (dft_control%qs_control%semi_empirical) THEN
1477 4990 : ALLOCATE (dispersion_env)
1478 : ! set general defaults
1479 : dispersion_env%doabc = .FALSE.
1480 : dispersion_env%c9cnst = .FALSE.
1481 : dispersion_env%lrc = .FALSE.
1482 : dispersion_env%srb = .FALSE.
1483 : dispersion_env%verbose = .FALSE.
1484 : NULLIFY (dispersion_env%c6ab, dispersion_env%maxci, dispersion_env%r0ab, dispersion_env%rcov, &
1485 : dispersion_env%r2r4, dispersion_env%cn, dispersion_env%cnkind, dispersion_env%cnlist, &
1486 : dispersion_env%d3_exclude_pair)
1487 : NULLIFY (dispersion_env%q_mesh, dispersion_env%kernel, dispersion_env%d2phi_dk2, &
1488 : dispersion_env%d2y_dx2, dispersion_env%dftd_section)
1489 : NULLIFY (dispersion_env%sab_vdw, dispersion_env%sab_cn)
1490 998 : IF (se_control%dispersion) THEN
1491 6 : dispersion_env%type = xc_vdw_fun_pairpot
1492 6 : dispersion_env%pp_type = vdw_pairpot_dftd3
1493 6 : dispersion_env%eps_cn = se_control%epscn
1494 6 : dispersion_env%s6 = se_control%sd3(1)
1495 6 : dispersion_env%sr6 = se_control%sd3(2)
1496 6 : dispersion_env%s8 = se_control%sd3(3)
1497 : dispersion_env%domol = .FALSE.
1498 6 : dispersion_env%kgc8 = 0._dp
1499 6 : dispersion_env%rc_disp = se_control%rcdisp
1500 6 : dispersion_env%exp_pre = 0._dp
1501 6 : dispersion_env%scaling = 0._dp
1502 6 : dispersion_env%nd3_exclude_pair = 0
1503 6 : dispersion_env%parameter_file_name = se_control%dispersion_parameter_file
1504 6 : CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, para_env=para_env)
1505 : ELSE
1506 992 : dispersion_env%type = xc_vdw_fun_none
1507 : END IF
1508 998 : CALL set_qs_env(qs_env, dispersion_env=dispersion_env)
1509 : END IF
1510 :
1511 : ! Initialize possible geomertical counterpoise correction potential
1512 : IF (dft_control%qs_control%method_id == do_method_gpw .OR. &
1513 : dft_control%qs_control%method_id == do_method_gapw .OR. &
1514 : dft_control%qs_control%method_id == do_method_gapw_xc .OR. &
1515 : dft_control%qs_control%method_id == do_method_lrigpw .OR. &
1516 7404 : dft_control%qs_control%method_id == do_method_rigpw .OR. &
1517 : dft_control%qs_control%method_id == do_method_ofgpw) THEN
1518 5200 : ALLOCATE (gcp_env)
1519 5200 : NULLIFY (xc_section)
1520 5200 : xc_section => section_vals_get_subs_vals(dft_section, "XC")
1521 5200 : CALL qs_gcp_env_set(gcp_env, xc_section)
1522 5200 : CALL qs_gcp_init(qs_env, gcp_env)
1523 5200 : CALL set_qs_env(qs_env, gcp_env=gcp_env)
1524 : END IF
1525 :
1526 : ! *** Allocate the MO data types ***
1527 7404 : CALL get_qs_kind_set(qs_kind_set, nsgf=n_ao, nelectron=nelectron)
1528 :
1529 : ! the total number of electrons
1530 7404 : nelectron = nelectron - dft_control%charge
1531 :
1532 7404 : IF (dft_control%multiplicity == 0) THEN
1533 6192 : IF (MODULO(nelectron, 2) == 0) THEN
1534 5717 : dft_control%multiplicity = 1
1535 : ELSE
1536 475 : dft_control%multiplicity = 2
1537 : END IF
1538 : END IF
1539 :
1540 7404 : multiplicity = dft_control%multiplicity
1541 :
1542 7404 : IF ((dft_control%nspins < 1) .OR. (dft_control%nspins > 2)) THEN
1543 0 : CPABORT("nspins should be 1 or 2 for the time being ...")
1544 : END IF
1545 :
1546 7404 : IF ((MODULO(nelectron, 2) /= 0) .AND. (dft_control%nspins == 1)) THEN
1547 12 : IF (.NOT. dft_control%qs_control%ofgpw .AND. .NOT. dft_control%smear) THEN
1548 0 : CPABORT("Use the LSD option for an odd number of electrons")
1549 : END IF
1550 : END IF
1551 :
1552 : ! The transition potential method to calculate XAS needs LSD
1553 7404 : IF (dft_control%do_xas_calculation) THEN
1554 42 : IF (dft_control%nspins == 1) THEN
1555 0 : CPABORT("Use the LSD option for XAS with transition potential")
1556 : END IF
1557 : END IF
1558 :
1559 : ! assigning the number of states per spin initial version, not yet very
1560 : ! general. Should work for an even number of electrons and a single
1561 : ! additional electron this set of options that requires full matrices,
1562 : ! however, makes things a bit ugly right now.... we try to make a
1563 : ! distinction between the number of electrons per spin and the number of
1564 : ! MOs per spin this should allow the use of fractional occupations later
1565 : ! on
1566 7404 : IF (dft_control%qs_control%ofgpw) THEN
1567 :
1568 0 : IF (dft_control%nspins == 1) THEN
1569 0 : maxocc = nelectron
1570 0 : nelectron_spin(1) = nelectron
1571 0 : nelectron_spin(2) = 0
1572 0 : n_mo(1) = 1
1573 0 : n_mo(2) = 0
1574 : ELSE
1575 0 : IF (MODULO(nelectron + multiplicity - 1, 2) /= 0) THEN
1576 0 : CPABORT("LSD: try to use a different multiplicity")
1577 : END IF
1578 0 : nelectron_spin(1) = (nelectron + multiplicity - 1)/2
1579 0 : nelectron_spin(2) = (nelectron - multiplicity + 1)/2
1580 0 : IF (nelectron_spin(1) < 0) THEN
1581 0 : CPABORT("LSD: too few electrons for this multiplicity")
1582 : END IF
1583 0 : maxocc = MAXVAL(nelectron_spin)
1584 0 : n_mo(1) = MIN(nelectron_spin(1), 1)
1585 0 : n_mo(2) = MIN(nelectron_spin(2), 1)
1586 : END IF
1587 :
1588 : ELSE
1589 :
1590 7404 : IF (dft_control%nspins == 1) THEN
1591 5799 : maxocc = 2.0_dp
1592 5799 : nelectron_spin(1) = nelectron
1593 5799 : nelectron_spin(2) = 0
1594 5799 : IF (MODULO(nelectron, 2) == 0) THEN
1595 5787 : n_mo(1) = nelectron/2
1596 : ELSE
1597 12 : n_mo(1) = INT(nelectron/2._dp) + 1
1598 : END IF
1599 5799 : n_mo(2) = 0
1600 : ELSE
1601 1605 : maxocc = 1.0_dp
1602 :
1603 : ! The simplist spin distribution is written here. Special cases will
1604 : ! need additional user input
1605 1605 : IF (MODULO(nelectron + multiplicity - 1, 2) /= 0) THEN
1606 0 : CPABORT("LSD: try to use a different multiplicity")
1607 : END IF
1608 :
1609 1605 : nelectron_spin(1) = (nelectron + multiplicity - 1)/2
1610 1605 : nelectron_spin(2) = (nelectron - multiplicity + 1)/2
1611 :
1612 1605 : IF (nelectron_spin(2) < 0) THEN
1613 0 : CPABORT("LSD: too few electrons for this multiplicity")
1614 : END IF
1615 :
1616 1605 : n_mo(1) = nelectron_spin(1)
1617 1605 : n_mo(2) = nelectron_spin(2)
1618 :
1619 : END IF
1620 :
1621 : END IF
1622 :
1623 : ! Read the total_zeff_corr here [SGh]
1624 7404 : CALL get_qs_kind_set(qs_kind_set, total_zeff_corr=total_zeff_corr)
1625 : ! store it in qs_env
1626 7404 : qs_env%total_zeff_corr = total_zeff_corr
1627 :
1628 : ! store the number of electrons once an for all
1629 : CALL qs_subsys_set(subsys, &
1630 : nelectron_total=nelectron, &
1631 7404 : nelectron_spin=nelectron_spin)
1632 :
1633 : ! Check and set number of added (unoccupied) MOs
1634 7404 : IF (dft_control%nspins == 2) THEN
1635 1605 : IF (scf_control%added_mos(2) < 0) THEN
1636 128 : n_mo_add = n_ao - n_mo(2) ! use all available MOs
1637 1477 : ELSEIF (scf_control%added_mos(2) > 0) THEN
1638 : n_mo_add = scf_control%added_mos(2)
1639 : ELSE
1640 1331 : n_mo_add = scf_control%added_mos(1)
1641 : END IF
1642 1605 : IF (n_mo_add > n_ao - n_mo(2)) THEN
1643 18 : CPWARN("More ADDED_MOs requested for beta spin than available.")
1644 : END IF
1645 1605 : scf_control%added_mos(2) = MIN(n_mo_add, n_ao - n_mo(2))
1646 1605 : n_mo(2) = n_mo(2) + scf_control%added_mos(2)
1647 : END IF
1648 :
1649 : ! proceed alpha orbitals after the beta orbitals; this is essential to avoid
1650 : ! reduction in the number of available unoccupied molecular orbitals.
1651 : ! E.g. n_ao = 10, nelectrons = 10, multiplicity = 3 implies n_mo(1) = 6, n_mo(2) = 4;
1652 : ! added_mos(1:2) = (6,undef) should increase the number of molecular orbitals as
1653 : ! n_mo(1) = min(n_ao, n_mo(1) + added_mos(1)) = 10, n_mo(2) = 10.
1654 : ! However, if we try to proceed alpha orbitals first, this leads us n_mo(1:2) = (10,8)
1655 : ! due to the following assignment instruction above:
1656 : ! IF (scf_control%added_mos(2) > 0) THEN ... ELSE; n_mo_add = scf_control%added_mos(1); END IF
1657 7404 : IF (scf_control%added_mos(1) < 0) THEN
1658 682 : scf_control%added_mos(1) = n_ao - n_mo(1) ! use all available MOs
1659 6722 : ELSEIF (scf_control%added_mos(1) > n_ao - n_mo(1)) THEN
1660 : CALL cp_warn(__LOCATION__, &
1661 : "More added MOs requested than available. "// &
1662 : "The full set of unoccupied MOs will be used. "// &
1663 : "Use 'ADDED_MOS -1' to always use all available MOs "// &
1664 98 : "and to get rid of this warning.")
1665 : END IF
1666 7404 : scf_control%added_mos(1) = MIN(scf_control%added_mos(1), n_ao - n_mo(1))
1667 7404 : n_mo(1) = n_mo(1) + scf_control%added_mos(1)
1668 :
1669 7404 : IF (dft_control%nspins == 2) THEN
1670 1605 : IF (n_mo(2) > n_mo(1)) &
1671 : CALL cp_warn(__LOCATION__, &
1672 : "More beta than alpha MOs requested. "// &
1673 0 : "The number of beta MOs will be reduced to the number alpha MOs.")
1674 1605 : n_mo(2) = MIN(n_mo(1), n_mo(2))
1675 1605 : CPASSERT(n_mo(1) >= nelectron_spin(1))
1676 1605 : CPASSERT(n_mo(2) >= nelectron_spin(2))
1677 : END IF
1678 :
1679 : ! kpoints
1680 7404 : CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints)
1681 7404 : IF (do_kpoints .AND. dft_control%nspins == 2) THEN
1682 : ! we need equal number of calculated states
1683 22 : IF (n_mo(2) /= n_mo(1)) &
1684 : CALL cp_warn(__LOCATION__, &
1685 : "Kpoints: Different number of MOs requested. "// &
1686 2 : "The number of beta MOs will be set to the number alpha MOs.")
1687 22 : n_mo(2) = n_mo(1)
1688 22 : CPASSERT(n_mo(1) >= nelectron_spin(1))
1689 22 : CPASSERT(n_mo(2) >= nelectron_spin(2))
1690 : END IF
1691 :
1692 : ! Compatibility checks for smearing
1693 7404 : IF (scf_control%smear%do_smear) THEN
1694 914 : IF (scf_control%added_mos(1) == 0) THEN
1695 0 : CPABORT("Extra MOs (ADDED_MOS) are required for smearing")
1696 : END IF
1697 : END IF
1698 :
1699 : ! *** Some options require that all MOs are computed ... ***
1700 : IF (BTEST(cp_print_key_should_output(logger%iter_info, dft_section, &
1701 : "PRINT%MO/CARTESIAN"), &
1702 : cp_p_file) .OR. &
1703 : (scf_control%level_shift /= 0.0_dp) .OR. &
1704 7404 : (scf_control%diagonalization%eps_jacobi /= 0.0_dp) .OR. &
1705 : (dft_control%roks .AND. (.NOT. scf_control%use_ot))) THEN
1706 7548 : n_mo(:) = n_ao
1707 : END IF
1708 :
1709 : ! Compatibility checks for ROKS
1710 7404 : IF (dft_control%roks .AND. (.NOT. scf_control%use_ot)) THEN
1711 36 : IF (scf_control%roks_scheme == general_roks) THEN
1712 0 : CPWARN("General ROKS scheme is not yet tested!")
1713 : END IF
1714 36 : IF (scf_control%smear%do_smear) THEN
1715 : CALL cp_abort(__LOCATION__, &
1716 : "The options ROKS and SMEAR are not compatible. "// &
1717 0 : "Try UKS instead of ROKS")
1718 : END IF
1719 : END IF
1720 7404 : IF (dft_control%low_spin_roks) THEN
1721 8 : SELECT CASE (dft_control%qs_control%method_id)
1722 : CASE DEFAULT
1723 : CASE (do_method_xtb, do_method_dftb)
1724 : CALL cp_abort(__LOCATION__, &
1725 0 : "xTB/DFTB methods are not compatible with low spin ROKS.")
1726 : CASE (do_method_rm1, do_method_am1, do_method_mndo, do_method_pm3, &
1727 : do_method_pm6, do_method_pm6fm, do_method_mndod, do_method_pnnl)
1728 : CALL cp_abort(__LOCATION__, &
1729 8 : "SE methods are not compatible with low spin ROKS.")
1730 : END SELECT
1731 : END IF
1732 :
1733 : ! in principle the restricted calculation could be performed
1734 : ! using just one set of MOs and special casing most of the code
1735 : ! right now we'll just take care of what is effectively an additional constraint
1736 : ! at as few places as possible, just duplicating the beta orbitals
1737 7404 : IF (dft_control%restricted .AND. (output_unit > 0)) THEN
1738 : ! it is really not yet tested till the end ! Joost
1739 23 : WRITE (output_unit, *) ""
1740 23 : WRITE (output_unit, *) " **************************************"
1741 23 : WRITE (output_unit, *) " restricted calculation cutting corners"
1742 23 : WRITE (output_unit, *) " experimental feature, check code "
1743 23 : WRITE (output_unit, *) " **************************************"
1744 : END IF
1745 :
1746 : ! no point in allocating these things here ?
1747 7404 : IF (dft_control%qs_control%do_ls_scf) THEN
1748 338 : NULLIFY (mos)
1749 : ELSE
1750 29855 : ALLOCATE (mos(dft_control%nspins))
1751 15723 : DO ispin = 1, dft_control%nspins
1752 : CALL allocate_mo_set(mo_set=mos(ispin), &
1753 : nao=n_ao, &
1754 : nmo=n_mo(ispin), &
1755 : nelectron=nelectron_spin(ispin), &
1756 : n_el_f=REAL(nelectron_spin(ispin), dp), &
1757 : maxocc=maxocc, &
1758 15723 : flexible_electron_count=dft_control%relax_multiplicity)
1759 : END DO
1760 : END IF
1761 :
1762 7404 : CALL set_qs_env(qs_env, mos=mos)
1763 :
1764 : ! allocate mos when switch_surf_dip is triggered [SGh]
1765 7404 : IF (dft_control%switch_surf_dip) THEN
1766 8 : ALLOCATE (mos_last_converged(dft_control%nspins))
1767 4 : DO ispin = 1, dft_control%nspins
1768 : CALL allocate_mo_set(mo_set=mos_last_converged(ispin), &
1769 : nao=n_ao, &
1770 : nmo=n_mo(ispin), &
1771 : nelectron=nelectron_spin(ispin), &
1772 : n_el_f=REAL(nelectron_spin(ispin), dp), &
1773 : maxocc=maxocc, &
1774 4 : flexible_electron_count=dft_control%relax_multiplicity)
1775 : END DO
1776 2 : CALL set_qs_env(qs_env, mos_last_converged=mos_last_converged)
1777 : END IF
1778 :
1779 7404 : IF (.NOT. be_silent) THEN
1780 : ! Print the DFT control parameters
1781 7398 : CALL write_dft_control(dft_control, dft_section)
1782 :
1783 : ! Print the vdW control parameters
1784 : IF (dft_control%qs_control%method_id == do_method_gpw .OR. &
1785 : dft_control%qs_control%method_id == do_method_gapw .OR. &
1786 : dft_control%qs_control%method_id == do_method_gapw_xc .OR. &
1787 : dft_control%qs_control%method_id == do_method_lrigpw .OR. &
1788 : dft_control%qs_control%method_id == do_method_rigpw .OR. &
1789 : dft_control%qs_control%method_id == do_method_dftb .OR. &
1790 : (dft_control%qs_control%method_id == do_method_xtb .AND. &
1791 7398 : (.NOT. dft_control%qs_control%xtb_control%do_tblite)) .OR. &
1792 : dft_control%qs_control%method_id == do_method_ofgpw) THEN
1793 6360 : CALL get_qs_env(qs_env, dispersion_env=dispersion_env)
1794 6360 : CALL qs_write_dispersion(qs_env, dispersion_env)
1795 : END IF
1796 :
1797 : ! Print the Quickstep control parameters
1798 7398 : CALL write_qs_control(dft_control%qs_control, dft_section)
1799 :
1800 : ! Print the ADMM control parameters
1801 7398 : IF (dft_control%do_admm) THEN
1802 450 : CALL write_admm_control(dft_control%admm_control, dft_section)
1803 : END IF
1804 :
1805 : ! Print XES/XAS control parameters
1806 7398 : IF (dft_control%do_xas_calculation) THEN
1807 42 : CALL cite_reference(Iannuzzi2007)
1808 : !CALL write_xas_control(dft_control%xas_control,dft_section)
1809 : END IF
1810 :
1811 : ! Print the unnormalized basis set information (input data)
1812 7398 : CALL write_gto_basis_sets(qs_kind_set, subsys_section)
1813 :
1814 : ! Print the atomic kind set
1815 7398 : CALL write_qs_kind_set(qs_kind_set, subsys_section)
1816 :
1817 : ! Print the molecule kind set
1818 7398 : CALL write_molecule_kind_set(molecule_kind_set, subsys_section)
1819 :
1820 : ! Print the total number of kinds, atoms, basis functions etc.
1821 7398 : CALL write_total_numbers(qs_kind_set, particle_set, qs_env%input)
1822 :
1823 : ! Print the atomic coordinates
1824 7398 : CALL write_qs_particle_coordinates(particle_set, qs_kind_set, subsys_section, label="QUICKSTEP")
1825 :
1826 : ! Print the interatomic distances
1827 7398 : CALL write_particle_distances(particle_set, cell, subsys_section)
1828 :
1829 : ! Print the requested structure data
1830 7398 : CALL write_structure_data(particle_set, cell, subsys_section)
1831 :
1832 : ! Print symmetry information
1833 7398 : CALL write_symmetry(particle_set, cell, subsys_section)
1834 :
1835 : ! Print the SCF parameters
1836 7398 : IF ((.NOT. dft_control%qs_control%do_ls_scf) .AND. &
1837 : (.NOT. dft_control%qs_control%do_almo_scf)) THEN
1838 6994 : CALL scf_c_write_parameters(scf_control, dft_section)
1839 : END IF
1840 : END IF
1841 :
1842 : ! Sets up pw_env, qs_charges, mpools ...
1843 7404 : CALL qs_env_setup(qs_env)
1844 :
1845 : ! Allocate and initialise rho0 soft on the global grid
1846 7404 : IF (dft_control%qs_control%method_id == do_method_gapw) THEN
1847 828 : CALL get_qs_env(qs_env=qs_env, pw_env=pw_env, rho0_mpole=rho0_mpole)
1848 828 : CALL rho0_s_grid_create(pw_env, rho0_mpole)
1849 : END IF
1850 :
1851 7404 : IF (output_unit > 0) CALL m_flush(output_unit)
1852 7404 : CALL timestop(handle)
1853 :
1854 74040 : END SUBROUTINE qs_init_subsys
1855 :
1856 : ! **************************************************************************************************
1857 : !> \brief Write the total number of kinds, atoms, etc. to the logical unit
1858 : !> number lunit.
1859 : !> \param qs_kind_set ...
1860 : !> \param particle_set ...
1861 : !> \param force_env_section ...
1862 : !> \author Creation (06.10.2000)
1863 : ! **************************************************************************************************
1864 7398 : SUBROUTINE write_total_numbers(qs_kind_set, particle_set, force_env_section)
1865 :
1866 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
1867 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1868 : TYPE(section_vals_type), POINTER :: force_env_section
1869 :
1870 : INTEGER :: maxlgto, maxlppl, maxlppnl, natom, ncgf, &
1871 : nkind, npgf, nset, nsgf, nshell, &
1872 : output_unit
1873 : TYPE(cp_logger_type), POINTER :: logger
1874 :
1875 7398 : NULLIFY (logger)
1876 7398 : logger => cp_get_default_logger()
1877 : output_unit = cp_print_key_unit_nr(logger, force_env_section, "PRINT%TOTAL_NUMBERS", &
1878 7398 : extension=".Log")
1879 :
1880 7398 : IF (output_unit > 0) THEN
1881 3723 : natom = SIZE(particle_set)
1882 3723 : nkind = SIZE(qs_kind_set)
1883 :
1884 : CALL get_qs_kind_set(qs_kind_set, &
1885 : maxlgto=maxlgto, &
1886 : ncgf=ncgf, &
1887 : npgf=npgf, &
1888 : nset=nset, &
1889 : nsgf=nsgf, &
1890 : nshell=nshell, &
1891 : maxlppl=maxlppl, &
1892 3723 : maxlppnl=maxlppnl)
1893 :
1894 : WRITE (UNIT=output_unit, FMT="(/,/,T2,A)") &
1895 3723 : "TOTAL NUMBERS AND MAXIMUM NUMBERS"
1896 :
1897 3723 : IF (nset + npgf + ncgf > 0) THEN
1898 : WRITE (UNIT=output_unit, FMT="(/,T3,A,(T30,A,T71,I10))") &
1899 3723 : "Total number of", &
1900 3723 : "- Atomic kinds: ", nkind, &
1901 3723 : "- Atoms: ", natom, &
1902 3723 : "- Shell sets: ", nset, &
1903 3723 : "- Shells: ", nshell, &
1904 3723 : "- Primitive Cartesian functions: ", npgf, &
1905 3723 : "- Cartesian basis functions: ", ncgf, &
1906 7446 : "- Spherical basis functions: ", nsgf
1907 0 : ELSE IF (nshell + nsgf > 0) THEN
1908 : WRITE (UNIT=output_unit, FMT="(/,T3,A,(T30,A,T71,I10))") &
1909 0 : "Total number of", &
1910 0 : "- Atomic kinds: ", nkind, &
1911 0 : "- Atoms: ", natom, &
1912 0 : "- Shells: ", nshell, &
1913 0 : "- Spherical basis functions: ", nsgf
1914 : ELSE
1915 : WRITE (UNIT=output_unit, FMT="(/,T3,A,(T30,A,T71,I10))") &
1916 0 : "Total number of", &
1917 0 : "- Atomic kinds: ", nkind, &
1918 0 : "- Atoms: ", natom
1919 : END IF
1920 :
1921 3723 : IF ((maxlppl > -1) .AND. (maxlppnl > -1)) THEN
1922 : WRITE (UNIT=output_unit, FMT="(/,T3,A,(T30,A,T75,I6))") &
1923 1888 : "Maximum angular momentum of the", &
1924 1888 : "- Orbital basis functions: ", maxlgto, &
1925 1888 : "- Local part of the GTH pseudopotential: ", maxlppl, &
1926 3776 : "- Non-local part of the GTH pseudopotential: ", maxlppnl
1927 1835 : ELSEIF (maxlppl > -1) THEN
1928 : WRITE (UNIT=output_unit, FMT="(/,T3,A,(T30,A,T75,I6))") &
1929 454 : "Maximum angular momentum of the", &
1930 454 : "- Orbital basis functions: ", maxlgto, &
1931 908 : "- Local part of the GTH pseudopotential: ", maxlppl
1932 : ELSE
1933 : WRITE (UNIT=output_unit, FMT="(/,T3,A,T75,I6)") &
1934 1381 : "Maximum angular momentum of the orbital basis functions: ", maxlgto
1935 : END IF
1936 :
1937 : ! LRI_AUX BASIS
1938 : CALL get_qs_kind_set(qs_kind_set, &
1939 : maxlgto=maxlgto, &
1940 : ncgf=ncgf, &
1941 : npgf=npgf, &
1942 : nset=nset, &
1943 : nsgf=nsgf, &
1944 : nshell=nshell, &
1945 3723 : basis_type="LRI_AUX")
1946 3723 : IF (nset + npgf + ncgf > 0) THEN
1947 : WRITE (UNIT=output_unit, FMT="(/,T3,A,/,T3,A,(T30,A,T71,I10))") &
1948 135 : "LRI_AUX Basis: ", &
1949 135 : "Total number of", &
1950 135 : "- Shell sets: ", nset, &
1951 135 : "- Shells: ", nshell, &
1952 135 : "- Primitive Cartesian functions: ", npgf, &
1953 135 : "- Cartesian basis functions: ", ncgf, &
1954 270 : "- Spherical basis functions: ", nsgf
1955 : WRITE (UNIT=output_unit, FMT="(T30,A,T75,I6)") &
1956 135 : " Maximum angular momentum ", maxlgto
1957 : END IF
1958 :
1959 : ! RI_HXC BASIS
1960 : CALL get_qs_kind_set(qs_kind_set, &
1961 : maxlgto=maxlgto, &
1962 : ncgf=ncgf, &
1963 : npgf=npgf, &
1964 : nset=nset, &
1965 : nsgf=nsgf, &
1966 : nshell=nshell, &
1967 3723 : basis_type="RI_HXC")
1968 3723 : IF (nset + npgf + ncgf > 0) THEN
1969 : WRITE (UNIT=output_unit, FMT="(/,T3,A,/,T3,A,(T30,A,T71,I10))") &
1970 111 : "RI_HXC Basis: ", &
1971 111 : "Total number of", &
1972 111 : "- Shell sets: ", nset, &
1973 111 : "- Shells: ", nshell, &
1974 111 : "- Primitive Cartesian functions: ", npgf, &
1975 111 : "- Cartesian basis functions: ", ncgf, &
1976 222 : "- Spherical basis functions: ", nsgf
1977 : WRITE (UNIT=output_unit, FMT="(T30,A,T75,I6)") &
1978 111 : " Maximum angular momentum ", maxlgto
1979 : END IF
1980 :
1981 : ! AUX_FIT BASIS
1982 : CALL get_qs_kind_set(qs_kind_set, &
1983 : maxlgto=maxlgto, &
1984 : ncgf=ncgf, &
1985 : npgf=npgf, &
1986 : nset=nset, &
1987 : nsgf=nsgf, &
1988 : nshell=nshell, &
1989 3723 : basis_type="AUX_FIT")
1990 3723 : IF (nset + npgf + ncgf > 0) THEN
1991 : WRITE (UNIT=output_unit, FMT="(/,T3,A,/,T3,A,(T30,A,T71,I10))") &
1992 336 : "AUX_FIT ADMM-Basis: ", &
1993 336 : "Total number of", &
1994 336 : "- Shell sets: ", nset, &
1995 336 : "- Shells: ", nshell, &
1996 336 : "- Primitive Cartesian functions: ", npgf, &
1997 336 : "- Cartesian basis functions: ", ncgf, &
1998 672 : "- Spherical basis functions: ", nsgf
1999 : WRITE (UNIT=output_unit, FMT="(T30,A,T75,I6)") &
2000 336 : " Maximum angular momentum ", maxlgto
2001 : END IF
2002 :
2003 : END IF
2004 : CALL cp_print_key_finished_output(output_unit, logger, force_env_section, &
2005 7398 : "PRINT%TOTAL_NUMBERS")
2006 :
2007 7398 : END SUBROUTINE write_total_numbers
2008 :
2009 : END MODULE qs_environment
|