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