(git:ed6f26b)
Loading...
Searching...
No Matches
qs_neighbor_lists.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Generate the atomic neighbor lists.
10!> \par History
11!> - List rebuild for sab_orb neighbor list (10.09.2002,MK)
12!> - List rebuild for all lists (25.09.2002,MK)
13!> - Row-wise parallelized version (16.06.2003,MK)
14!> - Row- and column-wise parallelized version (19.07.2003,MK)
15!> - bug fix for non-periodic case (23.02.06,MK)
16!> - major refactoring (25.07.10,jhu)
17!> \author Matthias Krack (08.10.1999,26.03.2002,16.06.2003)
18! **************************************************************************************************
27 USE cell_types, ONLY: cell_type,&
28 get_cell,&
29 pbc,&
36 USE cp_output_handling, ONLY: cp_p_file,&
49 USE input_constants, ONLY: &
56 USE kinds, ONLY: default_string_length,&
57 dp,&
58 int_8
59 USE kpoint_types, ONLY: kpoint_type
61 USE mathlib, ONLY: erfc_cutoff
67 USE periodic_table, ONLY: ptable
68 USE physcon, ONLY: bohr
74 USE qs_gcp_types, ONLY: qs_gcp_type
75 USE qs_kind_types, ONLY: get_qs_kind,&
78 USE qs_ks_types, ONLY: get_ks_env,&
81 USE qs_neighbor_list_types, ONLY: &
86 USE string_utilities, ONLY: compress,&
92 USE util, ONLY: locate,&
93 sort
94 USE xtb_types, ONLY: get_xtb_atom_param,&
96#include "./base/base_uses.f90"
97
98 IMPLICIT NONE
99
100 PRIVATE
101
102! **************************************************************************************************
104 INTEGER, DIMENSION(:), POINTER :: list => null(), &
105 list_local_a_index => null(), &
106 list_local_b_index => null(), &
107 list_1d => null(), &
108 list_a_mol => null(), &
109 list_b_mol => null()
110 END TYPE local_atoms_type
111! **************************************************************************************************
112
113 CHARACTER(len=*), PARAMETER, PRIVATE :: modulen = 'qs_neighbor_lists'
114
115 ! private counter, used to version qs neighbor lists
116 INTEGER, SAVE, PRIVATE :: last_qs_neighbor_list_id_nr = 0
117
118 ! Public subroutines
122CONTAINS
123
124! **************************************************************************************************
125!> \brief free the internals of atom2d
126!> \param atom2d ...
127!> \param
128! **************************************************************************************************
129 SUBROUTINE atom2d_cleanup(atom2d)
130 TYPE(local_atoms_type), DIMENSION(:) :: atom2d
131
132 CHARACTER(len=*), PARAMETER :: routinen = 'atom2d_cleanup'
133
134 INTEGER :: handle, ikind
135
136 CALL timeset(routinen, handle)
137 DO ikind = 1, SIZE(atom2d)
138 NULLIFY (atom2d(ikind)%list)
139 IF (ASSOCIATED(atom2d(ikind)%list_local_a_index)) THEN
140 DEALLOCATE (atom2d(ikind)%list_local_a_index)
141 END IF
142 IF (ASSOCIATED(atom2d(ikind)%list_local_b_index)) THEN
143 DEALLOCATE (atom2d(ikind)%list_local_b_index)
144 END IF
145 IF (ASSOCIATED(atom2d(ikind)%list_a_mol)) THEN
146 DEALLOCATE (atom2d(ikind)%list_a_mol)
147 END IF
148 IF (ASSOCIATED(atom2d(ikind)%list_b_mol)) THEN
149 DEALLOCATE (atom2d(ikind)%list_b_mol)
150 END IF
151 IF (ASSOCIATED(atom2d(ikind)%list_1d)) THEN
152 DEALLOCATE (atom2d(ikind)%list_1d)
153 END IF
154 END DO
155 CALL timestop(handle)
156
157 END SUBROUTINE atom2d_cleanup
158
159! **************************************************************************************************
160!> \brief Build some distribution structure of atoms, refactored from build_qs_neighbor_lists
161!> \param atom2d output
162!> \param distribution_1d ...
163!> \param distribution_2d ...
164!> \param atomic_kind_set ...
165!> \param molecule_set ...
166!> \param molecule_only ...
167!> \param particle_set ...
168!> \author JH
169! **************************************************************************************************
170 SUBROUTINE atom2d_build(atom2d, distribution_1d, distribution_2d, &
171 atomic_kind_set, molecule_set, molecule_only, particle_set)
172 TYPE(local_atoms_type), DIMENSION(:) :: atom2d
173 TYPE(distribution_1d_type), POINTER :: distribution_1d
174 TYPE(distribution_2d_type), POINTER :: distribution_2d
175 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
176 TYPE(molecule_type), DIMENSION(:), POINTER :: molecule_set
177 LOGICAL :: molecule_only
178 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
179
180 CHARACTER(len=*), PARAMETER :: routinen = 'atom2d_build'
181
182 INTEGER :: atom_a, handle, ia, iat, iatom, &
183 iatom_local, ikind, imol, natom, &
184 natom_a, natom_local_a, natom_local_b, &
185 nel, nkind
186 INTEGER, ALLOCATABLE, DIMENSION(:) :: atom2mol, atom_of_kind, listindex, &
187 listsort
188 INTEGER, DIMENSION(:), POINTER :: local_cols_array, local_rows_array
189
190 CALL timeset(routinen, handle)
191
192 nkind = SIZE(atomic_kind_set)
193 natom = SIZE(particle_set)
194 CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, atom_of_kind=atom_of_kind)
195
196 IF (molecule_only) THEN
197 ALLOCATE (atom2mol(natom))
198 DO imol = 1, SIZE(molecule_set)
199 DO iat = molecule_set(imol)%first_atom, molecule_set(imol)%last_atom
200 atom2mol(iat) = imol
201 END DO
202 END DO
203 END IF
204
205 DO ikind = 1, nkind
206 NULLIFY (atom2d(ikind)%list)
207 NULLIFY (atom2d(ikind)%list_local_a_index)
208 NULLIFY (atom2d(ikind)%list_local_b_index)
209 NULLIFY (atom2d(ikind)%list_1d)
210 NULLIFY (atom2d(ikind)%list_a_mol)
211 NULLIFY (atom2d(ikind)%list_b_mol)
212
213 CALL get_atomic_kind(atomic_kind_set(ikind), atom_list=atom2d(ikind)%list)
214
215 natom_a = SIZE(atom2d(ikind)%list)
216
217 natom_local_a = distribution_2d%n_local_rows(ikind)
218 natom_local_b = distribution_2d%n_local_cols(ikind)
219 local_rows_array => distribution_2d%local_rows(ikind)%array
220 local_cols_array => distribution_2d%local_cols(ikind)%array
221
222 nel = distribution_1d%n_el(ikind)
223 ALLOCATE (atom2d(ikind)%list_1d(nel))
224 DO iat = 1, nel
225 ia = distribution_1d%list(ikind)%array(iat)
226 atom2d(ikind)%list_1d(iat) = atom_of_kind(ia)
227 END DO
228
229 ALLOCATE (listsort(natom_a), listindex(natom_a))
230 listsort(1:natom_a) = atom2d(ikind)%list(1:natom_a)
231 CALL sort(listsort, natom_a, listindex)
232 ! Block rows
233 IF (natom_local_a > 0) THEN
234 ALLOCATE (atom2d(ikind)%list_local_a_index(natom_local_a))
235 ALLOCATE (atom2d(ikind)%list_a_mol(natom_local_a))
236 atom2d(ikind)%list_a_mol(:) = 0
237
238 ! Build index vector for mapping
239 DO iatom_local = 1, natom_local_a
240 atom_a = local_rows_array(iatom_local)
241 iatom = locate(listsort, atom_a)
242 atom2d(ikind)%list_local_a_index(iatom_local) = listindex(iatom)
243 IF (molecule_only) atom2d(ikind)%list_a_mol(iatom_local) = atom2mol(atom_a)
244 END DO
245
246 END IF
247
248 ! Block columns
249 IF (natom_local_b > 0) THEN
250
251 ALLOCATE (atom2d(ikind)%list_local_b_index(natom_local_b))
252 ALLOCATE (atom2d(ikind)%list_b_mol(natom_local_b))
253 atom2d(ikind)%list_b_mol(:) = 0
254
255 ! Build index vector for mapping
256 DO iatom_local = 1, natom_local_b
257 atom_a = local_cols_array(iatom_local)
258 iatom = locate(listsort, atom_a)
259 atom2d(ikind)%list_local_b_index(iatom_local) = listindex(iatom)
260 IF (molecule_only) atom2d(ikind)%list_b_mol(iatom_local) = atom2mol(atom_a)
261 END DO
262
263 END IF
264
265 DEALLOCATE (listsort, listindex)
266
267 END DO
268
269 CALL timestop(handle)
270
271 END SUBROUTINE atom2d_build
272
273! **************************************************************************************************
274!> \brief Build all the required neighbor lists for Quickstep.
275!> \param qs_env ...
276!> \param para_env ...
277!> \param molecular ...
278!> \param force_env_section ...
279!> \date 28.08.2000
280!> \par History
281!> - Major refactoring (25.07.2010,jhu)
282!> \author MK
283!> \version 1.0
284! **************************************************************************************************
285 SUBROUTINE build_qs_neighbor_lists(qs_env, para_env, molecular, force_env_section)
286 TYPE(qs_environment_type), POINTER :: qs_env
287 TYPE(mp_para_env_type), POINTER :: para_env
288 LOGICAL, OPTIONAL :: molecular
289 TYPE(section_vals_type), POINTER :: force_env_section
290
291 CHARACTER(len=*), PARAMETER :: routinen = 'build_qs_neighbor_lists'
292
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
308 TYPE(all_potential_type), POINTER :: all_potential
309 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
310 TYPE(cell_type), POINTER :: cell
311 TYPE(cp_logger_type), POINTER :: logger
312 TYPE(dft_control_type), POINTER :: dft_control
313 TYPE(distribution_1d_type), POINTER :: distribution_1d
314 TYPE(distribution_2d_type), POINTER :: distribution_2d
315 TYPE(ewald_environment_type), POINTER :: ewald_env
316 TYPE(gth_potential_type), POINTER :: gth_potential
317 TYPE(gto_basis_set_type), POINTER :: aux_basis_set, aux_fit_basis_set, &
318 orb_basis_set, ri_basis_set
319 TYPE(kpoint_type), POINTER :: kpoints
320 TYPE(local_atoms_type), ALLOCATABLE, DIMENSION(:) :: atom2d
321 TYPE(molecule_type), DIMENSION(:), POINTER :: molecule_set
322 TYPE(neighbor_list_set_p_type), DIMENSION(:), POINTER :: saa_list, sab_all, sab_almo, &
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
326 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
327 TYPE(paw_proj_set_type), POINTER :: paw_proj
328 TYPE(qs_dftb_atom_type), POINTER :: dftb_atom
329 TYPE(qs_dispersion_type), POINTER :: dispersion_env
330 TYPE(qs_gcp_type), POINTER :: gcp_env
331 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
332 TYPE(qs_ks_env_type), POINTER :: ks_env
333 TYPE(section_vals_type), POINTER :: hfx_sections, neighbor_list_section
334 TYPE(sgp_potential_type), POINTER :: sgp_potential
335 TYPE(xtb_atom_type), POINTER :: xtb_atom
336
337 CALL timeset(routinen, handle)
338 NULLIFY (logger)
339 logger => cp_get_default_logger()
340
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)
344
345 NULLIFY (sab_orb)
346 NULLIFY (sac_ae)
347 NULLIFY (sac_ppl)
348 NULLIFY (sac_lri)
349 NULLIFY (sap_ppnl)
350 NULLIFY (sap_oce)
351 NULLIFY (sab_se)
352 NULLIFY (sab_lrc)
353 NULLIFY (sab_tbe)
354 NULLIFY (sab_xtbe)
355 NULLIFY (sab_core)
356 NULLIFY (sab_xb)
357 NULLIFY (sab_xtb_pp)
358 NULLIFY (sab_xtb_nonbond)
359 NULLIFY (sab_all)
360 NULLIFY (sab_vdw)
361 NULLIFY (sab_cn)
362 NULLIFY (soo_list)
363 NULLIFY (sab_scp)
364 NULLIFY (sab_almo)
365 NULLIFY (sab_kp)
366 NULLIFY (sab_kp_nosym)
367
368 CALL get_qs_env(qs_env, &
369 ks_env=ks_env, &
370 atomic_kind_set=atomic_kind_set, &
371 qs_kind_set=qs_kind_set, &
372 cell=cell, &
373 kpoints=kpoints, &
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)
379
380 neighbor_list_section => section_vals_get_subs_vals(force_env_section, "DFT%PRINT%NEIGHBOR_LISTS")
381
382 ! This sets the id number of the qs neighbor lists, new lists, means new version
383 ! new version implies new sparsity of the matrices
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)
386
387 CALL get_ks_env(ks_env=ks_env, &
388 sab_orb=sab_orb, &
389 sac_ae=sac_ae, &
390 sac_ppl=sac_ppl, &
391 sac_lri=sac_lri, &
392 sab_vdw=sab_vdw, &
393 sap_ppnl=sap_ppnl, &
394 sap_oce=sap_oce, &
395 sab_se=sab_se, &
396 sab_lrc=sab_lrc, &
397 sab_tbe=sab_tbe, &
398 sab_xtbe=sab_xtbe, &
399 sab_core=sab_core, &
400 sab_xb=sab_xb, &
401 sab_xtb_pp=sab_xtb_pp, &
402 sab_xtb_nonbond=sab_xtb_nonbond, &
403 sab_scp=sab_scp, &
404 sab_all=sab_all, &
405 sab_almo=sab_almo, &
406 sab_kp=sab_kp, &
407 sab_kp_nosym=sab_kp_nosym)
408
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
414 lrigpw = (dft_control%qs_control%method_id == do_method_lrigpw)
415 rigpw = (dft_control%qs_control%method_id == do_method_rigpw)
416 lri_optbas = dft_control%qs_control%lri_optbas
417
418 ! molecular lists
419 molecule_only = .false.
420 IF (PRESENT(molecular)) molecule_only = molecular
421 ! minimum image convention (MIC)
422 mic = molecule_only
423 IF (dokp) THEN
424 ! no MIC for kpoints
425 mic = .false.
426 ELSEIF (nddo) THEN
427 ! enforce MIC for interaction lists in SE
428 mic = .true.
429 END IF
430 pdist = dft_control%qs_control%pairlist_radius
431
432 hfx_sections => section_vals_get_subs_vals(qs_env%input, "DFT%XC%HF")
433 CALL section_vals_get(hfx_sections, explicit=do_hfx)
434
435 CALL get_atomic_kind_set(atomic_kind_set, maxatom=maxatom)
436 CALL get_qs_kind_set(qs_kind_set, paw_atom_present=paw_atom_present, &
437 gth_potential_present=gth_potential_present, &
438 sgp_potential_present=sgp_potential_present, &
439 all_potential_present=all_potential_present)
440
441 CALL section_vals_val_get(qs_env%input, "DFT%SUBCELLS", r_val=subcells)
442
443 ! Allocate work storage
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
451 c_radius(:) = 0.0_dp
452 core_radius(:) = 0.0_dp
453 calpha(:) = 0.0_dp
454 zeff(:) = 0.0_dp
455
456 ALLOCATE (pair_radius(nkind, nkind))
457 IF (gth_potential_present .OR. sgp_potential_present) THEN
458 ALLOCATE (ppl_present(nkind), ppl_radius(nkind))
459 ppl_radius = 0.0_dp
460 ALLOCATE (ppnl_present(nkind), ppnl_radius(nkind))
461 ppnl_radius = 0.0_dp
462 END IF
463 IF (paw_atom_present) THEN
464 ALLOCATE (oce_present(nkind), oce_radius(nkind))
465 oce_radius = 0.0_dp
466 END IF
467 IF (all_potential_present .OR. sgp_potential_present) THEN
468 ALLOCATE (all_present(nkind), all_pot_rad(nkind))
469 all_pot_rad = 0.0_dp
470 END IF
471
472 ! Initialize the local data structures
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)
476
477 DO ikind = 1, nkind
478
479 CALL get_atomic_kind(atomic_kind_set(ikind), atom_list=atom2d(ikind)%list)
480
481 CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set, basis_type="ORB")
482 CALL get_qs_kind(qs_kind_set(ikind), basis_set=aux_basis_set, basis_type="AUX")
483 CALL get_qs_kind(qs_kind_set(ikind), basis_set=aux_fit_basis_set, basis_type="AUX_FIT")
484
485 CALL get_qs_kind(qs_kind_set(ikind), &
486 paw_proj_set=paw_proj, &
487 paw_atom=paw_atom, &
488 all_potential=all_potential, &
489 gth_potential=gth_potential, &
490 sgp_potential=sgp_potential)
491
492 IF (dftb) THEN
493 ! Set the interaction radius for the neighbor lists (DFTB case)
494 ! This includes all interactions (orbitals and short range pair potential) except vdW
495 CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_atom)
496 CALL get_dftb_atom_param(dftb_parameter=dftb_atom, &
497 cutoff=orb_radius(ikind), &
498 defined=orb_present(ikind))
499 ELSE
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))
503 ELSE
504 orb_present(ikind) = .false.
505 END IF
506 END IF
507
508 IF (ASSOCIATED(aux_basis_set)) THEN
509 aux_present(ikind) = .true.
510 ELSE
511 aux_present(ikind) = .false.
512 END IF
513
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))
517 ELSE
518 aux_fit_present(ikind) = .false.
519 END IF
520
521 ! core overlap
522 CALL get_qs_kind(qs_kind_set(ikind), &
523 alpha_core_charge=calpha(ikind), &
524 core_charge_radius=core_radius(ikind), &
525 zeff=zeff(ikind))
526 IF (zeff(ikind) /= 0._dp .AND. calpha(ikind) /= 0._dp) THEN
527 core_present(ikind) = .true.
528 ELSE
529 core_present(ikind) = .false.
530 END IF
531
532 ! Pseudopotentials
533 IF (gth_potential_present .OR. sgp_potential_present) THEN
534 IF (ASSOCIATED(gth_potential)) THEN
535 CALL get_potential(potential=gth_potential, &
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
541 CALL get_potential(potential=sgp_potential, &
542 ppl_present=ppl_present(ikind), &
543 ppl_radius=ppl_radius(ikind), &
544 ppnl_present=ppnl_present(ikind), &
545 ppnl_radius=ppnl_radius(ikind))
546 ELSE
547 ppl_present(ikind) = .false.
548 ppnl_present(ikind) = .false.
549 END IF
550 END IF
551
552 ! GAPW
553 IF (paw_atom_present) THEN
554 IF (paw_atom) THEN
555 oce_present(ikind) = .true.
556 CALL get_paw_proj_set(paw_proj_set=paw_proj, rcprj=oce_radius(ikind))
557 ELSE
558 oce_present(ikind) = .false.
559 END IF
560 END IF
561
562 ! Check the presence of an all electron potential or ERFC potential
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))
573 END IF
574 END IF
575 END IF
576
577 END DO
578
579 ! Build the orbital-orbital overlap neighbor lists
580 IF (pdist < 0.0_dp) THEN
581 pdist = max(plane_distance(1, 0, 0, cell), &
582 plane_distance(0, 1, 0, cell), &
583 plane_distance(0, 0, 1, cell))
584 END IF
585 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius, pdist)
586 CALL build_neighbor_lists(sab_orb, particle_set, atom2d, cell, pair_radius, &
587 mic=mic, subcells=subcells, molecular=molecule_only, nlname="sab_orb")
588 CALL set_ks_env(ks_env=ks_env, sab_orb=sab_orb)
589 CALL write_neighbor_lists(sab_orb, particle_set, cell, para_env, neighbor_list_section, &
590 "/SAB_ORB", "sab_orb", "ORBITAL ORBITAL")
591
592 ! Build orbital-orbital list containing all the pairs, to be used with
593 ! non-symmetric operators. Beware: the cutoff of the orbital-orbital overlap
594 ! might not be optimal. It should be verified for each operator.
595 IF (.NOT. (nddo .OR. dftb .OR. xtb)) THEN
596 CALL build_neighbor_lists(sab_all, particle_set, atom2d, cell, pair_radius, &
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)
599 END IF
600
601 ! Build the core-core overlap neighbor lists
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)
607 CALL write_neighbor_lists(sab_core, particle_set, cell, para_env, neighbor_list_section, &
608 "/SAB_CORE", "sab_core", "CORE CORE")
609 END IF
610
611 IF (dokp) THEN
612 ! We try to guess an integration radius for K-points
613 ! For non-HFX calculations we use the overlap list
614 ! For HFX we use the interaction radius of kinds (ORB or ADMM basis)
615 ! plus a range for the operator
616 IF (do_hfx) THEN
617
618 !case study on the HFX potential: TC, SR or Overlap?
619 CALL section_vals_val_get(hfx_sections, "INTERACTION_POTENTIAL%POTENTIAL_TYPE", i_val=hfx_pot)
620
621 SELECT CASE (hfx_pot)
622 CASE (do_potential_id)
623 roperator = 0.0_dp
625 CALL section_vals_val_get(hfx_sections, "INTERACTION_POTENTIAL%CUTOFF_RADIUS", r_val=roperator)
626 CASE (do_potential_short)
627 CALL section_vals_val_get(hfx_sections, "INTERACTION_POTENTIAL%OMEGA", r_val=omega)
628 CALL section_vals_val_get(hfx_sections, "SCREENING%EPS_SCHWARZ", r_val=eps_schwarz)
629 CALL erfc_cutoff(eps_schwarz, omega, roperator)
630 CASE DEFAULT
631 cpabort("HFX potential not available for K-points (NYI)")
632 END SELECT
633
634 IF (dft_control%do_admm) THEN
635 CALL pair_radius_setup(aux_fit_present, aux_fit_present, aux_fit_radius, aux_fit_radius, &
636 pair_radius)
637
638 !We cannot accept a pair radius smaller than the ORB overlap, for sanity reasons
639 ALLOCATE (pair_radius_lb(nkind, nkind))
640 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius_lb)
641 DO jkind = 1, nkind
642 DO ikind = 1, nkind
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
645 END DO
646 END DO
647 ELSE
648 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius)
649 END IF
650 pair_radius = pair_radius + cutoff_screen_factor*roperator
651 ELSE
652 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius)
653 END IF
654 CALL build_neighbor_lists(sab_kp, particle_set, atom2d, cell, pair_radius, &
655 subcells=subcells, nlname="sab_kp")
656 CALL set_ks_env(ks_env=ks_env, sab_kp=sab_kp)
657
658 IF (do_hfx) THEN
659 CALL build_neighbor_lists(sab_kp_nosym, particle_set, atom2d, cell, pair_radius, &
660 subcells=subcells, nlname="sab_kp_nosym", symmetric=.false.)
661 CALL set_ks_env(ks_env=ks_env, sab_kp_nosym=sab_kp_nosym)
662 END IF
663 END IF
664
665 ! Build orbital GTH-PPL operator overlap list
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)
669 CALL build_neighbor_lists(sac_ppl, particle_set, atom2d, cell, pair_radius, &
670 subcells=subcells, operator_type="ABC", nlname="sac_ppl")
671 CALL set_ks_env(ks_env=ks_env, sac_ppl=sac_ppl)
672 CALL write_neighbor_lists(sac_ppl, particle_set, cell, para_env, neighbor_list_section, &
673 "/SAC_PPL", "sac_ppl", "ORBITAL GTH-PPL")
674 IF (lrigpw) THEN
675 IF (qs_env%lri_env%ppl_ri) THEN
676 CALL build_neighbor_lists(sac_lri, particle_set, atom2d, cell, pair_radius, &
677 subcells=subcells, symmetric=.false., operator_type="PP", nlname="sac_lri")
678 CALL set_ks_env(ks_env=ks_env, sac_lri=sac_lri)
679 END IF
680 END IF
681 END IF
682
683 IF (any(ppnl_present)) THEN
684 CALL pair_radius_setup(orb_present, ppnl_present, orb_radius, ppnl_radius, pair_radius)
685 CALL build_neighbor_lists(sap_ppnl, particle_set, atom2d, cell, pair_radius, &
686 subcells=subcells, operator_type="ABBA", nlname="sap_ppnl")
687 CALL set_ks_env(ks_env=ks_env, sap_ppnl=sap_ppnl)
688 CALL write_neighbor_lists(sap_ppnl, particle_set, cell, para_env, neighbor_list_section, &
689 "/SAP_PPNL", "sap_ppnl", "ORBITAL GTH-PPNL")
690 END IF
691 END IF
692
693 IF (paw_atom_present) THEN
694 ! Build orbital-GAPW projector overlap list
695 IF (any(oce_present)) THEN
696 CALL pair_radius_setup(orb_present, oce_present, orb_radius, oce_radius, pair_radius)
697 CALL build_neighbor_lists(sap_oce, particle_set, atom2d, cell, pair_radius, &
698 subcells=subcells, operator_type="ABBA", nlname="sap_oce")
699 CALL set_ks_env(ks_env=ks_env, sap_oce=sap_oce)
700 CALL write_neighbor_lists(sap_oce, particle_set, cell, para_env, neighbor_list_section, &
701 "/SAP_OCE", "sap_oce", "ORBITAL(A) PAW-PRJ")
702 END IF
703 END IF
704
705 ! Build orbital-ERFC potential list
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)
709 CALL build_neighbor_lists(sac_ae, particle_set, atom2d, cell, pair_radius, &
710 subcells=subcells, operator_type="ABC", nlname="sac_ae")
711 CALL set_ks_env(ks_env=ks_env, sac_ae=sac_ae)
712 CALL write_neighbor_lists(sac_ae, particle_set, cell, para_env, neighbor_list_section, &
713 "/SAC_AE", "sac_ae", "ORBITAL ERFC POTENTIAL")
714 END IF
715 END IF
716
717 IF (nddo) THEN
718 ! Semi-empirical neighbor lists
719 default_present = .true.
720 c_radius = dft_control%qs_control%se_control%cutoff_cou
721 ! Build the neighbor lists for the Hartree terms
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
724 ! Use MIC for the periodic code of GKS
725 CALL build_neighbor_lists(sab_se, particle_set, atom2d, cell, pair_radius, mic=mic, &
726 subcells=subcells, nlname="sab_se")
727 ELSE
728 CALL build_neighbor_lists(sab_se, particle_set, atom2d, cell, pair_radius, &
729 subcells=subcells, nlname="sab_se")
730 END IF
731 CALL set_ks_env(ks_env=ks_env, sab_se=sab_se)
732 CALL write_neighbor_lists(sab_se, particle_set, cell, para_env, neighbor_list_section, &
733 "/SAB_SE", "sab_se", "HARTREE INTERACTIONS")
734
735 ! If requested build the SE long-range correction neighbor list
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)
740 CALL build_neighbor_lists(sab_lrc, particle_set, atom2d, cell, pair_radius, &
741 subcells=subcells, nlname="sab_lrc")
742 CALL set_ks_env(ks_env=ks_env, sab_lrc=sab_lrc)
743 CALL write_neighbor_lists(sab_lrc, particle_set, cell, para_env, neighbor_list_section, &
744 "/SAB_LRC", "sab_lrc", "SE LONG-RANGE CORRECTION")
745 END IF
746 END IF
747
748 IF (dftb) THEN
749 ! Build the neighbor lists for the DFTB Ewald methods
750 IF (dft_control%qs_control%dftb_control%do_ewald) THEN
751 CALL get_qs_env(qs_env=qs_env, ewald_env=ewald_env)
752 CALL ewald_env_get(ewald_env, rcut=rcut)
753 c_radius = rcut
754 CALL pair_radius_setup(orb_present, orb_present, c_radius, c_radius, pair_radius)
755 CALL build_neighbor_lists(sab_tbe, particle_set, atom2d, cell, pair_radius, mic=mic, &
756 subcells=subcells, nlname="sab_tbe")
757 CALL set_ks_env(ks_env=ks_env, sab_tbe=sab_tbe)
758 END IF
759
760 ! Build the neighbor lists for the DFTB vdW pair potential
761 IF (dft_control%qs_control%dftb_control%dispersion) THEN
762 IF (dft_control%qs_control%dftb_control%dispersion_type == dispersion_uff) THEN
763 DO ikind = 1, nkind
764 CALL get_qs_kind(qs_kind_set(ikind), dftb_parameter=dftb_atom)
765 CALL get_dftb_atom_param(dftb_parameter=dftb_atom, rcdisp=c_radius(ikind))
766 END DO
767 default_present = .true.
768 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
769 CALL build_neighbor_lists(sab_vdw, particle_set, atom2d, cell, pair_radius, &
770 subcells=subcells, nlname="sab_vdw")
771 CALL set_ks_env(ks_env=ks_env, sab_vdw=sab_vdw)
772 END IF
773 END IF
774 END IF
775
776 IF (xtb) THEN
777 ! Build the neighbor lists for the xTB Ewald method
778 IF (dft_control%qs_control%xtb_control%do_ewald) THEN
779 CALL get_qs_env(qs_env=qs_env, ewald_env=ewald_env)
780 CALL ewald_env_get(ewald_env, rcut=rcut)
781 c_radius = rcut
782 CALL pair_radius_setup(orb_present, orb_present, c_radius, c_radius, pair_radius)
783 CALL build_neighbor_lists(sab_tbe, particle_set, atom2d, cell, pair_radius, mic=mic, &
784 subcells=subcells, nlname="sab_tbe")
785 CALL set_ks_env(ks_env=ks_env, sab_tbe=sab_tbe)
786 END IF
787 ! Repulsive Potential
788 pair_radius(1:nkind, 1:nkind) = dft_control%qs_control%xtb_control%rcpair(1:nkind, 1:nkind)
789 default_present = .true.
790 CALL build_neighbor_lists(sab_xtb_pp, particle_set, atom2d, cell, pair_radius, &
791 subcells=subcells, nlname="sab_xtb_pp")
792 CALL set_ks_env(ks_env=ks_env, sab_xtb_pp=sab_xtb_pp)
793 ! SR part of Coulomb interaction
794 DO ikind = 1, nkind
795 CALL get_qs_kind(qs_kind_set(ikind), xtb_parameter=xtb_atom)
796 CALL get_xtb_atom_param(xtb_parameter=xtb_atom, rcut=c_radius(ikind))
797 END DO
798 default_present = .true.
799 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
800 CALL build_neighbor_lists(sab_xtbe, particle_set, atom2d, cell, pair_radius, &
801 subcells=subcells, nlname="sab_xtbe")
802 CALL set_ks_env(ks_env=ks_env, sab_xtbe=sab_xtbe)
803 ! XB list
804 ALLOCATE (xb1_atom(nkind), xb2_atom(nkind))
805 c_radius = 0.5_dp*dft_control%qs_control%xtb_control%xb_radius
806 DO ikind = 1, nkind
807 CALL get_atomic_kind(atomic_kind_set(ikind), z=zat)
808 IF (zat == 17 .OR. zat == 35 .OR. zat == 53 .OR. zat == 85) THEN
809 xb1_atom(ikind) = .true.
810 ELSE
811 xb1_atom(ikind) = .false.
812 END IF
813 IF (zat == 7 .OR. zat == 8 .OR. zat == 15 .OR. zat == 16) THEN
814 xb2_atom(ikind) = .true.
815 ELSE
816 xb2_atom(ikind) = .false.
817 END IF
818 END DO
819 CALL pair_radius_setup(xb1_atom, xb2_atom, c_radius, c_radius, pair_radius)
820 CALL build_neighbor_lists(sab_xb, particle_set, atom2d, cell, pair_radius, &
821 symmetric=.false., subcells=subcells, operator_type="PP", nlname="sab_xb")
822 CALL set_ks_env(ks_env=ks_env, sab_xb=sab_xb)
823 CALL write_neighbor_lists(sab_xb, particle_set, cell, para_env, neighbor_list_section, &
824 "/SAB_XB", "sab_xb", "XB bonding")
825
826 ! nonbonded interactions list
827 IF (dft_control%qs_control%xtb_control%do_nonbonded) THEN
828 ngp = SIZE(dft_control%qs_control%xtb_control%nonbonded%pot)
829 ALLOCATE (nonbond1_atom(nkind), nonbond2_atom(nkind))
830 nonbond1_atom = .false.
831 nonbond2_atom = .false.
832 DO ingp = 1, ngp
833 DO ikind = 1, nkind
834 rcut = sqrt(dft_control%qs_control%xtb_control%nonbonded%pot(ingp)%pot%rcutsq)
835 c_radius = rcut
836 CALL get_atomic_kind(atomic_kind_set(ikind), element_symbol=element_symbol)
837 CALL uppercase(element_symbol)
838 IF (trim(dft_control%qs_control%xtb_control%nonbonded%pot(ingp)%pot%at1) == trim(element_symbol)) THEN
839 nonbond1_atom(ikind) = .true.
840 DO jkind = 1, nkind
841 CALL get_atomic_kind(atomic_kind_set(jkind), element_symbol=element_symbol2)
842 CALL uppercase(element_symbol2)
843 IF (trim(dft_control%qs_control%xtb_control%nonbonded%pot(ingp)%pot%at2) == trim(element_symbol2)) THEN
844 nonbond2_atom(jkind) = .true.
845 END IF
846 END DO
847 END IF
848 END DO
849 CALL pair_radius_setup(nonbond1_atom, nonbond2_atom, c_radius, c_radius, pair_radius)
850 CALL build_neighbor_lists(sab_xtb_nonbond, particle_set, atom2d, cell, pair_radius, &
851 symmetric=.false., subcells=subcells, operator_type="PP", nlname="sab_xtb_nonbond")
852 CALL set_ks_env(ks_env=ks_env, sab_xtb_nonbond=sab_xtb_nonbond)
853 CALL write_neighbor_lists(sab_xtb_nonbond, particle_set, cell, para_env, neighbor_list_section, &
854 "/SAB_XTB_NONBOND", "sab_xtb_nonbond", "XTB NONBONDED INTERACTIONS")
855 END DO
856 END IF
857 END IF
858
859 ! Build the neighbor lists for the vdW pair potential
860 CALL get_qs_env(qs_env=qs_env, dispersion_env=dispersion_env)
861 sab_vdw => dispersion_env%sab_vdw
862 sab_cn => dispersion_env%sab_cn
863 IF (dispersion_env%type == xc_vdw_fun_pairpot .OR. xtb) THEN
864 IF (dispersion_env%pp_type == vdw_pairpot_dftd4) THEN
865 c_radius(:) = dispersion_env%rc_d4
866 ELSE
867 c_radius(:) = dispersion_env%rc_disp
868 END IF
869 default_present = .true. !include all atoms in vdW (even without basis)
870 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
871 CALL build_neighbor_lists(sab_vdw, particle_set, atom2d, cell, pair_radius, &
872 subcells=subcells, operator_type="PP", nlname="sab_vdw")
873 dispersion_env%sab_vdw => sab_vdw
874
875 ! Build the neighbor lists for coordination numbers as needed by the DFT-D3/D4 method
876 ! This is also needed for the xTB Hamiltonian
877 DO ikind = 1, nkind
878 CALL get_atomic_kind(atomic_kind_set(ikind), z=zat)
879 c_radius(ikind) = 4._dp*ptable(zat)%covalent_radius*bohr
880 END DO
881 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
882 CALL build_neighbor_lists(sab_cn, particle_set, atom2d, cell, pair_radius, &
883 subcells=subcells, operator_type="PP", nlname="sab_cn")
884 dispersion_env%sab_cn => sab_cn
885 END IF
886
887 ! Build the neighbor lists for the gCP pair potential
888 NULLIFY (gcp_env)
889 CALL get_qs_env(qs_env=qs_env, gcp_env=gcp_env)
890 IF (ASSOCIATED(gcp_env)) THEN
891 IF (gcp_env%do_gcp) THEN
892 sab_gcp => gcp_env%sab_gcp
893 DO ikind = 1, nkind
894 c_radius(ikind) = gcp_env%gcp_kind(ikind)%rcsto
895 END DO
896 CALL pair_radius_setup(orb_present, orb_present, c_radius, c_radius, pair_radius)
897 CALL build_neighbor_lists(sab_gcp, particle_set, atom2d, cell, pair_radius, &
898 subcells=subcells, operator_type="PP", nlname="sab_gcp")
899 gcp_env%sab_gcp => sab_gcp
900 ELSE
901 NULLIFY (gcp_env%sab_gcp)
902 END IF
903 END IF
904
905 IF (lrigpw .OR. lri_optbas) THEN
906 ! set neighborlists in lri_env environment
907 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius)
908 soo_list => qs_env%lri_env%soo_list
909 CALL build_neighbor_lists(soo_list, particle_set, atom2d, cell, pair_radius, &
910 mic=mic, molecular=molecule_only, subcells=subcells, nlname="soo_list")
911 qs_env%lri_env%soo_list => soo_list
912 CALL write_neighbor_lists(soo_list, particle_set, cell, para_env, neighbor_list_section, &
913 "/SOO_LIST", "soo_list", "ORBITAL ORBITAL (RI)")
914 ELSEIF (rigpw) THEN
915 ALLOCATE (ri_present(nkind), ri_radius(nkind))
916 ri_present = .false.
917 ri_radius = 0.0_dp
918 DO ikind = 1, nkind
919 CALL get_qs_kind(qs_kind_set(ikind), basis_set=ri_basis_set, basis_type="RI_HXC")
920 IF (ASSOCIATED(ri_basis_set)) THEN
921 ri_present(ikind) = .true.
922 CALL get_gto_basis_set(gto_basis_set=ri_basis_set, kind_radius=ri_radius(ikind))
923 ELSE
924 ri_present(ikind) = .false.
925 END IF
926 END DO
927 ! set neighborlists in lri_env environment
928 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius)
929 soo_list => qs_env%lri_env%soo_list
930 CALL build_neighbor_lists(soo_list, particle_set, atom2d, cell, pair_radius, &
931 mic=mic, molecular=molecule_only, subcells=subcells, nlname="soo_list")
932 qs_env%lri_env%soo_list => soo_list
933 !
934 CALL pair_radius_setup(ri_present, ri_present, ri_radius, ri_radius, pair_radius)
935 saa_list => qs_env%lri_env%saa_list
936 CALL build_neighbor_lists(saa_list, particle_set, atom2d, cell, pair_radius, &
937 mic=mic, molecular=molecule_only, subcells=subcells, nlname="saa_list")
938 qs_env%lri_env%saa_list => saa_list
939 !
940 CALL pair_radius_setup(ri_present, orb_present, ri_radius, orb_radius, pair_radius)
941 soa_list => qs_env%lri_env%soa_list
942 CALL build_neighbor_lists(soa_list, particle_set, atom2d, cell, pair_radius, &
943 mic=mic, symmetric=.false., molecular=molecule_only, &
944 subcells=subcells, operator_type="ABC", nlname="saa_list")
945 qs_env%lri_env%soa_list => soa_list
946 END IF
947
948 ! Build the neighbor lists for the ALMO delocalization
949 IF (almo) THEN
950 DO ikind = 1, nkind
951 CALL get_atomic_kind(atomic_kind_set(ikind), rcov=almo_rcov, rvdw=almo_rvdw)
952 ! multiply the radius by some hard-coded number
953 c_radius(ikind) = max(almo_rcov, almo_rvdw)*bohr* &
955 END DO
956 default_present = .true. !include all atoms (even without basis)
957 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
958 CALL build_neighbor_lists(sab_almo, particle_set, atom2d, cell, pair_radius, &
959 subcells=subcells, operator_type="PP", nlname="sab_almo")
960 CALL set_ks_env(ks_env=ks_env, sab_almo=sab_almo)
961 END IF
962
963 ! Print particle distribution
964 print_key_path = "PRINT%DISTRIBUTION"
965 IF (btest(cp_print_key_should_output(logger%iter_info, force_env_section, &
966 print_key_path), &
967 cp_p_file)) THEN
968 iw = cp_print_key_unit_nr(logger=logger, &
969 basis_section=force_env_section, &
970 print_key_path=print_key_path, &
971 extension=".out")
972 CALL write_neighbor_distribution(sab_orb, qs_kind_set, iw, para_env)
973 CALL cp_print_key_finished_output(unit_nr=iw, &
974 logger=logger, &
975 basis_section=force_env_section, &
976 print_key_path=print_key_path)
977 END IF
978
979 ! Release work storage
980 CALL atom2d_cleanup(atom2d)
981
982 DEALLOCATE (atom2d)
983 DEALLOCATE (orb_present, default_present, core_present)
984 DEALLOCATE (orb_radius, aux_fit_radius, c_radius, core_radius)
985 DEALLOCATE (calpha, zeff)
986 DEALLOCATE (pair_radius)
987 IF (gth_potential_present .OR. sgp_potential_present) THEN
988 DEALLOCATE (ppl_present, ppl_radius)
989 DEALLOCATE (ppnl_present, ppnl_radius)
990 END IF
991 IF (paw_atom_present) THEN
992 DEALLOCATE (oce_present, oce_radius)
993 END IF
994 IF (all_potential_present .OR. sgp_potential_present) THEN
995 DEALLOCATE (all_present, all_pot_rad)
996 END IF
997
998 CALL timestop(handle)
999
1000 END SUBROUTINE build_qs_neighbor_lists
1001
1002! **************************************************************************************************
1003!> \brief Build simple pair neighbor lists.
1004!> \param ab_list ...
1005!> \param particle_set ...
1006!> \param atom ...
1007!> \param cell ...
1008!> \param pair_radius ...
1009!> \param subcells ...
1010!> \param mic ...
1011!> \param symmetric ...
1012!> \param molecular ...
1013!> \param subset_of_mol ...
1014!> \param current_subset ...
1015!> \param operator_type ...
1016!> \param nlname ...
1017!> \param atomb_to_keep the list of atom indices to keep for pairs from the atom2d%b_list
1018!> \date 20.03.2002
1019!> \par History
1020!> - Major refactoring (25.07.2010,jhu)
1021!> - Added option to filter out atoms from list_b (08.2018, A. Bussy)
1022!> \author MK
1023!> \version 2.0
1024! **************************************************************************************************
1025 SUBROUTINE build_neighbor_lists(ab_list, particle_set, atom, cell, pair_radius, subcells, &
1026 mic, symmetric, molecular, subset_of_mol, current_subset, &
1027 operator_type, nlname, atomb_to_keep)
1028
1029 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1030 POINTER :: ab_list
1031 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1032 TYPE(local_atoms_type), DIMENSION(:), INTENT(IN) :: atom
1033 TYPE(cell_type), POINTER :: cell
1034 REAL(dp), DIMENSION(:, :), INTENT(IN) :: pair_radius
1035 REAL(dp), INTENT(IN) :: subcells
1036 LOGICAL, INTENT(IN), OPTIONAL :: mic, symmetric, molecular
1037 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: subset_of_mol
1038 INTEGER, OPTIONAL :: current_subset
1039 CHARACTER(LEN=*), INTENT(IN), OPTIONAL :: operator_type
1040 CHARACTER(LEN=*), INTENT(IN) :: nlname
1041 INTEGER, DIMENSION(:), INTENT(IN), OPTIONAL :: atomb_to_keep
1042
1043 CHARACTER(len=*), PARAMETER :: routinen = 'build_neighbor_lists'
1044
1045 TYPE local_lists
1046 INTEGER, DIMENSION(:), POINTER :: list
1047 END TYPE local_lists
1048
1049 INTEGER :: atom_a, atom_b, handle, i, iab, iatom, iatom_local, &
1050 iatom_subcell, icell, ikind, j, jatom, jatom_local, jcell, jkind, k, &
1051 kcell, maxat, mol_a, mol_b, nkind, otype, natom, inode, nnode, nentry
1052 INTEGER, DIMENSION(3) :: cell_b, ncell, nsubcell, periodic
1053 INTEGER, DIMENSION(:), POINTER :: index_list
1054 LOGICAL :: include_ab, my_mic, &
1055 my_molecular, my_symmetric, my_sort_atomb
1056 LOGICAL, ALLOCATABLE, DIMENSION(:) :: pres_a, pres_b
1057 REAL(dp) :: rab2, rab2_max, rab_max, rabm, deth, subcell_scale
1058 REAL(dp), DIMENSION(3) :: r, rab, ra, rb, sab_max, sb, &
1059 sb_pbc, sb_min, sb_max, rab_pbc, pd, sab_max_guard
1060 INTEGER, ALLOCATABLE, DIMENSION(:) :: nlista, nlistb
1061 TYPE(local_lists), DIMENSION(:), POINTER :: lista, listb
1062 TYPE(neighbor_list_p_type), &
1063 ALLOCATABLE, DIMENSION(:) :: kind_a
1064 TYPE(neighbor_list_set_type), POINTER :: neighbor_list_set
1065 TYPE(subcell_type), DIMENSION(:, :, :), &
1066 POINTER :: subcell
1067 REAL(kind=dp), DIMENSION(:, :), ALLOCATABLE :: r_pbc
1069 DIMENSION(:), POINTER :: nl_iterator
1070
1071 CALL timeset(routinen//"_"//trim(nlname), handle)
1072
1073 ! input options
1074 my_mic = .false.
1075 IF (PRESENT(mic)) my_mic = mic
1076 my_symmetric = .true.
1077 IF (PRESENT(symmetric)) my_symmetric = symmetric
1078 my_molecular = .false.
1079 ! if we have a molecular NL, MIC has to be used
1080 IF (PRESENT(molecular)) my_molecular = molecular
1081 ! check for operator types
1082 IF (PRESENT(operator_type)) THEN
1083 SELECT CASE (operator_type)
1084 CASE ("AB")
1085 otype = 1 ! simple overlap
1086 CASE ("ABC")
1087 otype = 2 ! for three center operators
1088 cpassert(.NOT. my_molecular)
1089 my_symmetric = .false.
1090 CASE ("ABBA")
1091 otype = 3 ! for separable nonlocal operators
1092 my_symmetric = .false.
1093 CASE ("PP")
1094 otype = 4 ! simple atomic pair potential list
1095 CASE default
1096 cpabort("")
1097 END SELECT
1098 ELSE
1099 ! default is a simple AB neighbor list
1100 otype = 1
1101 END IF
1102 my_sort_atomb = .false.
1103 IF (PRESENT(atomb_to_keep)) THEN
1104 my_sort_atomb = .true.
1105 END IF
1106
1107 nkind = SIZE(atom)
1108 ! Deallocate the old neighbor list structure
1109 CALL release_neighbor_list_sets(ab_list)
1110 ! Allocate and initialize the new neighbor list structure
1111 ALLOCATE (ab_list(nkind*nkind))
1112 DO iab = 1, SIZE(ab_list)
1113 NULLIFY (ab_list(iab)%neighbor_list_set)
1114 ab_list(iab)%nl_size = -1
1115 ab_list(iab)%nl_start = -1
1116 ab_list(iab)%nl_end = -1
1117 NULLIFY (ab_list(iab)%nlist_task)
1118 END DO
1119
1120 ! Allocate and initialize the kind availability
1121 ALLOCATE (pres_a(nkind), pres_b(nkind))
1122 DO ikind = 1, nkind
1123 pres_a(ikind) = any(pair_radius(ikind, :) > 0._dp)
1124 pres_b(ikind) = any(pair_radius(:, ikind) > 0._dp)
1125 END DO
1126
1127 ! create a copy of the pbc'ed coordinates
1128 natom = SIZE(particle_set)
1129 ALLOCATE (r_pbc(3, natom))
1130 DO i = 1, natom
1131 r_pbc(1:3, i) = pbc(particle_set(i)%r(1:3), cell)
1132 END DO
1133
1134 ! setup the local lists of atoms
1135 maxat = 0
1136 DO ikind = 1, nkind
1137 maxat = max(maxat, SIZE(atom(ikind)%list))
1138 END DO
1139 ALLOCATE (index_list(maxat))
1140 DO i = 1, maxat
1141 index_list(i) = i
1142 END DO
1143 ALLOCATE (lista(nkind), listb(nkind), nlista(nkind), nlistb(nkind))
1144 nlista = 0
1145 nlistb = 0
1146 DO ikind = 1, nkind
1147 NULLIFY (lista(ikind)%list, listb(ikind)%list)
1148 SELECT CASE (otype)
1149 CASE (1)
1150 IF (ASSOCIATED(atom(ikind)%list_local_a_index)) THEN
1151 lista(ikind)%list => atom(ikind)%list_local_a_index
1152 nlista(ikind) = SIZE(lista(ikind)%list)
1153 END IF
1154 IF (ASSOCIATED(atom(ikind)%list_local_b_index)) THEN
1155 listb(ikind)%list => atom(ikind)%list_local_b_index
1156 nlistb(ikind) = SIZE(listb(ikind)%list)
1157 END IF
1158 CASE (2)
1159 IF (ASSOCIATED(atom(ikind)%list_local_a_index)) THEN
1160 lista(ikind)%list => atom(ikind)%list_local_a_index
1161 nlista(ikind) = SIZE(lista(ikind)%list)
1162 END IF
1163 nlistb(ikind) = SIZE(atom(ikind)%list)
1164 listb(ikind)%list => index_list
1165 CASE (3)
1166 CALL combine_lists(lista(ikind)%list, nlista(ikind), ikind, atom)
1167 nlistb(ikind) = SIZE(atom(ikind)%list)
1168 listb(ikind)%list => index_list
1169 CASE (4)
1170 nlista(ikind) = SIZE(atom(ikind)%list_1d)
1171 lista(ikind)%list => atom(ikind)%list_1d
1172 nlistb(ikind) = SIZE(atom(ikind)%list)
1173 listb(ikind)%list => index_list
1174 CASE default
1175 cpabort("")
1176 END SELECT
1177 END DO
1178
1179 ! Determine max. number of local atoms
1180 maxat = 0
1181 DO ikind = 1, nkind
1182 maxat = max(maxat, nlista(ikind), nlistb(ikind))
1183 END DO
1184 ALLOCATE (kind_a(2*maxat))
1185
1186 ! Load informations about the simulation cell
1187 CALL get_cell(cell=cell, periodic=periodic, deth=deth)
1188
1189 ! Loop over all atomic kind pairs
1190 DO ikind = 1, nkind
1191 IF (.NOT. pres_a(ikind)) cycle
1192
1193 DO jkind = 1, nkind
1194 IF (.NOT. pres_b(jkind)) cycle
1195
1196 iab = ikind + nkind*(jkind - 1)
1197
1198 ! Calculate the square of the maximum interaction distance
1199 IF (pair_radius(ikind, jkind) <= 0._dp) cycle
1200 rab_max = pair_radius(ikind, jkind)
1201 IF (otype == 3) THEN
1202 ! Calculate the square of the maximum interaction distance
1203 ! for sac_max / ncell this must be the maximum over all kinds
1204 ! to be correct for three center terms involving different kinds
1205 rabm = maxval(pair_radius(:, jkind))
1206 ELSE
1207 rabm = rab_max
1208 END IF
1209 rab2_max = rabm*rabm
1210
1211 pd(1) = plane_distance(1, 0, 0, cell)
1212 pd(2) = plane_distance(0, 1, 0, cell)
1213 pd(3) = plane_distance(0, 0, 1, cell)
1214
1215 sab_max = rabm/pd
1216 sab_max_guard = 15.0_dp/pd
1217
1218 ! It makes sense to have fewer subcells for larger systems
1219 subcell_scale = ((125.0_dp**3)/deth)**(1.0_dp/6.0_dp)
1220
1221 ! guess the number of subcells for optimal performance,
1222 ! guard against crazy stuff triggered by very small rabm
1223 nsubcell(:) = int(max(1.0_dp, min(0.5_dp*subcells*subcell_scale/sab_max(:), &
1224 0.5_dp*subcells*subcell_scale/sab_max_guard(:))))
1225
1226 ! number of image cells to be considered
1227 ncell(:) = (int(sab_max(:)) + 1)*periodic(:)
1228
1229 CALL allocate_neighbor_list_set(neighbor_list_set=ab_list(iab)%neighbor_list_set, &
1230 symmetric=my_symmetric)
1231 neighbor_list_set => ab_list(iab)%neighbor_list_set
1232
1233 DO iatom_local = 1, nlista(ikind)
1234 iatom = lista(ikind)%list(iatom_local)
1235 atom_a = atom(ikind)%list(iatom)
1236 CALL add_neighbor_list(neighbor_list_set=neighbor_list_set, &
1237 atom=atom_a, &
1238 neighbor_list=kind_a(iatom_local)%neighbor_list)
1239 END DO
1240
1241 CALL allocate_subcell(subcell, nsubcell)
1242 DO iatom_local = 1, nlista(ikind)
1243 iatom = lista(ikind)%list(iatom_local)
1244 atom_a = atom(ikind)%list(iatom)
1245 r = r_pbc(:, atom_a)
1246 CALL give_ijk_subcell(r, i, j, k, cell, nsubcell)
1247 subcell(i, j, k)%natom = subcell(i, j, k)%natom + 1
1248 END DO
1249 DO k = 1, nsubcell(3)
1250 DO j = 1, nsubcell(2)
1251 DO i = 1, nsubcell(1)
1252 maxat = subcell(i, j, k)%natom + subcell(i, j, k)%natom/10
1253 ALLOCATE (subcell(i, j, k)%atom_list(maxat))
1254 subcell(i, j, k)%natom = 0
1255 END DO
1256 END DO
1257 END DO
1258 DO iatom_local = 1, nlista(ikind)
1259 iatom = lista(ikind)%list(iatom_local)
1260 atom_a = atom(ikind)%list(iatom)
1261 r = r_pbc(:, atom_a)
1262 CALL give_ijk_subcell(r, i, j, k, cell, nsubcell)
1263 subcell(i, j, k)%natom = subcell(i, j, k)%natom + 1
1264 subcell(i, j, k)%atom_list(subcell(i, j, k)%natom) = iatom_local
1265 END DO
1266
1267 DO jatom_local = 1, nlistb(jkind)
1268 jatom = listb(jkind)%list(jatom_local)
1269 atom_b = atom(jkind)%list(jatom)
1270 IF (my_sort_atomb .AND. .NOT. my_symmetric) THEN
1271 IF (.NOT. any(atomb_to_keep == atom_b)) cycle
1272 END IF
1273 IF (my_molecular) THEN
1274 mol_b = atom(jkind)%list_b_mol(jatom_local)
1275 IF (PRESENT(subset_of_mol)) THEN
1276 IF (subset_of_mol(mol_b) .NE. current_subset) cycle
1277 END IF
1278 END IF
1279 r = r_pbc(:, atom_b)
1280 CALL real_to_scaled(sb_pbc(:), r(:), cell)
1281
1282 loop2_kcell: DO kcell = -ncell(3), ncell(3)
1283 sb(3) = sb_pbc(3) + real(kcell, dp)
1284 sb_min(3) = sb(3) - sab_max(3)
1285 sb_max(3) = sb(3) + sab_max(3)
1286 IF (periodic(3) /= 0) THEN
1287 IF (sb_min(3) >= 0.5_dp) EXIT loop2_kcell
1288 IF (sb_max(3) < -0.5_dp) cycle loop2_kcell
1289 END IF
1290 cell_b(3) = kcell
1291
1292 loop2_jcell: DO jcell = -ncell(2), ncell(2)
1293 sb(2) = sb_pbc(2) + real(jcell, dp)
1294 sb_min(2) = sb(2) - sab_max(2)
1295 sb_max(2) = sb(2) + sab_max(2)
1296 IF (periodic(2) /= 0) THEN
1297 IF (sb_min(2) >= 0.5_dp) EXIT loop2_jcell
1298 IF (sb_max(2) < -0.5_dp) cycle loop2_jcell
1299 END IF
1300 cell_b(2) = jcell
1301
1302 loop2_icell: DO icell = -ncell(1), ncell(1)
1303 sb(1) = sb_pbc(1) + real(icell, dp)
1304 sb_min(1) = sb(1) - sab_max(1)
1305 sb_max(1) = sb(1) + sab_max(1)
1306 IF (periodic(1) /= 0) THEN
1307 IF (sb_min(1) >= 0.5_dp) EXIT loop2_icell
1308 IF (sb_max(1) < -0.5_dp) cycle loop2_icell
1309 END IF
1310 cell_b(1) = icell
1311
1312 CALL scaled_to_real(rb, sb, cell)
1313
1314 loop_k: DO k = 1, nsubcell(3)
1315 loop_j: DO j = 1, nsubcell(2)
1316 loop_i: DO i = 1, nsubcell(1)
1317
1318 ! FIXME for non-periodic systems, the whole subcell trick is skipped
1319 ! yielding a Natom**2 pair list build.
1320 IF (periodic(3) /= 0) THEN
1321 IF (sb_max(3) < subcell(i, j, k)%s_min(3)) EXIT loop_k
1322 IF (sb_min(3) >= subcell(i, j, k)%s_max(3)) cycle loop_k
1323 END IF
1324
1325 IF (periodic(2) /= 0) THEN
1326 IF (sb_max(2) < subcell(i, j, k)%s_min(2)) EXIT loop_j
1327 IF (sb_min(2) >= subcell(i, j, k)%s_max(2)) cycle loop_j
1328 END IF
1329
1330 IF (periodic(1) /= 0) THEN
1331 IF (sb_max(1) < subcell(i, j, k)%s_min(1)) EXIT loop_i
1332 IF (sb_min(1) >= subcell(i, j, k)%s_max(1)) cycle loop_i
1333 END IF
1334
1335 IF (subcell(i, j, k)%natom == 0) cycle
1336
1337 DO iatom_subcell = 1, subcell(i, j, k)%natom
1338 iatom_local = subcell(i, j, k)%atom_list(iatom_subcell)
1339 iatom = lista(ikind)%list(iatom_local)
1340 atom_a = atom(ikind)%list(iatom)
1341 IF (my_molecular) THEN
1342 mol_a = atom(ikind)%list_a_mol(iatom_local)
1343 IF (mol_a /= mol_b) cycle
1344 END IF
1345 IF (my_symmetric) THEN
1346 IF (atom_a > atom_b) THEN
1347 include_ab = (modulo(atom_a + atom_b, 2) /= 0)
1348 ELSE
1349 include_ab = (modulo(atom_a + atom_b, 2) == 0)
1350 END IF
1351 IF (my_sort_atomb) THEN
1352 IF ((.NOT. any(atomb_to_keep == atom_b)) .AND. &
1353 (.NOT. any(atomb_to_keep == atom_a))) THEN
1354 include_ab = .false.
1355 END IF
1356 END IF
1357 ELSE
1358 include_ab = .true.
1359 END IF
1360 IF (include_ab) THEN
1361 ra(:) = r_pbc(:, atom_a)
1362 rab(:) = rb(:) - ra(:)
1363 rab2 = rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3)
1364 IF (rab2 < rab2_max) THEN
1365 include_ab = .true.
1366 IF (my_mic) THEN
1367 ! only if rab is minimum image the pair will be included
1368 ! ideally the range of the pair list is < L/2 so
1369 ! that this never triggers
1370 rab_pbc(:) = pbc(rab(:), cell)
1371 IF (sum((rab_pbc - rab)**2) > epsilon(1.0_dp)) THEN
1372 include_ab = .false.
1373 END IF
1374 END IF
1375 IF (include_ab) THEN
1376 CALL add_neighbor_node( &
1377 neighbor_list=kind_a(iatom_local)%neighbor_list, &
1378 neighbor=atom_b, &
1379 cell=cell_b, &
1380 r=rab, &
1381 nkind=nkind)
1382 END IF
1383 END IF
1384 END IF
1385 END DO
1386
1387 END DO loop_i
1388 END DO loop_j
1389 END DO loop_k
1390
1391 END DO loop2_icell
1392 END DO loop2_jcell
1393 END DO loop2_kcell
1394
1395 END DO
1396
1397 CALL deallocate_subcell(subcell)
1398
1399 END DO
1400 END DO
1401
1402 SELECT CASE (otype)
1403 CASE (1:2, 4)
1404 CASE (3)
1405 DO ikind = 1, nkind
1406 DEALLOCATE (lista(ikind)%list)
1407 END DO
1408 CASE default
1409 cpabort("")
1410 END SELECT
1411 DEALLOCATE (kind_a, pres_a, pres_b, lista, listb, nlista, nlistb)
1412 DEALLOCATE (index_list)
1413 DEALLOCATE (r_pbc)
1414
1415 nentry = 0
1416 CALL neighbor_list_iterator_create(nl_iterator, ab_list)
1417 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1418 CALL get_iterator_info(nl_iterator, inode=inode, nnode=nnode)
1419 IF (inode == 1) nentry = nentry + nnode
1420 END DO
1421 CALL neighbor_list_iterator_release(nl_iterator)
1422 !
1423 ALLOCATE (ab_list(1)%nlist_task(nentry))
1424 ab_list(1)%nl_size = nentry
1425 DO iab = 2, SIZE(ab_list)
1426 ab_list(iab)%nl_size = nentry
1427 ab_list(iab)%nlist_task => ab_list(1)%nlist_task
1428 END DO
1429 !
1430 nentry = 0
1431 CALL neighbor_list_iterator_create(nl_iterator, ab_list)
1432 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1433 nentry = nentry + 1
1434 CALL get_iterator_task(nl_iterator, ab_list(1)%nlist_task(nentry))
1435 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, nkind=nkind)
1436 iab = (ikind - 1)*nkind + jkind
1437 IF (ab_list(iab)%nl_start < 0) ab_list(iab)%nl_start = nentry
1438 IF (ab_list(iab)%nl_end < 0) THEN
1439 ab_list(iab)%nl_end = nentry
1440 ELSE
1441 cpassert(ab_list(iab)%nl_end + 1 == nentry)
1442 ab_list(iab)%nl_end = nentry
1443 END IF
1444 END DO
1445 CALL neighbor_list_iterator_release(nl_iterator)
1446
1447 CALL timestop(handle)
1448
1449 END SUBROUTINE build_neighbor_lists
1450
1451! **************************************************************************************************
1452!> \brief Build a neighborlist
1453!> \param ab_list ...
1454!> \param basis_set_a ...
1455!> \param basis_set_b ...
1456!> \param qs_env ...
1457!> \param mic ...
1458!> \param symmetric ...
1459!> \param molecular ...
1460!> \param operator_type ...
1461!> \date 14.03.2016
1462!> \author JGH
1463! **************************************************************************************************
1464 SUBROUTINE setup_neighbor_list(ab_list, basis_set_a, basis_set_b, qs_env, &
1465 mic, symmetric, molecular, operator_type)
1466
1467 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1468 POINTER :: ab_list
1469 TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER :: basis_set_a
1470 TYPE(gto_basis_set_p_type), DIMENSION(:), &
1471 OPTIONAL, POINTER :: basis_set_b
1472 TYPE(qs_environment_type), POINTER :: qs_env
1473 LOGICAL, INTENT(IN), OPTIONAL :: mic, symmetric, molecular
1474 CHARACTER(LEN=*), INTENT(IN), OPTIONAL :: operator_type
1475
1476 CHARACTER(LEN=4) :: otype
1477 INTEGER :: ikind, nkind
1478 LOGICAL :: my_mic, my_molecular, my_symmetric
1479 LOGICAL, ALLOCATABLE, DIMENSION(:) :: a_present, b_present
1480 REAL(dp), ALLOCATABLE, DIMENSION(:) :: a_radius, b_radius
1481 REAL(dp), ALLOCATABLE, DIMENSION(:, :) :: pair_radius
1482 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
1483 TYPE(cell_type), POINTER :: cell
1484 TYPE(distribution_1d_type), POINTER :: distribution_1d
1485 TYPE(distribution_2d_type), POINTER :: distribution_2d
1486 TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER :: basis_a, basis_b
1487 TYPE(gto_basis_set_type), POINTER :: abas, bbas
1488 TYPE(local_atoms_type), ALLOCATABLE, DIMENSION(:) :: atom2d
1489 TYPE(molecule_type), DIMENSION(:), POINTER :: molecule_set
1490 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1491
1492 basis_a => basis_set_a
1493 IF (PRESENT(basis_set_b)) THEN
1494 basis_b => basis_set_b
1495 my_symmetric = .false.
1496 ELSE
1497 basis_b => basis_set_a
1498 my_symmetric = .true.
1499 END IF
1500 IF (PRESENT(symmetric)) my_symmetric = symmetric
1501
1502 IF (PRESENT(mic)) THEN
1503 my_mic = mic
1504 ELSE
1505 my_mic = .false.
1506 END IF
1507
1508 IF (PRESENT(molecular)) THEN
1509 my_molecular = molecular
1510 ELSE
1511 my_molecular = .false.
1512 END IF
1513
1514 IF (PRESENT(operator_type)) THEN
1515 otype = operator_type
1516 ELSE
1517 ! default is a simple AB neighbor list
1518 otype = "AB"
1519 END IF
1520
1521 nkind = SIZE(basis_a)
1522 ALLOCATE (a_present(nkind), b_present(nkind))
1523 a_present = .false.
1524 b_present = .false.
1525 ALLOCATE (a_radius(nkind), b_radius(nkind))
1526 a_radius = 0.0_dp
1527 b_radius = 0.0_dp
1528 DO ikind = 1, nkind
1529 IF (ASSOCIATED(basis_a(ikind)%gto_basis_set)) THEN
1530 a_present(ikind) = .true.
1531 abas => basis_a(ikind)%gto_basis_set
1532 CALL get_gto_basis_set(gto_basis_set=abas, kind_radius=a_radius(ikind))
1533 END IF
1534 IF (ASSOCIATED(basis_b(ikind)%gto_basis_set)) THEN
1535 b_present(ikind) = .true.
1536 bbas => basis_b(ikind)%gto_basis_set
1537 CALL get_gto_basis_set(gto_basis_set=bbas, kind_radius=b_radius(ikind))
1538 END IF
1539 END DO
1540
1541 ALLOCATE (pair_radius(nkind, nkind))
1542 pair_radius = 0.0_dp
1543 CALL pair_radius_setup(a_present, b_present, a_radius, b_radius, pair_radius)
1544
1545 CALL get_qs_env(qs_env, &
1546 atomic_kind_set=atomic_kind_set, &
1547 cell=cell, &
1548 distribution_2d=distribution_2d, &
1549 local_particles=distribution_1d, &
1550 particle_set=particle_set, &
1551 molecule_set=molecule_set)
1552
1553 ALLOCATE (atom2d(nkind))
1554 CALL atom2d_build(atom2d, distribution_1d, distribution_2d, atomic_kind_set, &
1555 molecule_set, my_molecular, particle_set=particle_set)
1556 CALL build_neighbor_lists(ab_list, particle_set, atom2d, cell, pair_radius, &
1557 mic=my_mic, symmetric=my_symmetric, molecular=my_molecular, &
1558 subcells=2.0_dp, nlname="AUX_NL")
1559
1560 CALL atom2d_cleanup(atom2d)
1561
1562 DEALLOCATE (a_present, b_present, a_radius, b_radius, pair_radius, atom2d)
1563
1564 END SUBROUTINE setup_neighbor_list
1565
1566! **************************************************************************************************
1567!> \brief ...
1568!> \param list ...
1569!> \param n ...
1570!> \param ikind ...
1571!> \param atom ...
1572! **************************************************************************************************
1573 SUBROUTINE combine_lists(list, n, ikind, atom)
1574 INTEGER, DIMENSION(:), POINTER :: list
1575 INTEGER, INTENT(OUT) :: n
1576 INTEGER, INTENT(IN) :: ikind
1577 TYPE(local_atoms_type), DIMENSION(:), INTENT(IN) :: atom
1578
1579 INTEGER :: i, ib, na, nb
1580 INTEGER, DIMENSION(:), POINTER :: lista, listb
1581
1582 cpassert(.NOT. ASSOCIATED(list))
1583
1584 lista => atom(ikind)%list_local_a_index
1585 listb => atom(ikind)%list_local_b_index
1586
1587 IF (ASSOCIATED(lista)) THEN
1588 na = SIZE(lista)
1589 ELSE
1590 na = 0
1591 END IF
1592
1593 IF (ASSOCIATED(listb)) THEN
1594 nb = SIZE(listb)
1595 ELSE
1596 nb = 0
1597 END IF
1598
1599 ALLOCATE (list(na + nb))
1600
1601 n = na
1602 IF (na .GT. 0) list(1:na) = lista(1:na)
1603 IF (nb .GT. 0) THEN
1604 loopb: DO ib = 1, nb
1605 DO i = 1, na
1606 IF (listb(ib) == list(i)) cycle loopb
1607 END DO
1608 n = n + 1
1609 list(n) = listb(ib)
1610 END DO loopb
1611 END IF
1612 END SUBROUTINE combine_lists
1613
1614! **************************************************************************************************
1615
1616! **************************************************************************************************
1617!> \brief ...
1618!> \param present_a ...
1619!> \param present_b ...
1620!> \param radius_a ...
1621!> \param radius_b ...
1622!> \param pair_radius ...
1623!> \param prmin ...
1624! **************************************************************************************************
1625 SUBROUTINE pair_radius_setup(present_a, present_b, radius_a, radius_b, pair_radius, prmin)
1626 LOGICAL, DIMENSION(:), INTENT(IN) :: present_a, present_b
1627 REAL(dp), DIMENSION(:), INTENT(IN) :: radius_a, radius_b
1628 REAL(dp), DIMENSION(:, :), INTENT(OUT) :: pair_radius
1629 REAL(dp), INTENT(IN), OPTIONAL :: prmin
1630
1631 INTEGER :: i, j, nkind
1632 REAL(dp) :: rrmin
1633
1634 nkind = SIZE(present_a)
1635
1636 pair_radius = 0._dp
1637
1638 rrmin = 0.0_dp
1639 IF (PRESENT(prmin)) rrmin = prmin
1640
1641 DO i = 1, nkind
1642 IF (.NOT. present_a(i)) cycle
1643 DO j = 1, nkind
1644 IF (.NOT. present_b(j)) cycle
1645 pair_radius(i, j) = radius_a(i) + radius_b(j)
1646 pair_radius(i, j) = max(pair_radius(i, j), rrmin)
1647 END DO
1648 END DO
1649
1650 END SUBROUTINE pair_radius_setup
1651
1652! **************************************************************************************************
1653!> \brief Print the distribution of the simple pair neighbor list.
1654!> \param ab ...
1655!> \param qs_kind_set ...
1656!> \param output_unit ...
1657!> \param para_env ...
1658!> \date 19.06.2003
1659!> \author MK
1660!> \version 1.0
1661! **************************************************************************************************
1662 SUBROUTINE write_neighbor_distribution(ab, qs_kind_set, output_unit, para_env)
1663 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1664 POINTER :: ab
1665 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
1666 INTEGER, INTENT(in) :: output_unit
1667 TYPE(mp_para_env_type), POINTER :: para_env
1668
1669 CHARACTER(len=*), PARAMETER :: routinen = 'write_neighbor_distribution'
1670 LOGICAL, PARAMETER :: full_output = .false.
1671
1672 INTEGER :: handle, ikind, inode, ipe, jkind, n, &
1673 nkind, nnode
1674 INTEGER(int_8) :: nblock_max, nblock_sum, nelement_max, &
1675 nelement_sum, tmp(2)
1676 INTEGER, ALLOCATABLE, DIMENSION(:) :: nblock, nelement, nnsgf
1677 TYPE(gto_basis_set_type), POINTER :: orb_basis_set
1679 DIMENSION(:), POINTER :: nl_iterator
1680
1681 CALL timeset(routinen, handle)
1682 associate(mype => para_env%mepos + 1, npe => para_env%num_pe)
1683
1684 ! Allocate work storage
1685 ALLOCATE (nblock(npe), nelement(npe))
1686 nblock(:) = 0
1687 nelement(:) = 0
1688 nkind = SIZE(qs_kind_set)
1689 ALLOCATE (nnsgf(nkind))
1690 nnsgf = 1
1691 DO ikind = 1, nkind
1692 CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set)
1693 IF (ASSOCIATED(orb_basis_set)) THEN
1694 CALL get_gto_basis_set(gto_basis_set=orb_basis_set, nsgf=nnsgf(ikind))
1695 END IF
1696 END DO
1697
1698 CALL neighbor_list_iterator_create(nl_iterator, ab)
1699 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1700 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, inode=inode, nnode=nnode)
1701 IF (inode == 1) THEN
1702 n = nnsgf(ikind)*nnsgf(jkind)
1703 nblock(mype) = nblock(mype) + nnode
1704 nelement(mype) = nelement(mype) + n*nnode
1705 END IF
1706 END DO
1707 CALL neighbor_list_iterator_release(nl_iterator)
1708
1709 IF (full_output) THEN
1710 ! XXXXXXXX should gather/scatter this on ionode
1711 CALL para_env%sum(nblock)
1712 CALL para_env%sum(nelement)
1713
1714 nblock_sum = sum(int(nblock, kind=int_8))
1715 nelement_sum = sum(int(nelement, kind=int_8))
1716 ELSE
1717 nblock_sum = nblock(mype)
1718 nblock_max = nblock(mype)
1719 nelement_sum = nelement(mype)
1720 nelement_max = nelement(mype)
1721 tmp = (/nblock_sum, nelement_sum/)
1722 CALL para_env%sum(tmp)
1723 nblock_sum = tmp(1); nelement_sum = tmp(2)
1724 tmp = (/nblock_max, nelement_max/)
1725 CALL para_env%max(tmp)
1726 nblock_max = tmp(1); nelement_max = tmp(2)
1727 END IF
1728
1729 IF (output_unit > 0) THEN
1730 IF (full_output) THEN
1731 WRITE (unit=output_unit, &
1732 fmt="(/,/,T2,A,/,/,T3,A,/,/,(T4,I6,T27,I10,T55,I10))") &
1733 "DISTRIBUTION OF THE NEIGHBOR LISTS", &
1734 "Process Number of particle pairs Number of matrix elements", &
1735 (ipe - 1, nblock(ipe), nelement(ipe), ipe=1, npe)
1736 WRITE (unit=output_unit, fmt="(/,T7,A3,T27,I10,T55,I10)") &
1737 "Sum", sum(nblock), sum(nelement)
1738 ELSE
1739 WRITE (unit=output_unit, fmt="(/,T2,A)") "DISTRIBUTION OF THE NEIGHBOR LISTS"
1740 WRITE (unit=output_unit, fmt="(T15,A,T68,I13)") "Total number of particle pairs:", nblock_sum
1741 WRITE (unit=output_unit, fmt="(T15,A,T68,I13)") "Total number of matrix elements:", nelement_sum
1742 WRITE (unit=output_unit, fmt="(T15,A,T68,I13)") "Average number of particle pairs:", (nblock_sum + npe - 1)/npe
1743 WRITE (unit=output_unit, fmt="(T15,A,T68,I13)") "Maximum number of particle pairs:", nblock_max
1744 WRITE (unit=output_unit, fmt="(T15,A,T68,I13)") "Average number of matrix element:", (nelement_sum + npe - 1)/npe
1745 WRITE (unit=output_unit, fmt="(T15,A,T68,I13)") "Maximum number of matrix elements:", nelement_max
1746 END IF
1747 END IF
1748 END associate
1749
1750 ! Release work storage
1751
1752 DEALLOCATE (nblock, nelement, nnsgf)
1753
1754 CALL timestop(handle)
1755
1756 END SUBROUTINE write_neighbor_distribution
1757
1758! **************************************************************************************************
1759!> \brief Write a set of neighbor lists to the output unit.
1760!> \param ab ...
1761!> \param particle_set ...
1762!> \param cell ...
1763!> \param para_env ...
1764!> \param neighbor_list_section ...
1765!> \param nl_type ...
1766!> \param middle_name ...
1767!> \param nlname ...
1768!> \date 04.03.2002
1769!> \par History
1770!> - Adapted to the new parallelized neighbor list version
1771!> (26.06.2003,MK)
1772!> \author MK
1773!> \version 1.0
1774! **************************************************************************************************
1775 SUBROUTINE write_neighbor_lists(ab, particle_set, cell, para_env, neighbor_list_section, &
1776 nl_type, middle_name, nlname)
1777
1778 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1779 POINTER :: ab
1780 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1781 TYPE(cell_type), POINTER :: cell
1782 TYPE(mp_para_env_type), POINTER :: para_env
1783 TYPE(section_vals_type), POINTER :: neighbor_list_section
1784 CHARACTER(LEN=*), INTENT(IN) :: nl_type, middle_name, nlname
1785
1786 CHARACTER(LEN=default_string_length) :: string, unit_str
1787 INTEGER :: iatom, inode, iw, jatom, nneighbor, nnode
1788 INTEGER, DIMENSION(3) :: cell_b
1789 REAL(dp) :: dab, unit_conv
1790 REAL(dp), DIMENSION(3) :: ra, rab, rb
1791 TYPE(cp_logger_type), POINTER :: logger
1792 TYPE(neighbor_list_iterator_p_type), &
1793 DIMENSION(:), POINTER :: nl_iterator
1794
1795 NULLIFY (logger)
1796 logger => cp_get_default_logger()
1797 IF (btest(cp_print_key_should_output(logger%iter_info, neighbor_list_section, &
1798 trim(nl_type)), &
1799 cp_p_file)) THEN
1800 iw = cp_print_key_unit_nr(logger=logger, &
1801 basis_section=neighbor_list_section, &
1802 print_key_path=trim(nl_type), &
1803 extension=".out", &
1804 middle_name=trim(middle_name), &
1805 local=.true., &
1806 log_filename=.false., &
1807 file_position="REWIND")
1808 associate(mype => para_env%mepos)
1809 CALL section_vals_val_get(neighbor_list_section, "UNIT", c_val=unit_str)
1810 unit_conv = cp_unit_from_cp2k(1.0_dp, trim(unit_str))
1811
1812 ! Print headline
1813 string = ""
1814 WRITE (unit=string, fmt="(A,I5,A)") &
1815 trim(nlname)//" IN "//trim(unit_str)//" (PROCESS", mype, ")"
1816 CALL compress(string)
1817 IF (iw > 0) WRITE (unit=iw, fmt="(/,/,T2,A)") trim(string)
1818
1819 nneighbor = 0
1820
1821 CALL neighbor_list_iterator_create(nl_iterator, ab)
1822 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1823 CALL get_iterator_info(nl_iterator, inode=inode, nnode=nnode, &
1824 iatom=iatom, jatom=jatom, cell=cell_b, r=rab)
1825 nneighbor = nneighbor + 1
1826 ra(:) = pbc(particle_set(iatom)%r, cell)
1827 rb(:) = ra(:) + rab(:)
1828 dab = sqrt(rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3))
1829 IF (iw > 0) THEN
1830 IF (inode == 1) THEN
1831 WRITE (unit=iw, fmt="(/,T2,I5,3X,I6,3X,3F12.6)") &
1832 iatom, nnode, ra(1:3)*unit_conv
1833 END IF
1834 WRITE (unit=iw, fmt="(T10,I6,3X,3I4,3F12.6,2X,F12.6)") &
1835 jatom, cell_b(1:3), rb(1:3)*unit_conv, dab*unit_conv
1836 END IF
1837 END DO
1838 CALL neighbor_list_iterator_release(nl_iterator)
1839
1840 string = ""
1841 WRITE (unit=string, fmt="(A,I12,A,I12)") &
1842 "Total number of neighbor interactions for process", mype, ":", &
1843 nneighbor
1844 CALL compress(string)
1845 IF (iw > 0) WRITE (unit=iw, fmt="(/,T2,A)") trim(string)
1846 CALL cp_print_key_finished_output(unit_nr=iw, &
1847 logger=logger, &
1848 basis_section=neighbor_list_section, &
1849 print_key_path=trim(nl_type), &
1850 local=.true.)
1851 END associate
1852 END IF
1853
1854 END SUBROUTINE write_neighbor_lists
1855
1856END MODULE qs_neighbor_lists
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
Types for all ALMO-based methods.
real(kind=dp), parameter, public almo_max_cutoff_multiplier
Definition atom.F:9
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind_set(atomic_kind_set, atom_of_kind, kind_of, natom_of_kind, maxatom, natom, nshell, fist_potential_present, shell_present, shell_adiabatic, shell_check_distance, damping_present)
Get attributes of an atomic kind set.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
subroutine, public get_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, lmin, lx, ly, lz, m, ncgf_set, npgf, nsgf_set, nshell, cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, last_cgf, last_sgf, n, gcc, maxco, maxl, maxpgf, maxsgf_set, maxshell, maxso, nco_sum, npgf_sum, nshell_sum, maxder, short_kind_radius, npgf_seg_sum)
...
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public scaled_to_real(r, s, cell)
Transform scaled cell coordinates real coordinates. r=h*s.
Definition cell_types.F:516
subroutine, public real_to_scaled(s, r, cell)
Transform real to scaled cell coordinates. s=h_inv*r.
Definition cell_types.F:486
subroutine, public get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, h, h_inv, symmetry_id, tag)
Get informations about a simulation cell.
Definition cell_types.F:195
real(kind=dp) function, public plane_distance(h, k, l, cell)
Calculate the distance between two lattice planes as defined by a triple of Miller indices (hkl).
Definition cell_types.F:252
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
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...
unit conversion facility
Definition cp_units.F:30
real(kind=dp) function, public cp_unit_from_cp2k(value, unit_str, defaults, power)
converts from the internal cp2k units to the given unit
Definition cp_units.F:1179
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
stores a mapping of 2D info (e.g. matrix) on a 2D processor distribution (i.e. blacs grid) where cpus...
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.
Definition of the atomic potential types.
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_method_rigpw
integer, parameter, public dispersion_uff
integer, parameter, public vdw_pairpot_dftd4
integer, parameter, public do_potential_truncated
integer, parameter, public do_potential_id
integer, parameter, public do_method_lrigpw
integer, parameter, public do_potential_short
integer, parameter, public do_se_is_slater
integer, parameter, public xc_vdw_fun_pairpot
objects that represent the structure of input sections and the data contained in an input section
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
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public int_8
Definition kinds.F:54
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
Types and basic routines needed for a kpoint calculation.
2- and 3-center electron repulsion integral routines based on libint2 Currently available operators: ...
real(kind=dp), parameter, public cutoff_screen_factor
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
Definition list.F:24
Collection of simple mathematical functions and subroutines.
Definition mathlib.F:15
subroutine, public erfc_cutoff(eps, omg, r_cutoff)
compute a truncation radius for the shortrange operator
Definition mathlib.F:1686
Interface to the message passing library MPI.
Define the data structure for the molecule information.
Define the data structure for the particle information.
subroutine, public get_paw_proj_set(paw_proj_set, csprj, chprj, first_prj, first_prjs, last_prj, local_oce_sphi_h, local_oce_sphi_s, maxl, ncgauprj, nsgauprj, nsatbas, nsotot, nprj, o2nindex, n2oindex, rcprj, rzetprj, zisomin, zetprj)
Get informations about a paw projectors set.
Periodic Table related data definitions.
type(atom), dimension(0:nelem), public ptable
Definition of physical constants:
Definition physcon.F:68
real(kind=dp), parameter, public bohr
Definition physcon.F:147
Definition of the DFTB parameter types.
Working with the DFTB parameter types.
subroutine, public get_dftb_atom_param(dftb_parameter, name, typ, defined, z, zeff, natorb, lmax, skself, occupation, eta, energy, cutoff, xi, di, rcdisp, dudq)
...
Definition of disperson types for DFT calculations.
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)
Get the QUICKSTEP environment.
Definition of gCP types for DFT calculations.
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, 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, 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)
Get attributes of an atomic kind set.
subroutine, public set_ks_env(ks_env, v_hartree_rspace, s_mstruct_changed, rho_changed, 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, 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, task_list, task_list_soft, subsys, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env)
...
subroutine, public get_ks_env(ks_env, v_hartree_rspace, s_mstruct_changed, rho_changed, 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, rho, rho_xc, vppl, rho_core, rho_nlcc, rho_nlcc_g, vee, neighbor_list_id, 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, task_list, task_list_soft, kpoints, do_kpoints, atomic_kind_set, qs_kind_set, cell, cell_ref, use_ref_cell, particle_set, energy, force, local_particles, local_molecules, molecule_kind_set, molecule_set, subsys, cp_subsys, virial, results, atprop, nkind, natom, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env, nelectron_total, nelectron_spin)
...
Define the neighbor list data types and the corresponding functionality.
subroutine, public get_iterator_task(iterator_set, task, mepos)
Captures the current state of the iterator in a neighbor_list_task_type.
subroutine, public release_neighbor_list_sets(nlists)
releases an array of neighbor_list_sets
subroutine, public add_neighbor_node(neighbor_list, neighbor, cell, r, exclusion_list, nkind)
Add a new neighbor list node to a neighbor list.
subroutine, public add_neighbor_list(neighbor_list_set, atom, neighbor_list)
Add a new neighbor list to a neighbor list set.
subroutine, public neighbor_list_iterator_create(iterator_set, nl, search, nthread)
Neighbor list iterator functions.
subroutine, public neighbor_list_iterator_release(iterator_set)
...
integer function, public neighbor_list_iterate(iterator_set, mepos)
...
subroutine, public allocate_neighbor_list_set(neighbor_list_set, symmetric)
Allocate and initialize a set of neighbor lists.
subroutine, public get_iterator_info(iterator_set, mepos, ikind, jkind, nkind, ilist, nlist, inode, nnode, iatom, jatom, r, cell)
...
Generate the atomic neighbor lists.
subroutine, public atom2d_cleanup(atom2d)
free the internals of atom2d
subroutine, public build_qs_neighbor_lists(qs_env, para_env, molecular, force_env_section)
Build all the required neighbor lists for Quickstep.
subroutine, public pair_radius_setup(present_a, present_b, radius_a, radius_b, pair_radius, prmin)
...
subroutine, public build_neighbor_lists(ab_list, particle_set, atom, cell, pair_radius, subcells, mic, symmetric, molecular, subset_of_mol, current_subset, operator_type, nlname, atomb_to_keep)
Build simple pair neighbor lists.
subroutine, public setup_neighbor_list(ab_list, basis_set_a, basis_set_b, qs_env, mic, symmetric, molecular, operator_type)
Build a neighborlist.
subroutine, public write_neighbor_lists(ab, particle_set, cell, para_env, neighbor_list_section, nl_type, middle_name, nlname)
Write a set of neighbor lists to the output unit.
subroutine, public atom2d_build(atom2d, distribution_1d, distribution_2d, atomic_kind_set, molecule_set, molecule_only, particle_set)
Build some distribution structure of atoms, refactored from build_qs_neighbor_lists.
Utilities for string manipulations.
subroutine, public compress(string, full)
Eliminate multiple space characters in a string. If full is .TRUE., then all spaces are eliminated.
elemental subroutine, public uppercase(string)
Convert all lower case characters in a string to upper case.
subcell types and allocation routines
subroutine, public deallocate_subcell(subcell)
Deallocate a subcell grid structure.
subroutine, public give_ijk_subcell(r, i, j, k, cell, nsubcell)
...
subroutine, public allocate_subcell(subcell, nsubcell, maxatom, cell)
Allocate and initialize a subcell grid structure for the atomic neighbor search.
All kind of helpful little routines.
Definition util.F:14
pure integer function, public locate(array, x)
Purpose: Given an array array(1:n), and given a value x, a value x_index is returned which is the ind...
Definition util.F:61
Definition of the xTB parameter types.
Definition xtb_types.F:20
subroutine, public get_xtb_atom_param(xtb_parameter, symbol, 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, occupation, electronegativity, chmax, en, kqat2, kcn, kq)
...
Definition xtb_types.F:199
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
type of a logger, at the moment it contains just a print level starting at which level it should be l...
structure to store local (to a processor) ordered lists of integers.
distributes pairs on a 2d grid of processors
Contains information about kpoints.
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...