(git:374b731)
Loading...
Searching...
No Matches
cp_control_utils.F
Go to the documentation of this file.
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!> \brief Utilities to set up the control types
10! **************************************************************************************************
12 USE bibliography, ONLY: &
17 USE cp_control_types, ONLY: &
22 USE cp_files, ONLY: close_file,&
28 USE cp_units, ONLY: cp_unit_from_cp2k,&
31 USE input_constants, ONLY: &
60 USE input_section_types, ONLY: &
63 USE kinds, ONLY: default_path_length,&
65 dp
66 USE mathconstants, ONLY: fourpi
71 USE util, ONLY: sort
75 USE xc_write_output, ONLY: xc_write
76#include "./base/base_uses.f90"
77
78 IMPLICIT NONE
79
80 PRIVATE
81
82 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'cp_control_utils'
83
84 PUBLIC :: read_dft_control, &
93CONTAINS
94
95! **************************************************************************************************
96!> \brief ...
97!> \param dft_control ...
98!> \param dft_section ...
99! **************************************************************************************************
100 SUBROUTINE read_dft_control(dft_control, dft_section)
101 TYPE(dft_control_type), POINTER :: dft_control
102 TYPE(section_vals_type), POINTER :: dft_section
103
104 CHARACTER(len=default_path_length) :: basis_set_file_name, potential_file_name
105 CHARACTER(LEN=default_string_length), &
106 DIMENSION(:), POINTER :: tmpstringlist
107 INTEGER :: admmtype, excitations, irep, isize, &
108 method_id, nrep, xc_deriv_method_id
109 LOGICAL :: do_ot, do_rtp, exopt1, exopt2, exopt3, &
110 explicit, is_present, l_param, not_se, &
111 was_present
112 REAL(kind=dp) :: density_cut, gradient_cut, tau_cut
113 REAL(kind=dp), DIMENSION(:), POINTER :: pol
114 TYPE(cp_logger_type), POINTER :: logger
115 TYPE(section_vals_type), POINTER :: maxwell_section, sccs_section, &
116 scf_section, tmp_section, &
117 xc_fun_section, xc_section
118
119 was_present = .false.
120
121 logger => cp_get_default_logger()
122
123 NULLIFY (tmp_section, xc_fun_section, xc_section)
124 ALLOCATE (dft_control)
125 CALL dft_control_create(dft_control)
126 ! determine wheather this is a semiempirical or DFTB run
127 ! --> (no XC section needs to be provided)
128 not_se = .true.
129 CALL section_vals_val_get(dft_section, "QS%METHOD", i_val=method_id)
130 SELECT CASE (method_id)
133 not_se = .false.
134 END SELECT
135 ! Check for XC section and XC_FUNCTIONAL section
136 xc_section => section_vals_get_subs_vals(dft_section, "XC")
137 CALL section_vals_get(xc_section, explicit=is_present)
138 IF (.NOT. is_present .AND. not_se) THEN
139 cpabort("XC section missing.")
140 END IF
141 IF (is_present) THEN
142 CALL section_vals_val_get(xc_section, "density_cutoff", r_val=density_cut)
143 CALL section_vals_val_get(xc_section, "gradient_cutoff", r_val=gradient_cut)
144 CALL section_vals_val_get(xc_section, "tau_cutoff", r_val=tau_cut)
145 ! Perform numerical stability checks and possibly correct the issues
146 IF (density_cut <= epsilon(0.0_dp)*100.0_dp) &
147 CALL cp_warn(__location__, &
148 "DENSITY_CUTOFF lower than 100*EPSILON, where EPSILON is the machine precision. "// &
149 "This may lead to numerical problems. Setting up shake_tol to 100*EPSILON! ")
150 density_cut = max(epsilon(0.0_dp)*100.0_dp, density_cut)
151 IF (gradient_cut <= epsilon(0.0_dp)*100.0_dp) &
152 CALL cp_warn(__location__, &
153 "GRADIENT_CUTOFF lower than 100*EPSILON, where EPSILON is the machine precision. "// &
154 "This may lead to numerical problems. Setting up shake_tol to 100*EPSILON! ")
155 gradient_cut = max(epsilon(0.0_dp)*100.0_dp, gradient_cut)
156 IF (tau_cut <= epsilon(0.0_dp)*100.0_dp) &
157 CALL cp_warn(__location__, &
158 "TAU_CUTOFF lower than 100*EPSILON, where EPSILON is the machine precision. "// &
159 "This may lead to numerical problems. Setting up shake_tol to 100*EPSILON! ")
160 tau_cut = max(epsilon(0.0_dp)*100.0_dp, tau_cut)
161 CALL section_vals_val_set(xc_section, "density_cutoff", r_val=density_cut)
162 CALL section_vals_val_set(xc_section, "gradient_cutoff", r_val=gradient_cut)
163 CALL section_vals_val_set(xc_section, "tau_cutoff", r_val=tau_cut)
164 END IF
165 xc_fun_section => section_vals_get_subs_vals(xc_section, "XC_FUNCTIONAL")
166 CALL section_vals_get(xc_fun_section, explicit=is_present)
167 IF (.NOT. is_present .AND. not_se) THEN
168 cpabort("XC_FUNCTIONAL section missing.")
169 END IF
170 scf_section => section_vals_get_subs_vals(dft_section, "SCF")
171 CALL section_vals_val_get(dft_section, "UKS", l_val=dft_control%uks)
172 CALL section_vals_val_get(dft_section, "ROKS", l_val=dft_control%roks)
173 IF (dft_control%uks .OR. dft_control%roks) THEN
174 dft_control%nspins = 2
175 ELSE
176 dft_control%nspins = 1
177 END IF
178
179 dft_control%lsd = (dft_control%nspins > 1)
180 dft_control%use_kinetic_energy_density = xc_uses_kinetic_energy_density(xc_fun_section, dft_control%lsd)
181
182 xc_deriv_method_id = section_get_ival(xc_section, "XC_GRID%XC_DERIV")
183 dft_control%drho_by_collocation = (xc_uses_norm_drho(xc_fun_section, dft_control%lsd) &
184 .AND. (xc_deriv_method_id == xc_deriv_collocate))
185 IF (dft_control%drho_by_collocation) THEN
186 cpabort("derivatives by collocation not implemented")
187 END IF
188
189 ! Automatic auxiliary basis set generation
190 CALL section_vals_val_get(dft_section, "AUTO_BASIS", n_rep_val=nrep)
191 DO irep = 1, nrep
192 CALL section_vals_val_get(dft_section, "AUTO_BASIS", i_rep_val=irep, c_vals=tmpstringlist)
193 IF (SIZE(tmpstringlist) == 2) THEN
194 CALL uppercase(tmpstringlist(2))
195 SELECT CASE (tmpstringlist(2))
196 CASE ("X")
197 isize = -1
198 CASE ("SMALL")
199 isize = 0
200 CASE ("MEDIUM")
201 isize = 1
202 CASE ("LARGE")
203 isize = 2
204 CASE ("HUGE")
205 isize = 3
206 CASE DEFAULT
207 cpwarn("Unknown basis size in AUTO_BASIS keyword:"//trim(tmpstringlist(1)))
208 END SELECT
209 !
210 SELECT CASE (tmpstringlist(1))
211 CASE ("X")
212 CASE ("RI_AUX")
213 dft_control%auto_basis_ri_aux = isize
214 CASE ("AUX_FIT")
215 dft_control%auto_basis_aux_fit = isize
216 CASE ("LRI_AUX")
217 dft_control%auto_basis_lri_aux = isize
218 CASE ("P_LRI_AUX")
219 dft_control%auto_basis_p_lri_aux = isize
220 CASE ("RI_HXC")
221 dft_control%auto_basis_ri_hxc = isize
222 CASE ("RI_XAS")
223 dft_control%auto_basis_ri_xas = isize
224 CASE ("RI_HFX")
225 dft_control%auto_basis_ri_hfx = isize
226 CASE DEFAULT
227 cpwarn("Unknown basis type in AUTO_BASIS keyword:"//trim(tmpstringlist(1)))
228 END SELECT
229 ELSE
230 CALL cp_abort(__location__, &
231 "AUTO_BASIS keyword in &DFT section has a wrong number of arguments.")
232 END IF
233 END DO
234
235 !! check if we do wavefunction fitting
236 tmp_section => section_vals_get_subs_vals(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD")
237 CALL section_vals_get(tmp_section, explicit=is_present)
238 dft_control%do_admm = is_present
239 dft_control%do_admm_mo = .false.
240 dft_control%do_admm_dm = .false.
241 IF (is_present) THEN
242 do_ot = .false.
243 CALL section_vals_val_get(scf_section, "OT%_SECTION_PARAMETERS_", l_val=do_ot)
244 CALL admm_control_create(dft_control%admm_control)
245
246 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%ADMM_TYPE", i_val=admmtype)
247 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%ADMM_PURIFICATION_METHOD", explicit=exopt1)
248 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%METHOD", explicit=exopt2)
249 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%EXCH_SCALING_MODEL", explicit=exopt3)
250 dft_control%admm_control%admm_type = admmtype
251 SELECT CASE (admmtype)
252 CASE (no_admm_type)
253 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%ADMM_PURIFICATION_METHOD", i_val=method_id)
254 dft_control%admm_control%purification_method = method_id
255 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%METHOD", i_val=method_id)
256 dft_control%admm_control%method = method_id
257 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%EXCH_SCALING_MODEL", i_val=method_id)
258 dft_control%admm_control%scaling_model = method_id
259 CASE (admm1_type)
260 ! METHOD BASIS_PROJECTION
261 ! ADMM_PURIFICATION_METHOD choose
262 ! EXCH_SCALING_MODEL NONE
263 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%ADMM_PURIFICATION_METHOD", i_val=method_id)
264 dft_control%admm_control%purification_method = method_id
265 dft_control%admm_control%method = do_admm_basis_projection
266 dft_control%admm_control%scaling_model = do_admm_exch_scaling_none
267 CASE (admm2_type)
268 ! METHOD BASIS_PROJECTION
269 ! ADMM_PURIFICATION_METHOD NONE
270 ! EXCH_SCALING_MODEL NONE
271 dft_control%admm_control%purification_method = do_admm_purify_none
272 dft_control%admm_control%method = do_admm_basis_projection
273 dft_control%admm_control%scaling_model = do_admm_exch_scaling_none
274 CASE (admms_type)
275 ! ADMM_PURIFICATION_METHOD NONE
276 ! METHOD CHARGE_CONSTRAINED_PROJECTION
277 ! EXCH_SCALING_MODEL MERLOT
278 dft_control%admm_control%purification_method = do_admm_purify_none
279 dft_control%admm_control%method = do_admm_charge_constrained_projection
280 dft_control%admm_control%scaling_model = do_admm_exch_scaling_merlot
281 CASE (admmp_type)
282 ! ADMM_PURIFICATION_METHOD NONE
283 ! METHOD BASIS_PROJECTION
284 ! EXCH_SCALING_MODEL MERLOT
285 dft_control%admm_control%purification_method = do_admm_purify_none
286 dft_control%admm_control%method = do_admm_basis_projection
287 dft_control%admm_control%scaling_model = do_admm_exch_scaling_merlot
288 CASE (admmq_type)
289 ! ADMM_PURIFICATION_METHOD NONE
290 ! METHOD CHARGE_CONSTRAINED_PROJECTION
291 ! EXCH_SCALING_MODEL NONE
292 dft_control%admm_control%purification_method = do_admm_purify_none
293 dft_control%admm_control%method = do_admm_charge_constrained_projection
294 dft_control%admm_control%scaling_model = do_admm_exch_scaling_none
295 CASE DEFAULT
296 CALL cp_abort(__location__, &
297 "ADMM_TYPE keyword in &AUXILIARY_DENSITY_MATRIX_METHOD section has a wrong value.")
298 END SELECT
299
300 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%EPS_FILTER", &
301 r_val=dft_control%admm_control%eps_filter)
302
303 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%EXCH_CORRECTION_FUNC", i_val=method_id)
304 dft_control%admm_control%aux_exch_func = method_id
305
306 ! parameters for X functional
307 dft_control%admm_control%aux_exch_func_param = .false.
308 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%OPTX_A1", explicit=explicit, &
309 r_val=dft_control%admm_control%aux_x_param(1))
310 IF (explicit) dft_control%admm_control%aux_exch_func_param = .true.
311 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%OPTX_A2", explicit=explicit, &
312 r_val=dft_control%admm_control%aux_x_param(2))
313 IF (explicit) dft_control%admm_control%aux_exch_func_param = .true.
314 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%OPTX_GAMMA", explicit=explicit, &
315 r_val=dft_control%admm_control%aux_x_param(3))
316 IF (explicit) dft_control%admm_control%aux_exch_func_param = .true.
317
318 CALL read_admm_block_list(dft_control%admm_control, dft_section)
319
320 ! check for double assignments
321 SELECT CASE (admmtype)
322 CASE (admm2_type)
323 IF (exopt2) CALL cp_warn(__location__, &
324 "Value of ADMM_PURIFICATION_METHOD keyword will be overwritten with ADMM_TYPE selections.")
325 IF (exopt3) CALL cp_warn(__location__, &
326 "Value of EXCH_SCALING_MODEL keyword will be overwritten with ADMM_TYPE selections.")
328 IF (exopt1) CALL cp_warn(__location__, &
329 "Value of METHOD keyword will be overwritten with ADMM_TYPE selections.")
330 IF (exopt2) CALL cp_warn(__location__, &
331 "Value of METHOD keyword will be overwritten with ADMM_TYPE selections.")
332 IF (exopt3) CALL cp_warn(__location__, &
333 "Value of EXCH_SCALING_MODEL keyword will be overwritten with ADMM_TYPE selections.")
334 END SELECT
335
336 ! In the case of charge-constrained projection (e.g. according to Merlot),
337 ! there is no purification needed and hence, do_admm_purify_none has to be set.
338
339 IF ((dft_control%admm_control%method == do_admm_blocking_purify_full .OR. &
340 dft_control%admm_control%method == do_admm_blocked_projection) &
341 .AND. dft_control%admm_control%scaling_model == do_admm_exch_scaling_merlot) THEN
342 cpabort("ADMM: Blocking and Merlot scaling are mutually exclusive.")
343 END IF
344
345 IF (dft_control%admm_control%method == do_admm_charge_constrained_projection .AND. &
346 dft_control%admm_control%purification_method /= do_admm_purify_none) THEN
347 CALL cp_abort(__location__, &
348 "ADMM: In the case of METHOD=CHARGE_CONSTRAINED_PROJECTION, "// &
349 "ADMM_PURIFICATION_METHOD=NONE has to be set.")
350 END IF
351
352 IF (dft_control%admm_control%purification_method == do_admm_purify_mo_diag .OR. &
353 dft_control%admm_control%purification_method == do_admm_purify_mo_no_diag) THEN
354 IF (dft_control%admm_control%method /= do_admm_basis_projection) &
355 cpabort("ADMM: Chosen purification requires BASIS_PROJECTION")
356
357 IF (.NOT. do_ot) cpabort("ADMM: MO-based purification requires OT.")
358 END IF
359
360 IF (dft_control%admm_control%purification_method == do_admm_purify_none_dm .OR. &
361 dft_control%admm_control%purification_method == do_admm_purify_mcweeny) THEN
362 dft_control%do_admm_dm = .true.
363 ELSE
364 dft_control%do_admm_mo = .true.
365 END IF
366 END IF
367
368 ! Set restricted to true, if both OT and ROKS are requested
369 !MK in principle dft_control%restricted could be dropped completely like the
370 !MK input key by using only dft_control%roks now
371 CALL section_vals_val_get(scf_section, "OT%_SECTION_PARAMETERS_", l_val=l_param)
372 dft_control%restricted = (dft_control%roks .AND. l_param)
373
374 CALL section_vals_val_get(dft_section, "CHARGE", i_val=dft_control%charge)
375 CALL section_vals_val_get(dft_section, "MULTIPLICITY", i_val=dft_control%multiplicity)
376 CALL section_vals_val_get(dft_section, "RELAX_MULTIPLICITY", r_val=dft_control%relax_multiplicity)
377 IF (dft_control%relax_multiplicity > 0.0_dp) THEN
378 IF (.NOT. dft_control%uks) &
379 CALL cp_abort(__location__, "The option RELAX_MULTIPLICITY is only valid for "// &
380 "unrestricted Kohn-Sham (UKS) calculations")
381 END IF
382
383 ! check for the presence of the low spin roks section
384 tmp_section => section_vals_get_subs_vals(dft_section, "LOW_SPIN_ROKS")
385 CALL section_vals_get(tmp_section, explicit=dft_control%low_spin_roks)
386
387 dft_control%sic_method_id = sic_none
388 dft_control%sic_scaling_a = 1.0_dp
389 dft_control%sic_scaling_b = 1.0_dp
390
391 ! DFT+U
392 dft_control%dft_plus_u = .false.
393 CALL section_vals_val_get(dft_section, "PLUS_U_METHOD", i_val=method_id)
394 dft_control%plus_u_method_id = method_id
395
396 ! Smearing in use
397 dft_control%smear = .false.
398
399 ! Surface dipole correction
400 dft_control%correct_surf_dip = .false.
401 CALL section_vals_val_get(dft_section, "SURFACE_DIPOLE_CORRECTION", l_val=dft_control%correct_surf_dip)
402 CALL section_vals_val_get(dft_section, "SURF_DIP_DIR", i_val=dft_control%dir_surf_dip)
403 dft_control%pos_dir_surf_dip = -1.0_dp
404 CALL section_vals_val_get(dft_section, "SURF_DIP_POS", r_val=dft_control%pos_dir_surf_dip)
405 ! another logical variable, surf_dip_correct_switch, is introduced for
406 ! implementation of "SURF_DIP_SWITCH" [SGh]
407 dft_control%switch_surf_dip = .false.
408 dft_control%surf_dip_correct_switch = dft_control%correct_surf_dip
409 CALL section_vals_val_get(dft_section, "SURF_DIP_SWITCH", l_val=dft_control%switch_surf_dip)
410 dft_control%correct_el_density_dip = .false.
411 CALL section_vals_val_get(dft_section, "CORE_CORR_DIP", l_val=dft_control%correct_el_density_dip)
412 IF (dft_control%correct_el_density_dip) THEN
413 IF (dft_control%correct_surf_dip) THEN
414 ! Do nothing, move on
415 ELSE
416 dft_control%correct_el_density_dip = .false.
417 cpwarn("CORE_CORR_DIP keyword is activated only if SURFACE_DIPOLE_CORRECTION is TRUE")
418 END IF
419 END IF
420
421 CALL section_vals_val_get(dft_section, "BASIS_SET_FILE_NAME", &
422 c_val=basis_set_file_name)
423 CALL section_vals_val_get(dft_section, "POTENTIAL_FILE_NAME", &
424 c_val=potential_file_name)
425
426 ! Read the input section
427 tmp_section => section_vals_get_subs_vals(dft_section, "sic")
428 CALL section_vals_val_get(tmp_section, "SIC_METHOD", &
429 i_val=dft_control%sic_method_id)
430 CALL section_vals_val_get(tmp_section, "ORBITAL_SET", &
431 i_val=dft_control%sic_list_id)
432 CALL section_vals_val_get(tmp_section, "SIC_SCALING_A", &
433 r_val=dft_control%sic_scaling_a)
434 CALL section_vals_val_get(tmp_section, "SIC_SCALING_B", &
435 r_val=dft_control%sic_scaling_b)
436
437 ! Determine if this is a TDDFPT run
438 CALL section_vals_val_get(dft_section, "EXCITATIONS", i_val=excitations)
439 dft_control%do_tddfpt_calculation = (excitations == tddfpt_excitations)
440 IF (dft_control%do_tddfpt_calculation) THEN
441 CALL tddfpt_control_create(dft_control%tddfpt_control)
442 END IF
443
444 do_rtp = .false.
445 tmp_section => section_vals_get_subs_vals(dft_section, "REAL_TIME_PROPAGATION")
446 CALL section_vals_get(tmp_section, explicit=is_present)
447 IF (is_present) THEN
448 CALL read_rtp_section(dft_control, tmp_section)
449 do_rtp = .true.
450 END IF
451
452 ! Read the input section
453 tmp_section => section_vals_get_subs_vals(dft_section, "XAS")
454 CALL section_vals_get(tmp_section, explicit=dft_control%do_xas_calculation)
455 IF (dft_control%do_xas_calculation) THEN
456 ! Override with section parameter
457 CALL section_vals_val_get(tmp_section, "_SECTION_PARAMETERS_", &
458 l_val=dft_control%do_xas_calculation)
459 END IF
460
461 tmp_section => section_vals_get_subs_vals(dft_section, "XAS_TDP")
462 CALL section_vals_get(tmp_section, explicit=dft_control%do_xas_tdp_calculation)
463 IF (dft_control%do_xas_tdp_calculation) THEN
464 ! Override with section parameter
465 CALL section_vals_val_get(tmp_section, "_SECTION_PARAMETERS_", &
466 l_val=dft_control%do_xas_tdp_calculation)
467 END IF
468
469 ! Read the finite field input section
470 dft_control%apply_efield = .false.
471 dft_control%apply_efield_field = .false. !this is for RTP
472 dft_control%apply_vector_potential = .false. !this is for RTP
473 tmp_section => section_vals_get_subs_vals(dft_section, "EFIELD")
474 CALL section_vals_get(tmp_section, n_repetition=nrep, explicit=is_present)
475 IF (is_present) THEN
476 ALLOCATE (dft_control%efield_fields(nrep))
477 CALL read_efield_sections(dft_control, tmp_section)
478 IF (do_rtp) THEN
479 IF (.NOT. dft_control%rtp_control%velocity_gauge) THEN
480 dft_control%apply_efield_field = .true.
481 ELSE
482 dft_control%apply_vector_potential = .true.
483 ! Use this input value of vector potential to (re)start RTP
484 dft_control%rtp_control%vec_pot = dft_control%efield_fields(1)%efield%vec_pot_initial
485 END IF
486 ELSE
487 dft_control%apply_efield = .true.
488 cpassert(nrep == 1)
489 END IF
490 END IF
491
492 ! Read the finite field input section for periodic fields
493 tmp_section => section_vals_get_subs_vals(dft_section, "PERIODIC_EFIELD")
494 CALL section_vals_get(tmp_section, explicit=dft_control%apply_period_efield)
495 IF (dft_control%apply_period_efield) THEN
496 ALLOCATE (dft_control%period_efield)
497 CALL section_vals_val_get(tmp_section, "POLARISATION", r_vals=pol)
498 dft_control%period_efield%polarisation(1:3) = pol(1:3)
499 CALL section_vals_val_get(tmp_section, "D_FILTER", r_vals=pol)
500 dft_control%period_efield%d_filter(1:3) = pol(1:3)
501 CALL section_vals_val_get(tmp_section, "INTENSITY", &
502 r_val=dft_control%period_efield%strength)
503 dft_control%period_efield%displacement_field = .false.
504 CALL section_vals_val_get(tmp_section, "DISPLACEMENT_FIELD", &
505 l_val=dft_control%period_efield%displacement_field)
506 ! periodic fields don't work with RTP
507 cpassert(.NOT. do_rtp)
508 IF (dft_control%period_efield%displacement_field) THEN
509 CALL cite_reference(stengel2009)
510 ELSE
511 CALL cite_reference(souza2002)
512 CALL cite_reference(umari2002)
513 END IF
514 END IF
515
516 ! Read the external potential input section
517 tmp_section => section_vals_get_subs_vals(dft_section, "EXTERNAL_POTENTIAL")
518 CALL section_vals_get(tmp_section, explicit=dft_control%apply_external_potential)
519 IF (dft_control%apply_external_potential) THEN
520 CALL expot_control_create(dft_control%expot_control)
521 CALL section_vals_val_get(tmp_section, "READ_FROM_CUBE", &
522 l_val=dft_control%expot_control%read_from_cube)
523 CALL section_vals_val_get(tmp_section, "STATIC", &
524 l_val=dft_control%expot_control%static)
525 CALL section_vals_val_get(tmp_section, "SCALING_FACTOR", &
526 r_val=dft_control%expot_control%scaling_factor)
527 ! External potential using Maxwell equation
528 maxwell_section => section_vals_get_subs_vals(tmp_section, "MAXWELL")
529 CALL section_vals_get(maxwell_section, explicit=is_present)
530 IF (is_present) THEN
531 dft_control%expot_control%maxwell_solver = .true.
532 CALL maxwell_control_create(dft_control%maxwell_control)
533 ! read the input values from Maxwell section
534 CALL section_vals_val_get(maxwell_section, "TEST_REAL", &
535 r_val=dft_control%maxwell_control%real_test)
536 CALL section_vals_val_get(maxwell_section, "TEST_INTEGER", &
537 i_val=dft_control%maxwell_control%int_test)
538 CALL section_vals_val_get(maxwell_section, "TEST_LOGICAL", &
539 l_val=dft_control%maxwell_control%log_test)
540 ELSE
541 dft_control%expot_control%maxwell_solver = .false.
542 END IF
543 END IF
544
545 ! Read the SCCS input section if present
546 sccs_section => section_vals_get_subs_vals(dft_section, "SCCS")
547 CALL section_vals_get(sccs_section, explicit=is_present)
548 IF (is_present) THEN
549 ! Check section parameter if SCCS is activated
550 CALL section_vals_val_get(sccs_section, "_SECTION_PARAMETERS_", &
551 l_val=dft_control%do_sccs)
552 IF (dft_control%do_sccs) THEN
553 ALLOCATE (dft_control%sccs_control)
554 CALL section_vals_val_get(sccs_section, "RELATIVE_PERMITTIVITY", &
555 r_val=dft_control%sccs_control%epsilon_solvent)
556 CALL section_vals_val_get(sccs_section, "ALPHA", &
557 r_val=dft_control%sccs_control%alpha_solvent)
558 CALL section_vals_val_get(sccs_section, "BETA", &
559 r_val=dft_control%sccs_control%beta_solvent)
560 CALL section_vals_val_get(sccs_section, "DELTA_RHO", &
561 r_val=dft_control%sccs_control%delta_rho)
562 CALL section_vals_val_get(sccs_section, "DERIVATIVE_METHOD", &
563 i_val=dft_control%sccs_control%derivative_method)
564 CALL section_vals_val_get(sccs_section, "METHOD", &
565 i_val=dft_control%sccs_control%method_id)
566 CALL section_vals_val_get(sccs_section, "EPS_SCCS", &
567 r_val=dft_control%sccs_control%eps_sccs)
568 CALL section_vals_val_get(sccs_section, "EPS_SCF", &
569 r_val=dft_control%sccs_control%eps_scf)
570 CALL section_vals_val_get(sccs_section, "GAMMA", &
571 r_val=dft_control%sccs_control%gamma_solvent)
572 CALL section_vals_val_get(sccs_section, "MAX_ITER", &
573 i_val=dft_control%sccs_control%max_iter)
574 CALL section_vals_val_get(sccs_section, "MIXING", &
575 r_val=dft_control%sccs_control%mixing)
576 SELECT CASE (dft_control%sccs_control%method_id)
577 CASE (sccs_andreussi)
578 tmp_section => section_vals_get_subs_vals(sccs_section, "ANDREUSSI")
579 CALL section_vals_val_get(tmp_section, "RHO_MAX", &
580 r_val=dft_control%sccs_control%rho_max)
581 CALL section_vals_val_get(tmp_section, "RHO_MIN", &
582 r_val=dft_control%sccs_control%rho_min)
583 IF (dft_control%sccs_control%rho_max < dft_control%sccs_control%rho_min) THEN
584 CALL cp_abort(__location__, &
585 "The SCCS parameter RHO_MAX is smaller than RHO_MIN. "// &
586 "Please, check your input!")
587 END IF
588 CALL cite_reference(andreussi2012)
590 tmp_section => section_vals_get_subs_vals(sccs_section, "FATTEBERT-GYGI")
591 CALL section_vals_val_get(tmp_section, "BETA", &
592 r_val=dft_control%sccs_control%beta)
593 IF (dft_control%sccs_control%beta < 0.5_dp) THEN
594 CALL cp_abort(__location__, &
595 "A value smaller than 0.5 for the SCCS parameter beta "// &
596 "causes numerical problems. Please, check your input!")
597 END IF
598 CALL section_vals_val_get(tmp_section, "RHO_ZERO", &
599 r_val=dft_control%sccs_control%rho_zero)
600 CALL cite_reference(fattebert2002)
601 CASE DEFAULT
602 cpabort("Invalid SCCS model specified. Please, check your input!")
603 END SELECT
604 CALL cite_reference(yin2017)
605 END IF
606 END IF
607
608 ! ZMP added input sections
609 ! Read the external density input section
610 tmp_section => section_vals_get_subs_vals(dft_section, "EXTERNAL_DENSITY")
611 CALL section_vals_get(tmp_section, explicit=dft_control%apply_external_density)
612
613 ! Read the external vxc input section
614 tmp_section => section_vals_get_subs_vals(dft_section, "EXTERNAL_VXC")
615 CALL section_vals_get(tmp_section, explicit=dft_control%apply_external_vxc)
616
617 END SUBROUTINE read_dft_control
618
619! **************************************************************************************************
620!> \brief ...
621!> \param qs_control ...
622!> \param dft_section ...
623! **************************************************************************************************
624 SUBROUTINE read_mgrid_section(qs_control, dft_section)
625
626 TYPE(qs_control_type), INTENT(INOUT) :: qs_control
627 TYPE(section_vals_type), POINTER :: dft_section
628
629 CHARACTER(len=*), PARAMETER :: routinen = 'read_mgrid_section'
630
631 INTEGER :: handle, igrid_level, ngrid_level
632 LOGICAL :: explicit, multigrid_set
633 REAL(dp) :: cutoff
634 REAL(dp), DIMENSION(:), POINTER :: cutofflist
635 TYPE(section_vals_type), POINTER :: mgrid_section
636
637 CALL timeset(routinen, handle)
638
639 NULLIFY (mgrid_section, cutofflist)
640 mgrid_section => section_vals_get_subs_vals(dft_section, "MGRID")
641
642 CALL section_vals_val_get(mgrid_section, "NGRIDS", i_val=ngrid_level)
643 CALL section_vals_val_get(mgrid_section, "MULTIGRID_SET", l_val=multigrid_set)
644 CALL section_vals_val_get(mgrid_section, "CUTOFF", r_val=cutoff)
645 CALL section_vals_val_get(mgrid_section, "PROGRESSION_FACTOR", r_val=qs_control%progression_factor)
646 CALL section_vals_val_get(mgrid_section, "COMMENSURATE", l_val=qs_control%commensurate_mgrids)
647 CALL section_vals_val_get(mgrid_section, "REALSPACE", l_val=qs_control%realspace_mgrids)
648 CALL section_vals_val_get(mgrid_section, "REL_CUTOFF", r_val=qs_control%relative_cutoff)
649 CALL section_vals_val_get(mgrid_section, "SKIP_LOAD_BALANCE_DISTRIBUTED", &
650 l_val=qs_control%skip_load_balance_distributed)
651
652 ! For SE and DFTB possibly override with new defaults
653 IF (qs_control%semi_empirical .OR. qs_control%dftb .OR. qs_control%xtb) THEN
654 ngrid_level = 1
655 multigrid_set = .false.
656 ! Override default cutoff value unless user specified an explicit argument..
657 CALL section_vals_val_get(mgrid_section, "CUTOFF", explicit=explicit, r_val=cutoff)
658 IF (.NOT. explicit) cutoff = 1.0_dp
659 END IF
660
661 ALLOCATE (qs_control%e_cutoff(ngrid_level))
662 qs_control%cutoff = cutoff
663
664 IF (multigrid_set) THEN
665 ! Read the values from input
666 IF (qs_control%commensurate_mgrids) THEN
667 cpabort("Do not specify cutoffs for the commensurate grids (NYI)")
668 END IF
669
670 CALL section_vals_val_get(mgrid_section, "MULTIGRID_CUTOFF", r_vals=cutofflist)
671 IF (ASSOCIATED(cutofflist)) THEN
672 IF (SIZE(cutofflist, 1) /= ngrid_level) THEN
673 cpabort("Number of multi-grids requested and number of cutoff values do not match")
674 END IF
675 DO igrid_level = 1, ngrid_level
676 qs_control%e_cutoff(igrid_level) = cutofflist(igrid_level)
677 END DO
678 END IF
679 ! set cutoff to smallest value in multgrid available with >= cutoff
680 DO igrid_level = ngrid_level, 1, -1
681 IF (qs_control%cutoff <= qs_control%e_cutoff(igrid_level)) THEN
682 qs_control%cutoff = qs_control%e_cutoff(igrid_level)
683 EXIT
684 END IF
685 ! set largest grid value to cutoff
686 IF (igrid_level == 1) THEN
687 qs_control%cutoff = qs_control%e_cutoff(1)
688 END IF
689 END DO
690 ELSE
691 IF (qs_control%commensurate_mgrids) qs_control%progression_factor = 4.0_dp
692 qs_control%e_cutoff(1) = qs_control%cutoff
693 DO igrid_level = 2, ngrid_level
694 qs_control%e_cutoff(igrid_level) = qs_control%e_cutoff(igrid_level - 1)/ &
695 qs_control%progression_factor
696 END DO
697 END IF
698 ! check that multigrids are ordered
699 DO igrid_level = 2, ngrid_level
700 IF (qs_control%e_cutoff(igrid_level) > qs_control%e_cutoff(igrid_level - 1)) THEN
701 cpabort("The cutoff values for the multi-grids are not ordered from large to small")
702 ELSE IF (qs_control%e_cutoff(igrid_level) == qs_control%e_cutoff(igrid_level - 1)) THEN
703 cpabort("The same cutoff value was specified for two multi-grids")
704 END IF
705 END DO
706 CALL timestop(handle)
707 END SUBROUTINE read_mgrid_section
708
709! **************************************************************************************************
710!> \brief ...
711!> \param qs_control ...
712!> \param qs_section ...
713! **************************************************************************************************
714 SUBROUTINE read_qs_section(qs_control, qs_section)
715
716 TYPE(qs_control_type), INTENT(INOUT) :: qs_control
717 TYPE(section_vals_type), POINTER :: qs_section
718
719 CHARACTER(len=*), PARAMETER :: routinen = 'read_qs_section'
720
721 CHARACTER(LEN=default_string_length), &
722 DIMENSION(:), POINTER :: clist
723 INTEGER :: handle, itmp, j, jj, k, n_rep, n_var, &
724 ngauss, ngp, nrep
725 INTEGER, DIMENSION(:), POINTER :: tmplist
726 LOGICAL :: explicit, was_present
727 REAL(dp) :: tmp, tmpsqrt, value
728 REAL(dp), POINTER :: scal(:)
729 TYPE(section_vals_type), POINTER :: cdft_control_section, ddapc_restraint_section, &
730 dftb_parameter, dftb_section, genpot_section, lri_optbas_section, mull_section, &
731 nonbonded_section, s2_restraint_section, se_section, xtb_parameter, xtb_section
732
733 CALL timeset(routinen, handle)
734
735 was_present = .false.
736 NULLIFY (mull_section, ddapc_restraint_section, s2_restraint_section, &
737 se_section, dftb_section, xtb_section, dftb_parameter, xtb_parameter, lri_optbas_section, &
738 cdft_control_section, genpot_section)
739
740 mull_section => section_vals_get_subs_vals(qs_section, "MULLIKEN_RESTRAINT")
741 ddapc_restraint_section => section_vals_get_subs_vals(qs_section, "DDAPC_RESTRAINT")
742 s2_restraint_section => section_vals_get_subs_vals(qs_section, "S2_RESTRAINT")
743 se_section => section_vals_get_subs_vals(qs_section, "SE")
744 dftb_section => section_vals_get_subs_vals(qs_section, "DFTB")
745 xtb_section => section_vals_get_subs_vals(qs_section, "xTB")
746 dftb_parameter => section_vals_get_subs_vals(dftb_section, "PARAMETER")
747 xtb_parameter => section_vals_get_subs_vals(xtb_section, "PARAMETER")
748 lri_optbas_section => section_vals_get_subs_vals(qs_section, "OPTIMIZE_LRI_BASIS")
749 cdft_control_section => section_vals_get_subs_vals(qs_section, "CDFT")
750 nonbonded_section => section_vals_get_subs_vals(xtb_section, "NONBONDED")
751 genpot_section => section_vals_get_subs_vals(nonbonded_section, "GENPOT")
752
753 ! Setup all defaults values and overwrite input parameters
754 ! EPS_DEFAULT should set the target accuracy in the total energy (~per electron) or a closely related value
755 CALL section_vals_val_get(qs_section, "EPS_DEFAULT", r_val=value)
756 tmpsqrt = sqrt(value) ! a trick to work around a NAG 5.1 optimizer bug
757
758 ! random choice ?
759 qs_control%eps_core_charge = value/100.0_dp
760 ! correct if all Gaussians would have the same radius (overlap will be smaller than eps_pgf_orb**2).
761 ! Can be significantly in error if not... requires fully new screening/pairlist procedures
762 qs_control%eps_pgf_orb = tmpsqrt
763 qs_control%eps_kg_orb = qs_control%eps_pgf_orb
764 ! consistent since also a kind of overlap
765 qs_control%eps_ppnl = qs_control%eps_pgf_orb/100.0_dp
766 ! accuracy is basically set by the overlap, this sets an empirical shift
767 qs_control%eps_ppl = 1.0e-2_dp
768 !
769 qs_control%gapw_control%eps_cpc = value
770 ! expexted error in the density
771 qs_control%eps_rho_gspace = value
772 qs_control%eps_rho_rspace = value
773 ! error in the gradient, can be the sqrt of the error in the energy, ignored if map_consistent
774 qs_control%eps_gvg_rspace = tmpsqrt
775 !
776 CALL section_vals_val_get(qs_section, "EPS_CORE_CHARGE", n_rep_val=n_rep)
777 IF (n_rep /= 0) THEN
778 CALL section_vals_val_get(qs_section, "EPS_CORE_CHARGE", r_val=qs_control%eps_core_charge)
779 END IF
780 CALL section_vals_val_get(qs_section, "EPS_GVG_RSPACE", n_rep_val=n_rep)
781 IF (n_rep /= 0) THEN
782 CALL section_vals_val_get(qs_section, "EPS_GVG_RSPACE", r_val=qs_control%eps_gvg_rspace)
783 END IF
784 CALL section_vals_val_get(qs_section, "EPS_PGF_ORB", n_rep_val=n_rep)
785 IF (n_rep /= 0) THEN
786 CALL section_vals_val_get(qs_section, "EPS_PGF_ORB", r_val=qs_control%eps_pgf_orb)
787 END IF
788 CALL section_vals_val_get(qs_section, "EPS_KG_ORB", n_rep_val=n_rep)
789 IF (n_rep /= 0) THEN
790 CALL section_vals_val_get(qs_section, "EPS_KG_ORB", r_val=tmp)
791 qs_control%eps_kg_orb = sqrt(tmp)
792 END IF
793 CALL section_vals_val_get(qs_section, "EPS_PPL", n_rep_val=n_rep)
794 IF (n_rep /= 0) THEN
795 CALL section_vals_val_get(qs_section, "EPS_PPL", r_val=qs_control%eps_ppl)
796 END IF
797 CALL section_vals_val_get(qs_section, "EPS_PPNL", n_rep_val=n_rep)
798 IF (n_rep /= 0) THEN
799 CALL section_vals_val_get(qs_section, "EPS_PPNL", r_val=qs_control%eps_ppnl)
800 END IF
801 CALL section_vals_val_get(qs_section, "EPS_RHO", n_rep_val=n_rep)
802 IF (n_rep /= 0) THEN
803 CALL section_vals_val_get(qs_section, "EPS_RHO", r_val=qs_control%eps_rho_gspace)
804 qs_control%eps_rho_rspace = qs_control%eps_rho_gspace
805 END IF
806 CALL section_vals_val_get(qs_section, "EPS_RHO_RSPACE", n_rep_val=n_rep)
807 IF (n_rep /= 0) THEN
808 CALL section_vals_val_get(qs_section, "EPS_RHO_RSPACE", r_val=qs_control%eps_rho_rspace)
809 END IF
810 CALL section_vals_val_get(qs_section, "EPS_RHO_GSPACE", n_rep_val=n_rep)
811 IF (n_rep /= 0) THEN
812 CALL section_vals_val_get(qs_section, "EPS_RHO_GSPACE", r_val=qs_control%eps_rho_gspace)
813 END IF
814 CALL section_vals_val_get(qs_section, "EPS_FILTER_MATRIX", n_rep_val=n_rep)
815 IF (n_rep /= 0) THEN
816 CALL section_vals_val_get(qs_section, "EPS_FILTER_MATRIX", r_val=qs_control%eps_filter_matrix)
817 END IF
818 CALL section_vals_val_get(qs_section, "EPS_CPC", n_rep_val=n_rep)
819 IF (n_rep /= 0) THEN
820 CALL section_vals_val_get(qs_section, "EPS_CPC", r_val=qs_control%gapw_control%eps_cpc)
821 END IF
822
823 CALL section_vals_val_get(qs_section, "EPSFIT", r_val=qs_control%gapw_control%eps_fit)
824 CALL section_vals_val_get(qs_section, "EPSISO", r_val=qs_control%gapw_control%eps_iso)
825 CALL section_vals_val_get(qs_section, "EPSSVD", r_val=qs_control%gapw_control%eps_svd)
826 CALL section_vals_val_get(qs_section, "EPSRHO0", r_val=qs_control%gapw_control%eps_Vrho0)
827 CALL section_vals_val_get(qs_section, "ALPHA0_HARD", r_val=qs_control%gapw_control%alpha0_hard)
828 qs_control%gapw_control%lrho1_eq_lrho0 = .false.
829 qs_control%gapw_control%alpha0_hard_from_input = .false.
830 IF (qs_control%gapw_control%alpha0_hard /= 0.0_dp) qs_control%gapw_control%alpha0_hard_from_input = .true.
831 CALL section_vals_val_get(qs_section, "FORCE_PAW", l_val=qs_control%gapw_control%force_paw)
832 CALL section_vals_val_get(qs_section, "MAX_RAD_LOCAL", r_val=qs_control%gapw_control%max_rad_local)
833
834 CALL section_vals_val_get(qs_section, "MIN_PAIR_LIST_RADIUS", r_val=qs_control%pairlist_radius)
835
836 CALL section_vals_val_get(qs_section, "LS_SCF", l_val=qs_control%do_ls_scf)
837 CALL section_vals_val_get(qs_section, "ALMO_SCF", l_val=qs_control%do_almo_scf)
838 CALL section_vals_val_get(qs_section, "KG_METHOD", l_val=qs_control%do_kg)
839
840 ! Logicals
841 CALL section_vals_val_get(qs_section, "REF_EMBED_SUBSYS", l_val=qs_control%ref_embed_subsys)
842 CALL section_vals_val_get(qs_section, "CLUSTER_EMBED_SUBSYS", l_val=qs_control%cluster_embed_subsys)
843 CALL section_vals_val_get(qs_section, "HIGH_LEVEL_EMBED_SUBSYS", l_val=qs_control%high_level_embed_subsys)
844 CALL section_vals_val_get(qs_section, "DFET_EMBEDDED", l_val=qs_control%dfet_embedded)
845 CALL section_vals_val_get(qs_section, "DMFET_EMBEDDED", l_val=qs_control%dmfet_embedded)
846
847 ! Integers gapw
848 CALL section_vals_val_get(qs_section, "LMAXN1", i_val=qs_control%gapw_control%lmax_sphere)
849 CALL section_vals_val_get(qs_section, "LMAXN0", i_val=qs_control%gapw_control%lmax_rho0)
850 CALL section_vals_val_get(qs_section, "LADDN0", i_val=qs_control%gapw_control%ladd_rho0)
851 CALL section_vals_val_get(qs_section, "QUADRATURE", i_val=qs_control%gapw_control%quadrature)
852 ! GAPW 1c basis
853 CALL section_vals_val_get(qs_section, "GAPW_1C_BASIS", i_val=qs_control%gapw_control%basis_1c)
854 IF (qs_control%gapw_control%basis_1c /= gapw_1c_orb) THEN
855 qs_control%gapw_control%eps_svd = max(qs_control%gapw_control%eps_svd, 1.e-12_dp)
856 END IF
857
858 ! Integers grids
859 CALL section_vals_val_get(qs_section, "PW_GRID", i_val=itmp)
860 SELECT CASE (itmp)
862 qs_control%pw_grid_opt%spherical = .true.
863 qs_control%pw_grid_opt%fullspace = .false.
865 qs_control%pw_grid_opt%spherical = .false.
866 qs_control%pw_grid_opt%fullspace = .true.
868 qs_control%pw_grid_opt%spherical = .false.
869 qs_control%pw_grid_opt%fullspace = .false.
870 END SELECT
871
872 ! Method for PPL calculation
873 CALL section_vals_val_get(qs_section, "CORE_PPL", i_val=itmp)
874 qs_control%do_ppl_method = itmp
875
876 CALL section_vals_val_get(qs_section, "PW_GRID_LAYOUT", i_vals=tmplist)
877 qs_control%pw_grid_opt%distribution_layout = tmplist
878 CALL section_vals_val_get(qs_section, "PW_GRID_BLOCKED", i_val=qs_control%pw_grid_opt%blocked)
879
880 !Integers extrapolation
881 CALL section_vals_val_get(qs_section, "EXTRAPOLATION", i_val=qs_control%wf_interpolation_method_nr)
882 CALL section_vals_val_get(qs_section, "EXTRAPOLATION_ORDER", i_val=qs_control%wf_extrapolation_order)
883
884 !Method
885 CALL section_vals_val_get(qs_section, "METHOD", i_val=qs_control%method_id)
886 qs_control%gapw = .false.
887 qs_control%gapw_xc = .false.
888 qs_control%gpw = .false.
889 qs_control%pao = .false.
890 qs_control%dftb = .false.
891 qs_control%xtb = .false.
892 qs_control%semi_empirical = .false.
893 qs_control%ofgpw = .false.
894 qs_control%lrigpw = .false.
895 qs_control%rigpw = .false.
896 SELECT CASE (qs_control%method_id)
897 CASE (do_method_gapw)
898 CALL cite_reference(lippert1999)
899 CALL cite_reference(krack2000)
900 qs_control%gapw = .true.
901 CASE (do_method_gapw_xc)
902 qs_control%gapw_xc = .true.
903 CASE (do_method_gpw)
904 CALL cite_reference(lippert1997)
905 CALL cite_reference(vandevondele2005a)
906 qs_control%gpw = .true.
907 CASE (do_method_ofgpw)
908 qs_control%ofgpw = .true.
909 CASE (do_method_lrigpw)
910 qs_control%lrigpw = .true.
911 CASE (do_method_rigpw)
912 qs_control%rigpw = .true.
913 CASE (do_method_dftb)
914 qs_control%dftb = .true.
915 CALL cite_reference(porezag1995)
916 CALL cite_reference(seifert1996)
917 CASE (do_method_xtb)
918 qs_control%xtb = .true.
919 CALL cite_reference(grimme2017)
920 CASE (do_method_mndo)
921 CALL cite_reference(dewar1977)
922 qs_control%semi_empirical = .true.
923 CASE (do_method_am1)
924 CALL cite_reference(dewar1985)
925 qs_control%semi_empirical = .true.
926 CASE (do_method_pm3)
927 CALL cite_reference(stewart1989)
928 qs_control%semi_empirical = .true.
929 CASE (do_method_pnnl)
930 CALL cite_reference(schenter2008)
931 qs_control%semi_empirical = .true.
932 CASE (do_method_pm6)
933 CALL cite_reference(stewart2007)
934 qs_control%semi_empirical = .true.
935 CASE (do_method_pm6fm)
936 CALL cite_reference(vanvoorhis2015)
937 qs_control%semi_empirical = .true.
938 CASE (do_method_pdg)
939 CALL cite_reference(repasky2002)
940 qs_control%semi_empirical = .true.
941 CASE (do_method_rm1)
942 CALL cite_reference(rocha2006)
943 qs_control%semi_empirical = .true.
944 CASE (do_method_mndod)
945 CALL cite_reference(dewar1977)
946 CALL cite_reference(thiel1992)
947 qs_control%semi_empirical = .true.
948 END SELECT
949
950 CALL section_vals_get(mull_section, explicit=qs_control%mulliken_restraint)
951
952 IF (qs_control%mulliken_restraint) THEN
953 CALL section_vals_val_get(mull_section, "STRENGTH", r_val=qs_control%mulliken_restraint_control%strength)
954 CALL section_vals_val_get(mull_section, "TARGET", r_val=qs_control%mulliken_restraint_control%target)
955 CALL section_vals_val_get(mull_section, "ATOMS", n_rep_val=n_rep)
956 jj = 0
957 DO k = 1, n_rep
958 CALL section_vals_val_get(mull_section, "ATOMS", i_rep_val=k, i_vals=tmplist)
959 jj = jj + SIZE(tmplist)
960 END DO
961 qs_control%mulliken_restraint_control%natoms = jj
962 IF (qs_control%mulliken_restraint_control%natoms < 1) &
963 cpabort("Need at least 1 atom to use mulliken constraints")
964 ALLOCATE (qs_control%mulliken_restraint_control%atoms(qs_control%mulliken_restraint_control%natoms))
965 jj = 0
966 DO k = 1, n_rep
967 CALL section_vals_val_get(mull_section, "ATOMS", i_rep_val=k, i_vals=tmplist)
968 DO j = 1, SIZE(tmplist)
969 jj = jj + 1
970 qs_control%mulliken_restraint_control%atoms(jj) = tmplist(j)
971 END DO
972 END DO
973 END IF
974 CALL section_vals_get(ddapc_restraint_section, n_repetition=nrep, explicit=qs_control%ddapc_restraint)
975 IF (qs_control%ddapc_restraint) THEN
976 ALLOCATE (qs_control%ddapc_restraint_control(nrep))
977 CALL read_ddapc_section(qs_control, qs_section=qs_section)
978 qs_control%ddapc_restraint_is_spin = .false.
979 qs_control%ddapc_explicit_potential = .false.
980 END IF
981
982 CALL section_vals_get(s2_restraint_section, explicit=qs_control%s2_restraint)
983 IF (qs_control%s2_restraint) THEN
984 CALL section_vals_val_get(s2_restraint_section, "STRENGTH", &
985 r_val=qs_control%s2_restraint_control%strength)
986 CALL section_vals_val_get(s2_restraint_section, "TARGET", &
987 r_val=qs_control%s2_restraint_control%target)
988 CALL section_vals_val_get(s2_restraint_section, "FUNCTIONAL_FORM", &
989 i_val=qs_control%s2_restraint_control%functional_form)
990 END IF
991
992 CALL section_vals_get(cdft_control_section, explicit=qs_control%cdft)
993 IF (qs_control%cdft) THEN
994 CALL read_cdft_control_section(qs_control, cdft_control_section)
995 END IF
996
997 ! Semi-empirical code
998 IF (qs_control%semi_empirical) THEN
999 CALL section_vals_val_get(se_section, "ORTHOGONAL_BASIS", &
1000 l_val=qs_control%se_control%orthogonal_basis)
1001 CALL section_vals_val_get(se_section, "DELTA", &
1002 r_val=qs_control%se_control%delta)
1003 CALL section_vals_val_get(se_section, "ANALYTICAL_GRADIENTS", &
1004 l_val=qs_control%se_control%analytical_gradients)
1005 CALL section_vals_val_get(se_section, "FORCE_KDSO-D_EXCHANGE", &
1006 l_val=qs_control%se_control%force_kdsod_EX)
1007 ! Integral Screening
1008 CALL section_vals_val_get(se_section, "INTEGRAL_SCREENING", &
1009 i_val=qs_control%se_control%integral_screening)
1010 IF (qs_control%method_id == do_method_pnnl) THEN
1011 IF (qs_control%se_control%integral_screening /= do_se_is_slater) &
1012 CALL cp_warn(__location__, &
1013 "PNNL semi-empirical parameterization supports only the Slater type "// &
1014 "integral scheme. Revert to Slater and continue the calculation.")
1015 qs_control%se_control%integral_screening = do_se_is_slater
1016 END IF
1017 ! Global Arrays variable
1018 CALL section_vals_val_get(se_section, "GA%NCELLS", &
1019 i_val=qs_control%se_control%ga_ncells)
1020 ! Long-Range correction
1021 CALL section_vals_val_get(se_section, "LR_CORRECTION%CUTOFF", &
1022 r_val=qs_control%se_control%cutoff_lrc)
1023 qs_control%se_control%taper_lrc = qs_control%se_control%cutoff_lrc
1024 CALL section_vals_val_get(se_section, "LR_CORRECTION%RC_TAPER", &
1025 explicit=explicit)
1026 IF (explicit) THEN
1027 CALL section_vals_val_get(se_section, "LR_CORRECTION%RC_TAPER", &
1028 r_val=qs_control%se_control%taper_lrc)
1029 END IF
1030 CALL section_vals_val_get(se_section, "LR_CORRECTION%RC_RANGE", &
1031 r_val=qs_control%se_control%range_lrc)
1032 ! Coulomb
1033 CALL section_vals_val_get(se_section, "COULOMB%CUTOFF", &
1034 r_val=qs_control%se_control%cutoff_cou)
1035 qs_control%se_control%taper_cou = qs_control%se_control%cutoff_cou
1036 CALL section_vals_val_get(se_section, "COULOMB%RC_TAPER", &
1037 explicit=explicit)
1038 IF (explicit) THEN
1039 CALL section_vals_val_get(se_section, "COULOMB%RC_TAPER", &
1040 r_val=qs_control%se_control%taper_cou)
1041 END IF
1042 CALL section_vals_val_get(se_section, "COULOMB%RC_RANGE", &
1043 r_val=qs_control%se_control%range_cou)
1044 ! Exchange
1045 CALL section_vals_val_get(se_section, "EXCHANGE%CUTOFF", &
1046 r_val=qs_control%se_control%cutoff_exc)
1047 qs_control%se_control%taper_exc = qs_control%se_control%cutoff_exc
1048 CALL section_vals_val_get(se_section, "EXCHANGE%RC_TAPER", &
1049 explicit=explicit)
1050 IF (explicit) THEN
1051 CALL section_vals_val_get(se_section, "EXCHANGE%RC_TAPER", &
1052 r_val=qs_control%se_control%taper_exc)
1053 END IF
1054 CALL section_vals_val_get(se_section, "EXCHANGE%RC_RANGE", &
1055 r_val=qs_control%se_control%range_exc)
1056 ! Screening (only if the integral scheme is of dumped type)
1057 IF (qs_control%se_control%integral_screening == do_se_is_kdso_d) THEN
1058 CALL section_vals_val_get(se_section, "SCREENING%RC_TAPER", &
1059 r_val=qs_control%se_control%taper_scr)
1060 CALL section_vals_val_get(se_section, "SCREENING%RC_RANGE", &
1061 r_val=qs_control%se_control%range_scr)
1062 END IF
1063 ! Periodic Type Calculation
1064 CALL section_vals_val_get(se_section, "PERIODIC", &
1065 i_val=qs_control%se_control%periodic_type)
1066 SELECT CASE (qs_control%se_control%periodic_type)
1067 CASE (do_se_lr_none)
1068 qs_control%se_control%do_ewald = .false.
1069 qs_control%se_control%do_ewald_r3 = .false.
1070 qs_control%se_control%do_ewald_gks = .false.
1071 CASE (do_se_lr_ewald)
1072 qs_control%se_control%do_ewald = .true.
1073 qs_control%se_control%do_ewald_r3 = .false.
1074 qs_control%se_control%do_ewald_gks = .false.
1075 CASE (do_se_lr_ewald_gks)
1076 qs_control%se_control%do_ewald = .false.
1077 qs_control%se_control%do_ewald_r3 = .false.
1078 qs_control%se_control%do_ewald_gks = .true.
1079 IF (qs_control%method_id /= do_method_pnnl) &
1080 CALL cp_abort(__location__, &
1081 "A periodic semi-empirical calculation was requested with a long-range "// &
1082 "summation on the single integral evaluation. This scheme is supported "// &
1083 "only by the PNNL parameterization.")
1084 CASE (do_se_lr_ewald_r3)
1085 qs_control%se_control%do_ewald = .true.
1086 qs_control%se_control%do_ewald_r3 = .true.
1087 qs_control%se_control%do_ewald_gks = .false.
1088 IF (qs_control%se_control%integral_screening /= do_se_is_kdso) &
1089 CALL cp_abort(__location__, &
1090 "A periodic semi-empirical calculation was requested with a long-range "// &
1091 "summation for the slowly convergent part 1/R^3, which is not congruent "// &
1092 "with the integral screening chosen. The only integral screening supported "// &
1093 "by this periodic type calculation is the standard Klopman-Dewar-Sabelli-Ohno.")
1094 END SELECT
1095
1096 ! dispersion pair potentials
1097 CALL section_vals_val_get(se_section, "DISPERSION", &
1098 l_val=qs_control%se_control%dispersion)
1099 CALL section_vals_val_get(se_section, "DISPERSION_RADIUS", &
1100 r_val=qs_control%se_control%rcdisp)
1101 CALL section_vals_val_get(se_section, "COORDINATION_CUTOFF", &
1102 r_val=qs_control%se_control%epscn)
1103 CALL section_vals_val_get(se_section, "D3_SCALING", r_vals=scal)
1104 qs_control%se_control%sd3(1) = scal(1)
1105 qs_control%se_control%sd3(2) = scal(2)
1106 qs_control%se_control%sd3(3) = scal(3)
1107 CALL section_vals_val_get(se_section, "DISPERSION_PARAMETER_FILE", &
1108 c_val=qs_control%se_control%dispersion_parameter_file)
1109
1110 ! Stop the execution for non-implemented features
1111 IF (qs_control%se_control%periodic_type == do_se_lr_ewald_r3) THEN
1112 cpabort("EWALD_R3 not implemented yet!")
1113 END IF
1114
1115 IF (qs_control%method_id == do_method_mndo .OR. &
1116 qs_control%method_id == do_method_am1 .OR. &
1117 qs_control%method_id == do_method_mndod .OR. &
1118 qs_control%method_id == do_method_pdg .OR. &
1119 qs_control%method_id == do_method_pm3 .OR. &
1120 qs_control%method_id == do_method_pm6 .OR. &
1121 qs_control%method_id == do_method_pm6fm .OR. &
1122 qs_control%method_id == do_method_pnnl .OR. &
1123 qs_control%method_id == do_method_rm1) THEN
1124 qs_control%se_control%orthogonal_basis = .true.
1125 END IF
1126 END IF
1127
1128 ! DFTB code
1129 IF (qs_control%dftb) THEN
1130 CALL section_vals_val_get(dftb_section, "ORTHOGONAL_BASIS", &
1131 l_val=qs_control%dftb_control%orthogonal_basis)
1132 CALL section_vals_val_get(dftb_section, "SELF_CONSISTENT", &
1133 l_val=qs_control%dftb_control%self_consistent)
1134 CALL section_vals_val_get(dftb_section, "DISPERSION", &
1135 l_val=qs_control%dftb_control%dispersion)
1136 CALL section_vals_val_get(dftb_section, "DIAGONAL_DFTB3", &
1137 l_val=qs_control%dftb_control%dftb3_diagonal)
1138 CALL section_vals_val_get(dftb_section, "HB_SR_GAMMA", &
1139 l_val=qs_control%dftb_control%hb_sr_damp)
1140 CALL section_vals_val_get(dftb_section, "EPS_DISP", &
1141 r_val=qs_control%dftb_control%eps_disp)
1142 CALL section_vals_val_get(dftb_section, "DO_EWALD", explicit=explicit)
1143 IF (explicit) THEN
1144 CALL section_vals_val_get(dftb_section, "DO_EWALD", &
1145 l_val=qs_control%dftb_control%do_ewald)
1146 ELSE
1147 qs_control%dftb_control%do_ewald = (qs_control%periodicity /= 0)
1148 END IF
1149 CALL section_vals_val_get(dftb_parameter, "PARAM_FILE_PATH", &
1150 c_val=qs_control%dftb_control%sk_file_path)
1151 CALL section_vals_val_get(dftb_parameter, "PARAM_FILE_NAME", &
1152 c_val=qs_control%dftb_control%sk_file_list)
1153 CALL section_vals_val_get(dftb_parameter, "HB_SR_PARAM", &
1154 r_val=qs_control%dftb_control%hb_sr_para)
1155 CALL section_vals_val_get(dftb_parameter, "SK_FILE", n_rep_val=n_var)
1156 ALLOCATE (qs_control%dftb_control%sk_pair_list(3, n_var))
1157 DO k = 1, n_var
1158 CALL section_vals_val_get(dftb_parameter, "SK_FILE", i_rep_val=k, &
1159 c_vals=clist)
1160 qs_control%dftb_control%sk_pair_list(1:3, k) = clist(1:3)
1161 END DO
1162 ! Dispersion type
1163 CALL section_vals_val_get(dftb_parameter, "DISPERSION_TYPE", &
1164 i_val=qs_control%dftb_control%dispersion_type)
1165 CALL section_vals_val_get(dftb_parameter, "UFF_FORCE_FIELD", &
1166 c_val=qs_control%dftb_control%uff_force_field)
1167 ! D3 Dispersion
1168 CALL section_vals_val_get(dftb_parameter, "DISPERSION_RADIUS", &
1169 r_val=qs_control%dftb_control%rcdisp)
1170 CALL section_vals_val_get(dftb_parameter, "COORDINATION_CUTOFF", &
1171 r_val=qs_control%dftb_control%epscn)
1172 CALL section_vals_val_get(dftb_parameter, "D2_EXP_PRE", &
1173 r_val=qs_control%dftb_control%exp_pre)
1174 CALL section_vals_val_get(dftb_parameter, "D2_SCALING", &
1175 r_val=qs_control%dftb_control%scaling)
1176 CALL section_vals_val_get(dftb_parameter, "D3_SCALING", r_vals=scal)
1177 qs_control%dftb_control%sd3(1) = scal(1)
1178 qs_control%dftb_control%sd3(2) = scal(2)
1179 qs_control%dftb_control%sd3(3) = scal(3)
1180 CALL section_vals_val_get(dftb_parameter, "D3BJ_SCALING", r_vals=scal)
1181 qs_control%dftb_control%sd3bj(1) = scal(1)
1182 qs_control%dftb_control%sd3bj(2) = scal(2)
1183 qs_control%dftb_control%sd3bj(3) = scal(3)
1184 qs_control%dftb_control%sd3bj(4) = scal(4)
1185 CALL section_vals_val_get(dftb_parameter, "DISPERSION_PARAMETER_FILE", &
1186 c_val=qs_control%dftb_control%dispersion_parameter_file)
1187
1188 IF (qs_control%dftb_control%dispersion) CALL cite_reference(zhechkov2005)
1189 IF (qs_control%dftb_control%self_consistent) CALL cite_reference(elstner1998)
1190 IF (qs_control%dftb_control%hb_sr_damp) CALL cite_reference(hu2007)
1191 END IF
1192
1193 ! xTB code
1194 IF (qs_control%xtb) THEN
1195 CALL section_vals_val_get(xtb_section, "DO_EWALD", explicit=explicit)
1196 IF (explicit) THEN
1197 CALL section_vals_val_get(xtb_section, "DO_EWALD", &
1198 l_val=qs_control%xtb_control%do_ewald)
1199 ELSE
1200 qs_control%xtb_control%do_ewald = (qs_control%periodicity /= 0)
1201 END IF
1202 CALL section_vals_val_get(xtb_section, "STO_NG", i_val=ngauss)
1203 qs_control%xtb_control%sto_ng = ngauss
1204 CALL section_vals_val_get(xtb_section, "HYDROGEN_STO_NG", i_val=ngauss)
1205 qs_control%xtb_control%h_sto_ng = ngauss
1206 CALL section_vals_val_get(xtb_parameter, "PARAM_FILE_PATH", &
1207 c_val=qs_control%xtb_control%parameter_file_path)
1208 CALL section_vals_val_get(xtb_parameter, "PARAM_FILE_NAME", &
1209 c_val=qs_control%xtb_control%parameter_file_name)
1210 ! D3 Dispersion
1211 CALL section_vals_val_get(xtb_parameter, "DISPERSION_RADIUS", &
1212 r_val=qs_control%xtb_control%rcdisp)
1213 CALL section_vals_val_get(xtb_parameter, "COORDINATION_CUTOFF", &
1214 r_val=qs_control%xtb_control%epscn)
1215 CALL section_vals_val_get(xtb_parameter, "D3BJ_SCALING", r_vals=scal)
1216 qs_control%xtb_control%s6 = scal(1)
1217 qs_control%xtb_control%s8 = scal(2)
1218 CALL section_vals_val_get(xtb_parameter, "D3BJ_PARAM", r_vals=scal)
1219 qs_control%xtb_control%a1 = scal(1)
1220 qs_control%xtb_control%a2 = scal(2)
1221 CALL section_vals_val_get(xtb_parameter, "DISPERSION_PARAMETER_FILE", &
1222 c_val=qs_control%xtb_control%dispersion_parameter_file)
1223 ! global parameters
1224 CALL section_vals_val_get(xtb_parameter, "HUCKEL_CONSTANTS", r_vals=scal)
1225 qs_control%xtb_control%ks = scal(1)
1226 qs_control%xtb_control%kp = scal(2)
1227 qs_control%xtb_control%kd = scal(3)
1228 qs_control%xtb_control%ksp = scal(4)
1229 qs_control%xtb_control%k2sh = scal(5)
1230 CALL section_vals_val_get(xtb_parameter, "COULOMB_CONSTANTS", r_vals=scal)
1231 qs_control%xtb_control%kg = scal(1)
1232 qs_control%xtb_control%kf = scal(2)
1233 CALL section_vals_val_get(xtb_parameter, "CN_CONSTANTS", r_vals=scal)
1234 qs_control%xtb_control%kcns = scal(1)
1235 qs_control%xtb_control%kcnp = scal(2)
1236 qs_control%xtb_control%kcnd = scal(3)
1237 CALL section_vals_val_get(xtb_parameter, "EN_CONSTANT", r_vals=scal)
1238 qs_control%xtb_control%ken = scal(1)
1239 ! XB
1240 CALL section_vals_val_get(xtb_section, "USE_HALOGEN_CORRECTION", &
1241 l_val=qs_control%xtb_control%xb_interaction)
1242 CALL section_vals_val_get(xtb_parameter, "HALOGEN_BINDING", r_vals=scal)
1243 qs_control%xtb_control%kxr = scal(1)
1244 qs_control%xtb_control%kx2 = scal(2)
1245 ! NONBONDED interactions
1246 CALL section_vals_val_get(xtb_section, "DO_NONBONDED", &
1247 l_val=qs_control%xtb_control%do_nonbonded)
1248 CALL section_vals_get(nonbonded_section, explicit=explicit)
1249 IF (explicit .AND. qs_control%xtb_control%do_nonbonded) THEN
1250 CALL section_vals_get(genpot_section, explicit=explicit, n_repetition=ngp)
1251 IF (explicit) THEN
1252 CALL pair_potential_reallocate(qs_control%xtb_control%nonbonded, 1, ngp, gp=.true.)
1253 CALL read_gp_section(qs_control%xtb_control%nonbonded, genpot_section, 0)
1254 END IF
1255 END IF !nonbonded
1256 ! SR Coulomb
1257 CALL section_vals_val_get(xtb_section, "OLD_COULOMB_DAMPING", &
1258 l_val=qs_control%xtb_control%old_coulomb_damping)
1259 CALL section_vals_val_get(xtb_parameter, "COULOMB_SR_CUT", r_vals=scal)
1260 qs_control%xtb_control%coulomb_sr_cut = scal(1)
1261 CALL section_vals_val_get(xtb_parameter, "COULOMB_SR_EPS", r_vals=scal)
1262 qs_control%xtb_control%coulomb_sr_eps = scal(1)
1263 ! XB_radius
1264 CALL section_vals_val_get(xtb_parameter, "XB_RADIUS", r_val=qs_control%xtb_control%xb_radius)
1265 ! Kab
1266 CALL section_vals_val_get(xtb_parameter, "KAB_PARAM", n_rep_val=n_rep)
1267 ! For debug purposes
1268 CALL section_vals_val_get(xtb_section, "COULOMB_INTERACTION", &
1269 l_val=qs_control%xtb_control%coulomb_interaction)
1270 CALL section_vals_val_get(xtb_section, "COULOMB_LR", &
1271 l_val=qs_control%xtb_control%coulomb_lr)
1272 CALL section_vals_val_get(xtb_section, "TB3_INTERACTION", &
1273 l_val=qs_control%xtb_control%tb3_interaction)
1274 ! Check for bad atomic charges
1275 CALL section_vals_val_get(xtb_section, "CHECK_ATOMIC_CHARGES", &
1276 l_val=qs_control%xtb_control%check_atomic_charges)
1277
1278 qs_control%xtb_control%kab_nval = n_rep
1279 IF (n_rep > 0) THEN
1280 ALLOCATE (qs_control%xtb_control%kab_param(3, n_rep))
1281 ALLOCATE (qs_control%xtb_control%kab_types(2, n_rep))
1282 ALLOCATE (qs_control%xtb_control%kab_vals(n_rep))
1283 DO j = 1, n_rep
1284 CALL section_vals_val_get(xtb_parameter, "KAB_PARAM", i_rep_val=j, c_vals=clist)
1285 qs_control%xtb_control%kab_param(1, j) = clist(1)
1286 CALL get_ptable_info(clist(1), &
1287 ielement=qs_control%xtb_control%kab_types(1, j))
1288 qs_control%xtb_control%kab_param(2, j) = clist(2)
1289 CALL get_ptable_info(clist(2), &
1290 ielement=qs_control%xtb_control%kab_types(2, j))
1291 qs_control%xtb_control%kab_param(3, j) = clist(3)
1292 READ (clist(3), '(F10.0)') qs_control%xtb_control%kab_vals(j)
1293 END DO
1294 END IF
1295 END IF
1296
1297 ! Optimize LRI basis set
1298 CALL section_vals_get(lri_optbas_section, explicit=qs_control%lri_optbas)
1299
1300 CALL timestop(handle)
1301 END SUBROUTINE read_qs_section
1302
1303! **************************************************************************************************
1304!> \brief ...
1305!> \param t_control ...
1306!> \param dft_section ...
1307! **************************************************************************************************
1308 SUBROUTINE read_tddfpt_control(t_control, dft_section)
1309 TYPE(tddfpt_control_type) :: t_control
1310 TYPE(section_vals_type), POINTER :: dft_section
1311
1312 LOGICAL :: kenergy_den
1313 TYPE(section_vals_type), POINTER :: sic_section, t_section
1314
1315 kenergy_den = .false.
1316 NULLIFY (sic_section, t_section)
1317 t_section => section_vals_get_subs_vals(dft_section, "TDDFPT")
1318
1319 CALL section_vals_val_get(t_section, "CONVERGENCE", r_val=t_control%tolerance)
1320 CALL section_vals_val_get(t_section, "NEV", i_val=t_control%n_ev)
1321 CALL section_vals_val_get(t_section, "MAX_KV", i_val=t_control%max_kv)
1322 CALL section_vals_val_get(t_section, "RESTARTS", i_val=t_control%n_restarts)
1323 CALL section_vals_val_get(t_section, "NREORTHO", i_val=t_control%n_reortho)
1324 CALL section_vals_val_get(t_section, "RES_ETYPE", i_val=t_control%res_etype)
1325 CALL section_vals_val_get(t_section, "DIAG_METHOD", i_val=t_control%diag_method)
1326 CALL section_vals_val_get(t_section, "KERNEL", l_val=t_control%do_kernel)
1327 CALL section_vals_val_get(t_section, "LSD_SINGLETS", l_val=t_control%lsd_singlets)
1328 CALL section_vals_val_get(t_section, "INVERT_S", l_val=t_control%invert_S)
1329 CALL section_vals_val_get(t_section, "PRECOND", l_val=t_control%precond)
1330 CALL section_vals_val_get(t_section, "OE_CORR", i_val=t_control%oe_corr)
1331
1332 t_control%use_kinetic_energy_density = .false.
1333 sic_section => section_vals_get_subs_vals(t_section, "SIC")
1334 CALL section_vals_val_get(sic_section, "SIC_METHOD", i_val=t_control%sic_method_id)
1335 CALL section_vals_val_get(sic_section, "ORBITAL_SET", i_val=t_control%sic_list_id)
1336 CALL section_vals_val_get(sic_section, "SIC_SCALING_A", r_val=t_control%sic_scaling_a)
1337 CALL section_vals_val_get(sic_section, "SIC_SCALING_B", r_val=t_control%sic_scaling_b)
1338
1339 END SUBROUTINE read_tddfpt_control
1340
1341! **************************************************************************************************
1342!> \brief Read TDDFPT-related input parameters.
1343!> \param t_control TDDFPT control parameters
1344!> \param t_section TDDFPT input section
1345!> \param qs_control Quickstep control parameters
1346! **************************************************************************************************
1347 SUBROUTINE read_tddfpt2_control(t_control, t_section, qs_control)
1348 TYPE(tddfpt2_control_type), POINTER :: t_control
1349 TYPE(section_vals_type), POINTER :: t_section
1350 TYPE(qs_control_type), POINTER :: qs_control
1351
1352 CHARACTER(LEN=*), PARAMETER :: routinen = 'read_tddfpt2_control'
1353
1354 CHARACTER(LEN=default_string_length), &
1355 DIMENSION(:), POINTER :: tmpstringlist
1356 INTEGER :: handle, irep, isize, nrep
1357 INTEGER, ALLOCATABLE, DIMENSION(:) :: inds
1358 LOGICAL :: do_ewald, do_exchange, expl, explicit, &
1359 multigrid_set
1360 REAL(kind=dp) :: filter, fval, hfx
1361 TYPE(section_vals_type), POINTER :: dipole_section, mgrid_section, &
1362 soc_section, stda_section, xc_func, &
1363 xc_section
1364
1365 CALL timeset(routinen, handle)
1366
1367 CALL section_vals_val_get(t_section, "_SECTION_PARAMETERS_", l_val=t_control%enabled)
1368
1369 CALL section_vals_val_get(t_section, "NSTATES", i_val=t_control%nstates)
1370 CALL section_vals_val_get(t_section, "MAX_ITER", i_val=t_control%niters)
1371 CALL section_vals_val_get(t_section, "MAX_KV", i_val=t_control%nkvs)
1372 CALL section_vals_val_get(t_section, "NLUMO", i_val=t_control%nlumo)
1373 CALL section_vals_val_get(t_section, "NPROC_STATE", i_val=t_control%nprocs)
1374 CALL section_vals_val_get(t_section, "KERNEL", i_val=t_control%kernel)
1375 CALL section_vals_val_get(t_section, "OE_CORR", i_val=t_control%oe_corr)
1376 CALL section_vals_val_get(t_section, "EV_SHIFT", r_val=t_control%ev_shift)
1377 CALL section_vals_val_get(t_section, "EOS_SHIFT", r_val=t_control%eos_shift)
1378
1379 CALL section_vals_val_get(t_section, "CONVERGENCE", r_val=t_control%conv)
1380 CALL section_vals_val_get(t_section, "MIN_AMPLITUDE", r_val=t_control%min_excitation_amplitude)
1381 CALL section_vals_val_get(t_section, "ORTHOGONAL_EPS", r_val=t_control%orthogonal_eps)
1382
1383 CALL section_vals_val_get(t_section, "RESTART", l_val=t_control%is_restart)
1384 CALL section_vals_val_get(t_section, "RKS_TRIPLETS", l_val=t_control%rks_triplets)
1385 CALL section_vals_val_get(t_section, "DO_LRIGPW", l_val=t_control%do_lrigpw)
1386 CALL section_vals_val_get(t_section, "ADMM_KERNEL_CORRECTION_SYMMETRIC", l_val=t_control%admm_symm)
1387 CALL section_vals_val_get(t_section, "ADMM_KERNEL_XC_CORRECTION", l_val=t_control%admm_xc_correction)
1388
1389 ! read automatically generated auxiliary basis for LRI
1390 CALL section_vals_val_get(t_section, "AUTO_BASIS", n_rep_val=nrep)
1391 DO irep = 1, nrep
1392 CALL section_vals_val_get(t_section, "AUTO_BASIS", i_rep_val=irep, c_vals=tmpstringlist)
1393 IF (SIZE(tmpstringlist) == 2) THEN
1394 CALL uppercase(tmpstringlist(2))
1395 SELECT CASE (tmpstringlist(2))
1396 CASE ("X")
1397 isize = -1
1398 CASE ("SMALL")
1399 isize = 0
1400 CASE ("MEDIUM")
1401 isize = 1
1402 CASE ("LARGE")
1403 isize = 2
1404 CASE ("HUGE")
1405 isize = 3
1406 CASE DEFAULT
1407 cpabort("Unknown basis size in AUTO_BASIS keyword:"//trim(tmpstringlist(1)))
1408 END SELECT
1409 !
1410 SELECT CASE (tmpstringlist(1))
1411 CASE ("X")
1412 CASE ("P_LRI_AUX")
1413 t_control%auto_basis_p_lri_aux = isize
1414 CASE DEFAULT
1415 cpabort("Unknown basis type in AUTO_BASIS keyword:"//trim(tmpstringlist(1)))
1416 END SELECT
1417 ELSE
1418 CALL cp_abort(__location__, &
1419 "AUTO_BASIS keyword in &PROPERTIES &TDDFT section has a wrong number of arguments.")
1420 END IF
1421 END DO
1422
1423 IF (t_control%conv < 0) &
1424 t_control%conv = abs(t_control%conv)
1425
1426 ! DIPOLE_MOMENTS subsection
1427 dipole_section => section_vals_get_subs_vals(t_section, "DIPOLE_MOMENTS")
1428 CALL section_vals_val_get(dipole_section, "DIPOLE_FORM", explicit=explicit)
1429 IF (explicit) THEN
1430 CALL section_vals_val_get(dipole_section, "DIPOLE_FORM", i_val=t_control%dipole_form)
1431 ELSE
1432 t_control%dipole_form = 0
1433 END IF
1434 CALL section_vals_val_get(dipole_section, "REFERENCE", i_val=t_control%dipole_reference)
1435 CALL section_vals_val_get(dipole_section, "REFERENCE_POINT", explicit=explicit)
1436 IF (explicit) THEN
1437 CALL section_vals_val_get(dipole_section, "REFERENCE_POINT", r_vals=t_control%dipole_ref_point)
1438 ELSE
1439 NULLIFY (t_control%dipole_ref_point)
1440 IF (t_control%dipole_form == tddfpt_dipole_length .AND. t_control%dipole_reference == use_mom_ref_user) THEN
1441 cpabort("User-defined reference point should be given explicitly")
1442 END IF
1443 END IF
1444
1445 !SOC subsection
1446 soc_section => section_vals_get_subs_vals(t_section, "SOC")
1447 CALL section_vals_get(soc_section, explicit=explicit)
1448 IF (explicit) THEN
1449 t_control%do_soc = .true.
1450 END IF
1451
1452 ! MGRID subsection
1453 mgrid_section => section_vals_get_subs_vals(t_section, "MGRID")
1454 CALL section_vals_get(mgrid_section, explicit=t_control%mgrid_is_explicit)
1455
1456 IF (t_control%mgrid_is_explicit) THEN
1457 CALL section_vals_val_get(mgrid_section, "NGRIDS", i_val=t_control%mgrid_ngrids, explicit=explicit)
1458 IF (.NOT. explicit) t_control%mgrid_ngrids = SIZE(qs_control%e_cutoff)
1459
1460 CALL section_vals_val_get(mgrid_section, "CUTOFF", r_val=t_control%mgrid_cutoff, explicit=explicit)
1461 IF (.NOT. explicit) t_control%mgrid_cutoff = qs_control%cutoff
1462
1463 CALL section_vals_val_get(mgrid_section, "PROGRESSION_FACTOR", &
1464 r_val=t_control%mgrid_progression_factor, explicit=explicit)
1465 IF (explicit) THEN
1466 IF (t_control%mgrid_progression_factor <= 1.0_dp) &
1467 CALL cp_abort(__location__, &
1468 "Progression factor should be greater then 1.0 to ensure multi-grid ordering")
1469 ELSE
1470 t_control%mgrid_progression_factor = qs_control%progression_factor
1471 END IF
1472
1473 CALL section_vals_val_get(mgrid_section, "COMMENSURATE", l_val=t_control%mgrid_commensurate_mgrids, explicit=explicit)
1474 IF (.NOT. explicit) t_control%mgrid_commensurate_mgrids = qs_control%commensurate_mgrids
1475 IF (t_control%mgrid_commensurate_mgrids) THEN
1476 IF (explicit) THEN
1477 t_control%mgrid_progression_factor = 4.0_dp
1478 ELSE
1479 t_control%mgrid_progression_factor = qs_control%progression_factor
1480 END IF
1481 END IF
1482
1483 CALL section_vals_val_get(mgrid_section, "REL_CUTOFF", r_val=t_control%mgrid_relative_cutoff, explicit=explicit)
1484 IF (.NOT. explicit) t_control%mgrid_relative_cutoff = qs_control%relative_cutoff
1485
1486 CALL section_vals_val_get(mgrid_section, "MULTIGRID_SET", l_val=multigrid_set, explicit=explicit)
1487 IF (.NOT. explicit) multigrid_set = .false.
1488 IF (multigrid_set) THEN
1489 CALL section_vals_val_get(mgrid_section, "MULTIGRID_CUTOFF", r_vals=t_control%mgrid_e_cutoff)
1490 ELSE
1491 NULLIFY (t_control%mgrid_e_cutoff)
1492 END IF
1493
1494 CALL section_vals_val_get(mgrid_section, "REALSPACE", l_val=t_control%mgrid_realspace_mgrids, explicit=explicit)
1495 IF (.NOT. explicit) t_control%mgrid_realspace_mgrids = qs_control%realspace_mgrids
1496
1497 CALL section_vals_val_get(mgrid_section, "SKIP_LOAD_BALANCE_DISTRIBUTED", &
1498 l_val=t_control%mgrid_skip_load_balance, explicit=explicit)
1499 IF (.NOT. explicit) t_control%mgrid_skip_load_balance = qs_control%skip_load_balance_distributed
1500
1501 IF (ASSOCIATED(t_control%mgrid_e_cutoff)) THEN
1502 IF (SIZE(t_control%mgrid_e_cutoff) /= t_control%mgrid_ngrids) &
1503 cpabort("Inconsistent values for number of multi-grids")
1504
1505 ! sort multi-grids in descending order according to their cutoff values
1506 t_control%mgrid_e_cutoff = -t_control%mgrid_e_cutoff
1507 ALLOCATE (inds(t_control%mgrid_ngrids))
1508 CALL sort(t_control%mgrid_e_cutoff, t_control%mgrid_ngrids, inds)
1509 DEALLOCATE (inds)
1510 t_control%mgrid_e_cutoff = -t_control%mgrid_e_cutoff
1511 END IF
1512 END IF
1513
1514 ! expand XC subsection (if given explicitly)
1515 xc_section => section_vals_get_subs_vals(t_section, "XC")
1516 xc_func => section_vals_get_subs_vals(xc_section, "XC_FUNCTIONAL")
1517 CALL section_vals_get(xc_func, explicit=explicit)
1518 IF (explicit) &
1519 CALL xc_functionals_expand(xc_func, xc_section)
1520
1521 ! sTDA subsection
1522 stda_section => section_vals_get_subs_vals(t_section, "STDA")
1523 IF (t_control%kernel == tddfpt_kernel_stda) THEN
1524 t_control%stda_control%hfx_fraction = 0.0_dp
1525 t_control%stda_control%do_exchange = .true.
1526 t_control%stda_control%eps_td_filter = 1.e-10_dp
1527 t_control%stda_control%mn_alpha = -99.0_dp
1528 t_control%stda_control%mn_beta = -99.0_dp
1529 ! set default for Ewald method (on/off) dependent on periodicity
1530 SELECT CASE (qs_control%periodicity)
1531 CASE (0)
1532 t_control%stda_control%do_ewald = .false.
1533 CASE (1)
1534 t_control%stda_control%do_ewald = .true.
1535 CASE (2)
1536 t_control%stda_control%do_ewald = .true.
1537 CASE (3)
1538 t_control%stda_control%do_ewald = .true.
1539 CASE DEFAULT
1540 cpabort("Illegal value for periodiciy")
1541 END SELECT
1542 CALL section_vals_get(stda_section, explicit=explicit)
1543 IF (explicit) THEN
1544 CALL section_vals_val_get(stda_section, "HFX_FRACTION", r_val=hfx, explicit=expl)
1545 IF (expl) t_control%stda_control%hfx_fraction = hfx
1546 CALL section_vals_val_get(stda_section, "EPS_TD_FILTER", r_val=filter, explicit=expl)
1547 IF (expl) t_control%stda_control%eps_td_filter = filter
1548 CALL section_vals_val_get(stda_section, "DO_EWALD", l_val=do_ewald, explicit=expl)
1549 IF (expl) t_control%stda_control%do_ewald = do_ewald
1550 CALL section_vals_val_get(stda_section, "DO_EXCHANGE", l_val=do_exchange, explicit=expl)
1551 IF (expl) t_control%stda_control%do_exchange = do_exchange
1552 CALL section_vals_val_get(stda_section, "MATAGA_NISHIMOTO_CEXP", r_val=fval)
1553 t_control%stda_control%mn_alpha = fval
1554 CALL section_vals_val_get(stda_section, "MATAGA_NISHIMOTO_XEXP", r_val=fval)
1555 t_control%stda_control%mn_beta = fval
1556 END IF
1557 CALL section_vals_val_get(stda_section, "COULOMB_SR_CUT", r_val=fval)
1558 t_control%stda_control%coulomb_sr_cut = fval
1559 CALL section_vals_val_get(stda_section, "COULOMB_SR_EPS", r_val=fval)
1560 t_control%stda_control%coulomb_sr_eps = fval
1561 END IF
1562
1563 CALL timestop(handle)
1564 END SUBROUTINE read_tddfpt2_control
1565
1566! **************************************************************************************************
1567!> \brief Write the DFT control parameters to the output unit.
1568!> \param dft_control ...
1569!> \param dft_section ...
1570! **************************************************************************************************
1571 SUBROUTINE write_dft_control(dft_control, dft_section)
1572 TYPE(dft_control_type), POINTER :: dft_control
1573 TYPE(section_vals_type), POINTER :: dft_section
1574
1575 CHARACTER(len=*), PARAMETER :: routinen = 'write_dft_control'
1576
1577 CHARACTER(LEN=20) :: tmpstr
1578 INTEGER :: handle, output_unit
1579 REAL(kind=dp) :: density_cut, density_smooth_cut_range, &
1580 gradient_cut, tau_cut
1581 TYPE(cp_logger_type), POINTER :: logger
1582 TYPE(enumeration_type), POINTER :: enum
1583 TYPE(keyword_type), POINTER :: keyword
1584 TYPE(section_type), POINTER :: section
1585 TYPE(section_vals_type), POINTER :: xc_section
1586
1587 IF (dft_control%qs_control%semi_empirical) RETURN
1588 IF (dft_control%qs_control%dftb) RETURN
1589 IF (dft_control%qs_control%xtb) THEN
1590 CALL write_xtb_control(dft_control%qs_control%xtb_control, dft_section)
1591 RETURN
1592 END IF
1593 CALL timeset(routinen, handle)
1594
1595 NULLIFY (logger)
1596 logger => cp_get_default_logger()
1597
1598 output_unit = cp_print_key_unit_nr(logger, dft_section, &
1599 "PRINT%DFT_CONTROL_PARAMETERS", extension=".Log")
1600
1601 IF (output_unit > 0) THEN
1602
1603 xc_section => section_vals_get_subs_vals(dft_section, "XC")
1604
1605 IF (dft_control%uks) THEN
1606 WRITE (unit=output_unit, fmt="(/,T2,A,T78,A)") &
1607 "DFT| Spin unrestricted (spin-polarized) Kohn-Sham calculation", "UKS"
1608 ELSE IF (dft_control%roks) THEN
1609 WRITE (unit=output_unit, fmt="(/,T2,A,T77,A)") &
1610 "DFT| Spin restricted open Kohn-Sham calculation", "ROKS"
1611 ELSE
1612 WRITE (unit=output_unit, fmt="(/,T2,A,T78,A)") &
1613 "DFT| Spin restricted Kohn-Sham (RKS) calculation", "RKS"
1614 END IF
1615
1616 WRITE (unit=output_unit, fmt="(T2,A,T76,I5)") &
1617 "DFT| Multiplicity", dft_control%multiplicity
1618 WRITE (unit=output_unit, fmt="(T2,A,T76,I5)") &
1619 "DFT| Number of spin states", dft_control%nspins
1620
1621 WRITE (unit=output_unit, fmt="(T2,A,T76,I5)") &
1622 "DFT| Charge", dft_control%charge
1623
1624 IF (dft_control%sic_method_id .NE. sic_none) CALL cite_reference(vandevondele2005b)
1625 SELECT CASE (dft_control%sic_method_id)
1626 CASE (sic_none)
1627 tmpstr = "NO"
1628 CASE (sic_mauri_spz)
1629 tmpstr = "SPZ/MAURI SIC"
1630 CASE (sic_mauri_us)
1631 tmpstr = "US/MAURI SIC"
1632 CASE (sic_ad)
1633 tmpstr = "AD SIC"
1634 CASE (sic_eo)
1635 tmpstr = "Explicit Orbital SIC"
1636 CASE DEFAULT
1637 ! fix throughout the cp2k for this option
1638 cpabort("SIC option unknown")
1639 END SELECT
1640
1641 WRITE (unit=output_unit, fmt="(T2,A,T61,A20)") &
1642 "DFT| Self-interaction correction (SIC)", adjustr(trim(tmpstr))
1643
1644 IF (dft_control%sic_method_id /= sic_none) THEN
1645 WRITE (unit=output_unit, fmt="(T2,A,T66,ES15.6)") &
1646 "DFT| SIC scaling parameter a", dft_control%sic_scaling_a, &
1647 "DFT| SIC scaling parameter b", dft_control%sic_scaling_b
1648 END IF
1649
1650 IF (dft_control%sic_method_id == sic_eo) THEN
1651 IF (dft_control%sic_list_id == sic_list_all) THEN
1652 WRITE (unit=output_unit, fmt="(T2,A,T66,A)") &
1653 "DFT| SIC orbitals", "ALL"
1654 END IF
1655 IF (dft_control%sic_list_id == sic_list_unpaired) THEN
1656 WRITE (unit=output_unit, fmt="(T2,A,T66,A)") &
1657 "DFT| SIC orbitals", "UNPAIRED"
1658 END IF
1659 END IF
1660
1661 CALL section_vals_val_get(xc_section, "density_cutoff", r_val=density_cut)
1662 CALL section_vals_val_get(xc_section, "gradient_cutoff", r_val=gradient_cut)
1663 CALL section_vals_val_get(xc_section, "tau_cutoff", r_val=tau_cut)
1664 CALL section_vals_val_get(xc_section, "density_smooth_cutoff_range", r_val=density_smooth_cut_range)
1665
1666 WRITE (unit=output_unit, fmt="(T2,A,T66,ES15.6)") &
1667 "DFT| Cutoffs: density ", density_cut, &
1668 "DFT| gradient", gradient_cut, &
1669 "DFT| tau ", tau_cut, &
1670 "DFT| cutoff_smoothing_range", density_smooth_cut_range
1671 CALL section_vals_val_get(xc_section, "XC_GRID%XC_SMOOTH_RHO", &
1672 c_val=tmpstr)
1673 WRITE (output_unit, '( A, T61, A )') &
1674 " DFT| XC density smoothing ", adjustr(tmpstr)
1675 CALL section_vals_val_get(xc_section, "XC_GRID%XC_DERIV", &
1676 c_val=tmpstr)
1677 WRITE (output_unit, '( A, T61, A )') &
1678 " DFT| XC derivatives ", adjustr(tmpstr)
1679 IF (dft_control%dft_plus_u) THEN
1680 NULLIFY (enum, keyword, section)
1681 CALL create_dft_section(section)
1682 keyword => section_get_keyword(section, "PLUS_U_METHOD")
1683 CALL keyword_get(keyword, enum=enum)
1684 WRITE (unit=output_unit, fmt="(/,T2,A,T41,A40)") &
1685 "DFT+U| Method", adjustr(trim(enum_i2c(enum, dft_control%plus_u_method_id)))
1686 WRITE (unit=output_unit, fmt="(T2,A)") &
1687 "DFT+U| Check atomic kind information for details"
1688 CALL section_release(section)
1689 END IF
1690
1691 WRITE (unit=output_unit, fmt="(A)") ""
1692 CALL xc_write(output_unit, xc_section, dft_control%lsd)
1693
1694 IF (dft_control%apply_period_efield) THEN
1695 WRITE (unit=output_unit, fmt="(A)") ""
1696 IF (dft_control%period_efield%displacement_field) THEN
1697 WRITE (unit=output_unit, fmt="(T2,A)") &
1698 "PERIODIC_EFIELD| Use displacement field formulation"
1699 WRITE (unit=output_unit, fmt="(T2,A,T66,1X,ES14.6)") &
1700 "PERIODIC_EFIELD| Displacement field filter: x", &
1701 dft_control%period_efield%d_filter(1), &
1702 "PERIODIC_EFIELD| y", &
1703 dft_control%period_efield%d_filter(2), &
1704 "PERIODIC_EFIELD| z", &
1705 dft_control%period_efield%d_filter(3)
1706 END IF
1707 WRITE (unit=output_unit, fmt="(T2,A,T66,1X,ES14.6)") &
1708 "PERIODIC_EFIELD| Polarisation vector: x", &
1709 dft_control%period_efield%polarisation(1), &
1710 "PERIODIC_EFIELD| y", &
1711 dft_control%period_efield%polarisation(2), &
1712 "PERIODIC_EFIELD| z", &
1713 dft_control%period_efield%polarisation(3), &
1714 "PERIODIC_EFIELD| Intensity [a.u.]:", &
1715 dft_control%period_efield%strength
1716 IF (sqrt(dot_product(dft_control%period_efield%polarisation, &
1717 dft_control%period_efield%polarisation)) < epsilon(0.0_dp)) THEN
1718 cpabort("Invalid (too small) polarisation vector specified for PERIODIC_EFIELD")
1719 END IF
1720 END IF
1721
1722 IF (dft_control%do_sccs) THEN
1723 WRITE (unit=output_unit, fmt="(/,T2,A)") &
1724 "SCCS| Self-consistent continuum solvation model"
1725 WRITE (unit=output_unit, fmt="(T2,A,T61,ES20.6)") &
1726 "SCCS| Relative permittivity of the solvent (medium)", &
1727 dft_control%sccs_control%epsilon_solvent, &
1728 "SCCS| Absolute permittivity [a.u.]", &
1729 dft_control%sccs_control%epsilon_solvent/fourpi
1730 SELECT CASE (dft_control%sccs_control%method_id)
1731 CASE (sccs_andreussi)
1732 WRITE (unit=output_unit, fmt="(T2,A,/,(T2,A,T61,ES20.6))") &
1733 "SCCS| Dielectric function proposed by Andreussi et al.", &
1734 "SCCS| rho_max", dft_control%sccs_control%rho_max, &
1735 "SCCS| rho_min", dft_control%sccs_control%rho_min
1736 CASE (sccs_fattebert_gygi)
1737 WRITE (unit=output_unit, fmt="(T2,A,/,(T2,A,T61,ES20.6))") &
1738 "SCCS| Dielectric function proposed by Fattebert and Gygi", &
1739 "SCCS| beta", dft_control%sccs_control%beta, &
1740 "SCCS| rho_zero", dft_control%sccs_control%rho_zero
1741 CASE DEFAULT
1742 cpabort("Invalid SCCS model specified. Please, check your input!")
1743 END SELECT
1744 SELECT CASE (dft_control%sccs_control%derivative_method)
1745 CASE (sccs_derivative_fft)
1746 WRITE (unit=output_unit, fmt="(T2,A,T46,A35)") &
1747 "SCCS| Numerical derivative calculation", &
1748 adjustr("FFT")
1749 CASE (sccs_derivative_cd3)
1750 WRITE (unit=output_unit, fmt="(T2,A,T46,A35)") &
1751 "SCCS| Numerical derivative calculation", &
1752 adjustr("3-point stencil central differences")
1753 CASE (sccs_derivative_cd5)
1754 WRITE (unit=output_unit, fmt="(T2,A,T46,A35)") &
1755 "SCCS| Numerical derivative calculation", &
1756 adjustr("5-point stencil central differences")
1757 CASE (sccs_derivative_cd7)
1758 WRITE (unit=output_unit, fmt="(T2,A,T46,A35)") &
1759 "SCCS| Numerical derivative calculation", &
1760 adjustr("7-point stencil central differences")
1761 CASE DEFAULT
1762 CALL cp_abort(__location__, &
1763 "Invalid derivative method specified for SCCS model. "// &
1764 "Please, check your input!")
1765 END SELECT
1766 WRITE (unit=output_unit, fmt="(T2,A,T61,ES20.6)") &
1767 "SCCS| Repulsion parameter alpha [mN/m] = [dyn/cm]", &
1768 cp_unit_from_cp2k(dft_control%sccs_control%alpha_solvent, "mN/m")
1769 WRITE (unit=output_unit, fmt="(T2,A,T61,ES20.6)") &
1770 "SCCS| Dispersion parameter beta [GPa]", &
1771 cp_unit_from_cp2k(dft_control%sccs_control%beta_solvent, "GPa")
1772 WRITE (unit=output_unit, fmt="(T2,A,T61,ES20.6)") &
1773 "SCCS| Surface tension gamma [mN/m] = [dyn/cm]", &
1774 cp_unit_from_cp2k(dft_control%sccs_control%gamma_solvent, "mN/m")
1775 WRITE (unit=output_unit, fmt="(T2,A,T61,ES20.6)") &
1776 "SCCS| Mixing parameter applied during the iteration cycle", &
1777 dft_control%sccs_control%mixing
1778 WRITE (unit=output_unit, fmt="(T2,A,T61,ES20.6)") &
1779 "SCCS| Tolerance for the convergence of the SCCS iteration cycle", &
1780 dft_control%sccs_control%eps_sccs
1781 WRITE (unit=output_unit, fmt="(T2,A,T61,I20)") &
1782 "SCCS| Maximum number of iteration steps", &
1783 dft_control%sccs_control%max_iter
1784 WRITE (unit=output_unit, fmt="(T2,A,T61,ES20.6)") &
1785 "SCCS| SCF convergence threshold for starting the SCCS iteration", &
1786 dft_control%sccs_control%eps_scf
1787 WRITE (unit=output_unit, fmt="(T2,A,T61,ES20.6)") &
1788 "SCCS| Numerical increment for the cavity surface calculation", &
1789 dft_control%sccs_control%delta_rho
1790 END IF
1791
1792 WRITE (unit=output_unit, fmt="(A)") ""
1793
1794 END IF
1795
1796 CALL cp_print_key_finished_output(output_unit, logger, dft_section, &
1797 "PRINT%DFT_CONTROL_PARAMETERS")
1798
1799 CALL timestop(handle)
1800
1801 END SUBROUTINE write_dft_control
1802
1803! **************************************************************************************************
1804!> \brief Write the ADMM control parameters to the output unit.
1805!> \param admm_control ...
1806!> \param dft_section ...
1807! **************************************************************************************************
1808 SUBROUTINE write_admm_control(admm_control, dft_section)
1809 TYPE(admm_control_type), POINTER :: admm_control
1810 TYPE(section_vals_type), POINTER :: dft_section
1811
1812 INTEGER :: iounit
1813 TYPE(cp_logger_type), POINTER :: logger
1814
1815 NULLIFY (logger)
1816 logger => cp_get_default_logger()
1817
1818 iounit = cp_print_key_unit_nr(logger, dft_section, &
1819 "PRINT%DFT_CONTROL_PARAMETERS", extension=".Log")
1820
1821 IF (iounit > 0) THEN
1822
1823 SELECT CASE (admm_control%admm_type)
1824 CASE (no_admm_type)
1825 WRITE (unit=iounit, fmt="(/,T2,A,T77,A)") "ADMM| Specific ADMM type specified", "NONE"
1826 CASE (admm1_type)
1827 WRITE (unit=iounit, fmt="(/,T2,A,T76,A)") "ADMM| Specific ADMM type specified", "ADMM1"
1828 CASE (admm2_type)
1829 WRITE (unit=iounit, fmt="(/,T2,A,T76,A)") "ADMM| Specific ADMM type specified", "ADMM2"
1830 CASE (admms_type)
1831 WRITE (unit=iounit, fmt="(/,T2,A,T76,A)") "ADMM| Specific ADMM type specified", "ADMMS"
1832 CASE (admmp_type)
1833 WRITE (unit=iounit, fmt="(/,T2,A,T76,A)") "ADMM| Specific ADMM type specified", "ADMMP"
1834 CASE (admmq_type)
1835 WRITE (unit=iounit, fmt="(/,T2,A,T76,A)") "ADMM| Specific ADMM type specified", "ADMMQ"
1836 CASE DEFAULT
1837 cpabort("admm_type")
1838 END SELECT
1839
1840 SELECT CASE (admm_control%purification_method)
1841 CASE (do_admm_purify_none)
1842 WRITE (unit=iounit, fmt="(T2,A,T77,A)") "ADMM| Density matrix purification method", "NONE"
1844 WRITE (unit=iounit, fmt="(T2,A,T75,A)") "ADMM| Density matrix purification method", "Cauchy"
1846 WRITE (unit=iounit, fmt="(T2,A,T66,A)") "ADMM| Density matrix purification method", "Cauchy subspace"
1848 WRITE (unit=iounit, fmt="(T2,A,T63,A)") "ADMM| Density matrix purification method", "MO diagonalization"
1850 WRITE (unit=iounit, fmt="(T2,A,T71,A)") "ADMM| Density matrix purification method", "MO no diag"
1852 WRITE (unit=iounit, fmt="(T2,A,T74,A)") "ADMM| Density matrix purification method", "McWeeny"
1854 WRITE (unit=iounit, fmt="(T2,A,T73,A)") "ADMM| Density matrix purification method", "NONE(DM)"
1855 CASE DEFAULT
1856 cpabort("admm_purification_method")
1857 END SELECT
1858
1859 SELECT CASE (admm_control%method)
1861 WRITE (unit=iounit, fmt="(T2,A)") "ADMM| Orbital projection on ADMM basis"
1863 WRITE (unit=iounit, fmt="(T2,A)") "ADMM| Blocked Fock matrix projection with full purification"
1865 WRITE (unit=iounit, fmt="(T2,A)") "ADMM| Blocked Fock matrix projection"
1867 WRITE (unit=iounit, fmt="(T2,A)") "ADMM| Orbital projection with charge constrain"
1868 CASE DEFAULT
1869 cpabort("admm method")
1870 END SELECT
1871
1872 SELECT CASE (admm_control%scaling_model)
1875 WRITE (unit=iounit, fmt="(T2,A)") "ADMM| Use Merlot (2014) scaling model"
1876 CASE DEFAULT
1877 cpabort("admm scaling_model")
1878 END SELECT
1879
1880 WRITE (unit=iounit, fmt="(T2,A,T61,G20.10)") "ADMM| eps_filter", admm_control%eps_filter
1881
1882 SELECT CASE (admm_control%aux_exch_func)
1884 WRITE (unit=iounit, fmt="(T2,A)") "ADMM| No exchange functional correction term used"
1886 WRITE (unit=iounit, fmt="(T2,A,T74,A)") "ADMM| Exchange functional in correction term", "(W)PBEX"
1888 WRITE (unit=iounit, fmt="(T2,A,T77,A)") "ADMM| Exchange functional in correction term", "PBEX"
1890 WRITE (unit=iounit, fmt="(T2,A,T77,A)") "ADMM| Exchange functional in correction term", "OPTX"
1892 WRITE (unit=iounit, fmt="(T2,A,T74,A)") "ADMM| Exchange functional in correction term", "Becke88"
1894 WRITE (unit=iounit, fmt="(T2,A,T74,A)") "ADMM| Exchange functional in correction term", "SlaterX"
1895 CASE DEFAULT
1896 cpabort("admm aux_exch_func")
1897 END SELECT
1898
1899 WRITE (unit=iounit, fmt="(A)") ""
1900
1901 END IF
1902
1903 CALL cp_print_key_finished_output(iounit, logger, dft_section, &
1904 "PRINT%DFT_CONTROL_PARAMETERS")
1905 END SUBROUTINE write_admm_control
1906
1907! **************************************************************************************************
1908!> \brief Write the xTB control parameters to the output unit.
1909!> \param xtb_control ...
1910!> \param dft_section ...
1911! **************************************************************************************************
1912 SUBROUTINE write_xtb_control(xtb_control, dft_section)
1913 TYPE(xtb_control_type), POINTER :: xtb_control
1914 TYPE(section_vals_type), POINTER :: dft_section
1915
1916 CHARACTER(len=*), PARAMETER :: routinen = 'write_xtb_control'
1917
1918 INTEGER :: handle, output_unit
1919 TYPE(cp_logger_type), POINTER :: logger
1920
1921 CALL timeset(routinen, handle)
1922 NULLIFY (logger)
1923 logger => cp_get_default_logger()
1924
1925 output_unit = cp_print_key_unit_nr(logger, dft_section, &
1926 "PRINT%DFT_CONTROL_PARAMETERS", extension=".Log")
1927
1928 IF (output_unit > 0) THEN
1929
1930 WRITE (unit=output_unit, fmt="(/,T2,A,T31,A50)") &
1931 "xTB| Parameter file", adjustr(trim(xtb_control%parameter_file_name))
1932 WRITE (unit=output_unit, fmt="(T2,A,T71,I10)") &
1933 "xTB| Basis expansion STO-NG", xtb_control%sto_ng
1934 WRITE (unit=output_unit, fmt="(T2,A,T71,I10)") &
1935 "xTB| Basis expansion STO-NG for Hydrogen", xtb_control%h_sto_ng
1936 WRITE (unit=output_unit, fmt="(T2,A,T71,L10)") &
1937 "xTB| Halogen interaction potential", xtb_control%xb_interaction
1938 WRITE (unit=output_unit, fmt="(T2,A,T71,F10.3)") &
1939 "xTB| Halogen interaction potential cutoff radius", xtb_control%xb_radius
1940 WRITE (unit=output_unit, fmt="(T2,A,T71,L10)") &
1941 "xTB| Nonbonded interactions", xtb_control%do_nonbonded
1942 WRITE (unit=output_unit, fmt="(T2,A,T31,A50)") &
1943 "xTB| D3 Dispersion: Parameter file", adjustr(trim(xtb_control%dispersion_parameter_file))
1944 WRITE (unit=output_unit, fmt="(T2,A,T51,3F10.3)") &
1945 "xTB| Huckel constants ks kp kd", xtb_control%ks, xtb_control%kp, xtb_control%kd
1946 WRITE (unit=output_unit, fmt="(T2,A,T61,2F10.3)") &
1947 "xTB| Huckel constants ksp k2sh", xtb_control%ksp, xtb_control%k2sh
1948 WRITE (unit=output_unit, fmt="(T2,A,T71,F10.3)") &
1949 "xTB| Mataga-Nishimoto exponent", xtb_control%kg
1950 WRITE (unit=output_unit, fmt="(T2,A,T71,F10.3)") &
1951 "xTB| Repulsion potential exponent", xtb_control%kf
1952 WRITE (unit=output_unit, fmt="(T2,A,T51,3F10.3)") &
1953 "xTB| Coordination number scaling kcn(s) kcn(p) kcn(d)", &
1954 xtb_control%kcns, xtb_control%kcnp, xtb_control%kcnd
1955 WRITE (unit=output_unit, fmt="(T2,A,T71,F10.3)") &
1956 "xTB| Electronegativity scaling", xtb_control%ken
1957 WRITE (unit=output_unit, fmt="(T2,A,T61,2F10.3)") &
1958 "xTB| Halogen potential scaling kxr kx2", xtb_control%kxr, xtb_control%kx2
1959 WRITE (unit=output_unit, fmt="(/)")
1960
1961 END IF
1962
1963 CALL cp_print_key_finished_output(output_unit, logger, dft_section, &
1964 "PRINT%DFT_CONTROL_PARAMETERS")
1965
1966 CALL timestop(handle)
1967
1968 END SUBROUTINE write_xtb_control
1969
1970! **************************************************************************************************
1971!> \brief Purpose: Write the QS control parameters to the output unit.
1972!> \param qs_control ...
1973!> \param dft_section ...
1974! **************************************************************************************************
1975 SUBROUTINE write_qs_control(qs_control, dft_section)
1976 TYPE(qs_control_type), INTENT(IN) :: qs_control
1977 TYPE(section_vals_type), POINTER :: dft_section
1978
1979 CHARACTER(len=*), PARAMETER :: routinen = 'write_qs_control'
1980
1981 CHARACTER(len=20) :: method, quadrature
1982 INTEGER :: handle, i, igrid_level, ngrid_level, &
1983 output_unit
1984 TYPE(cp_logger_type), POINTER :: logger
1985 TYPE(ddapc_restraint_type), POINTER :: ddapc_restraint_control
1986 TYPE(enumeration_type), POINTER :: enum
1987 TYPE(keyword_type), POINTER :: keyword
1988 TYPE(section_type), POINTER :: qs_section
1989 TYPE(section_vals_type), POINTER :: print_section_vals, qs_section_vals
1990
1991 IF (qs_control%semi_empirical) RETURN
1992 IF (qs_control%dftb) RETURN
1993 IF (qs_control%xtb) RETURN
1994 CALL timeset(routinen, handle)
1995 NULLIFY (logger, print_section_vals, qs_section, qs_section_vals)
1996 logger => cp_get_default_logger()
1997 print_section_vals => section_vals_get_subs_vals(dft_section, "PRINT")
1998 qs_section_vals => section_vals_get_subs_vals(dft_section, "QS")
1999 CALL section_vals_get(qs_section_vals, section=qs_section)
2000
2001 NULLIFY (enum, keyword)
2002 keyword => section_get_keyword(qs_section, "METHOD")
2003 CALL keyword_get(keyword, enum=enum)
2004 method = enum_i2c(enum, qs_control%method_id)
2005
2006 NULLIFY (enum, keyword)
2007 keyword => section_get_keyword(qs_section, "QUADRATURE")
2008 CALL keyword_get(keyword, enum=enum)
2009 quadrature = enum_i2c(enum, qs_control%gapw_control%quadrature)
2010
2011 output_unit = cp_print_key_unit_nr(logger, print_section_vals, &
2012 "DFT_CONTROL_PARAMETERS", extension=".Log")
2013 IF (output_unit > 0) THEN
2014 ngrid_level = SIZE(qs_control%e_cutoff)
2015 WRITE (unit=output_unit, fmt="(/,T2,A,T61,A20)") &
2016 "QS| Method:", adjustr(method)
2017 IF (qs_control%pw_grid_opt%spherical) THEN
2018 WRITE (unit=output_unit, fmt="(T2,A,T61,A)") &
2019 "QS| Density plane wave grid type", " SPHERICAL HALFSPACE"
2020 ELSE IF (qs_control%pw_grid_opt%fullspace) THEN
2021 WRITE (unit=output_unit, fmt="(T2,A,T57,A)") &
2022 "QS| Density plane wave grid type", " NON-SPHERICAL FULLSPACE"
2023 ELSE
2024 WRITE (unit=output_unit, fmt="(T2,A,T57,A)") &
2025 "QS| Density plane wave grid type", " NON-SPHERICAL HALFSPACE"
2026 END IF
2027 WRITE (unit=output_unit, fmt="(T2,A,T71,I10)") &
2028 "QS| Number of grid levels:", SIZE(qs_control%e_cutoff)
2029 IF (ngrid_level == 1) THEN
2030 WRITE (unit=output_unit, fmt="(T2,A,T71,F10.1)") &
2031 "QS| Density cutoff [a.u.]:", qs_control%e_cutoff(1)
2032 ELSE
2033 WRITE (unit=output_unit, fmt="(T2,A,T71,F10.1)") &
2034 "QS| Density cutoff [a.u.]:", qs_control%cutoff
2035 IF (qs_control%commensurate_mgrids) &
2036 WRITE (unit=output_unit, fmt="(T2,A)") "QS| Using commensurate multigrids"
2037 WRITE (unit=output_unit, fmt="(T2,A,T71,F10.1)") &
2038 "QS| Multi grid cutoff [a.u.]: 1) grid level", qs_control%e_cutoff(1)
2039 WRITE (unit=output_unit, fmt="(T2,A,I3,A,T71,F10.1)") &
2040 ("QS| ", igrid_level, ") grid level", &
2041 qs_control%e_cutoff(igrid_level), &
2042 igrid_level=2, SIZE(qs_control%e_cutoff))
2043 END IF
2044 IF (qs_control%pao) THEN
2045 WRITE (unit=output_unit, fmt="(T2,A)") "QS| PAO active"
2046 END IF
2047 WRITE (unit=output_unit, fmt="(T2,A,T71,F10.1)") &
2048 "QS| Grid level progression factor:", qs_control%progression_factor
2049 WRITE (unit=output_unit, fmt="(T2,A,T71,F10.1)") &
2050 "QS| Relative density cutoff [a.u.]:", qs_control%relative_cutoff
2051 WRITE (unit=output_unit, fmt="(T2,A,T73,ES8.1)") &
2052 "QS| Interaction thresholds: eps_pgf_orb:", &
2053 qs_control%eps_pgf_orb, &
2054 "QS| eps_filter_matrix:", &
2055 qs_control%eps_filter_matrix, &
2056 "QS| eps_core_charge:", &
2057 qs_control%eps_core_charge, &
2058 "QS| eps_rho_gspace:", &
2059 qs_control%eps_rho_gspace, &
2060 "QS| eps_rho_rspace:", &
2061 qs_control%eps_rho_rspace, &
2062 "QS| eps_gvg_rspace:", &
2063 qs_control%eps_gvg_rspace, &
2064 "QS| eps_ppl:", &
2065 qs_control%eps_ppl, &
2066 "QS| eps_ppnl:", &
2067 qs_control%eps_ppnl
2068 IF (qs_control%gapw) THEN
2069 SELECT CASE (qs_control%gapw_control%basis_1c)
2070 CASE (gapw_1c_orb)
2071 WRITE (unit=output_unit, fmt="(T2,A)") &
2072 "QS| GAPW| One center basis from orbital basis primitives"
2073 CASE (gapw_1c_small)
2074 WRITE (unit=output_unit, fmt="(T2,A)") &
2075 "QS| GAPW| One center basis extended with primitives (small:s)"
2076 CASE (gapw_1c_medium)
2077 WRITE (unit=output_unit, fmt="(T2,A)") &
2078 "QS| GAPW| One center basis extended with primitives (medium:sp)"
2079 CASE (gapw_1c_large)
2080 WRITE (unit=output_unit, fmt="(T2,A)") &
2081 "QS| GAPW| One center basis extended with primitives (large:spd)"
2082 CASE (gapw_1c_very_large)
2083 WRITE (unit=output_unit, fmt="(T2,A)") &
2084 "QS| GAPW| One center basis extended with primitives (very large:spdf)"
2085 CASE DEFAULT
2086 cpabort("basis_1c incorrect")
2087 END SELECT
2088 WRITE (unit=output_unit, fmt="(T2,A,T73,ES8.1)") &
2089 "QS| GAPW| eps_fit:", &
2090 qs_control%gapw_control%eps_fit, &
2091 "QS| GAPW| eps_iso:", &
2092 qs_control%gapw_control%eps_iso, &
2093 "QS| GAPW| eps_svd:", &
2094 qs_control%gapw_control%eps_svd, &
2095 "QS| GAPW| eps_cpc:", &
2096 qs_control%gapw_control%eps_cpc
2097 WRITE (unit=output_unit, fmt="(T2,A,T61,A20)") &
2098 "QS| GAPW| atom-r-grid: quadrature:", &
2099 adjustr(quadrature)
2100 WRITE (unit=output_unit, fmt="(T2,A,T71,I10)") &
2101 "QS| GAPW| atom-s-grid: max l :", &
2102 qs_control%gapw_control%lmax_sphere, &
2103 "QS| GAPW| max_l_rho0 :", &
2104 qs_control%gapw_control%lmax_rho0
2105 IF (qs_control%gapw_control%non_paw_atoms) THEN
2106 WRITE (unit=output_unit, fmt="(T2,A)") &
2107 "QS| GAPW| At least one kind is NOT PAW, i.e. it has only soft AO "
2108 END IF
2109 IF (qs_control%gapw_control%nopaw_as_gpw) THEN
2110 WRITE (unit=output_unit, fmt="(T2,A)") &
2111 "QS| GAPW| The NOT PAW atoms are treated fully GPW"
2112 END IF
2113 END IF
2114 IF (qs_control%gapw_xc) THEN
2115 SELECT CASE (qs_control%gapw_control%basis_1c)
2116 CASE (gapw_1c_orb)
2117 WRITE (unit=output_unit, fmt="(T2,A)") &
2118 "QS| GAPW_XC| One center basis from orbital basis primitives"
2119 CASE (gapw_1c_small)
2120 WRITE (unit=output_unit, fmt="(T2,A)") &
2121 "QS| GAPW_XC| One center basis extended with primitives (small:s)"
2122 CASE (gapw_1c_medium)
2123 WRITE (unit=output_unit, fmt="(T2,A)") &
2124 "QS| GAPW_XC| One center basis extended with primitives (medium:sp)"
2125 CASE (gapw_1c_large)
2126 WRITE (unit=output_unit, fmt="(T2,A)") &
2127 "QS| GAPW_XC| One center basis extended with primitives (large:spd)"
2128 CASE (gapw_1c_very_large)
2129 WRITE (unit=output_unit, fmt="(T2,A)") &
2130 "QS| GAPW_XC| One center basis extended with primitives (very large:spdf)"
2131 CASE DEFAULT
2132 cpabort("basis_1c incorrect")
2133 END SELECT
2134 WRITE (unit=output_unit, fmt="(T2,A,T73,ES8.1)") &
2135 "QS| GAPW_XC| eps_fit:", &
2136 qs_control%gapw_control%eps_fit, &
2137 "QS| GAPW_XC| eps_iso:", &
2138 qs_control%gapw_control%eps_iso, &
2139 "QS| GAPW_XC| eps_svd:", &
2140 qs_control%gapw_control%eps_svd
2141 WRITE (unit=output_unit, fmt="(T2,A,T55,A30)") &
2142 "QS| GAPW_XC|atom-r-grid: quadrature:", &
2143 enum_i2c(enum, qs_control%gapw_control%quadrature)
2144 WRITE (unit=output_unit, fmt="(T2,A,T71,I10)") &
2145 "QS| GAPW_XC| atom-s-grid: max l :", &
2146 qs_control%gapw_control%lmax_sphere
2147 END IF
2148 IF (qs_control%mulliken_restraint) THEN
2149 WRITE (unit=output_unit, fmt="(T2,A,T73,ES8.1)") &
2150 "QS| Mulliken restraint target", qs_control%mulliken_restraint_control%target
2151 WRITE (unit=output_unit, fmt="(T2,A,T73,ES8.1)") &
2152 "QS| Mulliken restraint strength", qs_control%mulliken_restraint_control%strength
2153 WRITE (unit=output_unit, fmt="(T2,A,T73,I8)") &
2154 "QS| Mulliken restraint atoms: ", qs_control%mulliken_restraint_control%natoms
2155 WRITE (unit=output_unit, fmt="(5I8)") qs_control%mulliken_restraint_control%atoms
2156 END IF
2157 IF (qs_control%ddapc_restraint) THEN
2158 DO i = 1, SIZE(qs_control%ddapc_restraint_control)
2159 ddapc_restraint_control => qs_control%ddapc_restraint_control(i)
2160 IF (SIZE(qs_control%ddapc_restraint_control) .GT. 1) &
2161 WRITE (unit=output_unit, fmt="(T2,A,T3,I8)") &
2162 "QS| parameters for DDAPC restraint number", i
2163 WRITE (unit=output_unit, fmt="(T2,A,T73,ES8.1)") &
2164 "QS| ddapc restraint target", ddapc_restraint_control%target
2165 WRITE (unit=output_unit, fmt="(T2,A,T73,ES8.1)") &
2166 "QS| ddapc restraint strength", ddapc_restraint_control%strength
2167 WRITE (unit=output_unit, fmt="(T2,A,T73,I8)") &
2168 "QS| ddapc restraint atoms: ", ddapc_restraint_control%natoms
2169 WRITE (unit=output_unit, fmt="(5I8)") ddapc_restraint_control%atoms
2170 WRITE (unit=output_unit, fmt="(T2,A)") "Coefficients:"
2171 WRITE (unit=output_unit, fmt="(5F6.2)") ddapc_restraint_control%coeff
2172 SELECT CASE (ddapc_restraint_control%functional_form)
2173 CASE (do_ddapc_restraint)
2174 WRITE (unit=output_unit, fmt="(T2,A,T61,A20)") &
2175 "QS| ddapc restraint functional form :", "RESTRAINT"
2176 CASE (do_ddapc_constraint)
2177 WRITE (unit=output_unit, fmt="(T2,A,T61,A20)") &
2178 "QS| ddapc restraint functional form :", "CONSTRAINT"
2179 CASE DEFAULT
2180 cpabort("Unknown ddapc restraint")
2181 END SELECT
2182 END DO
2183 END IF
2184 IF (qs_control%s2_restraint) THEN
2185 WRITE (unit=output_unit, fmt="(T2,A,T73,ES8.1)") &
2186 "QS| s2 restraint target", qs_control%s2_restraint_control%target
2187 WRITE (unit=output_unit, fmt="(T2,A,T73,ES8.1)") &
2188 "QS| s2 restraint strength", qs_control%s2_restraint_control%strength
2189 SELECT CASE (qs_control%s2_restraint_control%functional_form)
2190 CASE (do_s2_restraint)
2191 WRITE (unit=output_unit, fmt="(T2,A,T61,A20)") &
2192 "QS| s2 restraint functional form :", "RESTRAINT"
2193 cpabort("Not yet implemented")
2194 CASE (do_s2_constraint)
2195 WRITE (unit=output_unit, fmt="(T2,A,T61,A20)") &
2196 "QS| s2 restraint functional form :", "CONSTRAINT"
2197 CASE DEFAULT
2198 cpabort("Unknown ddapc restraint")
2199 END SELECT
2200 END IF
2201 END IF
2202 CALL cp_print_key_finished_output(output_unit, logger, print_section_vals, &
2203 "DFT_CONTROL_PARAMETERS")
2204
2205 CALL timestop(handle)
2206
2207 END SUBROUTINE write_qs_control
2208
2209! **************************************************************************************************
2210!> \brief reads the input parameters needed for ddapc.
2211!> \param qs_control ...
2212!> \param qs_section ...
2213!> \param ddapc_restraint_section ...
2214!> \author fschiff
2215!> \note
2216!> either reads DFT%QS%DDAPC_RESTRAINT or PROPERTIES%ET_coupling
2217!> if(qs_section is present the DFT part is read, if ddapc_restraint_section
2218!> is present ET_COUPLING is read. Avoid having both!!!
2219! **************************************************************************************************
2220 SUBROUTINE read_ddapc_section(qs_control, qs_section, ddapc_restraint_section)
2221
2222 TYPE(qs_control_type), INTENT(INOUT) :: qs_control
2223 TYPE(section_vals_type), OPTIONAL, POINTER :: qs_section, ddapc_restraint_section
2224
2225 INTEGER :: i, j, jj, k, n_rep
2226 INTEGER, DIMENSION(:), POINTER :: tmplist
2227 REAL(kind=dp), DIMENSION(:), POINTER :: rtmplist
2228 TYPE(ddapc_restraint_type), POINTER :: ddapc_restraint_control
2229 TYPE(section_vals_type), POINTER :: ddapc_section
2230
2231 IF (PRESENT(ddapc_restraint_section)) THEN
2232 IF (ASSOCIATED(qs_control%ddapc_restraint_control)) THEN
2233 IF (SIZE(qs_control%ddapc_restraint_control) .GE. 2) &
2234 cpabort("ET_COUPLING cannot be used in combination with a normal restraint")
2235 ELSE
2236 ddapc_section => ddapc_restraint_section
2237 ALLOCATE (qs_control%ddapc_restraint_control(1))
2238 END IF
2239 END IF
2240
2241 IF (PRESENT(qs_section)) THEN
2242 NULLIFY (ddapc_section)
2243 ddapc_section => section_vals_get_subs_vals(qs_section, &
2244 "DDAPC_RESTRAINT")
2245 END IF
2246
2247 DO i = 1, SIZE(qs_control%ddapc_restraint_control)
2248
2249 CALL ddapc_control_create(qs_control%ddapc_restraint_control(i))
2250 ddapc_restraint_control => qs_control%ddapc_restraint_control(i)
2251
2252 CALL section_vals_val_get(ddapc_section, "STRENGTH", i_rep_section=i, &
2253 r_val=ddapc_restraint_control%strength)
2254 CALL section_vals_val_get(ddapc_section, "TARGET", i_rep_section=i, &
2255 r_val=ddapc_restraint_control%target)
2256 CALL section_vals_val_get(ddapc_section, "FUNCTIONAL_FORM", i_rep_section=i, &
2257 i_val=ddapc_restraint_control%functional_form)
2258 CALL section_vals_val_get(ddapc_section, "ATOMS", i_rep_section=i, &
2259 n_rep_val=n_rep)
2260 CALL section_vals_val_get(ddapc_section, "TYPE_OF_DENSITY", i_rep_section=i, &
2261 i_val=ddapc_restraint_control%density_type)
2262
2263 jj = 0
2264 DO k = 1, n_rep
2265 CALL section_vals_val_get(ddapc_section, "ATOMS", i_rep_section=i, &
2266 i_rep_val=k, i_vals=tmplist)
2267 DO j = 1, SIZE(tmplist)
2268 jj = jj + 1
2269 END DO
2270 END DO
2271 IF (jj < 1) cpabort("Need at least 1 atom to use ddapc constraints")
2272 ddapc_restraint_control%natoms = jj
2273 IF (ASSOCIATED(ddapc_restraint_control%atoms)) &
2274 DEALLOCATE (ddapc_restraint_control%atoms)
2275 ALLOCATE (ddapc_restraint_control%atoms(ddapc_restraint_control%natoms))
2276 jj = 0
2277 DO k = 1, n_rep
2278 CALL section_vals_val_get(ddapc_section, "ATOMS", i_rep_section=i, &
2279 i_rep_val=k, i_vals=tmplist)
2280 DO j = 1, SIZE(tmplist)
2281 jj = jj + 1
2282 ddapc_restraint_control%atoms(jj) = tmplist(j)
2283 END DO
2284 END DO
2285
2286 IF (ASSOCIATED(ddapc_restraint_control%coeff)) &
2287 DEALLOCATE (ddapc_restraint_control%coeff)
2288 ALLOCATE (ddapc_restraint_control%coeff(ddapc_restraint_control%natoms))
2289 ddapc_restraint_control%coeff = 1.0_dp
2290
2291 CALL section_vals_val_get(ddapc_section, "COEFF", i_rep_section=i, &
2292 n_rep_val=n_rep)
2293 jj = 0
2294 DO k = 1, n_rep
2295 CALL section_vals_val_get(ddapc_section, "COEFF", i_rep_section=i, &
2296 i_rep_val=k, r_vals=rtmplist)
2297 DO j = 1, SIZE(rtmplist)
2298 jj = jj + 1
2299 IF (jj > ddapc_restraint_control%natoms) &
2300 cpabort("Need the same number of coeff as there are atoms ")
2301 ddapc_restraint_control%coeff(jj) = rtmplist(j)
2302 END DO
2303 END DO
2304 IF (jj < ddapc_restraint_control%natoms .AND. jj .NE. 0) &
2305 cpabort("Need no or the same number of coeff as there are atoms.")
2306 END DO
2307 k = 0
2308 DO i = 1, SIZE(qs_control%ddapc_restraint_control)
2309 IF (qs_control%ddapc_restraint_control(i)%functional_form == &
2310 do_ddapc_constraint) k = k + 1
2311 END DO
2312 IF (k == 2) CALL cp_abort(__location__, &
2313 "Only a single constraint possible yet, try to use restraints instead ")
2314
2315 END SUBROUTINE read_ddapc_section
2316
2317! **************************************************************************************************
2318!> \brief ...
2319!> \param dft_control ...
2320!> \param efield_section ...
2321! **************************************************************************************************
2322 SUBROUTINE read_efield_sections(dft_control, efield_section)
2323 TYPE(dft_control_type), POINTER :: dft_control
2324 TYPE(section_vals_type), POINTER :: efield_section
2325
2326 CHARACTER(len=default_path_length) :: file_name
2327 INTEGER :: i, io, j, n, unit_nr
2328 REAL(kind=dp), DIMENSION(:), POINTER :: tmp_vals
2329 TYPE(efield_type), POINTER :: efield
2330 TYPE(section_vals_type), POINTER :: tmp_section
2331
2332 DO i = 1, SIZE(dft_control%efield_fields)
2333 NULLIFY (dft_control%efield_fields(i)%efield)
2334 ALLOCATE (dft_control%efield_fields(i)%efield)
2335 efield => dft_control%efield_fields(i)%efield
2336 NULLIFY (efield%envelop_i_vars, efield%envelop_r_vars)
2337 CALL section_vals_val_get(efield_section, "INTENSITY", i_rep_section=i, &
2338 r_val=efield%strength)
2339
2340 CALL section_vals_val_get(efield_section, "POLARISATION", i_rep_section=i, &
2341 r_vals=tmp_vals)
2342 ALLOCATE (efield%polarisation(SIZE(tmp_vals)))
2343 efield%polarisation = tmp_vals
2344 CALL section_vals_val_get(efield_section, "PHASE", i_rep_section=i, &
2345 r_val=efield%phase_offset)
2346 CALL section_vals_val_get(efield_section, "ENVELOP", i_rep_section=i, &
2347 i_val=efield%envelop_id)
2348 CALL section_vals_val_get(efield_section, "WAVELENGTH", i_rep_section=i, &
2349 r_val=efield%wavelength)
2350 CALL section_vals_val_get(efield_section, "VEC_POT_INITIAL", i_rep_section=i, &
2351 r_vals=tmp_vals)
2352 efield%vec_pot_initial = tmp_vals
2353
2354 IF (efield%envelop_id == constant_env) THEN
2355 ALLOCATE (efield%envelop_i_vars(2))
2356 tmp_section => section_vals_get_subs_vals(efield_section, "CONSTANT_ENV", i_rep_section=i)
2357 CALL section_vals_val_get(tmp_section, "START_STEP", &
2358 i_val=efield%envelop_i_vars(1))
2359 CALL section_vals_val_get(tmp_section, "END_STEP", &
2360 i_val=efield%envelop_i_vars(2))
2361 ELSE IF (efield%envelop_id == gaussian_env) THEN
2362 ALLOCATE (efield%envelop_r_vars(2))
2363 tmp_section => section_vals_get_subs_vals(efield_section, "GAUSSIAN_ENV", i_rep_section=i)
2364 CALL section_vals_val_get(tmp_section, "T0", &
2365 r_val=efield%envelop_r_vars(1))
2366 CALL section_vals_val_get(tmp_section, "SIGMA", &
2367 r_val=efield%envelop_r_vars(2))
2368 ELSE IF (efield%envelop_id == ramp_env) THEN
2369 ALLOCATE (efield%envelop_i_vars(4))
2370 tmp_section => section_vals_get_subs_vals(efield_section, "RAMP_ENV", i_rep_section=i)
2371 CALL section_vals_val_get(tmp_section, "START_STEP_IN", &
2372 i_val=efield%envelop_i_vars(1))
2373 CALL section_vals_val_get(tmp_section, "END_STEP_IN", &
2374 i_val=efield%envelop_i_vars(2))
2375 CALL section_vals_val_get(tmp_section, "START_STEP_OUT", &
2376 i_val=efield%envelop_i_vars(3))
2377 CALL section_vals_val_get(tmp_section, "END_STEP_OUT", &
2378 i_val=efield%envelop_i_vars(4))
2379 ELSE IF (efield%envelop_id == custom_env) THEN
2380 tmp_section => section_vals_get_subs_vals(efield_section, "CUSTOM_ENV", i_rep_section=i)
2381 CALL section_vals_val_get(tmp_section, "EFIELD_FILE_NAME", c_val=file_name)
2382 CALL open_file(file_name=trim(file_name), file_action="READ", file_status="OLD", unit_number=unit_nr)
2383 !Determine the number of lines in file
2384 n = 0
2385 DO WHILE (.true.)
2386 READ (unit_nr, *, iostat=io)
2387 IF (io /= 0) EXIT
2388 n = n + 1
2389 END DO
2390 rewind(unit_nr)
2391 ALLOCATE (efield%envelop_r_vars(n + 1))
2392 !Store the timestep of the list in the first entry of the r_vars
2393 CALL section_vals_val_get(tmp_section, "TIMESTEP", r_val=efield%envelop_r_vars(1))
2394 !Read the file
2395 DO j = 2, n + 1
2396 READ (unit_nr, *) efield%envelop_r_vars(j)
2397 efield%envelop_r_vars(j) = cp_unit_to_cp2k(efield%envelop_r_vars(j), "volt/m")
2398 END DO
2399 CALL close_file(unit_nr)
2400 END IF
2401 END DO
2402 END SUBROUTINE read_efield_sections
2403
2404! **************************************************************************************************
2405!> \brief reads the input parameters needed real time propagation
2406!> \param dft_control ...
2407!> \param rtp_section ...
2408!> \author fschiff
2409! **************************************************************************************************
2410 SUBROUTINE read_rtp_section(dft_control, rtp_section)
2411
2412 TYPE(dft_control_type), INTENT(INOUT) :: dft_control
2413 TYPE(section_vals_type), POINTER :: rtp_section
2414
2415 INTEGER, DIMENSION(:), POINTER :: tmp
2416 LOGICAL :: is_present
2417 TYPE(section_vals_type), POINTER :: proj_mo_section
2418
2419 ALLOCATE (dft_control%rtp_control)
2420 CALL section_vals_val_get(rtp_section, "MAX_ITER", &
2421 i_val=dft_control%rtp_control%max_iter)
2422 CALL section_vals_val_get(rtp_section, "MAT_EXP", &
2423 i_val=dft_control%rtp_control%mat_exp)
2424 CALL section_vals_val_get(rtp_section, "ASPC_ORDER", &
2425 i_val=dft_control%rtp_control%aspc_order)
2426 CALL section_vals_val_get(rtp_section, "EXP_ACCURACY", &
2427 r_val=dft_control%rtp_control%eps_exp)
2428 CALL section_vals_val_get(rtp_section, "PROPAGATOR", &
2429 i_val=dft_control%rtp_control%propagator)
2430 CALL section_vals_val_get(rtp_section, "EPS_ITER", &
2431 r_val=dft_control%rtp_control%eps_ener)
2432 CALL section_vals_val_get(rtp_section, "INITIAL_WFN", &
2433 i_val=dft_control%rtp_control%initial_wfn)
2434 CALL section_vals_val_get(rtp_section, "HFX_BALANCE_IN_CORE", &
2435 l_val=dft_control%rtp_control%hfx_redistribute)
2436 CALL section_vals_val_get(rtp_section, "APPLY_WFN_MIX_INIT_RESTART", &
2437 l_val=dft_control%rtp_control%apply_wfn_mix_init_restart)
2438 CALL section_vals_val_get(rtp_section, "APPLY_DELTA_PULSE", &
2439 l_val=dft_control%rtp_control%apply_delta_pulse)
2440 CALL section_vals_val_get(rtp_section, "APPLY_DELTA_PULSE_MAG", &
2441 l_val=dft_control%rtp_control%apply_delta_pulse_mag)
2442 CALL section_vals_val_get(rtp_section, "VELOCITY_GAUGE", &
2443 l_val=dft_control%rtp_control%velocity_gauge)
2444 CALL section_vals_val_get(rtp_section, "VG_COM_NL", &
2445 l_val=dft_control%rtp_control%nl_gauge_transform)
2446 CALL section_vals_val_get(rtp_section, "PERIODIC", &
2447 l_val=dft_control%rtp_control%periodic)
2448 CALL section_vals_val_get(rtp_section, "DENSITY_PROPAGATION", &
2449 l_val=dft_control%rtp_control%linear_scaling)
2450 CALL section_vals_val_get(rtp_section, "MCWEENY_MAX_ITER", &
2451 i_val=dft_control%rtp_control%mcweeny_max_iter)
2452 CALL section_vals_val_get(rtp_section, "ACCURACY_REFINEMENT", &
2453 i_val=dft_control%rtp_control%acc_ref)
2454 CALL section_vals_val_get(rtp_section, "MCWEENY_EPS", &
2455 r_val=dft_control%rtp_control%mcweeny_eps)
2456 CALL section_vals_val_get(rtp_section, "DELTA_PULSE_SCALE", &
2457 r_val=dft_control%rtp_control%delta_pulse_scale)
2458 CALL section_vals_val_get(rtp_section, "DELTA_PULSE_DIRECTION", &
2459 i_vals=tmp)
2460 dft_control%rtp_control%delta_pulse_direction = tmp
2461 CALL section_vals_val_get(rtp_section, "SC_CHECK_START", &
2462 i_val=dft_control%rtp_control%sc_check_start)
2463 proj_mo_section => section_vals_get_subs_vals(rtp_section, "PRINT%PROJECTION_MO")
2464 CALL section_vals_get(proj_mo_section, explicit=is_present)
2465 IF (is_present) THEN
2466 IF (dft_control%rtp_control%linear_scaling) &
2467 CALL cp_abort(__location__, &
2468 "You have defined a time dependent projection of mos, but "// &
2469 "only the density matrix is propagated (DENSITY_PROPAGATION "// &
2470 ".TRUE.). Please either use MO-based real time DFT or do not "// &
2471 "define any PRINT%PROJECTION_MO section")
2472 dft_control%rtp_control%is_proj_mo = .true.
2473 ELSE
2474 dft_control%rtp_control%is_proj_mo = .false.
2475 END IF
2476
2477 END SUBROUTINE read_rtp_section
2478
2479! **************************************************************************************************
2480!> \brief Parses the BLOCK_LIST keywords from the ADMM section
2481!> \param admm_control ...
2482!> \param dft_section ...
2483! **************************************************************************************************
2484 SUBROUTINE read_admm_block_list(admm_control, dft_section)
2485 TYPE(admm_control_type), POINTER :: admm_control
2486 TYPE(section_vals_type), POINTER :: dft_section
2487
2488 INTEGER :: irep, list_size, n_rep
2489 INTEGER, DIMENSION(:), POINTER :: tmplist
2490
2491 NULLIFY (tmplist)
2492
2493 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%BLOCK_LIST", &
2494 n_rep_val=n_rep)
2495
2496 ALLOCATE (admm_control%blocks(n_rep))
2497
2498 DO irep = 1, n_rep
2499 CALL section_vals_val_get(dft_section, "AUXILIARY_DENSITY_MATRIX_METHOD%BLOCK_LIST", &
2500 i_rep_val=irep, i_vals=tmplist)
2501 list_size = SIZE(tmplist)
2502 ALLOCATE (admm_control%blocks(irep)%list(list_size))
2503 admm_control%blocks(irep)%list(:) = tmplist(:)
2504 END DO
2505
2506 END SUBROUTINE read_admm_block_list
2507
2508END MODULE cp_control_utils
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public vandevondele2005b
integer, save, public umari2002
integer, save, public yin2017
integer, save, public stewart2007
integer, save, public stengel2009
integer, save, public vandevondele2005a
integer, save, public grimme2017
integer, save, public lippert1999
integer, save, public dewar1977
integer, save, public elstner1998
integer, save, public vanvoorhis2015
integer, save, public repasky2002
integer, save, public hu2007
integer, save, public andreussi2012
integer, save, public rocha2006
integer, save, public lippert1997
integer, save, public fattebert2002
integer, save, public thiel1992
integer, save, public porezag1995
integer, save, public souza2002
integer, save, public schenter2008
integer, save, public krack2000
integer, save, public dewar1985
integer, save, public stewart1989
integer, save, public seifert1996
integer, save, public zhechkov2005
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dft_control_create(dft_control)
allocates and perform a very basic initialization
subroutine, public expot_control_create(expot_control)
...
subroutine, public maxwell_control_create(maxwell_control)
...
subroutine, public ddapc_control_create(ddapc_restraint_control)
create the ddapc_restraint_type
subroutine, public admm_control_create(admm_control)
...
subroutine, public tddfpt_control_create(tddfpt_control)
...
Utilities to set up the control types.
subroutine, public write_qs_control(qs_control, dft_section)
Purpose: Write the QS control parameters to the output unit.
subroutine, public read_tddfpt2_control(t_control, t_section, qs_control)
Read TDDFPT-related input parameters.
subroutine, public read_dft_control(dft_control, dft_section)
...
subroutine, public read_qs_section(qs_control, qs_section)
...
subroutine, public write_admm_control(admm_control, dft_section)
Write the ADMM control parameters to the output unit.
subroutine, public read_tddfpt_control(t_control, dft_section)
...
subroutine, public write_dft_control(dft_control, dft_section)
Write the DFT control parameters to the output unit.
subroutine, public read_ddapc_section(qs_control, qs_section, ddapc_restraint_section)
reads the input parameters needed for ddapc.
subroutine, public read_mgrid_section(qs_control, dft_section)
...
Utility routines to open and close files. Tracking of preconnections.
Definition cp_files.F:16
subroutine, public open_file(file_name, file_status, file_form, file_action, file_position, file_pad, unit_number, debug, skip_get_unit_number, file_access)
Opens the requested file using a free unit number.
Definition cp_files.F:308
subroutine, public close_file(unit_number, file_status, keep_preconnection)
Close an open file given by its logical unit number. Optionally, keep the file and unit preconnected.
Definition cp_files.F:119
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
unit conversion facility
Definition cp_units.F:30
real(kind=dp) function, public cp_unit_from_cp2k(value, unit_str, defaults, power)
converts from the internal cp2k units to the given unit
Definition cp_units.F:1179
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
Definition cp_units.F:1150
subroutine, public read_gp_section(nonbonded, section, start)
Reads the GENPOT - generic potential section.
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public sic_list_unpaired
integer, parameter, public sic_mauri_spz
integer, parameter, public do_method_ofgpw
integer, parameter, public do_admm_purify_mo_no_diag
integer, parameter, public do_se_lr_ewald_gks
integer, parameter, public do_admm_aux_exch_func_opt_libxc
integer, parameter, public do_s2_restraint
integer, parameter, public do_admm_purify_none
integer, parameter, public do_method_rigpw
integer, parameter, public do_s2_constraint
integer, parameter, public use_mom_ref_user
integer, parameter, public do_method_gpw
integer, parameter, public gapw_1c_large
integer, parameter, public do_method_pdg
integer, parameter, public do_admm_purify_none_dm
integer, parameter, public do_method_pnnl
integer, parameter, public do_ddapc_constraint
integer, parameter, public do_se_lr_none
integer, parameter, public do_admm_purify_mcweeny
integer, parameter, public do_se_lr_ewald
integer, parameter, public do_admm_blocking_purify_full
integer, parameter, public ramp_env
integer, parameter, public do_se_is_kdso_d
integer, parameter, public gapw_1c_medium
integer, parameter, public do_admm_aux_exch_func_sx_libxc
integer, parameter, public admm2_type
integer, parameter, public sic_list_all
integer, parameter, public constant_env
integer, parameter, public tddfpt_excitations
integer, parameter, public sic_eo
integer, parameter, public sccs_derivative_cd5
integer, parameter, public do_admm_aux_exch_func_bee
integer, parameter, public no_admm_type
integer, parameter, public do_admm_blocked_projection
integer, parameter, public do_admm_basis_projection
integer, parameter, public do_method_rm1
integer, parameter, public do_admm_aux_exch_func_default_libxc
integer, parameter, public do_admm_aux_exch_func_opt
integer, parameter, public gapw_1c_small
integer, parameter, public do_admm_aux_exch_func_none
integer, parameter, public do_admm_purify_cauchy_subspace
integer, parameter, public do_method_pm3
integer, parameter, public do_admm_aux_exch_func_bee_libxc
integer, parameter, public admm1_type
integer, parameter, public do_admm_aux_exch_func_pbex_libxc
integer, parameter, public do_method_mndo
integer, parameter, public gapw_1c_orb
integer, parameter, public do_admm_aux_exch_func_default
integer, parameter, public do_pwgrid_ns_fullspace
integer, parameter, public gapw_1c_very_large
integer, parameter, public do_method_gapw
integer, parameter, public admms_type
integer, parameter, public do_admm_charge_constrained_projection
integer, parameter, public do_admm_purify_cauchy
integer, parameter, public sccs_fattebert_gygi
integer, parameter, public gaussian_env
integer, parameter, public sccs_derivative_cd7
integer, parameter, public do_method_mndod
integer, parameter, public do_method_am1
integer, parameter, public do_method_dftb
integer, parameter, public tddfpt_dipole_length
integer, parameter, public sccs_derivative_fft
integer, parameter, public tddfpt_kernel_stda
integer, parameter, public do_pwgrid_spherical
integer, parameter, public do_se_lr_ewald_r3
integer, parameter, public do_se_is_kdso
integer, parameter, public do_admm_purify_mo_diag
integer, parameter, public do_method_lrigpw
integer, parameter, public sic_mauri_us
integer, parameter, public sic_none
integer, parameter, public do_se_is_slater
integer, parameter, public custom_env
integer, parameter, public do_method_xtb
integer, parameter, public do_ddapc_restraint
integer, parameter, public do_pwgrid_ns_halfspace
integer, parameter, public sccs_derivative_cd3
integer, parameter, public do_method_pm6fm
integer, parameter, public admmq_type
integer, parameter, public sccs_andreussi
integer, parameter, public sic_ad
integer, parameter, public do_admm_exch_scaling_none
integer, parameter, public admmp_type
integer, parameter, public do_method_gapw_xc
integer, parameter, public do_admm_exch_scaling_merlot
integer, parameter, public real_time_propagation
integer, parameter, public numerical
integer, parameter, public do_method_pm6
integer, parameter, public do_admm_aux_exch_func_pbex
integer, parameter, public slater
checks the input and perform some automatic "magic" on it
subroutine, public xc_functionals_expand(functionals, xc_section)
expand a shortcutted functional section
function that build the dft section of the input
subroutine, public create_dft_section(section)
creates the dft section
represents an enumeration, i.e. a mapping between integers and strings
character(len=default_string_length) function, public enum_i2c(enum, i)
maps an integer to a string
represents keywords in an input
subroutine, public keyword_get(keyword, names, usage, description, type_of_var, n_var, default_value, lone_keyword_value, repeats, enum, citations)
...
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_vals_val_set(section_vals, keyword_name, i_rep_section, i_rep_val, val, l_val, i_val, r_val, c_val, l_vals_ptr, i_vals_ptr, r_vals_ptr, c_vals_ptr)
sets the requested value
integer function, public section_get_ival(section_vals, keyword_name)
...
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
recursive subroutine, public section_release(section)
releases the given keyword list (see doc/ReferenceCounting.html)
recursive type(keyword_type) function, pointer, public section_get_keyword(section, keyword_name)
returns the requested keyword
subroutine, public section_vals_get(section_vals, ref_count, n_repetition, n_subs_vals_rep, section, explicit)
returns various attributes about the section_vals
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
integer, parameter, public default_path_length
Definition kinds.F:58
Definition of mathematical constants and functions.
real(kind=dp), parameter, public fourpi
subroutine, public pair_potential_reallocate(p, lb1_new, ub1_new, lj, lj_charmm, williams, goodwin, eam, quip, nequip, allegro, bmhft, bmhftd, ipbv, buck4r, buckmo, gp, tersoff, siepmann, gal, gal21, tab, deepmd)
Cleans the potential parameter type.
Periodic Table related data definitions.
subroutine, public get_ptable_info(symbol, number, amass, ielement, covalent_radius, metallic_radius, vdw_radius, found)
Pass information about the kind given the element symbol.
Utility subroutines for CDFT calculations.
subroutine, public read_cdft_control_section(qs_control, cdft_control_section)
reads the input parameters needed for CDFT with OT
Utilities for string manipulations.
elemental subroutine, public uppercase(string)
Convert all lower case characters in a string to upper case.
All kind of helpful little routines.
Definition util.F:14
input constants for xc
integer, parameter, public xc_deriv_collocate
Writes information on XC functionals to output.
subroutine, public xc_write(iounit, xc_section, lsd)
...
Exchange and Correlation functional calculations.
Definition xc.F:17
logical function, public xc_uses_norm_drho(xc_fun_section, lsd)
...
Definition xc.F:111
logical function, public xc_uses_kinetic_energy_density(xc_fun_section, lsd)
...
Definition xc.F:91
type of a logger, at the moment it contains just a print level starting at which level it should be l...
represent a keyword in the input
represent a section of the input file