288 LOGICAL,
OPTIONAL :: molecular
291 CHARACTER(len=*),
PARAMETER :: routinen =
'build_qs_neighbor_lists'
293 CHARACTER(LEN=2) :: element_symbol, element_symbol2
294 CHARACTER(LEN=default_string_length) :: print_key_path
295 INTEGER :: handle, hfx_pot, ikind, ingp, iw, jkind, &
296 maxatom, ngp, nkind, zat
297 LOGICAL :: all_potential_present, almo, dftb, do_hfx, dokp, gth_potential_present, &
298 lri_optbas, lrigpw, mic, molecule_only, nddo, paw_atom, paw_atom_present, rigpw, &
299 sgp_potential_present, xtb
300 LOGICAL,
ALLOCATABLE,
DIMENSION(:) :: all_present, aux_fit_present, aux_present, &
301 core_present, default_present, nonbond1_atom, nonbond2_atom, oce_present, orb_present, &
302 ppl_present, ppnl_present, ri_present, xb1_atom, xb2_atom
303 REAL(
dp) :: almo_rcov, almo_rvdw, eps_schwarz, &
304 omega, pdist, rcut, roperator, subcells
305 REAL(
dp),
ALLOCATABLE,
DIMENSION(:) :: all_pot_rad, aux_fit_radius, c_radius, calpha, &
306 core_radius, oce_radius, orb_radius, ppl_radius, ppnl_radius, ri_radius, zeff
307 REAL(
dp),
ALLOCATABLE,
DIMENSION(:, :) :: pair_radius, pair_radius_lb
318 orb_basis_set, ri_basis_set
323 sab_cn, sab_core, sab_gcp, sab_kp, sab_kp_nosym, sab_lrc, sab_orb, sab_scp, sab_se, &
324 sab_tbe, sab_vdw, sab_xb, sab_xtb_nonbond, sab_xtb_pp, sab_xtbe, sac_ae, sac_lri, &
325 sac_ppl, sap_oce, sap_ppnl, soa_list, soo_list
331 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
337 CALL timeset(routinen, handle)
341 NULLIFY (atomic_kind_set, qs_kind_set, cell, neighbor_list_section, &
342 distribution_1d, distribution_2d, gth_potential, sgp_potential, orb_basis_set, &
343 particle_set, molecule_set, dft_control, ks_env)
358 NULLIFY (sab_xtb_nonbond)
366 NULLIFY (sab_kp_nosym)
370 atomic_kind_set=atomic_kind_set, &
371 qs_kind_set=qs_kind_set, &
374 distribution_2d=distribution_2d, &
375 local_particles=distribution_1d, &
376 particle_set=particle_set, &
377 molecule_set=molecule_set, &
378 dft_control=dft_control)
384 last_qs_neighbor_list_id_nr = last_qs_neighbor_list_id_nr + 1
385 CALL set_ks_env(ks_env=ks_env, neighbor_list_id=last_qs_neighbor_list_id_nr)
401 sab_xtb_pp=sab_xtb_pp, &
402 sab_xtb_nonbond=sab_xtb_nonbond, &
407 sab_kp_nosym=sab_kp_nosym)
409 dokp = (kpoints%nkp > 0)
410 nddo = dft_control%qs_control%semi_empirical
411 dftb = dft_control%qs_control%dftb
412 xtb = dft_control%qs_control%xtb
413 almo = dft_control%qs_control%do_almo_scf
416 lri_optbas = dft_control%qs_control%lri_optbas
419 molecule_only = .false.
420 IF (
PRESENT(molecular)) molecule_only = molecular
430 pdist = dft_control%qs_control%pairlist_radius
437 gth_potential_present=gth_potential_present, &
438 sgp_potential_present=sgp_potential_present, &
439 all_potential_present=all_potential_present)
444 nkind =
SIZE(atomic_kind_set)
445 ALLOCATE (orb_present(nkind), aux_fit_present(nkind), aux_present(nkind), &
446 default_present(nkind), core_present(nkind))
447 ALLOCATE (orb_radius(nkind), aux_fit_radius(nkind), c_radius(nkind), &
448 core_radius(nkind), calpha(nkind), zeff(nkind))
449 orb_radius(:) = 0.0_dp
450 aux_fit_radius(:) = 0.0_dp
452 core_radius(:) = 0.0_dp
456 ALLOCATE (pair_radius(nkind, nkind))
457 IF (gth_potential_present .OR. sgp_potential_present)
THEN
458 ALLOCATE (ppl_present(nkind), ppl_radius(nkind))
460 ALLOCATE (ppnl_present(nkind), ppnl_radius(nkind))
463 IF (paw_atom_present)
THEN
464 ALLOCATE (oce_present(nkind), oce_radius(nkind))
467 IF (all_potential_present .OR. sgp_potential_present)
THEN
468 ALLOCATE (all_present(nkind), all_pot_rad(nkind))
473 ALLOCATE (atom2d(nkind))
474 CALL atom2d_build(atom2d, distribution_1d, distribution_2d, atomic_kind_set, &
475 molecule_set, molecule_only, particle_set=particle_set)
479 CALL get_atomic_kind(atomic_kind_set(ikind), atom_list=atom2d(ikind)%list)
481 CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set, basis_type=
"ORB")
483 CALL get_qs_kind(qs_kind_set(ikind), basis_set=aux_fit_basis_set, basis_type=
"AUX_FIT")
486 paw_proj_set=paw_proj, &
488 all_potential=all_potential, &
489 gth_potential=gth_potential, &
490 sgp_potential=sgp_potential)
495 CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_atom)
497 cutoff=orb_radius(ikind), &
498 defined=orb_present(ikind))
500 IF (
ASSOCIATED(orb_basis_set))
THEN
501 orb_present(ikind) = .true.
502 CALL get_gto_basis_set(gto_basis_set=orb_basis_set, kind_radius=orb_radius(ikind))
504 orb_present(ikind) = .false.
509 aux_present(ikind) = .true.
511 aux_present(ikind) = .false.
514 IF (
ASSOCIATED(aux_fit_basis_set))
THEN
515 aux_fit_present(ikind) = .true.
516 CALL get_gto_basis_set(gto_basis_set=aux_fit_basis_set, kind_radius=aux_fit_radius(ikind))
518 aux_fit_present(ikind) = .false.
523 alpha_core_charge=calpha(ikind), &
524 core_charge_radius=core_radius(ikind), &
526 IF (zeff(ikind) /= 0._dp .AND. calpha(ikind) /= 0._dp)
THEN
527 core_present(ikind) = .true.
529 core_present(ikind) = .false.
533 IF (gth_potential_present .OR. sgp_potential_present)
THEN
534 IF (
ASSOCIATED(gth_potential))
THEN
536 ppl_present=ppl_present(ikind), &
537 ppl_radius=ppl_radius(ikind), &
538 ppnl_present=ppnl_present(ikind), &
539 ppnl_radius=ppnl_radius(ikind))
540 ELSE IF (
ASSOCIATED(sgp_potential))
THEN
542 ppl_present=ppl_present(ikind), &
543 ppl_radius=ppl_radius(ikind), &
544 ppnl_present=ppnl_present(ikind), &
545 ppnl_radius=ppnl_radius(ikind))
547 ppl_present(ikind) = .false.
548 ppnl_present(ikind) = .false.
553 IF (paw_atom_present)
THEN
555 oce_present(ikind) = .true.
558 oce_present(ikind) = .false.
563 IF (all_potential_present .OR. sgp_potential_present)
THEN
564 all_present(ikind) = .false.
565 all_pot_rad(ikind) = 0.0_dp
566 IF (
ASSOCIATED(all_potential))
THEN
567 all_present(ikind) = .true.
568 CALL get_potential(potential=all_potential, core_charge_radius=all_pot_rad(ikind))
569 ELSE IF (
ASSOCIATED(sgp_potential))
THEN
570 IF (sgp_potential%ecp_local)
THEN
571 all_present(ikind) = .true.
572 CALL get_potential(potential=sgp_potential, core_charge_radius=all_pot_rad(ikind))
580 IF (pdist < 0.0_dp)
THEN
585 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius, pdist)
587 mic=mic, subcells=subcells, molecular=molecule_only, nlname=
"sab_orb")
588 CALL set_ks_env(ks_env=ks_env, sab_orb=sab_orb)
590 "/SAB_ORB",
"sab_orb",
"ORBITAL ORBITAL")
595 IF (.NOT. (nddo .OR. dftb .OR. xtb))
THEN
597 mic=mic, symmetric=.false., subcells=subcells, molecular=molecule_only, nlname=
"sab_all")
598 CALL set_ks_env(ks_env=ks_env, sab_all=sab_all)
602 IF (.NOT. (nddo .OR. dftb .OR. xtb))
THEN
603 CALL pair_radius_setup(core_present, core_present, core_radius, core_radius, pair_radius)
604 CALL build_neighbor_lists(sab_core, particle_set, atom2d, cell, pair_radius, subcells=subcells, &
605 operator_type=
"PP", nlname=
"sab_core")
606 CALL set_ks_env(ks_env=ks_env, sab_core=sab_core)
608 "/SAB_CORE",
"sab_core",
"CORE CORE")
621 SELECT CASE (hfx_pot)
631 cpabort(
"HFX potential not available for K-points (NYI)")
634 IF (dft_control%do_admm)
THEN
635 CALL pair_radius_setup(aux_fit_present, aux_fit_present, aux_fit_radius, aux_fit_radius, &
639 ALLOCATE (pair_radius_lb(nkind, nkind))
640 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius_lb)
643 IF (pair_radius(ikind, jkind) +
cutoff_screen_factor*roperator .LE. pair_radius_lb(ikind, jkind)) &
644 pair_radius(ikind, jkind) = pair_radius_lb(ikind, jkind) - roperator
648 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius)
652 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius)
655 subcells=subcells, nlname=
"sab_kp")
660 subcells=subcells, nlname=
"sab_kp_nosym", symmetric=.false.)
661 CALL set_ks_env(ks_env=ks_env, sab_kp_nosym=sab_kp_nosym)
666 IF (gth_potential_present .OR. sgp_potential_present)
THEN
667 IF (any(ppl_present))
THEN
668 CALL pair_radius_setup(orb_present, ppl_present, orb_radius, ppl_radius, pair_radius)
670 subcells=subcells, operator_type=
"ABC", nlname=
"sac_ppl")
671 CALL set_ks_env(ks_env=ks_env, sac_ppl=sac_ppl)
673 "/SAC_PPL",
"sac_ppl",
"ORBITAL GTH-PPL")
675 IF (qs_env%lri_env%ppl_ri)
THEN
677 subcells=subcells, symmetric=.false., operator_type=
"PP", nlname=
"sac_lri")
678 CALL set_ks_env(ks_env=ks_env, sac_lri=sac_lri)
683 IF (any(ppnl_present))
THEN
684 CALL pair_radius_setup(orb_present, ppnl_present, orb_radius, ppnl_radius, pair_radius)
686 subcells=subcells, operator_type=
"ABBA", nlname=
"sap_ppnl")
687 CALL set_ks_env(ks_env=ks_env, sap_ppnl=sap_ppnl)
689 "/SAP_PPNL",
"sap_ppnl",
"ORBITAL GTH-PPNL")
693 IF (paw_atom_present)
THEN
695 IF (any(oce_present))
THEN
696 CALL pair_radius_setup(orb_present, oce_present, orb_radius, oce_radius, pair_radius)
698 subcells=subcells, operator_type=
"ABBA", nlname=
"sap_oce")
699 CALL set_ks_env(ks_env=ks_env, sap_oce=sap_oce)
701 "/SAP_OCE",
"sap_oce",
"ORBITAL(A) PAW-PRJ")
706 IF (.NOT. (nddo .OR. dftb .OR. xtb))
THEN
707 IF (all_potential_present .OR. sgp_potential_present)
THEN
708 CALL pair_radius_setup(orb_present, all_present, orb_radius, all_pot_rad, pair_radius)
710 subcells=subcells, operator_type=
"ABC", nlname=
"sac_ae")
713 "/SAC_AE",
"sac_ae",
"ORBITAL ERFC POTENTIAL")
719 default_present = .true.
720 c_radius = dft_control%qs_control%se_control%cutoff_cou
722 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
723 IF (dft_control%qs_control%se_control%do_ewald_gks)
THEN
726 subcells=subcells, nlname=
"sab_se")
729 subcells=subcells, nlname=
"sab_se")
733 "/SAB_SE",
"sab_se",
"HARTREE INTERACTIONS")
736 IF ((dft_control%qs_control%se_control%do_ewald) .AND. &
737 (dft_control%qs_control%se_control%integral_screening /=
do_se_is_slater))
THEN
738 c_radius = dft_control%qs_control%se_control%cutoff_lrc
739 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
741 subcells=subcells, nlname=
"sab_lrc")
742 CALL set_ks_env(ks_env=ks_env, sab_lrc=sab_lrc)
744 "/SAB_LRC",
"sab_lrc",
"SE LONG-RANGE CORRECTION")
750 IF (dft_control%qs_control%dftb_control%do_ewald)
THEN
751 CALL get_qs_env(qs_env=qs_env, ewald_env=ewald_env)
754 CALL pair_radius_setup(orb_present, orb_present, c_radius, c_radius, pair_radius)
756 subcells=subcells, nlname=
"sab_tbe")
757 CALL set_ks_env(ks_env=ks_env, sab_tbe=sab_tbe)
761 IF (dft_control%qs_control%dftb_control%dispersion)
THEN
762 IF (dft_control%qs_control%dftb_control%dispersion_type ==
dispersion_uff)
THEN
764 CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_atom)
767 default_present = .true.
768 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
770 subcells=subcells, nlname=
"sab_vdw")
771 CALL set_ks_env(ks_env=ks_env, sab_vdw=sab_vdw)
776 IF (xtb .AND. (.NOT. dft_control%qs_control%xtb_control%do_tblite))
THEN
778 IF (dft_control%qs_control%xtb_control%do_ewald)
THEN
779 CALL get_qs_env(qs_env=qs_env, ewald_env=ewald_env)
782 CALL pair_radius_setup(orb_present, orb_present, c_radius, c_radius, pair_radius)
784 subcells=subcells, nlname=
"sab_tbe")
785 CALL set_ks_env(ks_env=ks_env, sab_tbe=sab_tbe)
788 pair_radius(1:nkind, 1:nkind) = dft_control%qs_control%xtb_control%rcpair(1:nkind, 1:nkind)
789 default_present = .true.
791 subcells=subcells, nlname=
"sab_xtb_pp")
792 CALL set_ks_env(ks_env=ks_env, sab_xtb_pp=sab_xtb_pp)
795 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom)
798 default_present = .true.
799 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
801 subcells=subcells, nlname=
"sab_xtbe")
802 CALL set_ks_env(ks_env=ks_env, sab_xtbe=sab_xtbe)
804 ALLOCATE (xb1_atom(nkind), xb2_atom(nkind))
805 c_radius = 0.5_dp*dft_control%qs_control%xtb_control%xb_radius
808 IF (zat == 17 .OR. zat == 35 .OR. zat == 53 .OR. zat == 85)
THEN
809 xb1_atom(ikind) = .true.
811 xb1_atom(ikind) = .false.
813 IF (zat == 7 .OR. zat == 8 .OR. zat == 15 .OR. zat == 16)
THEN
814 xb2_atom(ikind) = .true.
816 xb2_atom(ikind) = .false.
821 symmetric=.false., subcells=subcells, operator_type=
"PP", nlname=
"sab_xb")
824 "/SAB_XB",
"sab_xb",
"XB bonding")
827 IF (dft_control%qs_control%xtb_control%do_nonbonded &
828 .AND. (.NOT. dft_control%qs_control%xtb_control%do_tblite))
THEN
829 ngp =
SIZE(dft_control%qs_control%xtb_control%nonbonded%pot)
830 ALLOCATE (nonbond1_atom(nkind), nonbond2_atom(nkind))
831 nonbond1_atom = .false.
832 nonbond2_atom = .false.
835 rcut = sqrt(dft_control%qs_control%xtb_control%nonbonded%pot(ingp)%pot%rcutsq)
837 CALL get_atomic_kind(atomic_kind_set(ikind), element_symbol=element_symbol)
839 IF (trim(dft_control%qs_control%xtb_control%nonbonded%pot(ingp)%pot%at1) == trim(element_symbol))
THEN
840 nonbond1_atom(ikind) = .true.
842 CALL get_atomic_kind(atomic_kind_set(jkind), element_symbol=element_symbol2)
844 IF (trim(dft_control%qs_control%xtb_control%nonbonded%pot(ingp)%pot%at2) == trim(element_symbol2))
THEN
845 nonbond2_atom(jkind) = .true.
850 CALL pair_radius_setup(nonbond1_atom, nonbond2_atom, c_radius, c_radius, pair_radius)
852 symmetric=.false., subcells=subcells, operator_type=
"PP", nlname=
"sab_xtb_nonbond")
853 CALL set_ks_env(ks_env=ks_env, sab_xtb_nonbond=sab_xtb_nonbond)
854 CALL write_neighbor_lists(sab_xtb_nonbond, particle_set, cell, para_env, neighbor_list_section, &
855 "/SAB_XTB_NONBOND",
"sab_xtb_nonbond",
"XTB NONBONDED INTERACTIONS")
861 IF (.NOT. dft_control%qs_control%xtb_control%do_tblite)
THEN
862 CALL get_qs_env(qs_env=qs_env, dispersion_env=dispersion_env)
863 sab_vdw => dispersion_env%sab_vdw
864 sab_cn => dispersion_env%sab_cn
867 c_radius(:) = dispersion_env%rc_d4
869 c_radius(:) = dispersion_env%rc_disp
871 default_present = .true.
872 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
874 subcells=subcells, operator_type=
"PP", nlname=
"sab_vdw")
875 dispersion_env%sab_vdw => sab_vdw
881 c_radius(ikind) = 4._dp*
ptable(zat)%covalent_radius*
bohr
883 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
885 subcells=subcells, operator_type=
"PP", nlname=
"sab_cn")
886 dispersion_env%sab_cn => sab_cn
892 CALL get_qs_env(qs_env=qs_env, gcp_env=gcp_env)
893 IF (
ASSOCIATED(gcp_env))
THEN
894 IF (gcp_env%do_gcp)
THEN
895 sab_gcp => gcp_env%sab_gcp
897 c_radius(ikind) = gcp_env%gcp_kind(ikind)%rcsto
899 CALL pair_radius_setup(orb_present, orb_present, c_radius, c_radius, pair_radius)
901 subcells=subcells, operator_type=
"PP", nlname=
"sab_gcp")
902 gcp_env%sab_gcp => sab_gcp
904 NULLIFY (gcp_env%sab_gcp)
908 IF (lrigpw .OR. lri_optbas)
THEN
910 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius)
911 soo_list => qs_env%lri_env%soo_list
913 mic=mic, molecular=molecule_only, subcells=subcells, nlname=
"soo_list")
914 qs_env%lri_env%soo_list => soo_list
916 "/SOO_LIST",
"soo_list",
"ORBITAL ORBITAL (RI)")
918 ALLOCATE (ri_present(nkind), ri_radius(nkind))
922 CALL get_qs_kind(qs_kind_set(ikind), basis_set=ri_basis_set, basis_type=
"RI_HXC")
923 IF (
ASSOCIATED(ri_basis_set))
THEN
924 ri_present(ikind) = .true.
927 ri_present(ikind) = .false.
931 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius)
932 soo_list => qs_env%lri_env%soo_list
934 mic=mic, molecular=molecule_only, subcells=subcells, nlname=
"soo_list")
935 qs_env%lri_env%soo_list => soo_list
937 CALL pair_radius_setup(ri_present, ri_present, ri_radius, ri_radius, pair_radius)
938 saa_list => qs_env%lri_env%saa_list
940 mic=mic, molecular=molecule_only, subcells=subcells, nlname=
"saa_list")
941 qs_env%lri_env%saa_list => saa_list
943 CALL pair_radius_setup(ri_present, orb_present, ri_radius, orb_radius, pair_radius)
944 soa_list => qs_env%lri_env%soa_list
946 mic=mic, symmetric=.false., molecular=molecule_only, &
947 subcells=subcells, operator_type=
"ABC", nlname=
"saa_list")
948 qs_env%lri_env%soa_list => soa_list
954 CALL get_atomic_kind(atomic_kind_set(ikind), rcov=almo_rcov, rvdw=almo_rvdw)
956 c_radius(ikind) = max(almo_rcov, almo_rvdw)*
bohr* &
959 default_present = .true.
960 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
962 subcells=subcells, operator_type=
"PP", nlname=
"sab_almo")
963 CALL set_ks_env(ks_env=ks_env, sab_almo=sab_almo)
967 print_key_path =
"PRINT%DISTRIBUTION"
972 basis_section=force_env_section, &
973 print_key_path=print_key_path, &
975 CALL write_neighbor_distribution(sab_orb, qs_kind_set, iw, para_env)
978 basis_section=force_env_section, &
979 print_key_path=print_key_path)
986 DEALLOCATE (orb_present, default_present, core_present)
987 DEALLOCATE (orb_radius, aux_fit_radius, c_radius, core_radius)
988 DEALLOCATE (calpha, zeff)
989 DEALLOCATE (pair_radius)
990 IF (gth_potential_present .OR. sgp_potential_present)
THEN
991 DEALLOCATE (ppl_present, ppl_radius)
992 DEALLOCATE (ppnl_present, ppnl_radius)
994 IF (paw_atom_present)
THEN
995 DEALLOCATE (oce_present, oce_radius)
997 IF (all_potential_present .OR. sgp_potential_present)
THEN
998 DEALLOCATE (all_present, all_pot_rad)
1001 CALL timestop(handle)
1029 mic, symmetric, molecular, subset_of_mol, current_subset, &
1030 operator_type, nlname, atomb_to_keep)
1037 REAL(
dp),
DIMENSION(:, :),
INTENT(IN) :: pair_radius
1038 REAL(
dp),
INTENT(IN) :: subcells
1039 LOGICAL,
INTENT(IN),
OPTIONAL :: mic, symmetric, molecular
1040 INTEGER,
DIMENSION(:),
OPTIONAL,
POINTER :: subset_of_mol
1041 INTEGER,
OPTIONAL :: current_subset
1042 CHARACTER(LEN=*),
INTENT(IN),
OPTIONAL :: operator_type
1043 CHARACTER(LEN=*),
INTENT(IN) :: nlname
1044 INTEGER,
DIMENSION(:),
INTENT(IN),
OPTIONAL :: atomb_to_keep
1046 CHARACTER(len=*),
PARAMETER :: routinen =
'build_neighbor_lists'
1049 INTEGER,
DIMENSION(:),
POINTER ::
list
1050 END TYPE local_lists
1052 INTEGER :: atom_a, atom_b, handle, i, iab, iatom, iatom_local, &
1053 iatom_subcell, icell, ikind, j, jatom, jatom_local, jcell, jkind, k, &
1054 kcell, maxat, mol_a, mol_b, nkind, otype, natom, inode, nnode, nentry
1055 INTEGER,
DIMENSION(3) :: cell_b, ncell, nsubcell, periodic
1056 INTEGER,
DIMENSION(:),
POINTER :: index_list
1057 LOGICAL :: include_ab, my_mic, &
1058 my_molecular, my_symmetric, my_sort_atomb
1059 LOGICAL,
ALLOCATABLE,
DIMENSION(:) :: pres_a, pres_b
1060 REAL(
dp) :: rab2, rab2_max, rab_max, rabm, deth, subcell_scale
1061 REAL(
dp),
DIMENSION(3) :: r, rab, ra, rb, sab_max, sb, &
1062 sb_pbc, sb_min, sb_max, rab_pbc, pd, sab_max_guard
1063 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: nlista, nlistb
1064 TYPE(local_lists),
DIMENSION(:),
POINTER :: lista, listb
1066 ALLOCATABLE,
DIMENSION(:) :: kind_a
1070 REAL(kind=
dp),
DIMENSION(:, :),
ALLOCATABLE :: r_pbc
1072 DIMENSION(:),
POINTER :: nl_iterator
1074 CALL timeset(routinen//
"_"//trim(nlname), handle)
1078 IF (
PRESENT(mic)) my_mic = mic
1079 my_symmetric = .true.
1080 IF (
PRESENT(symmetric)) my_symmetric = symmetric
1081 my_molecular = .false.
1083 IF (
PRESENT(molecular)) my_molecular = molecular
1085 IF (
PRESENT(operator_type))
THEN
1086 SELECT CASE (operator_type)
1091 cpassert(.NOT. my_molecular)
1092 my_symmetric = .false.
1095 my_symmetric = .false.
1105 my_sort_atomb = .false.
1106 IF (
PRESENT(atomb_to_keep))
THEN
1107 my_sort_atomb = .true.
1114 ALLOCATE (ab_list(nkind*nkind))
1115 DO iab = 1,
SIZE(ab_list)
1116 NULLIFY (ab_list(iab)%neighbor_list_set)
1117 ab_list(iab)%nl_size = -1
1118 ab_list(iab)%nl_start = -1
1119 ab_list(iab)%nl_end = -1
1120 NULLIFY (ab_list(iab)%nlist_task)
1124 ALLOCATE (pres_a(nkind), pres_b(nkind))
1126 pres_a(ikind) = any(pair_radius(ikind, :) > 0._dp)
1127 pres_b(ikind) = any(pair_radius(:, ikind) > 0._dp)
1131 natom =
SIZE(particle_set)
1132 ALLOCATE (r_pbc(3, natom))
1134 r_pbc(1:3, i) =
pbc(particle_set(i)%r(1:3), cell)
1140 maxat = max(maxat,
SIZE(
atom(ikind)%list))
1142 ALLOCATE (index_list(maxat))
1146 ALLOCATE (lista(nkind), listb(nkind), nlista(nkind), nlistb(nkind))
1150 NULLIFY (lista(ikind)%list, listb(ikind)%list)
1153 IF (
ASSOCIATED(
atom(ikind)%list_local_a_index))
THEN
1154 lista(ikind)%list =>
atom(ikind)%list_local_a_index
1155 nlista(ikind) =
SIZE(lista(ikind)%list)
1157 IF (
ASSOCIATED(
atom(ikind)%list_local_b_index))
THEN
1158 listb(ikind)%list =>
atom(ikind)%list_local_b_index
1159 nlistb(ikind) =
SIZE(listb(ikind)%list)
1162 IF (
ASSOCIATED(
atom(ikind)%list_local_a_index))
THEN
1163 lista(ikind)%list =>
atom(ikind)%list_local_a_index
1164 nlista(ikind) =
SIZE(lista(ikind)%list)
1166 nlistb(ikind) =
SIZE(
atom(ikind)%list)
1167 listb(ikind)%list => index_list
1169 CALL combine_lists(lista(ikind)%list, nlista(ikind), ikind,
atom)
1170 nlistb(ikind) =
SIZE(
atom(ikind)%list)
1171 listb(ikind)%list => index_list
1173 nlista(ikind) =
SIZE(
atom(ikind)%list_1d)
1174 lista(ikind)%list =>
atom(ikind)%list_1d
1175 nlistb(ikind) =
SIZE(
atom(ikind)%list)
1176 listb(ikind)%list => index_list
1185 maxat = max(maxat, nlista(ikind), nlistb(ikind))
1187 ALLOCATE (kind_a(2*maxat))
1190 CALL get_cell(cell=cell, periodic=periodic, deth=deth)
1194 IF (.NOT. pres_a(ikind)) cycle
1197 IF (.NOT. pres_b(jkind)) cycle
1199 iab = ikind + nkind*(jkind - 1)
1202 IF (pair_radius(ikind, jkind) <= 0._dp) cycle
1203 rab_max = pair_radius(ikind, jkind)
1204 IF (otype == 3)
THEN
1208 rabm = maxval(pair_radius(:, jkind))
1212 rab2_max = rabm*rabm
1219 sab_max_guard = 15.0_dp/pd
1222 subcell_scale = ((125.0_dp**3)/deth)**(1.0_dp/6.0_dp)
1226 nsubcell(:) = int(max(1.0_dp, min(0.5_dp*subcells*subcell_scale/sab_max(:), &
1227 0.5_dp*subcells*subcell_scale/sab_max_guard(:))))
1230 ncell(:) = (int(sab_max(:)) + 1)*periodic(:)
1233 symmetric=my_symmetric)
1234 neighbor_list_set => ab_list(iab)%neighbor_list_set
1236 DO iatom_local = 1, nlista(ikind)
1237 iatom = lista(ikind)%list(iatom_local)
1238 atom_a =
atom(ikind)%list(iatom)
1241 neighbor_list=kind_a(iatom_local)%neighbor_list)
1245 DO iatom_local = 1, nlista(ikind)
1246 iatom = lista(ikind)%list(iatom_local)
1247 atom_a =
atom(ikind)%list(iatom)
1248 r = r_pbc(:, atom_a)
1250 subcell(i, j, k)%natom = subcell(i, j, k)%natom + 1
1252 DO k = 1, nsubcell(3)
1253 DO j = 1, nsubcell(2)
1254 DO i = 1, nsubcell(1)
1255 maxat = subcell(i, j, k)%natom + subcell(i, j, k)%natom/10
1256 ALLOCATE (subcell(i, j, k)%atom_list(maxat))
1257 subcell(i, j, k)%natom = 0
1261 DO iatom_local = 1, nlista(ikind)
1262 iatom = lista(ikind)%list(iatom_local)
1263 atom_a =
atom(ikind)%list(iatom)
1264 r = r_pbc(:, atom_a)
1266 subcell(i, j, k)%natom = subcell(i, j, k)%natom + 1
1267 subcell(i, j, k)%atom_list(subcell(i, j, k)%natom) = iatom_local
1270 DO jatom_local = 1, nlistb(jkind)
1271 jatom = listb(jkind)%list(jatom_local)
1272 atom_b =
atom(jkind)%list(jatom)
1273 IF (my_sort_atomb .AND. .NOT. my_symmetric)
THEN
1274 IF (.NOT. any(atomb_to_keep == atom_b)) cycle
1276 IF (my_molecular)
THEN
1277 mol_b =
atom(jkind)%list_b_mol(jatom_local)
1278 IF (
PRESENT(subset_of_mol))
THEN
1279 IF (subset_of_mol(mol_b) .NE. current_subset) cycle
1282 r = r_pbc(:, atom_b)
1285 loop2_kcell:
DO kcell = -ncell(3), ncell(3)
1286 sb(3) = sb_pbc(3) + real(kcell,
dp)
1287 sb_min(3) = sb(3) - sab_max(3)
1288 sb_max(3) = sb(3) + sab_max(3)
1289 IF (periodic(3) /= 0)
THEN
1290 IF (sb_min(3) >= 0.5_dp)
EXIT loop2_kcell
1291 IF (sb_max(3) < -0.5_dp) cycle loop2_kcell
1295 loop2_jcell:
DO jcell = -ncell(2), ncell(2)
1296 sb(2) = sb_pbc(2) + real(jcell,
dp)
1297 sb_min(2) = sb(2) - sab_max(2)
1298 sb_max(2) = sb(2) + sab_max(2)
1299 IF (periodic(2) /= 0)
THEN
1300 IF (sb_min(2) >= 0.5_dp)
EXIT loop2_jcell
1301 IF (sb_max(2) < -0.5_dp) cycle loop2_jcell
1305 loop2_icell:
DO icell = -ncell(1), ncell(1)
1306 sb(1) = sb_pbc(1) + real(icell,
dp)
1307 sb_min(1) = sb(1) - sab_max(1)
1308 sb_max(1) = sb(1) + sab_max(1)
1309 IF (periodic(1) /= 0)
THEN
1310 IF (sb_min(1) >= 0.5_dp)
EXIT loop2_icell
1311 IF (sb_max(1) < -0.5_dp) cycle loop2_icell
1317 loop_k:
DO k = 1, nsubcell(3)
1318 loop_j:
DO j = 1, nsubcell(2)
1319 loop_i:
DO i = 1, nsubcell(1)
1323 IF (periodic(3) /= 0)
THEN
1324 IF (sb_max(3) < subcell(i, j, k)%s_min(3))
EXIT loop_k
1325 IF (sb_min(3) >= subcell(i, j, k)%s_max(3)) cycle loop_k
1328 IF (periodic(2) /= 0)
THEN
1329 IF (sb_max(2) < subcell(i, j, k)%s_min(2))
EXIT loop_j
1330 IF (sb_min(2) >= subcell(i, j, k)%s_max(2)) cycle loop_j
1333 IF (periodic(1) /= 0)
THEN
1334 IF (sb_max(1) < subcell(i, j, k)%s_min(1))
EXIT loop_i
1335 IF (sb_min(1) >= subcell(i, j, k)%s_max(1)) cycle loop_i
1338 IF (subcell(i, j, k)%natom == 0) cycle
1340 DO iatom_subcell = 1, subcell(i, j, k)%natom
1341 iatom_local = subcell(i, j, k)%atom_list(iatom_subcell)
1342 iatom = lista(ikind)%list(iatom_local)
1343 atom_a =
atom(ikind)%list(iatom)
1344 IF (my_molecular)
THEN
1345 mol_a =
atom(ikind)%list_a_mol(iatom_local)
1346 IF (mol_a /= mol_b) cycle
1348 IF (my_symmetric)
THEN
1349 IF (atom_a > atom_b)
THEN
1350 include_ab = (
modulo(atom_a + atom_b, 2) /= 0)
1352 include_ab = (
modulo(atom_a + atom_b, 2) == 0)
1354 IF (my_sort_atomb)
THEN
1355 IF ((.NOT. any(atomb_to_keep == atom_b)) .AND. &
1356 (.NOT. any(atomb_to_keep == atom_a)))
THEN
1357 include_ab = .false.
1363 IF (include_ab)
THEN
1364 ra(:) = r_pbc(:, atom_a)
1365 rab(:) = rb(:) - ra(:)
1366 rab2 = rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3)
1367 IF (rab2 < rab2_max)
THEN
1373 rab_pbc(:) =
pbc(rab(:), cell)
1374 IF (sum((rab_pbc - rab)**2) > epsilon(1.0_dp))
THEN
1375 include_ab = .false.
1378 IF (include_ab)
THEN
1380 neighbor_list=kind_a(iatom_local)%neighbor_list, &
1409 DEALLOCATE (lista(ikind)%list)
1414 DEALLOCATE (kind_a, pres_a, pres_b, lista, listb, nlista, nlistb)
1415 DEALLOCATE (index_list)
1422 IF (inode == 1) nentry = nentry + nnode
1426 ALLOCATE (ab_list(1)%nlist_task(nentry))
1427 ab_list(1)%nl_size = nentry
1428 DO iab = 2,
SIZE(ab_list)
1429 ab_list(iab)%nl_size = nentry
1430 ab_list(iab)%nlist_task => ab_list(1)%nlist_task
1439 iab = (ikind - 1)*nkind + jkind
1440 IF (ab_list(iab)%nl_start < 0) ab_list(iab)%nl_start = nentry
1441 IF (ab_list(iab)%nl_end < 0)
THEN
1442 ab_list(iab)%nl_end = nentry
1444 cpassert(ab_list(iab)%nl_end + 1 == nentry)
1445 ab_list(iab)%nl_end = nentry
1450 CALL timestop(handle)
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, 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, 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, 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, ecoul_1c, rho0_s_rs, rho0_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, tb_tblite)
Get the QUICKSTEP environment.