(git:cccd2f3)
Loading...
Searching...
No Matches
qs_environment.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \par History
10!> - Merged with the Quickstep MODULE method_specification (17.01.2002,MK)
11!> - USE statements cleaned, added
12!> (25.09.2002,MK)
13!> - Added more LSD structure (01.2003,Joost VandeVondele)
14!> - New molecule data types introduced (Sep. 2003,MK)
15!> - Cleaning; getting rid of pnode (02.10.2003,MK)
16!> - Sub-system setup added (08.10.2003,MK)
17!> \author MK (18.05.2000)
18! **************************************************************************************************
30 USE bibliography, ONLY: iannuzzi2006,&
32 cite_reference,&
34 USE cell_types, ONLY: cell_type
44 USE cp_control_utils, ONLY: &
53 USE cp_output_handling, ONLY: cp_p_file,&
63 USE ec_environment, ONLY: ec_env_create,&
82 USE gamma, ONLY: init_md_ftable
85 USE header, ONLY: dftb_header,&
86 qs_header,&
87 se_header,&
92 USE input_constants, ONLY: &
109 USE kinds, ONLY: default_string_length,&
110 dp
114 USE kpoint_types, ONLY: get_kpoint_info,&
123 USE machine, ONLY: m_flush
124 USE mathconstants, ONLY: pi
129 USE mp2_setup, ONLY: read_mp2_section
130 USE mp2_types, ONLY: mp2_env_create,&
139 USE physcon, ONLY: kelvin
140 USE pw_env_types, ONLY: pw_env_type
157 USE qs_gcp_types, ONLY: qs_gcp_type
158 USE qs_gcp_utils, ONLY: qs_gcp_env_set,&
171 USE qs_kind_types, ONLY: &
175 USE qs_ks_types, ONLY: qs_ks_env_create,&
179 USE qs_mo_types, ONLY: allocate_mo_set,&
182 USE qs_rho0_methods, ONLY: init_rho0
187 USE qs_subsys_types, ONLY: qs_subsys_get,&
210 USE tblite_interface, ONLY: tb_get_basis,&
212 tb_init_wf,&
215 USE xtb_parameters, ONLY: init_xtb_basis,&
221#include "./base/base_uses.f90"
222
223 IMPLICIT NONE
224
225 PRIVATE
226
227 ! *** Global parameters ***
228 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_environment'
229
230 ! *** Public subroutines ***
231 PUBLIC :: qs_init
232
233CONTAINS
234
235! **************************************************************************************************
236!> \brief Read the input and the database files for the setup of the
237!> QUICKSTEP environment.
238!> \param qs_env ...
239!> \param para_env ...
240!> \param root_section ...
241!> \param globenv ...
242!> \param cp_subsys ...
243!> \param kpoint_env ...
244!> \param qmmm ...
245!> \param qmmm_env_qm ...
246!> \param force_env_section ...
247!> \param subsys_section ...
248!> \param use_motion_section ...
249!> \param silent ...
250!> \author Creation (22.05.2000,MK)
251! **************************************************************************************************
252 SUBROUTINE qs_init(qs_env, para_env, root_section, globenv, cp_subsys, kpoint_env, &
253 qmmm, qmmm_env_qm, force_env_section, subsys_section, &
254 use_motion_section, silent)
255
256 TYPE(qs_environment_type), POINTER :: qs_env
257 TYPE(mp_para_env_type), POINTER :: para_env
258 TYPE(section_vals_type), OPTIONAL, POINTER :: root_section
259 TYPE(global_environment_type), OPTIONAL, POINTER :: globenv
260 TYPE(cp_subsys_type), OPTIONAL, POINTER :: cp_subsys
261 TYPE(kpoint_type), OPTIONAL, POINTER :: kpoint_env
262 LOGICAL, INTENT(IN), OPTIONAL :: qmmm
263 TYPE(qmmm_env_qm_type), OPTIONAL, POINTER :: qmmm_env_qm
264 TYPE(section_vals_type), POINTER :: force_env_section, subsys_section
265 LOGICAL, INTENT(IN) :: use_motion_section
266 LOGICAL, INTENT(IN), OPTIONAL :: silent
267
268 CHARACTER(LEN=default_string_length) :: basis_type
269 INTEGER :: ikind, method_id, nelectron_total, &
270 nkind, nkp_grid(3)
271 LOGICAL :: do_active_space, do_admm_rpa, do_bse, do_debug_fdiff, do_debug_forces, &
272 do_debug_stress_tensor, do_ec_hfx, do_et, do_exx, do_gw, do_hfx, do_kpoints, &
273 do_linear_response, do_mp2, do_ri_mp2, do_ri_rpa, do_ri_sos_mp2, do_tddfpt, &
274 do_wfc_low_scaling, do_wfc_low_scaling_kpoints, do_xtb_tblite, is_identical, is_semi, &
275 mp2_present, my_qmmm, qmmm_decoupl, same_except_frac, use_ref_cell
276 REAL(kind=dp), DIMENSION(:, :), POINTER :: rtmat
277 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
278 TYPE(cell_type), POINTER :: my_cell, my_cell_ref
279 TYPE(cp_blacs_env_type), POINTER :: blacs_env
280 TYPE(dft_control_type), POINTER :: dft_control
281 TYPE(distribution_1d_type), POINTER :: local_particles
282 TYPE(energy_correction_type), POINTER :: ec_env
283 TYPE(excited_energy_type), POINTER :: exstate_env
284 TYPE(harris_type), POINTER :: harris_env
285 TYPE(kpoint_type), POINTER :: kpoints
286 TYPE(lri_environment_type), POINTER :: lri_env
287 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
288 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
289 TYPE(qs_ks_env_type), POINTER :: ks_env
290 TYPE(qs_subsys_type), POINTER :: subsys
291 TYPE(qs_wf_history_type), POINTER :: wf_history
292 TYPE(rel_control_type), POINTER :: rel_control
293 TYPE(scf_control_type), POINTER :: scf_control
294 TYPE(section_vals_type), POINTER :: active_space_section, dft_section, ec_hfx_section, &
295 ec_section, et_coupling_section, gw_section, hfx_section, kpoint_section, mp2_section, &
296 rpa_hfx_section, tddfpt_section, transport_section
297
298 NULLIFY (my_cell, my_cell_ref, atomic_kind_set, particle_set, &
299 qs_kind_set, kpoint_section, dft_section, ec_section, &
300 subsys, ks_env, dft_control, blacs_env)
301
302 CALL set_qs_env(qs_env, input=force_env_section)
303 IF (.NOT. ASSOCIATED(subsys_section)) THEN
304 subsys_section => section_vals_get_subs_vals(force_env_section, "SUBSYS")
305 END IF
306
307 ! QMMM
308 my_qmmm = .false.
309 IF (PRESENT(qmmm)) my_qmmm = qmmm
310 qmmm_decoupl = .false.
311 IF (PRESENT(qmmm_env_qm)) THEN
312 IF (qmmm_env_qm%qmmm_coupl_type == do_qmmm_gauss .OR. &
313 qmmm_env_qm%qmmm_coupl_type == do_qmmm_swave) THEN
314 ! For GAUSS/SWAVE methods there could be a DDAPC decoupling requested
315 qmmm_decoupl = my_qmmm .AND. qmmm_env_qm%periodic .AND. qmmm_env_qm%multipole
316 END IF
317 qs_env%qmmm_env_qm => qmmm_env_qm
318 END IF
319 CALL set_qs_env(qs_env=qs_env, qmmm=my_qmmm)
320
321 ! Possibly initialize arrays for SE
322 CALL section_vals_val_get(force_env_section, "DFT%QS%METHOD", i_val=method_id)
323 SELECT CASE (method_id)
326 CALL init_se_intd_array()
327 is_semi = .true.
329 is_semi = .true.
330 CASE DEFAULT
331 is_semi = .false.
332 END SELECT
333
334 ALLOCATE (subsys)
335 CALL qs_subsys_create(subsys, para_env, &
336 force_env_section=force_env_section, &
337 subsys_section=subsys_section, &
338 use_motion_section=use_motion_section, &
339 root_section=root_section, &
340 cp_subsys=cp_subsys, &
341 elkind=is_semi, silent=silent)
342
343 ALLOCATE (ks_env)
344 CALL qs_ks_env_create(ks_env)
345 CALL set_ks_env(ks_env, subsys=subsys)
346 CALL set_qs_env(qs_env, ks_env=ks_env)
347
348 CALL qs_subsys_get(subsys, &
349 cell=my_cell, &
350 cell_ref=my_cell_ref, &
351 use_ref_cell=use_ref_cell, &
352 atomic_kind_set=atomic_kind_set, &
353 qs_kind_set=qs_kind_set, &
354 particle_set=particle_set)
355
356 CALL set_ks_env(ks_env, para_env=para_env)
357 IF (PRESENT(globenv)) THEN
358 CALL cp_blacs_env_create(blacs_env, para_env, globenv%blacs_grid_layout, &
359 globenv%blacs_repeatable)
360 ELSE
361 CALL cp_blacs_env_create(blacs_env, para_env)
362 END IF
363 CALL set_ks_env(ks_env, blacs_env=blacs_env)
364 CALL cp_blacs_env_release(blacs_env)
365
366 ! *** Setup the grids for the G-space Interpolation if any
367 CALL cp_ddapc_ewald_create(qs_env%cp_ddapc_ewald, qmmm_decoupl, my_cell, &
368 force_env_section, subsys_section, para_env)
369
370 ! kpoints
371 IF (PRESENT(kpoint_env)) THEN
372 kpoints => kpoint_env
373 CALL set_qs_env(qs_env=qs_env, kpoints=kpoints)
374 CALL kpoint_initialize(kpoints, particle_set, my_cell)
375 ELSE
376 NULLIFY (kpoints)
377 CALL kpoint_create(kpoints)
378 CALL set_qs_env(qs_env=qs_env, kpoints=kpoints)
379 kpoint_section => section_vals_get_subs_vals(qs_env%input, "DFT%KPOINTS")
380 CALL read_kpoint_section(kpoints, kpoint_section, my_cell%hmat)
381 do_hfx = .false.
382 hfx_section => section_vals_get_subs_vals(qs_env%input, "DFT%XC%HF")
383 CALL section_vals_get(hfx_section, explicit=do_hfx)
384 do_exx = .false.
385 rpa_hfx_section => section_vals_get_subs_vals(qs_env%input, "DFT%XC%WF_CORRELATION%RI_RPA%HF")
386 CALL section_vals_get(rpa_hfx_section, explicit=do_exx)
387 do_gw = .false.
388 gw_section => section_vals_get_subs_vals(qs_env%input, "DFT%XC%WF_CORRELATION%RI_RPA%GW")
389 CALL section_vals_get(gw_section, explicit=do_gw)
390 IF (.NOT. do_gw) THEN
391 gw_section => section_vals_get_subs_vals(qs_env%input, "PROPERTIES%BANDSTRUCTURE%GW")
392 CALL section_vals_get(gw_section, explicit=do_gw)
393 END IF
394 do_tddfpt = .false.
395 do_bse = .false.
396 tddfpt_section => section_vals_get_subs_vals(qs_env%input, "PROPERTIES%TDDFPT")
397 CALL section_vals_get(tddfpt_section, explicit=do_tddfpt)
398 IF (do_tddfpt) THEN
399 CALL section_vals_val_get(tddfpt_section, "DO_BSE", l_val=do_bse)
400 IF (.NOT. do_bse) &
401 CALL section_vals_val_get(tddfpt_section, "DO_BSE_W_ONLY", l_val=do_bse)
402 IF (.NOT. do_bse) &
403 CALL section_vals_val_get(tddfpt_section, "DO_BSE_GW_ONLY", l_val=do_bse)
404 END IF
405 do_active_space = .false.
406 active_space_section => section_vals_get_subs_vals(qs_env%input, "DFT%ACTIVE_SPACE")
407 CALL section_vals_get(active_space_section, explicit=do_active_space)
408 do_xtb_tblite = .false.
409 IF (method_id == do_method_xtb) THEN
410 CALL section_vals_val_get(qs_env%input, "DFT%QS%XTB%TBLITE%_SECTION_PARAMETERS_", &
411 l_val=do_xtb_tblite)
412 END IF
413 do_linear_response = .false.
414 IF (PRESENT(globenv)) do_linear_response = globenv%run_type_id == linear_response_run
415 do_debug_fdiff = .false.
416 IF (PRESENT(globenv)) do_debug_fdiff = globenv%run_type_id == debug_run
417 IF (do_debug_fdiff .AND. PRESENT(root_section)) THEN
418 CALL section_vals_val_get(root_section, "DEBUG%DEBUG_FORCES", &
419 l_val=do_debug_forces)
420 CALL section_vals_val_get(root_section, "DEBUG%DEBUG_STRESS_TENSOR", &
421 l_val=do_debug_stress_tensor)
422 do_debug_fdiff = do_debug_forces .OR. do_debug_stress_tensor
423 END IF
424 do_mp2 = .false.
425 do_ri_mp2 = .false.
426 do_ri_sos_mp2 = .false.
427 do_ri_rpa = .false.
428 do_wfc_low_scaling = .false.
429 do_wfc_low_scaling_kpoints = .false.
430 mp2_section => section_vals_get_subs_vals(qs_env%input, "DFT%XC%WF_CORRELATION")
431 CALL section_vals_get(mp2_section, explicit=mp2_present)
432 IF (mp2_present) THEN
433 CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%MP2%_SECTION_PARAMETERS_", &
434 l_val=do_mp2)
435 CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%RI_MP2%_SECTION_PARAMETERS_", &
436 l_val=do_ri_mp2)
437 CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%RI_SOS_MP2%_SECTION_PARAMETERS_", &
438 l_val=do_ri_sos_mp2)
439 CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%RI_RPA%_SECTION_PARAMETERS_", &
440 l_val=do_ri_rpa)
441 CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%LOW_SCALING%_SECTION_PARAMETERS_", &
442 l_val=do_wfc_low_scaling)
443 CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%LOW_SCALING%DO_KPOINTS", &
444 l_val=do_wfc_low_scaling_kpoints)
445 IF (.NOT. do_bse) &
446 CALL section_vals_val_get(qs_env%input, &
447 "DFT%XC%WF_CORRELATION%RI_RPA%GW%BSE%_SECTION_PARAMETERS_", &
448 l_val=do_bse)
449 END IF
450 CALL restrict_unsupported_atomic_kpoint_symmetry(kpoints, method_id, do_hfx, do_exx, do_gw, &
451 do_tddfpt, do_active_space, do_linear_response, &
452 do_debug_fdiff, &
453 do_mp2 .OR. do_ri_mp2 .OR. do_ri_sos_mp2, &
454 do_ri_rpa .AND. .NOT. do_gw, do_bse, &
455 do_wfc_low_scaling, do_wfc_low_scaling_kpoints, &
456 do_xtb_tblite)
457 CALL kpoint_initialize(kpoints, particle_set, my_cell)
458 dft_section => section_vals_get_subs_vals(qs_env%input, "DFT")
459 CALL write_kpoint_info(kpoints, dft_section=dft_section)
460 END IF
461
462 CALL qs_init_subsys(qs_env, para_env, subsys, my_cell, my_cell_ref, use_ref_cell, &
463 subsys_section, silent=silent)
464
465 CALL get_qs_env(qs_env, dft_control=dft_control)
466 IF (method_id == do_method_lrigpw .OR. dft_control%qs_control%lri_optbas) THEN
467 CALL get_qs_env(qs_env=qs_env, lri_env=lri_env)
468 CALL lri_env_basis("LRI", qs_env, lri_env, qs_kind_set)
469 ELSE IF (method_id == do_method_rigpw) THEN
470 CALL cp_warn(__location__, "Experimental code: "// &
471 "RIGPW should only be used for testing.")
472 CALL get_qs_env(qs_env=qs_env, lri_env=lri_env)
473 CALL lri_env_basis("RI", qs_env, lri_env, qs_kind_set)
474 END IF
475
476 IF (my_qmmm .AND. PRESENT(qmmm_env_qm) .AND. .NOT. dft_control%qs_control%commensurate_mgrids) THEN
477 IF (qmmm_env_qm%qmmm_coupl_type == do_qmmm_gauss .OR. qmmm_env_qm%qmmm_coupl_type == do_qmmm_swave) THEN
478 CALL cp_abort(__location__, "QM/MM with coupling GAUSS or S-WAVE requires "// &
479 "keyword FORCE_EVAL/DFT/MGRID/COMMENSURATE to be enabled.")
480 END IF
481 END IF
482
483 ! more kpoint stuff
484 CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints, blacs_env=blacs_env)
485 IF (do_kpoints) THEN
486 CALL kpoint_env_initialize(kpoints, para_env, blacs_env, with_aux_fit=dft_control%do_admm)
487 CALL kpoint_initialize_mos(kpoints, qs_env%mos)
488 CALL get_qs_env(qs_env=qs_env, wf_history=wf_history)
489 CALL wfi_create_for_kp(wf_history)
490 END IF
491 ! basis set symmetry rotations
492 IF (do_kpoints) THEN
493 CALL qs_basis_rotation(qs_env, kpoints)
494 END IF
495
496 do_hfx = .false.
497 hfx_section => section_vals_get_subs_vals(qs_env%input, "DFT%XC%HF")
498 CALL section_vals_get(hfx_section, explicit=do_hfx)
499 CALL get_qs_env(qs_env, dft_control=dft_control, scf_control=scf_control, nelectron_total=nelectron_total)
500 IF (do_hfx) THEN
501 ! Retrieve particle_set and atomic_kind_set (needed for both kinds of initialization)
502 nkp_grid = 1
503 IF (do_kpoints) CALL get_kpoint_info(kpoints, nkp_grid=nkp_grid)
504 IF (dft_control%do_admm) THEN
505 basis_type = 'AUX_FIT'
506 ELSE
507 basis_type = 'ORB'
508 END IF
509 CALL hfx_create(qs_env%x_data, para_env, hfx_section, atomic_kind_set, &
510 qs_kind_set, particle_set, dft_control, my_cell, orb_basis=basis_type, &
511 nelectron_total=nelectron_total, nkp_grid=nkp_grid)
512 END IF
513
514 mp2_section => section_vals_get_subs_vals(qs_env%input, "DFT%XC%WF_CORRELATION")
515 CALL section_vals_get(mp2_section, explicit=mp2_present)
516 IF (mp2_present) THEN
517 cpassert(ASSOCIATED(qs_env%mp2_env))
518 CALL read_mp2_section(qs_env%input, qs_env%mp2_env)
519 ! create the EXX section if necessary
520 do_exx = .false.
521 rpa_hfx_section => section_vals_get_subs_vals(qs_env%input, "DFT%XC%WF_CORRELATION%RI_RPA%HF")
522 CALL section_vals_get(rpa_hfx_section, explicit=do_exx)
523 IF (do_exx) THEN
524
525 ! do_exx in call of hfx_create decides whether to go without ADMM (do_exx=.TRUE.) or with
526 ! ADMM (do_exx=.FALSE.)
527 CALL section_vals_val_get(mp2_section, "RI_RPA%ADMM", l_val=do_admm_rpa)
528
529 ! Reuse the HFX integrals from the qs_env if applicable
530 qs_env%mp2_env%ri_rpa%reuse_hfx = .true.
531 IF (.NOT. do_hfx) qs_env%mp2_env%ri_rpa%reuse_hfx = .false.
532 CALL compare_hfx_sections(hfx_section, rpa_hfx_section, is_identical, same_except_frac)
533 IF (.NOT. (is_identical .OR. same_except_frac)) qs_env%mp2_env%ri_rpa%reuse_hfx = .false.
534 IF (dft_control%do_admm .AND. .NOT. do_admm_rpa) qs_env%mp2_env%ri_rpa%reuse_hfx = .false.
535
536 IF (.NOT. qs_env%mp2_env%ri_rpa%reuse_hfx) THEN
537 IF (do_admm_rpa) THEN
538 basis_type = 'AUX_FIT'
539 ELSE
540 basis_type = 'ORB'
541 END IF
542 CALL hfx_create(qs_env%mp2_env%ri_rpa%x_data, para_env, rpa_hfx_section, atomic_kind_set, &
543 qs_kind_set, particle_set, dft_control, my_cell, orb_basis=basis_type, &
544 nelectron_total=nelectron_total)
545 ELSE
546 qs_env%mp2_env%ri_rpa%x_data => qs_env%x_data
547 END IF
548 END IF
549 END IF
550
551 IF (dft_control%qs_control%do_kg) THEN
552 CALL cite_reference(iannuzzi2006)
553 CALL kg_env_create(qs_env, qs_env%kg_env, qs_kind_set, qs_env%input)
554 END IF
555
556 dft_section => section_vals_get_subs_vals(qs_env%input, "DFT")
557 CALL section_vals_val_get(dft_section, "EXCITED_STATES%_SECTION_PARAMETERS_", &
558 l_val=qs_env%excited_state)
559 NULLIFY (exstate_env)
560 CALL exstate_create(exstate_env, qs_env%excited_state, dft_section)
561 CALL set_qs_env(qs_env, exstate_env=exstate_env)
562
563 et_coupling_section => section_vals_get_subs_vals(qs_env%input, &
564 "PROPERTIES%ET_COUPLING")
565 CALL section_vals_get(et_coupling_section, explicit=do_et)
566 IF (do_et) CALL et_coupling_create(qs_env%et_coupling)
567
568 transport_section => section_vals_get_subs_vals(qs_env%input, "DFT%TRANSPORT")
569 CALL section_vals_get(transport_section, explicit=qs_env%do_transport)
570 IF (qs_env%do_transport) THEN
571 CALL transport_env_create(qs_env)
572 END IF
573
574 CALL get_qs_env(qs_env, harris_env=harris_env)
575 IF (qs_env%harris_method) THEN
576 ! initialize the Harris input density and potential integrals
577 CALL get_qs_env(qs_env, local_particles=local_particles)
578 CALL harris_rhoin_init(harris_env%rhoin, "RHOIN", qs_kind_set, atomic_kind_set, &
579 local_particles, dft_control%nspins)
580 ! Print information of the HARRIS section
581 CALL harris_write_input(harris_env)
582 END IF
583
584 NULLIFY (ec_env)
585 dft_section => section_vals_get_subs_vals(qs_env%input, "DFT")
586 CALL section_vals_val_get(dft_section, "ENERGY_CORRECTION%_SECTION_PARAMETERS_", &
587 l_val=qs_env%energy_correction)
588 ec_section => section_vals_get_subs_vals(qs_env%input, "DFT%ENERGY_CORRECTION")
589 CALL ec_env_create(qs_env, ec_env, dft_section, ec_section)
590 CALL set_qs_env(qs_env, ec_env=ec_env)
591
592 IF (qs_env%energy_correction) THEN
593 ! Energy correction with Hartree-Fock exchange
594 ec_hfx_section => section_vals_get_subs_vals(ec_section, "XC%HF")
595 CALL section_vals_get(ec_hfx_section, explicit=do_ec_hfx)
596
597 IF (ec_env%do_ec_hfx) THEN
598
599 ! kpoints and HFX not yet compatible
600 IF (ec_env%do_kpoints) THEN
601 CALL cp_abort(__location__, &
602 "Energy correction methods with hybrid functionals "// &
603 "and kpoints is not yet available.")
604 END IF
605
606 ! Hybrid functionals require same basis
607 IF (ec_env%basis_inconsistent) THEN
608 CALL cp_abort(__location__, &
609 "Energy correction methods with hybrid functionals: "// &
610 "correction and ground state need to use the same basis. "// &
611 "Checked by comparing basis set names only.")
612 END IF
613
614 ! Similar to RPA_HFX we can check if HFX integrals from the qs_env can be reused
615 IF (ec_env%do_ec_admm .AND. .NOT. dft_control%do_admm) THEN
616 CALL cp_abort(__location__, "Need an ADMM input section for ADMM EC to work")
617 END IF
618
619 ec_env%reuse_hfx = .true.
620 IF (.NOT. do_hfx) ec_env%reuse_hfx = .false.
621 CALL compare_hfx_sections(hfx_section, ec_hfx_section, is_identical, same_except_frac)
622 IF (.NOT. (is_identical .OR. same_except_frac)) ec_env%reuse_hfx = .false.
623 IF (dft_control%do_admm .AND. .NOT. ec_env%do_ec_admm) ec_env%reuse_hfx = .false.
624
625 IF (.NOT. ec_env%reuse_hfx) THEN
626 IF (ec_env%do_ec_admm) THEN
627 basis_type = 'AUX_FIT'
628 ELSE
629 basis_type = 'ORB'
630 END IF
631 CALL hfx_create(ec_env%x_data, para_env, ec_hfx_section, atomic_kind_set, &
632 qs_kind_set, particle_set, dft_control, my_cell, orb_basis=basis_type, &
633 nelectron_total=nelectron_total)
634 ELSE
635 ec_env%x_data => qs_env%x_data
636 END IF
637 END IF
638
639 ! Print information of the EC section
640 CALL ec_write_input(ec_env)
641
642 END IF
643
644 IF (dft_control%qs_control%do_almo_scf) THEN
645 CALL almo_scf_env_create(qs_env)
646 END IF
647
648 ! see if we have atomic relativistic corrections
649 CALL get_qs_env(qs_env, rel_control=rel_control)
650 IF (rel_control%rel_method /= rel_none) THEN
651 IF (rel_control%rel_transformation == rel_trans_atom) THEN
652 nkind = SIZE(atomic_kind_set)
653 DO ikind = 1, nkind
654 NULLIFY (rtmat)
655 CALL calculate_atomic_relkin(atomic_kind_set(ikind), qs_kind_set(ikind), rel_control, rtmat)
656 IF (ASSOCIATED(rtmat)) CALL set_qs_kind(qs_kind_set(ikind), reltmat=rtmat)
657 END DO
658 END IF
659 END IF
660
661 END SUBROUTINE qs_init
662
663! **************************************************************************************************
664!> \brief Restrict atomic k-point symmetry for methods not supporting it yet
665!> \param kpoints ...
666!> \param method_id ...
667!> \param do_hfx ...
668!> \param do_exx ...
669!> \param do_gw ...
670!> \param do_tddfpt ...
671!> \param do_active_space ...
672!> \param do_linear_response ...
673!> \param do_debug_fdiff ...
674!> \param do_mp2 ...
675!> \param do_rpa ...
676!> \param do_bse ...
677!> \param do_wfc_low_scaling ...
678!> \param do_wfc_low_scaling_kpoints ...
679!> \param do_xtb_tblite ...
680! **************************************************************************************************
681 SUBROUTINE restrict_unsupported_atomic_kpoint_symmetry(kpoints, method_id, do_hfx, do_exx, do_gw, &
682 do_tddfpt, do_active_space, do_linear_response, &
683 do_debug_fdiff, &
684 do_mp2, do_rpa, do_bse, do_wfc_low_scaling, &
685 do_wfc_low_scaling_kpoints, do_xtb_tblite)
686 TYPE(kpoint_type), POINTER :: kpoints
687 INTEGER, INTENT(IN) :: method_id
688 LOGICAL, INTENT(IN) :: do_hfx, do_exx, do_gw, do_tddfpt, do_active_space, &
689 do_linear_response, do_debug_fdiff, do_mp2, do_rpa, do_bse, do_wfc_low_scaling, &
690 do_wfc_low_scaling_kpoints, do_xtb_tblite
691
692 CHARACTER(LEN=default_string_length) :: kp_scheme, reason
693 LOGICAL :: full_grid, inversion_symmetry_only, &
694 kpoint_symmetry
695
696 reason = unsupported_kpoint_method_reason(method_id, do_gw, do_tddfpt, do_linear_response, &
697 do_mp2, do_bse, do_xtb_tblite)
698 IF (len_trim(reason) > 0) THEN
699 CALL get_kpoint_info(kpoints, kp_scheme=kp_scheme)
700 IF (len_trim(kp_scheme) > 0 .AND. trim(kp_scheme) /= "NONE") THEN
701 IF (trim(reason) == "GW") THEN
702 CALL cp_abort(__location__, &
703 "DFT%KPOINTS are not supported with GW; use "// &
704 "WF_CORRELATION%LOW_SCALING%KPOINTS and RI_RPA%GW%KPOINTS_SELF_ENERGY "// &
705 "for GW k-point sampling.")
706 ELSE
707 CALL cp_abort(__location__, &
708 "DFT%KPOINTS are not supported with "//trim(reason)// &
709 "; remove DFT%KPOINTS for these calculations.")
710 END IF
711 END IF
712 END IF
713 IF (do_active_space) THEN
714 CALL get_kpoint_info(kpoints, kp_scheme=kp_scheme)
715 IF (len_trim(kp_scheme) > 0 .AND. trim(kp_scheme) /= "NONE" .AND. &
716 trim(kp_scheme) /= "GAMMA") THEN
717 CALL cp_abort(__location__, &
718 "Only Gamma-point DFT%KPOINTS are supported with ACTIVE_SPACE; "// &
719 "use SCHEME GAMMA, SCHEME NONE, or remove DFT%KPOINTS.")
720 END IF
721 END IF
722
723 CALL get_kpoint_info(kpoints, symmetry=kpoint_symmetry, full_grid=full_grid, &
724 inversion_symmetry_only=inversion_symmetry_only)
725 IF (.NOT. (kpoint_symmetry .AND. .NOT. full_grid .AND. .NOT. inversion_symmetry_only)) RETURN
726
727 reason = unsupported_atomic_kpoint_symmetry_reason(method_id, do_hfx, do_exx, do_gw, &
728 do_tddfpt, do_active_space, do_linear_response, &
729 do_debug_fdiff, &
730 do_mp2, do_rpa, do_bse, do_wfc_low_scaling, &
731 do_wfc_low_scaling_kpoints, do_xtb_tblite)
732 IF (len_trim(reason) == 0) RETURN
733
734 CALL cp_warn(__location__, &
735 "Atomic k-point symmetry is currently not implemented for "//trim(reason)// &
736 "; restricting to inversion/time-reversal symmetry.")
737 CALL set_kpoint_info(kpoints, inversion_symmetry_only=.true.)
738
739 END SUBROUTINE restrict_unsupported_atomic_kpoint_symmetry
740
741! **************************************************************************************************
742!> \brief Return the reason why k-points are not enabled for a method
743!> \param method_id ...
744!> \param do_gw ...
745!> \param do_tddfpt ...
746!> \param do_linear_response ...
747!> \param do_mp2 ...
748!> \param do_bse ...
749!> \param do_xtb_tblite ...
750!> \return reason
751! **************************************************************************************************
752 FUNCTION unsupported_kpoint_method_reason(method_id, do_gw, do_tddfpt, do_linear_response, &
753 do_mp2, do_bse, do_xtb_tblite) RESULT(reason)
754 INTEGER, INTENT(IN) :: method_id
755 LOGICAL, INTENT(IN) :: do_gw, do_tddfpt, do_linear_response, &
756 do_mp2, do_bse, do_xtb_tblite
757 CHARACTER(LEN=default_string_length) :: reason
758
759 reason = ""
760 mark_used(do_gw)
761 mark_used(do_mp2)
762 mark_used(do_xtb_tblite)
763
764 IF (do_bse) THEN
765 reason = "BSE"
766 RETURN
767 END IF
768 IF (do_tddfpt) THEN
769 reason = "TDDFPT/TDDFT"
770 RETURN
771 END IF
772 IF (do_linear_response) THEN
773 reason = "LINEAR_RESPONSE/DFPT"
774 RETURN
775 END IF
776 SELECT CASE (method_id)
777 CASE (do_method_rigpw)
778 reason = "RIGPW"
779 CASE (do_method_ofgpw)
780 reason = "OFGPW"
783 reason = "semiempirical methods"
784 CASE DEFAULT
785 reason = ""
786 END SELECT
787
788 END FUNCTION unsupported_kpoint_method_reason
789
790! **************************************************************************************************
791!> \brief Return the reason why atomic k-point symmetry is not enabled
792!> \param method_id ...
793!> \param do_hfx ...
794!> \param do_exx ...
795!> \param do_gw ...
796!> \param do_tddfpt ...
797!> \param do_active_space ...
798!> \param do_linear_response ...
799!> \param do_debug_fdiff ...
800!> \param do_mp2 ...
801!> \param do_rpa ...
802!> \param do_bse ...
803!> \param do_wfc_low_scaling ...
804!> \param do_wfc_low_scaling_kpoints ...
805!> \param do_xtb_tblite ...
806!> \return reason
807! **************************************************************************************************
808 FUNCTION unsupported_atomic_kpoint_symmetry_reason(method_id, do_hfx, do_exx, do_gw, do_tddfpt, &
809 do_active_space, do_linear_response, do_debug_fdiff, &
810 do_mp2, do_rpa, do_bse, do_wfc_low_scaling, &
811 do_wfc_low_scaling_kpoints, do_xtb_tblite) RESULT(reason)
812 INTEGER, INTENT(IN) :: method_id
813 LOGICAL, INTENT(IN) :: do_hfx, do_exx, do_gw, do_tddfpt, do_active_space, &
814 do_linear_response, do_debug_fdiff, do_mp2, do_rpa, do_bse, do_wfc_low_scaling, &
815 do_wfc_low_scaling_kpoints, do_xtb_tblite
816 CHARACTER(LEN=default_string_length) :: reason
817
818 reason = ""
819 mark_used(do_debug_fdiff)
820 mark_used(do_xtb_tblite)
821
822 SELECT CASE (method_id)
823 CASE (do_method_lrigpw)
824 reason = "LRIGPW"
825 CASE (do_method_rigpw)
826 reason = "RIGPW"
829 reason = "semiempirical methods"
830 CASE DEFAULT
831 reason = ""
832 END SELECT
833
834 IF (len_trim(reason) > 0) RETURN
835 IF (do_bse) THEN
836 reason = "BSE"
837 ELSE IF (do_gw) THEN
838 reason = "GW"
839 ELSE IF (do_hfx .OR. do_exx) THEN
840 reason = "HFX/HF"
841 ELSE IF (do_tddfpt) THEN
842 reason = "TDDFPT/TDDFT"
843 ELSE IF (do_active_space) THEN
844 reason = "ACTIVE_SPACE"
845 ELSE IF (do_linear_response) THEN
846 reason = "LINEAR_RESPONSE/DFPT"
847 ELSE IF (do_mp2) THEN
848 reason = "MP2"
849 ELSE IF (do_rpa .AND. do_wfc_low_scaling_kpoints) THEN
850 reason = "LOW_SCALING RPA"
851 ELSE IF (do_wfc_low_scaling) THEN
852 reason = "LOW_SCALING WF_CORRELATION"
853 ELSE IF (do_rpa) THEN
854 reason = "RPA"
855 END IF
856
857 END FUNCTION unsupported_atomic_kpoint_symmetry_reason
858
859! **************************************************************************************************
860!> \brief Initialize the qs environment (subsys)
861!> \param qs_env ...
862!> \param para_env ...
863!> \param subsys ...
864!> \param cell ...
865!> \param cell_ref ...
866!> \param use_ref_cell ...
867!> \param subsys_section ...
868!> \param silent ...
869!> \author Creation (22.05.2000,MK)
870! **************************************************************************************************
871 SUBROUTINE qs_init_subsys(qs_env, para_env, subsys, cell, cell_ref, use_ref_cell, subsys_section, &
872 silent)
873
874 TYPE(qs_environment_type), POINTER :: qs_env
875 TYPE(mp_para_env_type), POINTER :: para_env
876 TYPE(qs_subsys_type), POINTER :: subsys
877 TYPE(cell_type), POINTER :: cell, cell_ref
878 LOGICAL, INTENT(in) :: use_ref_cell
879 TYPE(section_vals_type), POINTER :: subsys_section
880 LOGICAL, INTENT(in), OPTIONAL :: silent
881
882 CHARACTER(len=*), PARAMETER :: routinen = 'qs_init_subsys'
883
884 CHARACTER(len=2) :: element_symbol
885 INTEGER :: gfn_type, handle, ikind, ispin, iw, lmax_sphere, maxl, maxlgto, maxlgto_lri, &
886 maxlgto_nuc, maxlppl, maxlppnl, method_id, multiplicity, my_ival, n_ao, n_mo_add, natom, &
887 nelectron, ngauss, nkind, output_unit, sort_basis, tnadd_method
888 INTEGER, DIMENSION(2) :: n_mo, nelectron_spin
889 INTEGER, DIMENSION(5) :: occ
890 INTEGER, DIMENSION(:), POINTER :: mo_index_range
891 LOGICAL :: all_potential_present, be_silent, cneo_potential_present, do_kpoints, do_ri_hfx, &
892 do_ri_mp2, do_ri_rpa, do_ri_sos_mp2, do_rpa_ri_exx, do_wfc_im_time, e1terms, &
893 has_unit_metric, lribas, mp2_present, orb_gradient, paw_atom
894 REAL(kind=dp) :: alpha, ccore, ewald_rcut, fxx, maxocc, &
895 rc, rcut, total_zeff_corr, &
896 verlet_skin, zeff_correction
897 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
898 TYPE(cp_logger_type), POINTER :: logger
899 TYPE(dft_control_type), POINTER :: dft_control
900 TYPE(dftb_control_type), POINTER :: dftb_control
901 TYPE(distribution_1d_type), POINTER :: local_molecules, local_particles
902 TYPE(ewald_environment_type), POINTER :: ewald_env
903 TYPE(ewald_pw_type), POINTER :: ewald_pw
904 TYPE(fist_nonbond_env_type), POINTER :: se_nonbond_env
905 TYPE(gapw_control_type), POINTER :: gapw_control
906 TYPE(gto_basis_set_type), POINTER :: aux_fit_basis, lri_aux_basis, &
907 rhoin_basis, ri_aux_basis_set, &
908 ri_hfx_basis, ri_xas_basis, &
909 tmp_basis_set
910 TYPE(harris_type), POINTER :: harris_env
911 TYPE(local_rho_type), POINTER :: local_rho_set
912 TYPE(lri_environment_type), POINTER :: lri_env
913 TYPE(mo_set_type), DIMENSION(:), POINTER :: mos, mos_last_converged
914 TYPE(molecule_kind_type), DIMENSION(:), POINTER :: molecule_kind_set
915 TYPE(molecule_type), DIMENSION(:), POINTER :: molecule_set
916 TYPE(mp2_type), POINTER :: mp2_env
917 TYPE(nddo_mpole_type), POINTER :: se_nddo_mpole
918 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
919 TYPE(pw_env_type), POINTER :: pw_env
920 TYPE(qs_control_type), POINTER :: qs_control
921 TYPE(qs_dftb_pairpot_type), DIMENSION(:, :), &
922 POINTER :: dftb_potential
923 TYPE(qs_dispersion_type), POINTER :: dispersion_env
924 TYPE(qs_energy_type), POINTER :: energy
925 TYPE(qs_force_type), DIMENSION(:), POINTER :: force
926 TYPE(qs_gcp_type), POINTER :: gcp_env
927 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
928 TYPE(qs_kind_type), POINTER :: qs_kind
929 TYPE(qs_ks_env_type), POINTER :: ks_env
930 TYPE(qs_wf_history_type), POINTER :: wf_history
931 TYPE(rho0_mpole_type), POINTER :: rho0_mpole
932 TYPE(rho_atom_type), DIMENSION(:), POINTER :: rho_atom_set
933 TYPE(scf_control_type), POINTER :: scf_control
934 TYPE(se_taper_type), POINTER :: se_taper
935 TYPE(section_vals_type), POINTER :: dft_section, et_coupling_section, et_ddapc_section, &
936 ewald_section, harris_section, lri_section, mp2_section, nl_section, poisson_section, &
937 pp_section, print_section, qs_section, rixs_section, se_section, tddfpt_section, &
938 xc_section
939 TYPE(semi_empirical_control_type), POINTER :: se_control
940 TYPE(semi_empirical_si_type), POINTER :: se_store_int_env
941 TYPE(xtb_control_type), POINTER :: xtb_control
942
943 CALL timeset(routinen, handle)
944 NULLIFY (logger)
945 logger => cp_get_default_logger()
946 output_unit = cp_logger_get_default_io_unit(logger)
947
948 be_silent = .false.
949 IF (PRESENT(silent)) be_silent = silent
950
951 CALL cite_reference(cp2kqs2020)
952
953 ! Initialise the Quickstep environment
954 NULLIFY (mos, se_taper)
955 NULLIFY (dft_control)
956 NULLIFY (energy)
957 NULLIFY (force)
958 NULLIFY (local_molecules)
959 NULLIFY (local_particles)
960 NULLIFY (scf_control)
961 NULLIFY (dft_section)
962 NULLIFY (et_coupling_section)
963 NULLIFY (ks_env)
964 NULLIFY (mos_last_converged)
965 dft_section => section_vals_get_subs_vals(qs_env%input, "DFT")
966 qs_section => section_vals_get_subs_vals(dft_section, "QS")
967 et_coupling_section => section_vals_get_subs_vals(qs_env%input, "PROPERTIES%ET_COUPLING")
968 ! reimplemented TDDFPT
969 tddfpt_section => section_vals_get_subs_vals(qs_env%input, "PROPERTIES%TDDFPT")
970 rixs_section => section_vals_get_subs_vals(qs_env%input, "PROPERTIES%RIXS")
971
972 CALL qs_subsys_get(subsys, particle_set=particle_set, &
973 qs_kind_set=qs_kind_set, &
974 atomic_kind_set=atomic_kind_set, &
975 molecule_set=molecule_set, &
976 molecule_kind_set=molecule_kind_set)
977
978 ! Read the input section with the DFT control parameters
979 CALL read_dft_control(dft_control, dft_section)
980
981 ! Set periodicity flag
982 dft_control%qs_control%periodicity = sum(cell%perd)
983
984 ! Read the input section with the Quickstep control parameters
985 CALL read_qs_section(dft_control%qs_control, qs_section)
986
987 ! Print the Quickstep program banner (copyright and version number)
988 IF (.NOT. be_silent) THEN
989 iw = cp_print_key_unit_nr(logger, dft_section, "PRINT%PROGRAM_BANNER", extension=".Log")
990 CALL section_vals_val_get(qs_section, "METHOD", i_val=method_id)
991 SELECT CASE (method_id)
992 CASE DEFAULT
993 CALL qs_header(iw)
996 CALL se_header(iw)
997 CASE (do_method_dftb)
998 CALL dftb_header(iw)
999 CASE (do_method_xtb)
1000 IF (dft_control%qs_control%xtb_control%do_tblite) THEN
1001 CALL tblite_header(iw, dft_control%qs_control%xtb_control%tblite_method)
1002 ELSE
1003 gfn_type = dft_control%qs_control%xtb_control%gfn_type
1004 CALL xtb_header(iw, gfn_type)
1005 END IF
1006 END SELECT
1007 CALL cp_print_key_finished_output(iw, logger, dft_section, &
1008 "PRINT%PROGRAM_BANNER")
1009 END IF
1010
1011 IF (dft_control%do_sccs .AND. dft_control%qs_control%gapw) THEN
1012 cpabort("SCCS is not yet implemented with GAPW")
1013 END IF
1014 CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints)
1015 IF (do_kpoints) THEN
1016 ! reset some of the settings for wfn extrapolation for kpoints
1017 SELECT CASE (dft_control%qs_control%wf_interpolation_method_nr)
1020 CALL cp_warn(__location__, "Linear WFN-based extrapolation methods are not "// &
1021 "implemented for k-points. Switching to USE_PREV_WF.")
1022 dft_control%qs_control%wf_interpolation_method_nr = wfi_use_prev_wf_method_nr
1023 END SELECT
1024 END IF
1025
1026 ! Check if any kind of electron transfer calculation has to be performed
1027 CALL section_vals_val_get(et_coupling_section, "TYPE_OF_CONSTRAINT", i_val=my_ival)
1028 dft_control%qs_control%et_coupling_calc = .false.
1029 IF (my_ival == do_et_ddapc) THEN
1030 et_ddapc_section => section_vals_get_subs_vals(et_coupling_section, "DDAPC_RESTRAINT_A")
1031 dft_control%qs_control%et_coupling_calc = .true.
1032 dft_control%qs_control%ddapc_restraint = .true.
1033 CALL read_ddapc_section(dft_control%qs_control, ddapc_restraint_section=et_ddapc_section)
1034 END IF
1035
1036 CALL read_mgrid_section(dft_control%qs_control, dft_section)
1037
1038 ! Reimplemented TDDFPT
1039 CALL read_tddfpt2_control(dft_control%tddfpt2_control, tddfpt_section, dft_control%qs_control)
1040
1041 ! RIXS
1042 CALL section_vals_get(rixs_section, explicit=qs_env%do_rixs)
1043 IF (qs_env%do_rixs) THEN
1044 CALL read_rixs_control(dft_control%rixs_control, rixs_section, dft_control%qs_control)
1045 END IF
1046
1047 ! Create relativistic control section
1048 block
1049 TYPE(rel_control_type), POINTER :: rel_control
1050 ALLOCATE (rel_control)
1051 CALL rel_c_create(rel_control)
1052 CALL rel_c_read_parameters(rel_control, dft_section)
1053 CALL set_qs_env(qs_env, rel_control=rel_control)
1054 END block
1055
1056 ! Read DFTB parameter files
1057 IF (dft_control%qs_control%method_id == do_method_dftb) THEN
1058 NULLIFY (ewald_env, ewald_pw, dftb_potential)
1059 dftb_control => dft_control%qs_control%dftb_control
1060 CALL qs_dftb_param_init(atomic_kind_set, qs_kind_set, dftb_control, dftb_potential, &
1061 subsys_section=subsys_section, para_env=para_env)
1062 CALL set_qs_env(qs_env, dftb_potential=dftb_potential)
1063 ! check for Ewald
1064 IF (dftb_control%do_ewald) THEN
1065 ALLOCATE (ewald_env)
1066 CALL ewald_env_create(ewald_env, para_env)
1067 poisson_section => section_vals_get_subs_vals(dft_section, "POISSON")
1068 CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
1069 ewald_section => section_vals_get_subs_vals(poisson_section, "EWALD")
1070 print_section => section_vals_get_subs_vals(qs_env%input, "PRINT%GRID_INFORMATION")
1071 CALL get_qs_kind_set(qs_kind_set, basis_rcut=ewald_rcut)
1072 CALL read_ewald_section_tb(ewald_env, ewald_section, cell_ref%hmat, &
1073 cell_periodic=cell%perd)
1074 ALLOCATE (ewald_pw)
1075 CALL ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section=print_section)
1076 CALL set_qs_env(qs_env, ewald_env=ewald_env, ewald_pw=ewald_pw)
1077 END IF
1078 ELSEIF (dft_control%qs_control%method_id == do_method_xtb) THEN
1079 ! Read xTB parameter file
1080 xtb_control => dft_control%qs_control%xtb_control
1081 CALL get_qs_env(qs_env, nkind=nkind)
1082 IF (xtb_control%do_tblite) THEN
1083 ! put geometry to tblite
1084 CALL tb_init_geometry(qs_env, qs_env%tb_tblite)
1085 ! select tblite method
1086 CALL tb_set_calculator(qs_env%tb_tblite, xtb_control%tblite_method)
1087 !set up wave function
1088 CALL tb_init_wf(qs_env%tb_tblite)
1089 !get basis set
1090 DO ikind = 1, nkind
1091 qs_kind => qs_kind_set(ikind)
1092 ! Setup proper xTB parameters
1093 cpassert(.NOT. ASSOCIATED(qs_kind%xtb_parameter))
1094 CALL allocate_xtb_atom_param(qs_kind%xtb_parameter)
1095 ! Set default parameters
1096 CALL get_qs_kind(qs_kind, element_symbol=element_symbol)
1097
1098 NULLIFY (tmp_basis_set)
1099 CALL tb_get_basis(qs_env%tb_tblite, tmp_basis_set, element_symbol, qs_kind%xtb_parameter, occ)
1100 CALL add_basis_set_to_container(qs_kind%basis_sets, tmp_basis_set, "ORB")
1101 CALL set_xtb_atom_param(qs_kind%xtb_parameter, occupation=occ)
1102
1103 !setting the potential for the computation
1104 zeff_correction = 0.0_dp
1105 CALL init_potential(qs_kind%all_potential, itype="BARE", &
1106 zeff=real(sum(occ), dp), zeff_correction=zeff_correction)
1107 END DO
1108 ELSE
1109 NULLIFY (ewald_env, ewald_pw)
1110 DO ikind = 1, nkind
1111 qs_kind => qs_kind_set(ikind)
1112 ! Setup proper xTB parameters
1113 cpassert(.NOT. ASSOCIATED(qs_kind%xtb_parameter))
1114 CALL allocate_xtb_atom_param(qs_kind%xtb_parameter)
1115 ! Set default parameters
1116 gfn_type = dft_control%qs_control%xtb_control%gfn_type
1117 CALL get_qs_kind(qs_kind, element_symbol=element_symbol)
1118 CALL xtb_parameters_init(qs_kind%xtb_parameter, gfn_type, element_symbol, &
1119 xtb_control%parameter_file_path, xtb_control%parameter_file_name, &
1120 para_env)
1121 ! set dependent parameters
1122 CALL xtb_parameters_set(qs_kind%xtb_parameter)
1123 ! Generate basis set
1124 NULLIFY (tmp_basis_set)
1125 IF (qs_kind%xtb_parameter%z == 1) THEN
1126 ! special case hydrogen
1127 ngauss = xtb_control%h_sto_ng
1128 ELSE
1129 ngauss = xtb_control%sto_ng
1130 END IF
1131 IF (qs_kind%xtb_parameter%defined) THEN
1132 CALL init_xtb_basis(qs_kind%xtb_parameter, tmp_basis_set, ngauss)
1133 CALL add_basis_set_to_container(qs_kind%basis_sets, tmp_basis_set, "ORB")
1134 ELSE
1135 CALL set_qs_kind(qs_kind, ghost=.true.)
1136 IF (ASSOCIATED(qs_kind%all_potential)) THEN
1137 DEALLOCATE (qs_kind%all_potential%elec_conf)
1138 DEALLOCATE (qs_kind%all_potential)
1139 END IF
1140 END IF
1141 ! potential
1142 IF (qs_kind%xtb_parameter%defined) THEN
1143 zeff_correction = 0.0_dp
1144 CALL init_potential(qs_kind%all_potential, itype="BARE", &
1145 zeff=qs_kind%xtb_parameter%zeff, zeff_correction=zeff_correction)
1146 CALL get_potential(qs_kind%all_potential, alpha_core_charge=alpha)
1147 ccore = qs_kind%xtb_parameter%zeff*sqrt((alpha/pi)**3)
1148 CALL set_potential(qs_kind%all_potential, ccore_charge=ccore)
1149 qs_kind%xtb_parameter%zeff = qs_kind%xtb_parameter%zeff - zeff_correction
1150 END IF
1151 END DO
1152 !
1153 ! set repulsive potential range
1154 !
1155 ALLOCATE (xtb_control%rcpair(nkind, nkind))
1156 CALL xtb_pp_radius(qs_kind_set, xtb_control%rcpair, xtb_control%eps_pair, xtb_control%kf)
1157 ! check for Ewald
1158 IF (xtb_control%do_ewald) THEN
1159 ALLOCATE (ewald_env)
1160 CALL ewald_env_create(ewald_env, para_env)
1161 poisson_section => section_vals_get_subs_vals(dft_section, "POISSON")
1162 CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
1163 ewald_section => section_vals_get_subs_vals(poisson_section, "EWALD")
1164 print_section => section_vals_get_subs_vals(qs_env%input, "PRINT%GRID_INFORMATION")
1165 IF (gfn_type == 0) THEN
1166 CALL read_ewald_section_tb(ewald_env, ewald_section, cell_ref%hmat, &
1167 silent=silent, pset="EEQ", cell_periodic=cell%perd)
1168 ELSE
1169 CALL read_ewald_section_tb(ewald_env, ewald_section, cell_ref%hmat, &
1170 silent=silent, cell_periodic=cell%perd)
1171 END IF
1172 ALLOCATE (ewald_pw)
1173 CALL ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section=print_section)
1174 CALL set_qs_env(qs_env, ewald_env=ewald_env, ewald_pw=ewald_pw)
1175 END IF
1176 END IF
1177 END IF
1178 ! lri or ri env initialization
1179 lri_section => section_vals_get_subs_vals(qs_section, "LRIGPW")
1180 IF (dft_control%qs_control%method_id == do_method_lrigpw .OR. &
1181 dft_control%qs_control%lri_optbas .OR. &
1182 dft_control%qs_control%method_id == do_method_rigpw) THEN
1183 CALL lri_env_init(lri_env, lri_section)
1184 CALL set_qs_env(qs_env, lri_env=lri_env)
1185 END IF
1186
1187 ! Check basis and fill in missing parts
1188 CALL check_qs_kind_set(qs_kind_set, dft_control, subsys_section=subsys_section)
1189
1190 ! Check that no all-electron potential is present if GPW or GAPW_XC
1191 CALL get_qs_kind_set(qs_kind_set, all_potential_present=all_potential_present)
1192 IF ((dft_control%qs_control%method_id == do_method_gpw) .OR. &
1193 (dft_control%qs_control%method_id == do_method_gapw_xc) .OR. &
1194 (dft_control%qs_control%method_id == do_method_ofgpw)) THEN
1195 IF (all_potential_present) THEN
1196 cpabort("All-electron calculations with GPW, GAPW_XC, and OFGPW are not implemented")
1197 END IF
1198 END IF
1199
1200 ! Check that no cneo potential is present if not GAPW
1201 CALL get_qs_kind_set(qs_kind_set, cneo_potential_present=cneo_potential_present)
1202 IF (cneo_potential_present .AND. &
1203 dft_control%qs_control%method_id /= do_method_gapw) THEN
1204 cpabort("CNEO calculations require GAPW method")
1205 END IF
1206
1207 ! DFT+U
1208 CALL get_qs_kind_set(qs_kind_set, dft_plus_u_atom_present=dft_control%dft_plus_u)
1209
1210 IF (dft_control%do_admm) THEN
1211 ! Check if ADMM basis is available
1212 CALL get_qs_env(qs_env, nkind=nkind)
1213 DO ikind = 1, nkind
1214 NULLIFY (aux_fit_basis)
1215 qs_kind => qs_kind_set(ikind)
1216 CALL get_qs_kind(qs_kind, basis_set=aux_fit_basis, basis_type="AUX_FIT")
1217 IF (.NOT. (ASSOCIATED(aux_fit_basis))) THEN
1218 ! AUX_FIT basis set is not available
1219 cpabort("AUX_FIT basis set is not defined. ")
1220 END IF
1221 END DO
1222 END IF
1223
1224 lribas = .false.
1225 e1terms = .false.
1226 IF (dft_control%qs_control%method_id == do_method_lrigpw) THEN
1227 lribas = .true.
1228 CALL get_qs_env(qs_env, lri_env=lri_env)
1229 e1terms = lri_env%exact_1c_terms
1230 END IF
1231 IF (dft_control%qs_control%do_kg) THEN
1232 CALL section_vals_val_get(dft_section, "KG_METHOD%TNADD_METHOD", i_val=tnadd_method)
1233 IF (tnadd_method == kg_tnadd_embed_ri) lribas = .true.
1234 END IF
1235 IF (lribas) THEN
1236 ! Check if LRI_AUX basis is available, auto-generate if needed
1237 CALL get_qs_env(qs_env, nkind=nkind)
1238 DO ikind = 1, nkind
1239 NULLIFY (lri_aux_basis)
1240 qs_kind => qs_kind_set(ikind)
1241 CALL get_qs_kind(qs_kind, basis_set=lri_aux_basis, basis_type="LRI_AUX")
1242 IF (.NOT. (ASSOCIATED(lri_aux_basis))) THEN
1243 ! LRI_AUX basis set is not yet loaded
1244 CALL cp_warn(__location__, "Automatic Generation of LRI_AUX basis. "// &
1245 "This is experimental code.")
1246 ! Generate a default basis
1247 CALL create_lri_aux_basis_set(lri_aux_basis, qs_kind, dft_control%auto_basis_lri_aux, e1terms)
1248 CALL add_basis_set_to_container(qs_kind%basis_sets, lri_aux_basis, "LRI_AUX")
1249 END IF
1250 END DO
1251 END IF
1252
1253 CALL section_vals_val_get(qs_env%input, "DFT%XC%HF%RI%_SECTION_PARAMETERS_", l_val=do_ri_hfx)
1254 CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%RI_RPA%HF%RI%_SECTION_PARAMETERS_", &
1255 l_val=do_rpa_ri_exx)
1256 IF (do_ri_hfx .OR. do_rpa_ri_exx) THEN
1257 CALL get_qs_env(qs_env, nkind=nkind)
1258 CALL section_vals_val_get(qs_env%input, "DFT%SORT_BASIS", i_val=sort_basis)
1259 DO ikind = 1, nkind
1260 NULLIFY (ri_hfx_basis)
1261 qs_kind => qs_kind_set(ikind)
1262 CALL get_qs_kind(qs_kind=qs_kind, basis_set=ri_hfx_basis, &
1263 basis_type="RI_HFX")
1264 IF (.NOT. (ASSOCIATED(ri_hfx_basis))) THEN
1265 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto)
1266 IF (dft_control%do_admm) THEN
1267 CALL create_ri_aux_basis_set(ri_hfx_basis, qs_kind, dft_control%auto_basis_ri_hfx, &
1268 basis_type="AUX_FIT", basis_sort=sort_basis)
1269 ELSE
1270 CALL create_ri_aux_basis_set(ri_hfx_basis, qs_kind, dft_control%auto_basis_ri_hfx, &
1271 basis_sort=sort_basis)
1272 END IF
1273 CALL add_basis_set_to_container(qs_kind%basis_sets, ri_hfx_basis, "RI_HFX")
1274 END IF
1275 END DO
1276 END IF
1277
1278 IF (dft_control%qs_control%method_id == do_method_rigpw) THEN
1279 ! Check if RI_HXC basis is available, auto-generate if needed
1280 CALL get_qs_env(qs_env, nkind=nkind)
1281 DO ikind = 1, nkind
1282 NULLIFY (ri_hfx_basis)
1283 qs_kind => qs_kind_set(ikind)
1284 CALL get_qs_kind(qs_kind, basis_set=ri_hfx_basis, basis_type="RI_HXC")
1285 IF (.NOT. (ASSOCIATED(ri_hfx_basis))) THEN
1286 ! Generate a default basis
1287 CALL create_ri_aux_basis_set(ri_hfx_basis, qs_kind, dft_control%auto_basis_ri_hxc)
1288 CALL add_basis_set_to_container(qs_kind%basis_sets, ri_hfx_basis, "RI_HXC")
1289 END IF
1290 END DO
1291 END IF
1292
1293 ! Harris method
1294 NULLIFY (harris_env)
1295 CALL section_vals_val_get(dft_section, "HARRIS_METHOD%_SECTION_PARAMETERS_", &
1296 l_val=qs_env%harris_method)
1297 harris_section => section_vals_get_subs_vals(dft_section, "HARRIS_METHOD")
1298 CALL harris_env_create(qs_env, harris_env, harris_section)
1299 CALL set_qs_env(qs_env, harris_env=harris_env)
1300 !
1301 IF (qs_env%harris_method) THEN
1302 CALL get_qs_env(qs_env, nkind=nkind)
1303 ! Check if RI_HXC basis is available, auto-generate if needed
1304 DO ikind = 1, nkind
1305 NULLIFY (tmp_basis_set)
1306 qs_kind => qs_kind_set(ikind)
1307 CALL get_qs_kind(qs_kind, basis_set=rhoin_basis, basis_type="RHOIN")
1308 IF (.NOT. (ASSOCIATED(rhoin_basis))) THEN
1309 ! Generate a default basis
1310 CALL create_ri_aux_basis_set(tmp_basis_set, qs_kind, dft_control%auto_basis_ri_hxc)
1311 IF (qs_env%harris_env%density_source == hden_atomic) THEN
1312 CALL create_primitive_basis_set(tmp_basis_set, rhoin_basis, lmax=0)
1313 CALL deallocate_gto_basis_set(tmp_basis_set)
1314 ELSE
1315 rhoin_basis => tmp_basis_set
1316 END IF
1317 CALL add_basis_set_to_container(qs_kind%basis_sets, rhoin_basis, "RHOIN")
1318 END IF
1319 END DO
1320 END IF
1321
1322 mp2_section => section_vals_get_subs_vals(qs_env%input, "DFT%XC%WF_CORRELATION")
1323 CALL section_vals_get(mp2_section, explicit=mp2_present)
1324 IF (mp2_present) THEN
1325
1326 ! basis should be sorted for imaginary time RPA/GW
1327 CALL section_vals_val_get(qs_env%input, "DFT%SORT_BASIS", i_val=sort_basis)
1328 CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%LOW_SCALING%_SECTION_PARAMETERS_", &
1329 l_val=do_wfc_im_time)
1330
1331 IF (do_wfc_im_time .AND. sort_basis /= basis_sort_zet) THEN
1332 CALL cp_warn(__location__, &
1333 "Low-scaling RPA requires SORT_BASIS EXP keyword (in DFT input section) for good performance")
1334 END IF
1335
1336 ! Check if RI_AUX basis (for MP2/RPA) is given, auto-generate if not
1337 CALL mp2_env_create(qs_env%mp2_env)
1338 CALL get_qs_env(qs_env, mp2_env=mp2_env, nkind=nkind)
1339 CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%RI_MP2%_SECTION_PARAMETERS_", l_val=do_ri_mp2)
1340 CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%RI_SOS_MP2%_SECTION_PARAMETERS_", l_val=do_ri_sos_mp2)
1341 CALL section_vals_val_get(qs_env%input, "DFT%XC%WF_CORRELATION%RI_RPA%_SECTION_PARAMETERS_", l_val=do_ri_rpa)
1342 IF (do_ri_mp2 .OR. do_ri_sos_mp2 .OR. do_ri_rpa) THEN
1343 DO ikind = 1, nkind
1344 NULLIFY (ri_aux_basis_set)
1345 qs_kind => qs_kind_set(ikind)
1346 CALL get_qs_kind(qs_kind=qs_kind, basis_set=ri_aux_basis_set, &
1347 basis_type="RI_AUX")
1348 IF (.NOT. (ASSOCIATED(ri_aux_basis_set))) THEN
1349 ! RI_AUX basis set is not yet loaded
1350 ! Generate a default basis
1351 CALL create_ri_aux_basis_set(ri_aux_basis_set, qs_kind, dft_control%auto_basis_ri_aux, basis_sort=sort_basis)
1352 CALL add_basis_set_to_container(qs_kind%basis_sets, ri_aux_basis_set, "RI_AUX")
1353 ! Add a flag, which allows to check if the basis was generated
1354 ! when applying ERI_METHOD OS to mp2, ri-rpa, gw etc
1355 qs_env%mp2_env%ri_aux_auto_generated = .true.
1356 END IF
1357 END DO
1358 END IF
1359
1360 END IF
1361
1362 IF (dft_control%do_xas_tdp_calculation .OR. qs_env%do_rixs) THEN
1363 ! Check if RI_XAS basis is given, auto-generate if not
1364 CALL get_qs_env(qs_env, nkind=nkind)
1365 DO ikind = 1, nkind
1366 NULLIFY (ri_xas_basis)
1367 qs_kind => qs_kind_set(ikind)
1368 CALL get_qs_kind(qs_kind, basis_set=ri_xas_basis, basis_type="RI_XAS")
1369 IF (.NOT. ASSOCIATED(ri_xas_basis)) THEN
1370 ! Generate a default basis
1371 CALL create_ri_aux_basis_set(ri_xas_basis, qs_kind, dft_control%auto_basis_ri_xas)
1372 CALL add_basis_set_to_container(qs_kind%basis_sets, ri_xas_basis, "RI_XAS")
1373 END IF
1374 END DO
1375 END IF
1376
1377 ! Initialize the spherical harmonics and the orbital transformation matrices
1378 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto, maxlppl=maxlppl, maxlppnl=maxlppnl)
1379
1380 ! CNEO nuclear basis contributes to GAPW rho0
1381 IF (cneo_potential_present) THEN
1382 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto_nuc, basis_type="NUC")
1383 maxlgto = max(maxlgto, maxlgto_nuc)
1384 END IF
1385 lmax_sphere = dft_control%qs_control%gapw_control%lmax_sphere
1386 IF (lmax_sphere < 0) THEN
1387 lmax_sphere = 2*maxlgto
1388 dft_control%qs_control%gapw_control%lmax_sphere = lmax_sphere
1389 END IF
1390 IF (dft_control%qs_control%method_id == do_method_lrigpw .OR. dft_control%qs_control%lri_optbas) THEN
1391 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto_lri, basis_type="LRI_AUX")
1392 !take maxlgto from lri basis if larger (usually)
1393 maxlgto = max(maxlgto, maxlgto_lri)
1394 ELSE IF (dft_control%qs_control%method_id == do_method_rigpw) THEN
1395 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto_lri, basis_type="RI_HXC")
1396 maxlgto = max(maxlgto, maxlgto_lri)
1397 END IF
1398 IF (dft_control%do_xas_tdp_calculation .OR. qs_env%do_rixs) THEN
1399 !done as a precaution
1400 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto_lri, basis_type="RI_XAS")
1401 maxlgto = max(maxlgto, maxlgto_lri)
1402 END IF
1403 maxl = max(2*maxlgto, maxlppl, maxlppnl, lmax_sphere) + 1
1404
1405 CALL init_orbital_pointers(maxl)
1406
1407 CALL init_spherical_harmonics(maxl, 0)
1408
1409 ! Initialise the qs_kind_set
1410 CALL init_qs_kind_set(qs_kind_set)
1411
1412 ! Initialise GAPW soft basis and projectors
1413 IF (dft_control%qs_control%method_id == do_method_gapw .OR. &
1414 dft_control%qs_control%method_id == do_method_gapw_xc) THEN
1415 qs_control => dft_control%qs_control
1416 CALL init_gapw_basis_set(qs_kind_set, qs_control, qs_env%input)
1417 END IF
1418
1419 ! Initialise CNEO nuclear soft basis
1420 IF (cneo_potential_present) THEN
1421 CALL init_cneo_basis_set(qs_kind_set, qs_control)
1422 END IF
1423
1424 ! Initialize the pretabulation for the calculation of the
1425 ! incomplete Gamma function F_n(t) after McMurchie-Davidson
1426 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto)
1427 maxl = max(3*maxlgto + 1, 0)
1428 CALL init_md_ftable(maxl)
1429
1430 ! Initialize the atomic interaction radii
1431 CALL init_interaction_radii(dft_control%qs_control, qs_kind_set)
1432 !
1433 IF (dft_control%qs_control%method_id == do_method_xtb) THEN
1434 IF (.NOT. dft_control%qs_control%xtb_control%do_tblite) THEN
1435 ! cutoff radius
1436 CALL get_qs_env(qs_env, nkind=nkind)
1437 DO ikind = 1, nkind
1438 qs_kind => qs_kind_set(ikind)
1439 IF (qs_kind%xtb_parameter%defined) THEN
1440 CALL get_qs_kind(qs_kind, basis_set=tmp_basis_set)
1441 rcut = xtb_control%coulomb_sr_cut
1442 fxx = 2.0_dp*xtb_control%coulomb_sr_eps*qs_kind%xtb_parameter%eta**2
1443 fxx = 0.80_dp*(1.0_dp/fxx)**0.3333_dp
1444 rcut = min(rcut, xtb_control%coulomb_sr_cut)
1445 qs_kind%xtb_parameter%rcut = min(rcut, fxx)
1446 ELSE
1447 qs_kind%xtb_parameter%rcut = 0.0_dp
1448 END IF
1449 END DO
1450 END IF
1451 END IF
1452
1453 IF (.NOT. be_silent) THEN
1454 CALL write_pgf_orb_radii("orb", atomic_kind_set, qs_kind_set, subsys_section)
1455 CALL write_pgf_orb_radii("aux", atomic_kind_set, qs_kind_set, subsys_section)
1456 CALL write_pgf_orb_radii("lri", atomic_kind_set, qs_kind_set, subsys_section)
1457 CALL write_pgf_orb_radii("nuc", atomic_kind_set, qs_kind_set, subsys_section)
1458 CALL write_core_charge_radii(atomic_kind_set, qs_kind_set, subsys_section)
1459 CALL write_ppl_radii(atomic_kind_set, qs_kind_set, subsys_section)
1460 CALL write_ppnl_radii(atomic_kind_set, qs_kind_set, subsys_section)
1461 CALL write_paw_radii(atomic_kind_set, qs_kind_set, subsys_section)
1462 END IF
1463
1464 ! Distribute molecules and atoms using the new data structures
1465 CALL distribute_molecules_1d(atomic_kind_set=atomic_kind_set, &
1466 particle_set=particle_set, &
1467 local_particles=local_particles, &
1468 molecule_kind_set=molecule_kind_set, &
1469 molecule_set=molecule_set, &
1470 local_molecules=local_molecules, &
1471 force_env_section=qs_env%input)
1472
1473 ! SCF parameters
1474 ALLOCATE (scf_control)
1475 ! set (non)-self consistency
1476 IF (dft_control%qs_control%dftb) THEN
1477 scf_control%non_selfconsistent = .NOT. dft_control%qs_control%dftb_control%self_consistent
1478 END IF
1479 IF (dft_control%qs_control%xtb) THEN
1480 IF (dft_control%qs_control%xtb_control%do_tblite) THEN
1481 scf_control%non_selfconsistent = .false.
1482 ELSE
1483 scf_control%non_selfconsistent = (dft_control%qs_control%xtb_control%gfn_type == 0)
1484 END IF
1485 END IF
1486 IF (qs_env%harris_method) THEN
1487 scf_control%non_selfconsistent = .true.
1488 END IF
1489 CALL scf_c_create(scf_control)
1490 CALL scf_c_read_parameters(scf_control, dft_section)
1491
1492 ! Allocate the data structure for Quickstep energies
1493 CALL allocate_qs_energy(energy)
1494
1495 ! Check for orthogonal basis
1496 has_unit_metric = .false.
1497 IF (dft_control%qs_control%semi_empirical) THEN
1498 IF (dft_control%qs_control%se_control%orthogonal_basis) has_unit_metric = .true.
1499 END IF
1500 IF (dft_control%qs_control%dftb) THEN
1501 IF (dft_control%qs_control%dftb_control%orthogonal_basis) has_unit_metric = .true.
1502 END IF
1503 CALL set_qs_env(qs_env, has_unit_metric=has_unit_metric)
1504
1505 ! Activate the interpolation
1506 CALL wfi_create(wf_history, &
1507 interpolation_method_nr= &
1508 dft_control%qs_control%wf_interpolation_method_nr, &
1509 extrapolation_order=dft_control%qs_control%wf_extrapolation_order, &
1510 has_unit_metric=has_unit_metric)
1511
1512 ! Set the current Quickstep environment
1513 CALL set_qs_env(qs_env=qs_env, &
1514 scf_control=scf_control, &
1515 wf_history=wf_history)
1516
1517 CALL qs_subsys_set(subsys, &
1518 cell_ref=cell_ref, &
1519 use_ref_cell=use_ref_cell, &
1520 energy=energy, &
1521 force=force)
1522
1523 CALL get_qs_env(qs_env, ks_env=ks_env)
1524 CALL set_ks_env(ks_env, dft_control=dft_control)
1525
1526 CALL qs_subsys_set(subsys, local_molecules=local_molecules, &
1527 local_particles=local_particles, cell=cell)
1528
1529 CALL distribution_1d_release(local_particles)
1530 CALL distribution_1d_release(local_molecules)
1531 CALL wfi_release(wf_history)
1532
1533 CALL get_qs_env(qs_env=qs_env, &
1534 atomic_kind_set=atomic_kind_set, &
1535 dft_control=dft_control, &
1536 scf_control=scf_control)
1537
1538 ! Decide what conditions need mo_derivs
1539 ! right now, this only appears to be OT
1540 IF (dft_control%qs_control%do_ls_scf .OR. &
1541 dft_control%qs_control%do_almo_scf) THEN
1542 CALL set_qs_env(qs_env=qs_env, requires_mo_derivs=.false.)
1543 ELSE
1544 IF (scf_control%use_ot) THEN
1545 CALL set_qs_env(qs_env=qs_env, requires_mo_derivs=.true.)
1546 ELSE
1547 CALL set_qs_env(qs_env=qs_env, requires_mo_derivs=.false.)
1548 END IF
1549 END IF
1550
1551 ! XXXXXXX this is backwards XXXXXXXX
1552 IF (dft_control%qs_control%xtb_control%do_tblite) THEN
1553 IF (.NOT. scf_control%smear%do_smear) THEN
1554 ! set tblite default smearing
1555 scf_control%smear%do_smear = .true.
1556 scf_control%smear%method = smear_fermi_dirac
1557 scf_control%smear%electronic_temperature = 300._dp/kelvin
1558 scf_control%smear%eps_fermi_dirac = 1.e-6_dp
1559 END IF
1560 END IF
1561 dft_control%smear = scf_control%smear%do_smear
1562
1563 ! Periodic efield needs equal occupation and orbital gradients
1564 IF (.NOT. (dft_control%qs_control%dftb .OR. dft_control%qs_control%xtb)) THEN
1565 IF (dft_control%apply_period_efield) THEN
1566 CALL get_qs_env(qs_env=qs_env, requires_mo_derivs=orb_gradient)
1567 IF (.NOT. orb_gradient) THEN
1568 CALL cp_abort(__location__, "Periodic Efield needs orbital gradient and direct optimization."// &
1569 " Use the OT optimization method.")
1570 END IF
1571 IF (dft_control%smear) THEN
1572 CALL cp_abort(__location__, "Periodic Efield needs equal occupation numbers."// &
1573 " Smearing option is not possible.")
1574 END IF
1575 END IF
1576 END IF
1577
1578 ! Initialize the GAPW local densities and potentials
1579 IF (dft_control%qs_control%method_id == do_method_gapw .OR. &
1580 dft_control%qs_control%method_id == do_method_gapw_xc) THEN
1581 ! Allocate and initialize the set of atomic densities
1582 NULLIFY (rho_atom_set)
1583 gapw_control => dft_control%qs_control%gapw_control
1584 CALL init_rho_atom(rho_atom_set, atomic_kind_set, qs_kind_set, dft_control, para_env)
1585 CALL set_qs_env(qs_env=qs_env, rho_atom_set=rho_atom_set)
1586 IF (dft_control%qs_control%method_id /= do_method_gapw_xc) THEN
1587 CALL get_qs_env(qs_env=qs_env, local_rho_set=local_rho_set, natom=natom)
1588 ! Allocate and initialize the compensation density rho0
1589 CALL init_rho0(local_rho_set, qs_env, gapw_control)
1590 ! Allocate and Initialize the local coulomb term
1591 CALL init_coulomb_local(qs_env%hartree_local, natom)
1592 END IF
1593 ! NLCC
1594 CALL init_gapw_nlcc(qs_kind_set)
1595 ! Accurate XC integration
1596 IF (gapw_control%accurate_xcint) THEN
1597 cpassert(.NOT. ASSOCIATED(gapw_control%aw))
1598 CALL get_qs_env(qs_env, nkind=nkind)
1599 ALLOCATE (gapw_control%aw(nkind))
1600 alpha = gapw_control%aweights
1601 DO ikind = 1, nkind
1602 qs_kind => qs_kind_set(ikind)
1603 CALL get_qs_kind(qs_kind, hard_radius=rc, paw_atom=paw_atom)
1604 IF (paw_atom) THEN
1605 gapw_control%aw(ikind) = alpha*(1.2_dp/rc)**2
1606 ELSE
1607 gapw_control%aw(ikind) = 0.0_dp
1608 END IF
1609 END DO
1610 END IF
1611 ELSE IF (dft_control%qs_control%method_id == do_method_lrigpw) THEN
1612 ! allocate local ri environment
1613 ! nothing to do here?
1614 ELSE IF (dft_control%qs_control%method_id == do_method_rigpw) THEN
1615 ! allocate ri environment
1616 ! nothing to do here?
1617 ELSE IF (dft_control%qs_control%semi_empirical) THEN
1618 NULLIFY (se_store_int_env, se_nddo_mpole, se_nonbond_env)
1619 natom = SIZE(particle_set)
1620 se_section => section_vals_get_subs_vals(qs_section, "SE")
1621 se_control => dft_control%qs_control%se_control
1622
1623 ! Make the cutoff radii choice a bit smarter
1624 CALL se_cutoff_compatible(se_control, se_section, cell, output_unit)
1625
1626 SELECT CASE (dft_control%qs_control%method_id)
1627 CASE DEFAULT
1630 ! Neighbor lists have to be MAX(interaction range, orbital range)
1631 ! set new kind radius
1632 CALL init_se_nlradius(se_control, atomic_kind_set, qs_kind_set, subsys_section)
1633 END SELECT
1634 ! Initialize to zero the max multipole to treat in the EWALD scheme..
1635 se_control%max_multipole = do_multipole_none
1636 ! check for Ewald
1637 IF (se_control%do_ewald .OR. se_control%do_ewald_gks) THEN
1638 ALLOCATE (ewald_env)
1639 CALL ewald_env_create(ewald_env, para_env)
1640 poisson_section => section_vals_get_subs_vals(dft_section, "POISSON")
1641 CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
1642 ewald_section => section_vals_get_subs_vals(poisson_section, "EWALD")
1643 print_section => section_vals_get_subs_vals(qs_env%input, &
1644 "PRINT%GRID_INFORMATION")
1645 CALL read_ewald_section(ewald_env, ewald_section)
1646 ! Create ewald grids
1647 ALLOCATE (ewald_pw)
1648 CALL ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, &
1649 print_section=print_section)
1650 ! Initialize ewald grids
1651 CALL ewald_pw_grid_update(ewald_pw, ewald_env, cell%hmat)
1652 ! Setup the nonbond environment (real space part of Ewald)
1653 CALL ewald_env_get(ewald_env, rcut=ewald_rcut)
1654 ! Setup the maximum level of multipoles to be treated in the periodic SE scheme
1655 IF (se_control%do_ewald) THEN
1656 CALL ewald_env_get(ewald_env, max_multipole=se_control%max_multipole)
1657 END IF
1658 CALL section_vals_val_get(se_section, "NEIGHBOR_LISTS%VERLET_SKIN", &
1659 r_val=verlet_skin)
1660 ALLOCATE (se_nonbond_env)
1661 CALL fist_nonbond_env_create(se_nonbond_env, atomic_kind_set, do_nonbonded=.true., &
1662 do_electrostatics=.true., verlet_skin=verlet_skin, ewald_rcut=ewald_rcut, &
1663 ei_scale14=0.0_dp, vdw_scale14=0.0_dp, shift_cutoff=.false.)
1664 ! Create and Setup NDDO multipole environment
1665 CALL nddo_mpole_setup(se_nddo_mpole, natom)
1666 CALL set_qs_env(qs_env, ewald_env=ewald_env, ewald_pw=ewald_pw, &
1667 se_nonbond_env=se_nonbond_env, se_nddo_mpole=se_nddo_mpole)
1668 ! Handle the residual integral part 1/R^3
1669 CALL semi_empirical_expns3_setup(qs_kind_set, se_control, &
1670 dft_control%qs_control%method_id)
1671 END IF
1672 ! Taper function
1673 CALL se_taper_create(se_taper, se_control%integral_screening, se_control%do_ewald, &
1674 se_control%taper_cou, se_control%range_cou, &
1675 se_control%taper_exc, se_control%range_exc, &
1676 se_control%taper_scr, se_control%range_scr, &
1677 se_control%taper_lrc, se_control%range_lrc)
1678 CALL set_qs_env(qs_env, se_taper=se_taper)
1679 ! Store integral environment
1680 CALL semi_empirical_si_create(se_store_int_env, se_section)
1681 CALL set_qs_env(qs_env, se_store_int_env=se_store_int_env)
1682 END IF
1683
1684 ! Initialize possible dispersion parameters
1685 IF (dft_control%qs_control%method_id == do_method_gpw .OR. &
1686 dft_control%qs_control%method_id == do_method_gapw .OR. &
1687 dft_control%qs_control%method_id == do_method_gapw_xc .OR. &
1688 dft_control%qs_control%method_id == do_method_lrigpw .OR. &
1689 dft_control%qs_control%method_id == do_method_rigpw .OR. &
1690 dft_control%qs_control%method_id == do_method_ofgpw) THEN
1691 ALLOCATE (dispersion_env)
1692 NULLIFY (xc_section)
1693 xc_section => section_vals_get_subs_vals(dft_section, "XC")
1694 CALL qs_dispersion_env_set(dispersion_env, xc_section)
1695 IF (dispersion_env%type == xc_vdw_fun_pairpot) THEN
1696 NULLIFY (pp_section)
1697 pp_section => section_vals_get_subs_vals(xc_section, "VDW_POTENTIAL%PAIR_POTENTIAL")
1698 CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, pp_section, para_env)
1699 ELSE IF (dispersion_env%type == xc_vdw_fun_nonloc) THEN
1700 NULLIFY (nl_section)
1701 nl_section => section_vals_get_subs_vals(xc_section, "VDW_POTENTIAL%NON_LOCAL")
1702 CALL qs_dispersion_nonloc_init(dispersion_env, para_env)
1703 END IF
1704 CALL set_qs_env(qs_env, dispersion_env=dispersion_env)
1705 ELSE IF (dft_control%qs_control%method_id == do_method_dftb) THEN
1706 ALLOCATE (dispersion_env)
1707 ! set general defaults
1708 dispersion_env%doabc = .false.
1709 dispersion_env%c9cnst = .false.
1710 dispersion_env%lrc = .false.
1711 dispersion_env%srb = .false.
1712 dispersion_env%verbose = .false.
1713 NULLIFY (dispersion_env%c6ab, dispersion_env%maxci, dispersion_env%r0ab, dispersion_env%rcov, &
1714 dispersion_env%r2r4, dispersion_env%cn, dispersion_env%cnkind, dispersion_env%cnlist, &
1715 dispersion_env%d3_exclude_pair)
1716 NULLIFY (dispersion_env%q_mesh, dispersion_env%kernel, dispersion_env%d2phi_dk2, &
1717 dispersion_env%d2y_dx2, dispersion_env%dftd_section)
1718 NULLIFY (dispersion_env%sab_vdw, dispersion_env%sab_cn)
1719 IF (dftb_control%dispersion .AND. dftb_control%dispersion_type == dispersion_d3) THEN
1720 dispersion_env%type = xc_vdw_fun_pairpot
1721 dispersion_env%pp_type = vdw_pairpot_dftd3
1722 dispersion_env%eps_cn = dftb_control%epscn
1723 dispersion_env%s6 = dftb_control%sd3(1)
1724 dispersion_env%sr6 = dftb_control%sd3(2)
1725 dispersion_env%s8 = dftb_control%sd3(3)
1726 dispersion_env%domol = .false.
1727 dispersion_env%kgc8 = 0._dp
1728 dispersion_env%rc_disp = dftb_control%rcdisp
1729 dispersion_env%exp_pre = 0._dp
1730 dispersion_env%scaling = 0._dp
1731 dispersion_env%nd3_exclude_pair = 0
1732 dispersion_env%parameter_file_name = dftb_control%dispersion_parameter_file
1733 CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, para_env=para_env)
1734 ELSEIF (dftb_control%dispersion .AND. dftb_control%dispersion_type == dispersion_d3bj) THEN
1735 dispersion_env%type = xc_vdw_fun_pairpot
1736 dispersion_env%pp_type = vdw_pairpot_dftd3bj
1737 dispersion_env%eps_cn = dftb_control%epscn
1738 dispersion_env%s6 = dftb_control%sd3bj(1)
1739 dispersion_env%a1 = dftb_control%sd3bj(2)
1740 dispersion_env%s8 = dftb_control%sd3bj(3)
1741 dispersion_env%a2 = dftb_control%sd3bj(4)
1742 dispersion_env%domol = .false.
1743 dispersion_env%kgc8 = 0._dp
1744 dispersion_env%rc_disp = dftb_control%rcdisp
1745 dispersion_env%exp_pre = 0._dp
1746 dispersion_env%scaling = 0._dp
1747 dispersion_env%nd3_exclude_pair = 0
1748 dispersion_env%parameter_file_name = dftb_control%dispersion_parameter_file
1749 CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, para_env=para_env)
1750 ELSEIF (dftb_control%dispersion .AND. dftb_control%dispersion_type == dispersion_d2) THEN
1751 dispersion_env%type = xc_vdw_fun_pairpot
1752 dispersion_env%pp_type = vdw_pairpot_dftd2
1753 dispersion_env%exp_pre = dftb_control%exp_pre
1754 dispersion_env%scaling = dftb_control%scaling
1755 dispersion_env%parameter_file_name = dftb_control%dispersion_parameter_file
1756 dispersion_env%rc_disp = dftb_control%rcdisp
1757 CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, para_env=para_env)
1758 ELSE
1759 dispersion_env%type = xc_vdw_fun_none
1760 END IF
1761 CALL set_qs_env(qs_env, dispersion_env=dispersion_env)
1762 ELSE IF (dft_control%qs_control%method_id == do_method_xtb) THEN
1763 IF (.NOT. (dft_control%qs_control%xtb_control%do_tblite)) THEN
1764 ALLOCATE (dispersion_env)
1765 ! set general defaults
1766 dispersion_env%doabc = .false.
1767 dispersion_env%c9cnst = .false.
1768 dispersion_env%lrc = .false.
1769 dispersion_env%srb = .false.
1770 dispersion_env%verbose = .false.
1771 NULLIFY (dispersion_env%c6ab, dispersion_env%maxci, &
1772 dispersion_env%r0ab, dispersion_env%rcov, &
1773 dispersion_env%r2r4, dispersion_env%cn, &
1774 dispersion_env%cnkind, dispersion_env%cnlist, &
1775 dispersion_env%d3_exclude_pair)
1776 NULLIFY (dispersion_env%q_mesh, dispersion_env%kernel, dispersion_env%d2phi_dk2, &
1777 dispersion_env%d2y_dx2, dispersion_env%dftd_section)
1778 NULLIFY (dispersion_env%sab_vdw, dispersion_env%sab_cn)
1779 dispersion_env%type = xc_vdw_fun_pairpot
1780 dispersion_env%eps_cn = xtb_control%epscn
1781 dispersion_env%s6 = xtb_control%s6
1782 dispersion_env%s8 = xtb_control%s8
1783 dispersion_env%a1 = xtb_control%a1
1784 dispersion_env%a2 = xtb_control%a2
1785 dispersion_env%domol = .false.
1786 dispersion_env%kgc8 = 0._dp
1787 dispersion_env%rc_disp = xtb_control%rcdisp
1788 dispersion_env%rc_d4 = xtb_control%rcdisp
1789 dispersion_env%exp_pre = 0._dp
1790 dispersion_env%scaling = 0._dp
1791 dispersion_env%nd3_exclude_pair = 0
1792 dispersion_env%parameter_file_name = xtb_control%dispersion_parameter_file
1793 !
1794 SELECT CASE (xtb_control%vdw_type)
1796 dispersion_env%pp_type = vdw_pairpot_dftd3bj
1797 CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, para_env=para_env)
1798 IF (xtb_control%vdw_type == xtb_vdw_type_none) dispersion_env%type = xc_vdw_fun_none
1799 CASE (xtb_vdw_type_d4)
1800 dispersion_env%pp_type = vdw_pairpot_dftd4
1801 dispersion_env%ref_functional = "none"
1802 CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, &
1803 dispersion_env, para_env=para_env)
1804 dispersion_env%cnfun = 2
1805 CASE DEFAULT
1806 cpabort("vdw type")
1807 END SELECT
1808 CALL set_qs_env(qs_env, dispersion_env=dispersion_env)
1809 END IF
1810 ELSE IF (dft_control%qs_control%semi_empirical) THEN
1811 ALLOCATE (dispersion_env)
1812 ! set general defaults
1813 dispersion_env%doabc = .false.
1814 dispersion_env%c9cnst = .false.
1815 dispersion_env%lrc = .false.
1816 dispersion_env%srb = .false.
1817 dispersion_env%verbose = .false.
1818 NULLIFY (dispersion_env%c6ab, dispersion_env%maxci, dispersion_env%r0ab, dispersion_env%rcov, &
1819 dispersion_env%r2r4, dispersion_env%cn, dispersion_env%cnkind, dispersion_env%cnlist, &
1820 dispersion_env%d3_exclude_pair)
1821 NULLIFY (dispersion_env%q_mesh, dispersion_env%kernel, dispersion_env%d2phi_dk2, &
1822 dispersion_env%d2y_dx2, dispersion_env%dftd_section)
1823 NULLIFY (dispersion_env%sab_vdw, dispersion_env%sab_cn)
1824 IF (se_control%dispersion) THEN
1825 dispersion_env%type = xc_vdw_fun_pairpot
1826 dispersion_env%pp_type = vdw_pairpot_dftd3
1827 dispersion_env%eps_cn = se_control%epscn
1828 dispersion_env%s6 = se_control%sd3(1)
1829 dispersion_env%sr6 = se_control%sd3(2)
1830 dispersion_env%s8 = se_control%sd3(3)
1831 dispersion_env%domol = .false.
1832 dispersion_env%kgc8 = 0._dp
1833 dispersion_env%rc_disp = se_control%rcdisp
1834 dispersion_env%exp_pre = 0._dp
1835 dispersion_env%scaling = 0._dp
1836 dispersion_env%nd3_exclude_pair = 0
1837 dispersion_env%parameter_file_name = se_control%dispersion_parameter_file
1838 CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, para_env=para_env)
1839 ELSE
1840 dispersion_env%type = xc_vdw_fun_none
1841 END IF
1842 CALL set_qs_env(qs_env, dispersion_env=dispersion_env)
1843 END IF
1844
1845 ! Initialize possible geomertical counterpoise correction potential
1846 IF (dft_control%qs_control%method_id == do_method_gpw .OR. &
1847 dft_control%qs_control%method_id == do_method_gapw .OR. &
1848 dft_control%qs_control%method_id == do_method_gapw_xc .OR. &
1849 dft_control%qs_control%method_id == do_method_lrigpw .OR. &
1850 dft_control%qs_control%method_id == do_method_rigpw .OR. &
1851 dft_control%qs_control%method_id == do_method_ofgpw) THEN
1852 ALLOCATE (gcp_env)
1853 NULLIFY (xc_section)
1854 xc_section => section_vals_get_subs_vals(dft_section, "XC")
1855 CALL qs_gcp_env_set(gcp_env, xc_section)
1856 CALL qs_gcp_init(qs_env, gcp_env)
1857 CALL set_qs_env(qs_env, gcp_env=gcp_env)
1858 END IF
1859
1860 ! Allocate the MO data types
1861 CALL get_qs_kind_set(qs_kind_set, nsgf=n_ao, nelectron=nelectron)
1862
1863 ! The total number of electrons
1864 nelectron = nelectron - dft_control%charge
1865
1866 IF (dft_control%multiplicity == 0) THEN
1867 IF (modulo(nelectron, 2) == 0) THEN
1868 dft_control%multiplicity = 1
1869 ELSE
1870 dft_control%multiplicity = 2
1871 END IF
1872 END IF
1873
1874 multiplicity = dft_control%multiplicity
1875
1876 IF ((dft_control%nspins < 1) .OR. (dft_control%nspins > 2)) THEN
1877 cpabort("nspins should be 1 or 2 for the time being ...")
1878 END IF
1879
1880 IF ((modulo(nelectron, 2) /= 0) .AND. (dft_control%nspins == 1)) THEN
1881 IF (.NOT. dft_control%qs_control%ofgpw .AND. .NOT. dft_control%smear) THEN
1882 cpabort("Use the LSD option for an odd number of electrons")
1883 END IF
1884 END IF
1885
1886 ! The transition potential method to calculate XAS needs LSD
1887 IF (dft_control%do_xas_calculation) THEN
1888 IF (dft_control%nspins == 1) THEN
1889 cpabort("Use the LSD option for XAS with transition potential")
1890 END IF
1891 END IF
1892
1893 ! assigning the number of states per spin initial version, not yet very
1894 ! general. Should work for an even number of electrons and a single
1895 ! additional electron this set of options that requires full matrices,
1896 ! however, makes things a bit ugly right now.... we try to make a
1897 ! distinction between the number of electrons per spin and the number of
1898 ! MOs per spin this should allow the use of fractional occupations later on
1899 IF (dft_control%qs_control%ofgpw) THEN
1900
1901 IF (dft_control%nspins == 1) THEN
1902 maxocc = nelectron
1903 nelectron_spin(1) = nelectron
1904 nelectron_spin(2) = 0
1905 n_mo(1) = 1
1906 n_mo(2) = 0
1907 ELSE
1908 IF (modulo(nelectron + multiplicity - 1, 2) /= 0) THEN
1909 cpabort("LSD: try to use a different multiplicity")
1910 END IF
1911 nelectron_spin(1) = (nelectron + multiplicity - 1)/2
1912 nelectron_spin(2) = (nelectron - multiplicity + 1)/2
1913 IF (nelectron_spin(1) < 0) THEN
1914 cpabort("LSD: too few electrons for this multiplicity")
1915 END IF
1916 maxocc = maxval(nelectron_spin)
1917 n_mo(1) = min(nelectron_spin(1), 1)
1918 n_mo(2) = min(nelectron_spin(2), 1)
1919 END IF
1920
1921 ELSE
1922
1923 IF (dft_control%nspins == 1) THEN
1924 maxocc = 2.0_dp
1925 nelectron_spin(1) = nelectron
1926 nelectron_spin(2) = 0
1927 IF (modulo(nelectron, 2) == 0) THEN
1928 n_mo(1) = nelectron/2
1929 ELSE
1930 n_mo(1) = int(nelectron/2._dp) + 1
1931 END IF
1932 n_mo(2) = 0
1933 ELSE
1934 maxocc = 1.0_dp
1935
1936 ! The simplist spin distribution is written here. Special cases will
1937 ! need additional user input
1938 IF (modulo(nelectron + multiplicity - 1, 2) /= 0) THEN
1939 cpabort("LSD: try to use a different multiplicity")
1940 END IF
1941
1942 nelectron_spin(1) = (nelectron + multiplicity - 1)/2
1943 nelectron_spin(2) = (nelectron - multiplicity + 1)/2
1944
1945 IF (nelectron_spin(2) < 0) THEN
1946 cpabort("LSD: too few electrons for this multiplicity")
1947 END IF
1948
1949 n_mo(1) = nelectron_spin(1)
1950 n_mo(2) = nelectron_spin(2)
1951
1952 END IF
1953
1954 END IF
1955
1956 ! Read the total_zeff_corr here [SGh]
1957 CALL get_qs_kind_set(qs_kind_set, total_zeff_corr=total_zeff_corr)
1958 ! store it in qs_env
1959 qs_env%total_zeff_corr = total_zeff_corr
1960
1961 ! Store the number of electrons once and for all
1962 CALL qs_subsys_set(subsys, &
1963 nelectron_total=nelectron, &
1964 nelectron_spin=nelectron_spin)
1965
1966 ! Ensure that all orbitals requested for printout are added even
1967 ! if the keyword ADDED_MOS was not specified or set properly
1968 mo_index_range => section_get_ivals(dft_section, "PRINT%MO%MO_INDEX_RANGE")
1969 cpassert(ASSOCIATED(mo_index_range))
1970 IF (all(mo_index_range > 0)) THEN
1971 IF (mo_index_range(1) > mo_index_range(2)) THEN
1972 CALL cp_abort(__location__, &
1973 "The upper orbital index ("// &
1974 trim(adjustl(cp_to_string(mo_index_range(2))))// &
1975 ") of the MO_INDEX_RANGE should be equal or larger "// &
1976 "than the lower orbital index ("// &
1977 trim(adjustl(cp_to_string(mo_index_range(1))))// &
1978 ") for printout.")
1979 END IF
1980 ! Adapt ADDED_MOS automatically if needed for printout
1981 IF (.NOT. scf_control%use_ot) THEN
1982 scf_control%added_mos(1) = min(max(scf_control%added_mos(1), &
1983 mo_index_range(2) - n_mo(1)), &
1984 n_ao - n_mo(1))
1985 IF (dft_control%nspins == 2) THEN
1986 scf_control%added_mos(2) = min(max(scf_control%added_mos(2), &
1987 mo_index_range(2) - n_mo(2)), &
1988 n_ao - n_mo(2))
1989 END IF
1990 END IF
1991 ELSE IF (mo_index_range(2) < 0) THEN
1992 IF (.NOT. scf_control%use_ot) THEN
1993 ! Add all available orbitals
1994 scf_control%added_mos(1) = n_ao - n_mo(1)
1995 IF (dft_control%nspins == 2) THEN
1996 ! Ensure the same number for the spin-down (beta) orbitals
1997 scf_control%added_mos(2) = n_ao - n_mo(2)
1998 END IF
1999 END IF
2000 END IF
2001
2002 IF (dft_control%nspins == 2) THEN
2003 ! Check and set number of added (unoccupied) orbitals for beta spin
2004 IF (scf_control%added_mos(2) < 0) THEN
2005 n_mo_add = n_ao - n_mo(2) ! use all available MOs
2006 ELSE IF (scf_control%added_mos(2) > 0) THEN
2007 n_mo_add = scf_control%added_mos(2)
2008 ELSE
2009 n_mo_add = scf_control%added_mos(1)
2010 END IF
2011 IF (n_mo_add > n_ao - n_mo(2)) THEN
2012 cpwarn("More ADDED_MOs requested for beta spin than available.")
2013 END IF
2014 scf_control%added_mos(2) = min(n_mo_add, n_ao - n_mo(2))
2015 n_mo(2) = n_mo(2) + scf_control%added_mos(2)
2016 END IF
2017
2018 ! proceed alpha orbitals after the beta orbitals; this is essential to avoid
2019 ! reduction in the number of available unoccupied molecular orbitals.
2020 ! E.g. n_ao = 10, nelectrons = 10, multiplicity = 3 implies n_mo(1) = 6, n_mo(2) = 4;
2021 ! added_mos(1:2) = (6,undef) should increase the number of molecular orbitals as
2022 ! n_mo(1) = min(n_ao, n_mo(1) + added_mos(1)) = 10, n_mo(2) = 10.
2023 ! However, if we try to proceed alpha orbitals first, this leads us n_mo(1:2) = (10,8)
2024 ! due to the following assignment instruction above:
2025 ! IF (scf_control%added_mos(2) > 0) THEN ... ELSE; n_mo_add = scf_control%added_mos(1); END IF
2026 IF (dft_control%qs_control%xtb_control%do_tblite) THEN
2027 scf_control%added_mos(1) = n_ao - n_mo(1) ! tblite needs all MO's
2028 ELSE IF (scf_control%added_mos(1) < 0) THEN
2029 scf_control%added_mos(1) = n_ao - n_mo(1) ! use all available MOs
2030 ELSE IF (scf_control%added_mos(1) > n_ao - n_mo(1)) THEN
2031 CALL cp_warn(__location__, &
2032 "More added MOs requested than available. "// &
2033 "The full set of unoccupied MOs will be used. "// &
2034 "Use 'ADDED_MOS -1' to always use all available MOs "// &
2035 "and to get rid of this warning.")
2036 END IF
2037 scf_control%added_mos(1) = min(scf_control%added_mos(1), n_ao - n_mo(1))
2038 n_mo(1) = n_mo(1) + scf_control%added_mos(1)
2039
2040 IF (dft_control%nspins == 2) THEN
2041 IF (n_mo(2) > n_mo(1)) &
2042 CALL cp_warn(__location__, &
2043 "More beta than alpha MOs requested. "// &
2044 "The number of beta MOs will be reduced to the number alpha MOs.")
2045 n_mo(2) = min(n_mo(1), n_mo(2))
2046 cpassert(n_mo(1) >= nelectron_spin(1))
2047 cpassert(n_mo(2) >= nelectron_spin(2))
2048 END IF
2049
2050 ! kpoints
2051 CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints)
2052 IF (do_kpoints .AND. dft_control%nspins == 2) THEN
2053 ! we need equal number of calculated states
2054 IF (n_mo(2) /= n_mo(1)) &
2055 CALL cp_warn(__location__, &
2056 "Kpoints: Different number of MOs requested. "// &
2057 "The number of beta MOs will be set to the number alpha MOs.")
2058 n_mo(2) = n_mo(1)
2059 cpassert(n_mo(1) >= nelectron_spin(1))
2060 cpassert(n_mo(2) >= nelectron_spin(2))
2061 END IF
2062
2063 ! Compatibility checks for smearing
2064 IF (scf_control%smear%do_smear) THEN
2065 IF (scf_control%added_mos(1) == 0) THEN
2066 cpabort("Extra MOs (ADDED_MOS) are required for smearing")
2067 END IF
2068 END IF
2069
2070 ! Some options require that all MOs are computed ...
2071 IF (btest(cp_print_key_should_output(logger%iter_info, dft_section, &
2072 "PRINT%MO/CARTESIAN"), &
2073 cp_p_file) .OR. &
2074 (scf_control%level_shift /= 0.0_dp) .OR. &
2075 (scf_control%diagonalization%eps_jacobi /= 0.0_dp) .OR. &
2076 (dft_control%roks .AND. (.NOT. scf_control%use_ot))) THEN
2077 n_mo(:) = n_ao
2078 END IF
2079
2080 ! Compatibility checks for ROKS
2081 IF (dft_control%roks .AND. (.NOT. scf_control%use_ot)) THEN
2082 IF (scf_control%roks_scheme == general_roks) THEN
2083 cpwarn("General ROKS scheme is not yet tested!")
2084 END IF
2085 IF (scf_control%smear%do_smear) THEN
2086 CALL cp_abort(__location__, &
2087 "The options ROKS and SMEAR are not compatible. "// &
2088 "Try UKS instead of ROKS")
2089 END IF
2090 END IF
2091 IF (dft_control%low_spin_roks) THEN
2092 SELECT CASE (dft_control%qs_control%method_id)
2093 CASE DEFAULT
2095 CALL cp_abort(__location__, &
2096 "xTB/DFTB methods are not compatible with low spin ROKS.")
2099 CALL cp_abort(__location__, &
2100 "SE methods are not compatible with low spin ROKS.")
2101 END SELECT
2102 END IF
2103
2104 ! in principle the restricted calculation could be performed
2105 ! using just one set of MOs and special casing most of the code
2106 ! right now we'll just take care of what is effectively an additional constraint
2107 ! at as few places as possible, just duplicating the beta orbitals
2108 IF (dft_control%restricted .AND. (output_unit > 0)) THEN
2109 ! it is really not yet tested till the end ! Joost
2110 WRITE (output_unit, *) ""
2111 WRITE (output_unit, *) " **************************************"
2112 WRITE (output_unit, *) " restricted calculation cutting corners"
2113 WRITE (output_unit, *) " experimental feature, check code "
2114 WRITE (output_unit, *) " **************************************"
2115 END IF
2116
2117 ! no point in allocating these things here ?
2118 IF (dft_control%qs_control%do_ls_scf) THEN
2119 NULLIFY (mos)
2120 ELSE
2121 ALLOCATE (mos(dft_control%nspins))
2122 DO ispin = 1, dft_control%nspins
2123 CALL allocate_mo_set(mo_set=mos(ispin), &
2124 nao=n_ao, &
2125 nmo=n_mo(ispin), &
2126 nelectron=nelectron_spin(ispin), &
2127 n_el_f=real(nelectron_spin(ispin), dp), &
2128 maxocc=maxocc, &
2129 flexible_electron_count=dft_control%relax_multiplicity)
2130 END DO
2131 END IF
2132
2133 CALL set_qs_env(qs_env, mos=mos)
2134
2135 ! allocate mos when switch_surf_dip is triggered [SGh]
2136 IF (dft_control%switch_surf_dip) THEN
2137 ALLOCATE (mos_last_converged(dft_control%nspins))
2138 DO ispin = 1, dft_control%nspins
2139 CALL allocate_mo_set(mo_set=mos_last_converged(ispin), &
2140 nao=n_ao, &
2141 nmo=n_mo(ispin), &
2142 nelectron=nelectron_spin(ispin), &
2143 n_el_f=real(nelectron_spin(ispin), dp), &
2144 maxocc=maxocc, &
2145 flexible_electron_count=dft_control%relax_multiplicity)
2146 END DO
2147 CALL set_qs_env(qs_env, mos_last_converged=mos_last_converged)
2148 END IF
2149
2150 IF (.NOT. be_silent) THEN
2151 ! Print the DFT control parameters
2152 CALL write_dft_control(dft_control, dft_section)
2153
2154 ! Print the vdW control parameters
2155 IF (dft_control%qs_control%method_id == do_method_gpw .OR. &
2156 dft_control%qs_control%method_id == do_method_gapw .OR. &
2157 dft_control%qs_control%method_id == do_method_gapw_xc .OR. &
2158 dft_control%qs_control%method_id == do_method_lrigpw .OR. &
2159 dft_control%qs_control%method_id == do_method_rigpw .OR. &
2160 dft_control%qs_control%method_id == do_method_dftb .OR. &
2161 (dft_control%qs_control%method_id == do_method_xtb .AND. &
2162 (.NOT. dft_control%qs_control%xtb_control%do_tblite)) .OR. &
2163 dft_control%qs_control%method_id == do_method_ofgpw) THEN
2164 CALL get_qs_env(qs_env, dispersion_env=dispersion_env)
2165 CALL qs_write_dispersion(qs_env, dispersion_env)
2166 END IF
2167
2168 ! Print the Quickstep control parameters
2169 CALL write_qs_control(dft_control%qs_control, dft_section)
2170
2171 ! Print the ADMM control parameters
2172 IF (dft_control%do_admm) THEN
2173 CALL write_admm_control(dft_control%admm_control, dft_section)
2174 END IF
2175
2176 ! Print XES/XAS control parameters
2177 IF (dft_control%do_xas_calculation) THEN
2178 CALL cite_reference(iannuzzi2007)
2179 !CALL write_xas_control(dft_control%xas_control,dft_section)
2180 END IF
2181
2182 ! Print the unnormalized basis set information (input data)
2183 CALL write_gto_basis_sets(qs_kind_set, subsys_section)
2184
2185 ! Print the atomic kind set
2186 CALL write_qs_kind_set(qs_kind_set, subsys_section)
2187
2188 ! Print the molecule kind set
2189 CALL write_molecule_kind_set(molecule_kind_set, subsys_section)
2190
2191 ! Print the total number of kinds, atoms, basis functions etc.
2192 CALL write_total_numbers(qs_kind_set, particle_set, qs_env%input)
2193
2194 ! Print the atomic coordinates
2195 CALL write_qs_particle_coordinates(particle_set, qs_kind_set, subsys_section, label="QUICKSTEP")
2196
2197 ! Print the interatomic distances
2198 CALL write_particle_distances(particle_set, cell, subsys_section)
2199
2200 ! Print the requested structure data
2201 CALL write_structure_data(particle_set, cell, subsys_section)
2202
2203 ! Print symmetry information
2204 CALL write_symmetry(particle_set, cell, subsys_section)
2205
2206 ! Print the SCF parameters
2207 IF ((.NOT. dft_control%qs_control%do_ls_scf) .AND. &
2208 (.NOT. dft_control%qs_control%do_almo_scf)) THEN
2209 CALL scf_c_write_parameters(scf_control, dft_section)
2210 END IF
2211 END IF
2212
2213 ! Sets up pw_env, qs_charges, mpools ...
2214 CALL qs_env_setup(qs_env)
2215
2216 ! Allocate and initialise rho0 soft on the global grid
2217 IF (dft_control%qs_control%method_id == do_method_gapw) THEN
2218 CALL get_qs_env(qs_env=qs_env, pw_env=pw_env, rho0_mpole=rho0_mpole)
2219 CALL rho0_s_grid_create(pw_env, rho0_mpole)
2220 END IF
2221
2222 IF (output_unit > 0) CALL m_flush(output_unit)
2223 CALL timestop(handle)
2224
2225 END SUBROUTINE qs_init_subsys
2226
2227! **************************************************************************************************
2228!> \brief Write the total number of kinds, atoms, etc. to the logical unit
2229!> number lunit.
2230!> \param qs_kind_set ...
2231!> \param particle_set ...
2232!> \param force_env_section ...
2233!> \author Creation (06.10.2000)
2234! **************************************************************************************************
2235 SUBROUTINE write_total_numbers(qs_kind_set, particle_set, force_env_section)
2236
2237 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
2238 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
2239 TYPE(section_vals_type), POINTER :: force_env_section
2240
2241 INTEGER :: maxlgto, maxlppl, maxlppnl, natom, &
2242 natom_q, ncgf, nkind, nkind_q, npgf, &
2243 nset, nsgf, nshell, output_unit
2244 TYPE(cp_logger_type), POINTER :: logger
2245
2246 NULLIFY (logger)
2247 logger => cp_get_default_logger()
2248 output_unit = cp_print_key_unit_nr(logger, force_env_section, "PRINT%TOTAL_NUMBERS", &
2249 extension=".Log")
2250
2251 IF (output_unit > 0) THEN
2252 natom = SIZE(particle_set)
2253 nkind = SIZE(qs_kind_set)
2254
2255 CALL get_qs_kind_set(qs_kind_set, &
2256 maxlgto=maxlgto, &
2257 ncgf=ncgf, &
2258 npgf=npgf, &
2259 nset=nset, &
2260 nsgf=nsgf, &
2261 nshell=nshell, &
2262 maxlppl=maxlppl, &
2263 maxlppnl=maxlppnl)
2264
2265 WRITE (unit=output_unit, fmt="(/,/,T2,A)") &
2266 "TOTAL NUMBERS AND MAXIMUM NUMBERS"
2267
2268 IF (nset + npgf + ncgf > 0) THEN
2269 WRITE (unit=output_unit, fmt="(/,T3,A,(T30,A,T71,I10))") &
2270 "Total number of", &
2271 "- Atomic kinds: ", nkind, &
2272 "- Atoms: ", natom, &
2273 "- Shell sets: ", nset, &
2274 "- Shells: ", nshell, &
2275 "- Primitive Cartesian functions: ", npgf, &
2276 "- Cartesian basis functions: ", ncgf, &
2277 "- Spherical basis functions: ", nsgf
2278 ELSE IF (nshell + nsgf > 0) THEN
2279 WRITE (unit=output_unit, fmt="(/,T3,A,(T30,A,T71,I10))") &
2280 "Total number of", &
2281 "- Atomic kinds: ", nkind, &
2282 "- Atoms: ", natom, &
2283 "- Shells: ", nshell, &
2284 "- Spherical basis functions: ", nsgf
2285 ELSE
2286 WRITE (unit=output_unit, fmt="(/,T3,A,(T30,A,T71,I10))") &
2287 "Total number of", &
2288 "- Atomic kinds: ", nkind, &
2289 "- Atoms: ", natom
2290 END IF
2291
2292 IF ((maxlppl > -1) .AND. (maxlppnl > -1)) THEN
2293 WRITE (unit=output_unit, fmt="(/,T3,A,(T30,A,T75,I6))") &
2294 "Maximum angular momentum of the", &
2295 "- Orbital basis functions: ", maxlgto, &
2296 "- Local part of the GTH pseudopotential: ", maxlppl, &
2297 "- Non-local part of the GTH pseudopotential: ", maxlppnl
2298 ELSEIF (maxlppl > -1) THEN
2299 WRITE (unit=output_unit, fmt="(/,T3,A,(T30,A,T75,I6))") &
2300 "Maximum angular momentum of the", &
2301 "- Orbital basis functions: ", maxlgto, &
2302 "- Local part of the GTH pseudopotential: ", maxlppl
2303 ELSE
2304 WRITE (unit=output_unit, fmt="(/,T3,A,T75,I6)") &
2305 "Maximum angular momentum of the orbital basis functions: ", maxlgto
2306 END IF
2307
2308 ! LRI_AUX BASIS
2309 CALL get_qs_kind_set(qs_kind_set, &
2310 maxlgto=maxlgto, &
2311 ncgf=ncgf, &
2312 npgf=npgf, &
2313 nset=nset, &
2314 nsgf=nsgf, &
2315 nshell=nshell, &
2316 basis_type="LRI_AUX")
2317 IF (nset + npgf + ncgf > 0) THEN
2318 WRITE (unit=output_unit, fmt="(/,T3,A,/,T3,A,(T30,A,T71,I10))") &
2319 "LRI_AUX Basis: ", &
2320 "Total number of", &
2321 "- Shell sets: ", nset, &
2322 "- Shells: ", nshell, &
2323 "- Primitive Cartesian functions: ", npgf, &
2324 "- Cartesian basis functions: ", ncgf, &
2325 "- Spherical basis functions: ", nsgf
2326 WRITE (unit=output_unit, fmt="(T30,A,T75,I6)") &
2327 " Maximum angular momentum ", maxlgto
2328 END IF
2329
2330 ! RI_HXC BASIS
2331 CALL get_qs_kind_set(qs_kind_set, &
2332 maxlgto=maxlgto, &
2333 ncgf=ncgf, &
2334 npgf=npgf, &
2335 nset=nset, &
2336 nsgf=nsgf, &
2337 nshell=nshell, &
2338 basis_type="RI_HXC")
2339 IF (nset + npgf + ncgf > 0) THEN
2340 WRITE (unit=output_unit, fmt="(/,T3,A,/,T3,A,(T30,A,T71,I10))") &
2341 "RI_HXC Basis: ", &
2342 "Total number of", &
2343 "- Shell sets: ", nset, &
2344 "- Shells: ", nshell, &
2345 "- Primitive Cartesian functions: ", npgf, &
2346 "- Cartesian basis functions: ", ncgf, &
2347 "- Spherical basis functions: ", nsgf
2348 WRITE (unit=output_unit, fmt="(T30,A,T75,I6)") &
2349 " Maximum angular momentum ", maxlgto
2350 END IF
2351
2352 ! AUX_FIT BASIS
2353 CALL get_qs_kind_set(qs_kind_set, &
2354 maxlgto=maxlgto, &
2355 ncgf=ncgf, &
2356 npgf=npgf, &
2357 nset=nset, &
2358 nsgf=nsgf, &
2359 nshell=nshell, &
2360 basis_type="AUX_FIT")
2361 IF (nset + npgf + ncgf > 0) THEN
2362 WRITE (unit=output_unit, fmt="(/,T3,A,/,T3,A,(T30,A,T71,I10))") &
2363 "AUX_FIT ADMM-Basis: ", &
2364 "Total number of", &
2365 "- Shell sets: ", nset, &
2366 "- Shells: ", nshell, &
2367 "- Primitive Cartesian functions: ", npgf, &
2368 "- Cartesian basis functions: ", ncgf, &
2369 "- Spherical basis functions: ", nsgf
2370 WRITE (unit=output_unit, fmt="(T30,A,T75,I6)") &
2371 " Maximum angular momentum ", maxlgto
2372 END IF
2373
2374 ! NUCLEAR BASIS
2375 CALL get_qs_kind_set(qs_kind_set, &
2376 nkind_q=nkind_q, &
2377 natom_q=natom_q, &
2378 maxlgto=maxlgto, &
2379 ncgf=ncgf, &
2380 npgf=npgf, &
2381 nset=nset, &
2382 nsgf=nsgf, &
2383 nshell=nshell, &
2384 basis_type="NUC")
2385 IF (nset + npgf + ncgf > 0) THEN
2386 WRITE (unit=output_unit, fmt="(/,T3,A,/,T3,A,(T30,A,T71,I10))") &
2387 "Nuclear Basis: ", &
2388 "Total number of", &
2389 "- Quantum atomic kinds: ", nkind_q, &
2390 "- Quantum atoms: ", natom_q, &
2391 "- Shell sets: ", nset, &
2392 "- Shells: ", nshell, &
2393 "- Primitive Cartesian functions: ", npgf, &
2394 "- Cartesian basis functions: ", ncgf, &
2395 "- Spherical basis functions: ", nsgf
2396 WRITE (unit=output_unit, fmt="(T30,A,T75,I6)") &
2397 " Maximum angular momentum ", maxlgto
2398 END IF
2399
2400 END IF
2401 CALL cp_print_key_finished_output(output_unit, logger, force_env_section, &
2402 "PRINT%TOTAL_NUMBERS")
2403
2404 END SUBROUTINE write_total_numbers
2405
2406END MODULE qs_environment
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
almo_scf_env methods
subroutine, public almo_scf_env_create(qs_env)
Creation and basic initialization of the almo environment.
calculate the orbitals for a given atomic kind type
subroutine, public calculate_atomic_relkin(atomic_kind, qs_kind, rel_control, rtmat)
...
Define the atomic kind types and their sub types.
Automatic generation of auxiliary basis sets of different kind.
Definition auto_basis.F:15
subroutine, public create_lri_aux_basis_set(lri_aux_basis_set, qs_kind, basis_cntrl, exact_1c_terms, tda_kernel)
Create a LRI_AUX basis set using some heuristics.
Definition auto_basis.F:233
subroutine, public create_ri_aux_basis_set(ri_aux_basis_set, qs_kind, basis_cntrl, basis_type, basis_sort)
Create a RI_AUX basis set using some heuristics.
Definition auto_basis.F:57
subroutine, public add_basis_set_to_container(container, basis_set, basis_set_type)
...
integer, parameter, public basis_sort_zet
subroutine, public deallocate_gto_basis_set(gto_basis_set)
...
subroutine, public create_primitive_basis_set(basis_set, pbasis, lmax)
...
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public cp2kqs2020
integer, save, public iannuzzi2007
integer, save, public iannuzzi2006
Handles all functions related to the CELL.
Definition cell_types.F:15
methods related to the blacs parallel environment
subroutine, public cp_blacs_env_release(blacs_env)
releases the given blacs_env
subroutine, public cp_blacs_env_create(blacs_env, para_env, blacs_grid_layout, blacs_repeatable, row_major, grid_2d)
allocates and initializes a type that represent a blacs context
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
Utilities to set up the control types.
subroutine, public write_qs_control(qs_control, dft_section)
Purpose: Write the QS control parameters to the output unit.
subroutine, public read_rixs_control(rixs_control, rixs_section, qs_control)
Reads the input and stores in the rixs_control_type.
subroutine, public read_tddfpt2_control(t_control, t_section, qs_control)
Read TDDFPT-related input parameters.
subroutine, public read_dft_control(dft_control, dft_section)
...
subroutine, public read_qs_section(qs_control, qs_section)
...
subroutine, public write_admm_control(admm_control, dft_section)
Write the ADMM control parameters to the output unit.
subroutine, public write_dft_control(dft_control, dft_section)
Write the DFT control parameters to the output unit.
subroutine, public read_ddapc_section(qs_control, qs_section, ddapc_restraint_section)
reads the input parameters needed for ddapc.
subroutine, public read_mgrid_section(qs_control, dft_section)
...
contains information regarding the decoupling/recoupling method of Bloechl
subroutine, public cp_ddapc_ewald_create(cp_ddapc_ewald, qmmm_decoupl, qm_cell, force_env_section, subsys_section, para_env)
...
various routines to log and control the output. The idea is that decisions about where to log should ...
integer function, public cp_logger_get_default_io_unit(logger)
returns the unit nr for the ionode (-1 on all other processors) skips as well checks if the procs cal...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
integer, parameter, public cp_p_file
integer function, public cp_print_key_should_output(iteration_info, basis_section, print_key_path, used_print_key, first_time)
returns what should be done with the given property if btest(res,cp_p_store) then the property should...
types that represent a subsys, i.e. a part of the system
Work with symmetry.
Definition cp_symmetry.F:13
subroutine, public write_symmetry(particle_set, cell, input_section)
Write symmetry information to output.
Definition cp_symmetry.F:61
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
subroutine, public distribution_1d_release(distribution_1d)
releases the given distribution_1d
Distribution methods for atoms, particles, or molecules.
subroutine, public distribute_molecules_1d(atomic_kind_set, particle_set, local_particles, molecule_kind_set, molecule_set, local_molecules, force_env_section, prev_molecule_kind_set, prev_local_molecules)
Distribute molecules and particles.
Types needed for a for a Energy Correction.
Energy correction environment setup and handling.
subroutine, public ec_write_input(ec_env)
Print out the energy correction input section.
subroutine, public ec_env_create(qs_env, ec_env, dft_section, ec_section)
Allocates and intitializes ec_env.
Definition and initialisation of the et_coupling data type.
subroutine, public et_coupling_create(et_coupling)
...
subroutine, public ewald_env_set(ewald_env, ewald_type, alpha, epsilon, eps_pol, gmax, ns_max, precs, o_spline, para_env, poisson_section, interaction_cutoffs, cell_hmat)
Purpose: Set the EWALD environment.
subroutine, public ewald_env_create(ewald_env, para_env)
allocates and intitializes a ewald_env
subroutine, public read_ewald_section(ewald_env, ewald_section)
Purpose: read the EWALD section.
subroutine, public read_ewald_section_tb(ewald_env, ewald_section, hmat, silent, pset, cell_periodic)
Purpose: read the EWALD section for TB methods.
subroutine, public ewald_env_get(ewald_env, ewald_type, alpha, eps_pol, epsilon, gmax, ns_max, o_spline, group, para_env, poisson_section, precs, rcut, do_multipoles, max_multipole, do_ipol, max_ipol_iter, interaction_cutoffs, cell_hmat)
Purpose: Get the EWALD environment.
subroutine, public ewald_pw_grid_update(ewald_pw, ewald_env, cell_hmat)
Rescales pw_grids for given box, if necessary.
subroutine, public ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section)
creates the structure ewald_pw_type
Types for excited states potential energies.
subroutine, public exstate_create(ex_env, excited_state, dft_section)
Allocates and intitializes exstate_env.
Definition of the atomic potential types.
subroutine, public fist_nonbond_env_create(fist_nonbond_env, atomic_kind_set, potparm14, potparm, do_nonbonded, do_electrostatics, verlet_skin, ewald_rcut, ei_scale14, vdw_scale14, shift_cutoff)
allocates and intitializes a fist_nonbond_env
Calculation of the incomplete Gamma function F_n(t) for multi-center integrals over Cartesian Gaussia...
Definition gamma.F:15
subroutine, public init_md_ftable(nmax)
Initialize a table of F_n(t) values in the range 0 <= t <= 12 with a stepsize of 0....
Definition gamma.F:540
Define type storing the global information of a run. Keep the amount of stored data small....
subroutine, public init_coulomb_local(hartree_local, natom)
...
subroutine, public dftb_header(iw)
...
Definition header.F:193
subroutine, public qs_header(iw)
...
Definition header.F:310
subroutine, public se_header(iw)
...
Definition header.F:286
subroutine, public xtb_header(iw, gfn_type)
...
Definition header.F:217
subroutine, public tblite_header(iw, tb_type)
...
Definition header.F:252
Types and set/get functions for HFX.
Definition hfx_types.F:15
subroutine, public hfx_create(x_data, para_env, hfx_section, atomic_kind_set, qs_kind_set, particle_set, dft_control, cell, orb_basis, ri_basis, nelectron_total, nkp_grid)
This routine allocates and initializes all types in hfx_data
Definition hfx_types.F:596
subroutine, public compare_hfx_sections(hfx_section1, hfx_section2, is_identical, same_except_frac)
Compares the non-technical parts of two HFX input section and check whether they are the same Ignore ...
Definition hfx_types.F:2982
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public smear_fermi_dirac
integer, parameter, public do_method_ofgpw
integer, parameter, public xc_vdw_fun_none
integer, parameter, public dispersion_d3
integer, parameter, public do_method_rigpw
integer, parameter, public do_method_gpw
integer, parameter, public vdw_pairpot_dftd3
integer, parameter, public do_method_pdg
integer, parameter, public wfi_linear_wf_method_nr
integer, parameter, public linear_response_run
integer, parameter, public wfi_linear_ps_method_nr
integer, parameter, public do_method_pnnl
integer, parameter, public xc_vdw_fun_nonloc
integer, parameter, public vdw_pairpot_dftd4
integer, parameter, public xtb_vdw_type_d3
integer, parameter, public kg_tnadd_embed_ri
integer, parameter, public dispersion_d3bj
integer, parameter, public debug_run
integer, parameter, public wfi_gext_proj_qtr_nr
integer, parameter, public do_method_rm1
integer, parameter, public do_qmmm_swave
integer, parameter, public do_method_pm3
integer, parameter, public do_method_mndo
integer, parameter, public xtb_vdw_type_d4
integer, parameter, public do_method_gapw
integer, parameter, public do_method_mndod
integer, parameter, public rel_trans_atom
integer, parameter, public xtb_vdw_type_none
integer, parameter, public do_method_am1
integer, parameter, public do_method_dftb
integer, parameter, public wfi_use_prev_wf_method_nr
integer, parameter, public do_qmmm_gauss
integer, parameter, public rel_none
integer, parameter, public general_roks
integer, parameter, public vdw_pairpot_dftd2
integer, parameter, public wfi_gext_proj_nr
integer, parameter, public do_method_lrigpw
integer, parameter, public dispersion_d2
integer, parameter, public do_method_xtb
integer, parameter, public do_method_pm6fm
integer, parameter, public xc_vdw_fun_pairpot
integer, parameter, public do_et_ddapc
integer, parameter, public vdw_pairpot_dftd3bj
integer, parameter, public do_method_gapw_xc
integer, parameter, public do_method_pm6
integer, parameter, public hden_atomic
objects that represent the structure of input sections and the data contained in an input section
integer function, dimension(:), pointer, public section_get_ivals(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
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
Routines for a Kim-Gordon-like partitioning into molecular subunits.
subroutine, public kg_env_create(qs_env, kg_env, qs_kind_set, input)
Allocates and intitializes kg_env.
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
Routines needed for kpoint calculation.
subroutine, public kpoint_initialize_mos(kpoint, mos, added_mos, for_aux_fit)
Initialize a set of MOs and density matrix for each kpoint (kpoint group)
subroutine, public kpoint_initialize(kpoint, particle_set, cell)
Generate the kpoints and initialize the kpoint environment.
subroutine, public kpoint_env_initialize(kpoint, para_env, blacs_env, with_aux_fit)
Initialize the kpoint environment.
Types and basic routines needed for a kpoint calculation.
subroutine, public read_kpoint_section(kpoint, kpoint_section, a_vec)
Read the kpoint input section.
subroutine, public write_kpoint_info(kpoint, iounit, dft_section)
Write information on the kpoints to output.
subroutine, public set_kpoint_info(kpoint, kp_scheme, nkp_grid, kp_shift, symmetry, verbose, full_grid, use_real_wfn, eps_geo, parallel_group_size, kp_range, nkp, xkp, wkp, para_env, blacs_env_all, para_env_kp, para_env_inter_kp, blacs_env, kp_env, kp_aux_env, mpools, iogrp, nkp_groups, kp_dist, cell_to_index, index_to_cell, sab_nl, sab_nl_nosym, inversion_symmetry_only, symmetry_backend, symmetry_reduction_method)
Set information in a kpoint environment.
subroutine, public kpoint_create(kpoint)
Create a kpoint environment.
subroutine, public get_kpoint_info(kpoint, kp_scheme, nkp_grid, kp_shift, symmetry, verbose, full_grid, use_real_wfn, eps_geo, parallel_group_size, kp_range, nkp, xkp, wkp, para_env, blacs_env_all, para_env_kp, para_env_inter_kp, blacs_env, kp_env, kp_aux_env, mpools, iogrp, nkp_groups, kp_dist, cell_to_index, index_to_cell, sab_nl, sab_nl_nosym, inversion_symmetry_only, symmetry_backend, symmetry_reduction_method)
Retrieve information from a kpoint environment.
initializes the environment for lri lri : local resolution of the identity
subroutine, public lri_env_init(lri_env, lri_section)
initializes the lri env
subroutine, public lri_env_basis(ri_type, qs_env, lri_env, qs_kind_set)
initializes the lri env
contains the types and subroutines for dealing with the lri_env lri : local resolution of the identit...
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:136
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
Interface to the message passing library MPI.
Define the molecule kind structure types and the corresponding functionality.
subroutine, public write_molecule_kind_set(molecule_kind_set, subsys_section)
Write a moleculeatomic kind set data set to the output unit.
Define the data structure for the molecule information.
Types needed for MP2 calculations.
Definition mp2_setup.F:14
subroutine, public read_mp2_section(input, mp2_env)
...
Definition mp2_setup.F:62
Types needed for MP2 calculations.
Definition mp2_types.F:14
subroutine, public mp2_env_create(mp2_env)
...
Definition mp2_types.F:471
Multipole structure: for multipole (fixed and induced) in FF based MD.
integer, parameter, public do_multipole_none
Provides Cartesian and spherical orbital pointers and indices.
subroutine, public init_orbital_pointers(maxl)
Initialize or update the orbital pointers.
Calculation of the spherical harmonics and the corresponding orbital transformation matrices.
subroutine, public init_spherical_harmonics(maxl, output_unit)
Initialize or update the orbital transformation matrices.
Define methods related to particle_type.
subroutine, public write_qs_particle_coordinates(particle_set, qs_kind_set, subsys_section, label)
Write the atomic coordinates to the output unit.
subroutine, public write_structure_data(particle_set, cell, input_section)
Write structure data requested by a separate structure data input section to the output unit....
subroutine, public write_particle_distances(particle_set, cell, subsys_section)
Write the matrix of the particle distances to the output unit.
Define the data structure for the particle information.
Definition of physical constants:
Definition physcon.F:68
real(kind=dp), parameter, public kelvin
Definition physcon.F:165
container for various plainwaves related things
subroutine, public qs_basis_rotation(qs_env, kpoints, basis_type)
Construct basis set rotation matrices.
subroutine, public qs_dftb_param_init(atomic_kind_set, qs_kind_set, dftb_control, dftb_potential, subsys_section, para_env)
...
Definition of the DFTB parameter types.
Calculation of non local dispersion functionals Some routines adapted from: Copyright (C) 2001-2009 Q...
subroutine, public qs_dispersion_nonloc_init(dispersion_env, para_env)
...
Calculation of dispersion using pair potentials.
subroutine, public qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, pp_section, para_env)
...
Definition of disperson types for DFT calculations.
Set disperson types for DFT calculations.
subroutine, public qs_dispersion_env_set(dispersion_env, xc_section)
...
subroutine, public qs_write_dispersion(qs_env, dispersion_env, ounit)
...
subroutine, public allocate_qs_energy(qs_energy)
Allocate and/or initialise a Quickstep energy data structure.
qs_environment methods that use many other modules
subroutine, public qs_env_setup(qs_env)
initializes various components of the qs_env, that need only atomic_kind_set, cell,...
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, mimic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
subroutine, public set_qs_env(qs_env, super_cell, mos, qmmm, qmmm_periodic, mimic, ewald_env, ewald_pw, mpools, rho_external, external_vxc, mask, scf_control, rel_control, qs_charges, ks_env, ks_qmmm_env, wf_history, scf_env, active_space, input, oce, rho_atom_set, rho0_atom_set, rho0_mpole, run_rtp, rtp, rhoz_set, rhoz_tot, ecoul_1c, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, efield, rhoz_cneo_set, linres_control, xas_env, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, ls_scf_env, do_transport, transport_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, harris_env, gcp_env, mp2_env, bs_env, kg_env, force, kpoints, wanniercentres, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Set the QUICKSTEP environment.
subroutine, public qs_init(qs_env, para_env, root_section, globenv, cp_subsys, kpoint_env, qmmm, qmmm_env_qm, force_env_section, subsys_section, use_motion_section, silent)
Read the input and the database files for the setup of the QUICKSTEP environment.
Definition of gCP types for DFT calculations.
Set disperson types for DFT calculations.
subroutine, public qs_gcp_env_set(gcp_env, xc_section)
...
subroutine, public qs_gcp_init(qs_env, gcp_env)
...
Types needed for a for a Harris model calculation.
subroutine, public harris_rhoin_init(rhoin, basis_type, qs_kind_set, atomic_kind_set, local_particles, nspin)
...
Harris method environment setup and handling.
subroutine, public harris_write_input(harris_env)
Print out the Harris method input section.
subroutine, public harris_env_create(qs_env, harris_env, harris_section)
Allocates and intitializes harris_env.
Calculate the interaction radii for the operator matrix calculation.
subroutine, public write_pgf_orb_radii(basis, atomic_kind_set, qs_kind_set, subsys_section)
Write the orbital basis function radii to the output unit.
subroutine, public write_paw_radii(atomic_kind_set, qs_kind_set, subsys_section)
Write the radii of the one center projector.
subroutine, public write_ppnl_radii(atomic_kind_set, qs_kind_set, subsys_section)
Write the radii of the projector functions of the Goedecker pseudopotential (GTH, non-local part) to ...
subroutine, public init_se_nlradius(se_control, atomic_kind_set, qs_kind_set, subsys_section)
...
subroutine, public write_ppl_radii(atomic_kind_set, qs_kind_set, subsys_section)
Write the radii of the exponential functions of the Goedecker pseudopotential (GTH,...
subroutine, public init_interaction_radii(qs_control, qs_kind_set)
Initialize all the atomic kind radii for a given threshold value.
subroutine, public write_core_charge_radii(atomic_kind_set, qs_kind_set, subsys_section)
Write the radii of the core charge distributions to the output unit.
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, cneo_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zatom, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, monovalent, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
subroutine, public get_qs_kind_set(qs_kind_set, all_potential_present, tnadd_potential_present, gth_potential_present, sgp_potential_present, paw_atom_present, dft_plus_u_atom_present, maxcgf, maxsgf, maxco, maxco_proj, maxgtops, maxlgto, maxlprj, maxnset, maxsgf_set, ncgf, npgf, nset, nsgf, nshell, maxpol, maxlppl, maxlppnl, maxppnl, nelectron, maxder, max_ngrid_rad, max_sph_harm, maxg_iso_not0, lmax_rho0, basis_rcut, basis_type, total_zeff_corr, npgf_seg, cneo_potential_present, nkind_q, natom_q)
Get attributes of an atomic kind set.
subroutine, public init_gapw_nlcc(qs_kind_set)
...
subroutine, public init_cneo_basis_set(qs_kind_set, qs_control)
...
subroutine, public set_qs_kind(qs_kind, paw_atom, ghost, floating, hard_radius, hard0_radius, covalent_radius, vdw_radius, lmax_rho0, zeff, no_optimize, dispersion, u_minus_j, reltmat, dftb_parameter, xtb_parameter, elec_conf, pao_basis_size)
Set the components of an atomic kind data set.
subroutine, public write_qs_kind_set(qs_kind_set, subsys_section)
Write an atomic kind set data set to the output unit.
subroutine, public init_gapw_basis_set(qs_kind_set, qs_control, force_env_section, modify_qs_control)
...
subroutine, public init_qs_kind_set(qs_kind_set)
Initialise an atomic kind set data set.
subroutine, public check_qs_kind_set(qs_kind_set, dft_control, subsys_section)
...
subroutine, public write_gto_basis_sets(qs_kind_set, subsys_section)
Write all the GTO basis sets of an atomic kind set to the output unit (for the printing of the unnorm...
subroutine, public set_ks_env(ks_env, v_hartree_rspace, s_mstruct_changed, rho_changed, exc_accint, potential_changed, forces_up_to_date, complex_ks, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, kinetic, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_ks_im_kp, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, vee, neighbor_list_id, kpoints, sab_orb, sab_all, sac_ae, sac_ppl, sac_lri, sap_ppnl, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_vdw, sab_scp, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, task_list, task_list_soft, subsys, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env)
...
subroutine, public qs_ks_env_create(ks_env)
Allocates a new instance of ks_env.
Definition and initialisation of the mo data type.
Definition qs_mo_types.F:22
subroutine, public allocate_mo_set(mo_set, nao, nmo, nelectron, n_el_f, maxocc, flexible_electron_count)
Allocates a mo set and partially initializes it (nao,nmo,nelectron, and flexible_electron_count are v...
subroutine, public rho0_s_grid_create(pw_env, rho0_mpole)
...
subroutine, public init_rho0(local_rho_set, qs_env, gapw_control, zcore)
...
subroutine, public init_rho_atom(rho_atom_set, atomic_kind_set, qs_kind_set, dft_control, para_env)
...
Routines that work on qs_subsys_type.
subroutine, public qs_subsys_create(subsys, para_env, root_section, force_env_section, subsys_section, use_motion_section, cp_subsys, elkind, silent)
Creates a qs_subsys. Optionally an existsing cp_subsys is used.
types that represent a quickstep subsys
subroutine, public qs_subsys_get(subsys, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell, cell_ref, use_ref_cell, energy, force, qs_kind_set, cp_subsys, nelectron_total, nelectron_spin)
...
subroutine, public qs_subsys_set(subsys, cp_subsys, local_particles, local_molecules, cell, cell_ref, use_ref_cell, energy, force, qs_kind_set, nelectron_total, nelectron_spin)
...
Storage of past states of the qs_env. Methods to interpolate (or actually normally extrapolate) the n...
subroutine, public wfi_create_for_kp(wf_history)
Adapts wf_history storage flags for k-point calculations. For ASPC, switches from Gamma WFN storage t...
subroutine, public wfi_create(wf_history, interpolation_method_nr, extrapolation_order, has_unit_metric)
...
interpolate the wavefunctions to speed up the convergence when doing MD
subroutine, public wfi_release(wf_history)
releases a wf_history of a wavefunction (see doc/ReferenceCounting.html)
parameters that control a relativistic calculation
subroutine, public rel_c_create(rel_control)
allocates and initializes an rel control object with the default values
subroutine, public rel_c_read_parameters(rel_control, dft_section)
reads the parameters of the relativistic section into the given rel_control
parameters that control an scf iteration
subroutine, public scf_c_read_parameters(scf_control, inp_section)
reads the parameters of the scf section into the given scf_control
subroutine, public scf_c_write_parameters(scf_control, dft_section)
writes out the scf parameters
subroutine, public scf_c_create(scf_control)
allocates and initializes an scf control object with the default values
Methods for handling the 1/R^3 residual integral part.
subroutine, public semi_empirical_expns3_setup(qs_kind_set, se_control, method_id)
Setup the quantity necessary to handle the slowly convergent residual integral term 1/R^3.
Arrays of parameters used in the semi-empirical calculations \References Everywhere in this module TC...
subroutine, public init_se_intd_array()
Initialize all arrays used for the evaluation of the integrals.
Setup and Methods for semi-empirical multipole types.
subroutine, public nddo_mpole_setup(nddo_mpole, natom)
Setup NDDO multipole type.
Definition of the semi empirical multipole integral expansions types.
Type to store integrals for semi-empirical calculations.
subroutine, public semi_empirical_si_create(store_int_env, se_section, compression)
Allocate semi-empirical store integrals type.
Definition of the semi empirical parameter types.
subroutine, public se_taper_create(se_taper, integral_screening, do_ewald, taper_cou, range_cou, taper_exc, range_exc, taper_scr, range_scr, taper_lrc, range_lrc)
Creates the taper type used in SE calculations.
Working with the semi empirical parameter types.
subroutine, public se_cutoff_compatible(se_control, se_section, cell, output_unit)
Reset cutoffs trying to be somehow a bit smarter.
interface to tblite
subroutine, public tb_set_calculator(tb, typ)
...
subroutine, public tb_init_wf(tb)
initialize wavefunction ...
subroutine, public tb_init_geometry(qs_env, tb)
intialize geometry objects ...
subroutine, public tb_get_basis(tb, gto_basis_set, element_symbol, param, occ)
...
routines for DFT+NEGF calculations (coupling with the quantum transport code OMEN)
Definition transport.F:19
subroutine, public transport_env_create(qs_env)
creates the transport environment
Definition transport.F:124
Read xTB parameters.
subroutine, public xtb_parameters_set(param)
Read atom parameters for xTB Hamiltonian from input file.
subroutine, public xtb_parameters_init(param, gfn_type, element_symbol, parameter_file_path, parameter_file_name, para_env)
...
subroutine, public init_xtb_basis(param, gto_basis_set, ngauss)
...
xTB (repulsive) pair potentials Reference: Stefan Grimme, Christoph Bannwarth, Philip Shushkov JCTC 1...
subroutine, public xtb_pp_radius(qs_kind_set, ppradius, eps_pair, kfparam)
...
Definition of the xTB parameter types.
Definition xtb_types.F:20
subroutine, public allocate_xtb_atom_param(xtb_parameter)
...
Definition xtb_types.F:99
subroutine, public set_xtb_atom_param(xtb_parameter, aname, typ, defined, z, zeff, natorb, lmax, nao, lao, rcut, rcov, kx, eta, xgamma, alpha, zneff, nshell, nval, lval, kpoly, kappa, hen, zeta, xi, kappa0, alpg, electronegativity, occupation, chmax, en, kqat2, kcn, kq)
...
Definition xtb_types.F:299
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:60
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
type of a logger, at the moment it contains just a print level starting at which level it should be l...
represents a system: atoms, molecules, their pos,vel,...
structure to store local (to a processor) ordered lists of integers.
Contains information on the energy correction functional for KG.
Contains information on the excited states energy.
contains the initially parsed file and the initial parallel environment
Contains information about kpoints.
stores all the informations relevant to an mpi environment
contained for different pw related things
Contains information on the Harris method.
Provides all information about a quickstep kind.
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...
keeps track of the previous wavefunctions and can extrapolate them for the next step of md
contains the parameters needed by a relativistic calculation
Global Multipolar NDDO information type.
Taper type use in semi-empirical calculations.