(git:e5b1968)
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 .AND. (.NOT. dft_control%qs_control%xtb_control%do_tblite)) 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 &
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.
833 DO ingp = 1, ngp
834 DO ikind = 1, nkind
835 rcut = sqrt(dft_control%qs_control%xtb_control%nonbonded%pot(ingp)%pot%rcutsq)
836 c_radius = rcut
837 CALL get_atomic_kind(atomic_kind_set(ikind), element_symbol=element_symbol)
838 CALL uppercase(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.
841 DO jkind = 1, nkind
842 CALL get_atomic_kind(atomic_kind_set(jkind), element_symbol=element_symbol2)
843 CALL uppercase(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.
846 END IF
847 END DO
848 END IF
849 END DO
850 CALL pair_radius_setup(nonbond1_atom, nonbond2_atom, c_radius, c_radius, pair_radius)
851 CALL build_neighbor_lists(sab_xtb_nonbond, particle_set, atom2d, cell, 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")
856 END DO
857 END IF
858 END IF
859
860 ! Build the neighbor lists for the vdW pair potential
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
865 IF (dispersion_env%type == xc_vdw_fun_pairpot .OR. xtb) THEN
866 IF (dispersion_env%pp_type == vdw_pairpot_dftd4) THEN
867 c_radius(:) = dispersion_env%rc_d4
868 ELSE
869 c_radius(:) = dispersion_env%rc_disp
870 END IF
871 default_present = .true. !include all atoms in vdW (even without basis)
872 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
873 CALL build_neighbor_lists(sab_vdw, particle_set, atom2d, cell, pair_radius, &
874 subcells=subcells, operator_type="PP", nlname="sab_vdw")
875 dispersion_env%sab_vdw => sab_vdw
876
877 ! Build the neighbor lists for coordination numbers as needed by the DFT-D3/D4 method
878 ! This is also needed for the xTB Hamiltonian
879 DO ikind = 1, nkind
880 CALL get_atomic_kind(atomic_kind_set(ikind), z=zat)
881 c_radius(ikind) = 4._dp*ptable(zat)%covalent_radius*bohr
882 END DO
883 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
884 CALL build_neighbor_lists(sab_cn, particle_set, atom2d, cell, pair_radius, &
885 subcells=subcells, operator_type="PP", nlname="sab_cn")
886 dispersion_env%sab_cn => sab_cn
887 END IF
888 END IF
889
890 ! Build the neighbor lists for the gCP pair potential
891 NULLIFY (gcp_env)
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
896 DO ikind = 1, nkind
897 c_radius(ikind) = gcp_env%gcp_kind(ikind)%rcsto
898 END DO
899 CALL pair_radius_setup(orb_present, orb_present, c_radius, c_radius, pair_radius)
900 CALL build_neighbor_lists(sab_gcp, particle_set, atom2d, cell, pair_radius, &
901 subcells=subcells, operator_type="PP", nlname="sab_gcp")
902 gcp_env%sab_gcp => sab_gcp
903 ELSE
904 NULLIFY (gcp_env%sab_gcp)
905 END IF
906 END IF
907
908 IF (lrigpw .OR. lri_optbas) THEN
909 ! set neighborlists in lri_env environment
910 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius)
911 soo_list => qs_env%lri_env%soo_list
912 CALL build_neighbor_lists(soo_list, particle_set, atom2d, cell, pair_radius, &
913 mic=mic, molecular=molecule_only, subcells=subcells, nlname="soo_list")
914 qs_env%lri_env%soo_list => soo_list
915 CALL write_neighbor_lists(soo_list, particle_set, cell, para_env, neighbor_list_section, &
916 "/SOO_LIST", "soo_list", "ORBITAL ORBITAL (RI)")
917 ELSEIF (rigpw) THEN
918 ALLOCATE (ri_present(nkind), ri_radius(nkind))
919 ri_present = .false.
920 ri_radius = 0.0_dp
921 DO ikind = 1, 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.
925 CALL get_gto_basis_set(gto_basis_set=ri_basis_set, kind_radius=ri_radius(ikind))
926 ELSE
927 ri_present(ikind) = .false.
928 END IF
929 END DO
930 ! set neighborlists in lri_env environment
931 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius)
932 soo_list => qs_env%lri_env%soo_list
933 CALL build_neighbor_lists(soo_list, particle_set, atom2d, cell, pair_radius, &
934 mic=mic, molecular=molecule_only, subcells=subcells, nlname="soo_list")
935 qs_env%lri_env%soo_list => soo_list
936 !
937 CALL pair_radius_setup(ri_present, ri_present, ri_radius, ri_radius, pair_radius)
938 saa_list => qs_env%lri_env%saa_list
939 CALL build_neighbor_lists(saa_list, particle_set, atom2d, cell, pair_radius, &
940 mic=mic, molecular=molecule_only, subcells=subcells, nlname="saa_list")
941 qs_env%lri_env%saa_list => saa_list
942 !
943 CALL pair_radius_setup(ri_present, orb_present, ri_radius, orb_radius, pair_radius)
944 soa_list => qs_env%lri_env%soa_list
945 CALL build_neighbor_lists(soa_list, particle_set, atom2d, cell, pair_radius, &
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
949 END IF
950
951 ! Build the neighbor lists for the ALMO delocalization
952 IF (almo) THEN
953 DO ikind = 1, nkind
954 CALL get_atomic_kind(atomic_kind_set(ikind), rcov=almo_rcov, rvdw=almo_rvdw)
955 ! multiply the radius by some hard-coded number
956 c_radius(ikind) = max(almo_rcov, almo_rvdw)*bohr* &
958 END DO
959 default_present = .true. !include all atoms (even without basis)
960 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
961 CALL build_neighbor_lists(sab_almo, particle_set, atom2d, cell, pair_radius, &
962 subcells=subcells, operator_type="PP", nlname="sab_almo")
963 CALL set_ks_env(ks_env=ks_env, sab_almo=sab_almo)
964 END IF
965
966 ! Print particle distribution
967 print_key_path = "PRINT%DISTRIBUTION"
968 IF (btest(cp_print_key_should_output(logger%iter_info, force_env_section, &
969 print_key_path), &
970 cp_p_file)) THEN
971 iw = cp_print_key_unit_nr(logger=logger, &
972 basis_section=force_env_section, &
973 print_key_path=print_key_path, &
974 extension=".out")
975 CALL write_neighbor_distribution(sab_orb, qs_kind_set, iw, para_env)
976 CALL cp_print_key_finished_output(unit_nr=iw, &
977 logger=logger, &
978 basis_section=force_env_section, &
979 print_key_path=print_key_path)
980 END IF
981
982 ! Release work storage
983 CALL atom2d_cleanup(atom2d)
984
985 DEALLOCATE (atom2d)
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)
993 END IF
994 IF (paw_atom_present) THEN
995 DEALLOCATE (oce_present, oce_radius)
996 END IF
997 IF (all_potential_present .OR. sgp_potential_present) THEN
998 DEALLOCATE (all_present, all_pot_rad)
999 END IF
1000
1001 CALL timestop(handle)
1002
1003 END SUBROUTINE build_qs_neighbor_lists
1004
1005! **************************************************************************************************
1006!> \brief Build simple pair neighbor lists.
1007!> \param ab_list ...
1008!> \param particle_set ...
1009!> \param atom ...
1010!> \param cell ...
1011!> \param pair_radius ...
1012!> \param subcells ...
1013!> \param mic ...
1014!> \param symmetric ...
1015!> \param molecular ...
1016!> \param subset_of_mol ...
1017!> \param current_subset ...
1018!> \param operator_type ...
1019!> \param nlname ...
1020!> \param atomb_to_keep the list of atom indices to keep for pairs from the atom2d%b_list
1021!> \date 20.03.2002
1022!> \par History
1023!> - Major refactoring (25.07.2010,jhu)
1024!> - Added option to filter out atoms from list_b (08.2018, A. Bussy)
1025!> \author MK
1026!> \version 2.0
1027! **************************************************************************************************
1028 SUBROUTINE build_neighbor_lists(ab_list, particle_set, atom, cell, pair_radius, subcells, &
1029 mic, symmetric, molecular, subset_of_mol, current_subset, &
1030 operator_type, nlname, atomb_to_keep)
1031
1032 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1033 POINTER :: ab_list
1034 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1035 TYPE(local_atoms_type), DIMENSION(:), INTENT(IN) :: atom
1036 TYPE(cell_type), POINTER :: cell
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
1045
1046 CHARACTER(len=*), PARAMETER :: routinen = 'build_neighbor_lists'
1047
1048 TYPE local_lists
1049 INTEGER, DIMENSION(:), POINTER :: list
1050 END TYPE local_lists
1051
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
1065 TYPE(neighbor_list_p_type), &
1066 ALLOCATABLE, DIMENSION(:) :: kind_a
1067 TYPE(neighbor_list_set_type), POINTER :: neighbor_list_set
1068 TYPE(subcell_type), DIMENSION(:, :, :), &
1069 POINTER :: subcell
1070 REAL(kind=dp), DIMENSION(:, :), ALLOCATABLE :: r_pbc
1072 DIMENSION(:), POINTER :: nl_iterator
1073
1074 CALL timeset(routinen//"_"//trim(nlname), handle)
1075
1076 ! input options
1077 my_mic = .false.
1078 IF (PRESENT(mic)) my_mic = mic
1079 my_symmetric = .true.
1080 IF (PRESENT(symmetric)) my_symmetric = symmetric
1081 my_molecular = .false.
1082 ! if we have a molecular NL, MIC has to be used
1083 IF (PRESENT(molecular)) my_molecular = molecular
1084 ! check for operator types
1085 IF (PRESENT(operator_type)) THEN
1086 SELECT CASE (operator_type)
1087 CASE ("AB")
1088 otype = 1 ! simple overlap
1089 CASE ("ABC")
1090 otype = 2 ! for three center operators
1091 cpassert(.NOT. my_molecular)
1092 my_symmetric = .false.
1093 CASE ("ABBA")
1094 otype = 3 ! for separable nonlocal operators
1095 my_symmetric = .false.
1096 CASE ("PP")
1097 otype = 4 ! simple atomic pair potential list
1098 CASE default
1099 cpabort("")
1100 END SELECT
1101 ELSE
1102 ! default is a simple AB neighbor list
1103 otype = 1
1104 END IF
1105 my_sort_atomb = .false.
1106 IF (PRESENT(atomb_to_keep)) THEN
1107 my_sort_atomb = .true.
1108 END IF
1109
1110 nkind = SIZE(atom)
1111 ! Deallocate the old neighbor list structure
1112 CALL release_neighbor_list_sets(ab_list)
1113 ! Allocate and initialize the new neighbor list structure
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)
1121 END DO
1122
1123 ! Allocate and initialize the kind availability
1124 ALLOCATE (pres_a(nkind), pres_b(nkind))
1125 DO ikind = 1, nkind
1126 pres_a(ikind) = any(pair_radius(ikind, :) > 0._dp)
1127 pres_b(ikind) = any(pair_radius(:, ikind) > 0._dp)
1128 END DO
1129
1130 ! create a copy of the pbc'ed coordinates
1131 natom = SIZE(particle_set)
1132 ALLOCATE (r_pbc(3, natom))
1133 DO i = 1, natom
1134 r_pbc(1:3, i) = pbc(particle_set(i)%r(1:3), cell)
1135 END DO
1136
1137 ! setup the local lists of atoms
1138 maxat = 0
1139 DO ikind = 1, nkind
1140 maxat = max(maxat, SIZE(atom(ikind)%list))
1141 END DO
1142 ALLOCATE (index_list(maxat))
1143 DO i = 1, maxat
1144 index_list(i) = i
1145 END DO
1146 ALLOCATE (lista(nkind), listb(nkind), nlista(nkind), nlistb(nkind))
1147 nlista = 0
1148 nlistb = 0
1149 DO ikind = 1, nkind
1150 NULLIFY (lista(ikind)%list, listb(ikind)%list)
1151 SELECT CASE (otype)
1152 CASE (1)
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)
1156 END IF
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)
1160 END IF
1161 CASE (2)
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)
1165 END IF
1166 nlistb(ikind) = SIZE(atom(ikind)%list)
1167 listb(ikind)%list => index_list
1168 CASE (3)
1169 CALL combine_lists(lista(ikind)%list, nlista(ikind), ikind, atom)
1170 nlistb(ikind) = SIZE(atom(ikind)%list)
1171 listb(ikind)%list => index_list
1172 CASE (4)
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
1177 CASE default
1178 cpabort("")
1179 END SELECT
1180 END DO
1181
1182 ! Determine max. number of local atoms
1183 maxat = 0
1184 DO ikind = 1, nkind
1185 maxat = max(maxat, nlista(ikind), nlistb(ikind))
1186 END DO
1187 ALLOCATE (kind_a(2*maxat))
1188
1189 ! Load informations about the simulation cell
1190 CALL get_cell(cell=cell, periodic=periodic, deth=deth)
1191
1192 ! Loop over all atomic kind pairs
1193 DO ikind = 1, nkind
1194 IF (.NOT. pres_a(ikind)) cycle
1195
1196 DO jkind = 1, nkind
1197 IF (.NOT. pres_b(jkind)) cycle
1198
1199 iab = ikind + nkind*(jkind - 1)
1200
1201 ! Calculate the square of the maximum interaction distance
1202 IF (pair_radius(ikind, jkind) <= 0._dp) cycle
1203 rab_max = pair_radius(ikind, jkind)
1204 IF (otype == 3) THEN
1205 ! Calculate the square of the maximum interaction distance
1206 ! for sac_max / ncell this must be the maximum over all kinds
1207 ! to be correct for three center terms involving different kinds
1208 rabm = maxval(pair_radius(:, jkind))
1209 ELSE
1210 rabm = rab_max
1211 END IF
1212 rab2_max = rabm*rabm
1213
1214 pd(1) = plane_distance(1, 0, 0, cell)
1215 pd(2) = plane_distance(0, 1, 0, cell)
1216 pd(3) = plane_distance(0, 0, 1, cell)
1217
1218 sab_max = rabm/pd
1219 sab_max_guard = 15.0_dp/pd
1220
1221 ! It makes sense to have fewer subcells for larger systems
1222 subcell_scale = ((125.0_dp**3)/deth)**(1.0_dp/6.0_dp)
1223
1224 ! guess the number of subcells for optimal performance,
1225 ! guard against crazy stuff triggered by very small rabm
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(:))))
1228
1229 ! number of image cells to be considered
1230 ncell(:) = (int(sab_max(:)) + 1)*periodic(:)
1231
1232 CALL allocate_neighbor_list_set(neighbor_list_set=ab_list(iab)%neighbor_list_set, &
1233 symmetric=my_symmetric)
1234 neighbor_list_set => ab_list(iab)%neighbor_list_set
1235
1236 DO iatom_local = 1, nlista(ikind)
1237 iatom = lista(ikind)%list(iatom_local)
1238 atom_a = atom(ikind)%list(iatom)
1239 CALL add_neighbor_list(neighbor_list_set=neighbor_list_set, &
1240 atom=atom_a, &
1241 neighbor_list=kind_a(iatom_local)%neighbor_list)
1242 END DO
1243
1244 CALL allocate_subcell(subcell, nsubcell)
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)
1249 CALL give_ijk_subcell(r, i, j, k, cell, nsubcell)
1250 subcell(i, j, k)%natom = subcell(i, j, k)%natom + 1
1251 END DO
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
1258 END DO
1259 END DO
1260 END DO
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)
1265 CALL give_ijk_subcell(r, i, j, k, cell, nsubcell)
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
1268 END DO
1269
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
1275 END IF
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
1280 END IF
1281 END IF
1282 r = r_pbc(:, atom_b)
1283 CALL real_to_scaled(sb_pbc(:), r(:), cell)
1284
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
1292 END IF
1293 cell_b(3) = kcell
1294
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
1302 END IF
1303 cell_b(2) = jcell
1304
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
1312 END IF
1313 cell_b(1) = icell
1314
1315 CALL scaled_to_real(rb, sb, cell)
1316
1317 loop_k: DO k = 1, nsubcell(3)
1318 loop_j: DO j = 1, nsubcell(2)
1319 loop_i: DO i = 1, nsubcell(1)
1320
1321 ! FIXME for non-periodic systems, the whole subcell trick is skipped
1322 ! yielding a Natom**2 pair list build.
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
1326 END IF
1327
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
1331 END IF
1332
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
1336 END IF
1337
1338 IF (subcell(i, j, k)%natom == 0) cycle
1339
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
1347 END IF
1348 IF (my_symmetric) THEN
1349 IF (atom_a > atom_b) THEN
1350 include_ab = (modulo(atom_a + atom_b, 2) /= 0)
1351 ELSE
1352 include_ab = (modulo(atom_a + atom_b, 2) == 0)
1353 END IF
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.
1358 END IF
1359 END IF
1360 ELSE
1361 include_ab = .true.
1362 END IF
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
1368 include_ab = .true.
1369 IF (my_mic) THEN
1370 ! only if rab is minimum image the pair will be included
1371 ! ideally the range of the pair list is < L/2 so
1372 ! that this never triggers
1373 rab_pbc(:) = pbc(rab(:), cell)
1374 IF (sum((rab_pbc - rab)**2) > epsilon(1.0_dp)) THEN
1375 include_ab = .false.
1376 END IF
1377 END IF
1378 IF (include_ab) THEN
1379 CALL add_neighbor_node( &
1380 neighbor_list=kind_a(iatom_local)%neighbor_list, &
1381 neighbor=atom_b, &
1382 cell=cell_b, &
1383 r=rab, &
1384 nkind=nkind)
1385 END IF
1386 END IF
1387 END IF
1388 END DO
1389
1390 END DO loop_i
1391 END DO loop_j
1392 END DO loop_k
1393
1394 END DO loop2_icell
1395 END DO loop2_jcell
1396 END DO loop2_kcell
1397
1398 END DO
1399
1400 CALL deallocate_subcell(subcell)
1401
1402 END DO
1403 END DO
1404
1405 SELECT CASE (otype)
1406 CASE (1:2, 4)
1407 CASE (3)
1408 DO ikind = 1, nkind
1409 DEALLOCATE (lista(ikind)%list)
1410 END DO
1411 CASE default
1412 cpabort("")
1413 END SELECT
1414 DEALLOCATE (kind_a, pres_a, pres_b, lista, listb, nlista, nlistb)
1415 DEALLOCATE (index_list)
1416 DEALLOCATE (r_pbc)
1417
1418 nentry = 0
1419 CALL neighbor_list_iterator_create(nl_iterator, ab_list)
1420 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1421 CALL get_iterator_info(nl_iterator, inode=inode, nnode=nnode)
1422 IF (inode == 1) nentry = nentry + nnode
1423 END DO
1424 CALL neighbor_list_iterator_release(nl_iterator)
1425 !
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
1431 END DO
1432 !
1433 nentry = 0
1434 CALL neighbor_list_iterator_create(nl_iterator, ab_list)
1435 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1436 nentry = nentry + 1
1437 CALL get_iterator_task(nl_iterator, ab_list(1)%nlist_task(nentry))
1438 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, nkind=nkind)
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
1443 ELSE
1444 cpassert(ab_list(iab)%nl_end + 1 == nentry)
1445 ab_list(iab)%nl_end = nentry
1446 END IF
1447 END DO
1448 CALL neighbor_list_iterator_release(nl_iterator)
1449
1450 CALL timestop(handle)
1451
1452 END SUBROUTINE build_neighbor_lists
1453
1454! **************************************************************************************************
1455!> \brief Build a neighborlist
1456!> \param ab_list ...
1457!> \param basis_set_a ...
1458!> \param basis_set_b ...
1459!> \param qs_env ...
1460!> \param mic ...
1461!> \param symmetric ...
1462!> \param molecular ...
1463!> \param operator_type ...
1464!> \date 14.03.2016
1465!> \author JGH
1466! **************************************************************************************************
1467 SUBROUTINE setup_neighbor_list(ab_list, basis_set_a, basis_set_b, qs_env, &
1468 mic, symmetric, molecular, operator_type)
1469
1470 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1471 POINTER :: ab_list
1472 TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER :: basis_set_a
1473 TYPE(gto_basis_set_p_type), DIMENSION(:), &
1474 OPTIONAL, POINTER :: basis_set_b
1475 TYPE(qs_environment_type), POINTER :: qs_env
1476 LOGICAL, INTENT(IN), OPTIONAL :: mic, symmetric, molecular
1477 CHARACTER(LEN=*), INTENT(IN), OPTIONAL :: operator_type
1478
1479 CHARACTER(LEN=4) :: otype
1480 INTEGER :: ikind, nkind
1481 LOGICAL :: my_mic, my_molecular, my_symmetric
1482 LOGICAL, ALLOCATABLE, DIMENSION(:) :: a_present, b_present
1483 REAL(dp), ALLOCATABLE, DIMENSION(:) :: a_radius, b_radius
1484 REAL(dp), ALLOCATABLE, DIMENSION(:, :) :: pair_radius
1485 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
1486 TYPE(cell_type), POINTER :: cell
1487 TYPE(distribution_1d_type), POINTER :: distribution_1d
1488 TYPE(distribution_2d_type), POINTER :: distribution_2d
1489 TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER :: basis_a, basis_b
1490 TYPE(gto_basis_set_type), POINTER :: abas, bbas
1491 TYPE(local_atoms_type), ALLOCATABLE, DIMENSION(:) :: atom2d
1492 TYPE(molecule_type), DIMENSION(:), POINTER :: molecule_set
1493 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1494
1495 basis_a => basis_set_a
1496 IF (PRESENT(basis_set_b)) THEN
1497 basis_b => basis_set_b
1498 my_symmetric = .false.
1499 ELSE
1500 basis_b => basis_set_a
1501 my_symmetric = .true.
1502 END IF
1503 IF (PRESENT(symmetric)) my_symmetric = symmetric
1504
1505 IF (PRESENT(mic)) THEN
1506 my_mic = mic
1507 ELSE
1508 my_mic = .false.
1509 END IF
1510
1511 IF (PRESENT(molecular)) THEN
1512 my_molecular = molecular
1513 ELSE
1514 my_molecular = .false.
1515 END IF
1516
1517 IF (PRESENT(operator_type)) THEN
1518 otype = operator_type
1519 ELSE
1520 ! default is a simple AB neighbor list
1521 otype = "AB"
1522 END IF
1523
1524 nkind = SIZE(basis_a)
1525 ALLOCATE (a_present(nkind), b_present(nkind))
1526 a_present = .false.
1527 b_present = .false.
1528 ALLOCATE (a_radius(nkind), b_radius(nkind))
1529 a_radius = 0.0_dp
1530 b_radius = 0.0_dp
1531 DO ikind = 1, nkind
1532 IF (ASSOCIATED(basis_a(ikind)%gto_basis_set)) THEN
1533 a_present(ikind) = .true.
1534 abas => basis_a(ikind)%gto_basis_set
1535 CALL get_gto_basis_set(gto_basis_set=abas, kind_radius=a_radius(ikind))
1536 END IF
1537 IF (ASSOCIATED(basis_b(ikind)%gto_basis_set)) THEN
1538 b_present(ikind) = .true.
1539 bbas => basis_b(ikind)%gto_basis_set
1540 CALL get_gto_basis_set(gto_basis_set=bbas, kind_radius=b_radius(ikind))
1541 END IF
1542 END DO
1543
1544 ALLOCATE (pair_radius(nkind, nkind))
1545 pair_radius = 0.0_dp
1546 CALL pair_radius_setup(a_present, b_present, a_radius, b_radius, pair_radius)
1547
1548 CALL get_qs_env(qs_env, &
1549 atomic_kind_set=atomic_kind_set, &
1550 cell=cell, &
1551 distribution_2d=distribution_2d, &
1552 local_particles=distribution_1d, &
1553 particle_set=particle_set, &
1554 molecule_set=molecule_set)
1555
1556 ALLOCATE (atom2d(nkind))
1557 CALL atom2d_build(atom2d, distribution_1d, distribution_2d, atomic_kind_set, &
1558 molecule_set, my_molecular, particle_set=particle_set)
1559 CALL build_neighbor_lists(ab_list, particle_set, atom2d, cell, pair_radius, &
1560 mic=my_mic, symmetric=my_symmetric, molecular=my_molecular, &
1561 subcells=2.0_dp, nlname="AUX_NL")
1562
1563 CALL atom2d_cleanup(atom2d)
1564
1565 DEALLOCATE (a_present, b_present, a_radius, b_radius, pair_radius, atom2d)
1566
1567 END SUBROUTINE setup_neighbor_list
1568
1569! **************************************************************************************************
1570!> \brief ...
1571!> \param list ...
1572!> \param n ...
1573!> \param ikind ...
1574!> \param atom ...
1575! **************************************************************************************************
1576 SUBROUTINE combine_lists(list, n, ikind, atom)
1577 INTEGER, DIMENSION(:), POINTER :: list
1578 INTEGER, INTENT(OUT) :: n
1579 INTEGER, INTENT(IN) :: ikind
1580 TYPE(local_atoms_type), DIMENSION(:), INTENT(IN) :: atom
1581
1582 INTEGER :: i, ib, na, nb
1583 INTEGER, DIMENSION(:), POINTER :: lista, listb
1584
1585 cpassert(.NOT. ASSOCIATED(list))
1586
1587 lista => atom(ikind)%list_local_a_index
1588 listb => atom(ikind)%list_local_b_index
1589
1590 IF (ASSOCIATED(lista)) THEN
1591 na = SIZE(lista)
1592 ELSE
1593 na = 0
1594 END IF
1595
1596 IF (ASSOCIATED(listb)) THEN
1597 nb = SIZE(listb)
1598 ELSE
1599 nb = 0
1600 END IF
1601
1602 ALLOCATE (list(na + nb))
1603
1604 n = na
1605 IF (na .GT. 0) list(1:na) = lista(1:na)
1606 IF (nb .GT. 0) THEN
1607 loopb: DO ib = 1, nb
1608 DO i = 1, na
1609 IF (listb(ib) == list(i)) cycle loopb
1610 END DO
1611 n = n + 1
1612 list(n) = listb(ib)
1613 END DO loopb
1614 END IF
1615 END SUBROUTINE combine_lists
1616
1617! **************************************************************************************************
1618
1619! **************************************************************************************************
1620!> \brief ...
1621!> \param present_a ...
1622!> \param present_b ...
1623!> \param radius_a ...
1624!> \param radius_b ...
1625!> \param pair_radius ...
1626!> \param prmin ...
1627! **************************************************************************************************
1628 SUBROUTINE pair_radius_setup(present_a, present_b, radius_a, radius_b, pair_radius, prmin)
1629 LOGICAL, DIMENSION(:), INTENT(IN) :: present_a, present_b
1630 REAL(dp), DIMENSION(:), INTENT(IN) :: radius_a, radius_b
1631 REAL(dp), DIMENSION(:, :), INTENT(OUT) :: pair_radius
1632 REAL(dp), INTENT(IN), OPTIONAL :: prmin
1633
1634 INTEGER :: i, j, nkind
1635 REAL(dp) :: rrmin
1636
1637 nkind = SIZE(present_a)
1638
1639 pair_radius = 0._dp
1640
1641 rrmin = 0.0_dp
1642 IF (PRESENT(prmin)) rrmin = prmin
1643
1644 DO i = 1, nkind
1645 IF (.NOT. present_a(i)) cycle
1646 DO j = 1, nkind
1647 IF (.NOT. present_b(j)) cycle
1648 pair_radius(i, j) = radius_a(i) + radius_b(j)
1649 pair_radius(i, j) = max(pair_radius(i, j), rrmin)
1650 END DO
1651 END DO
1652
1653 END SUBROUTINE pair_radius_setup
1654
1655! **************************************************************************************************
1656!> \brief Print the distribution of the simple pair neighbor list.
1657!> \param ab ...
1658!> \param qs_kind_set ...
1659!> \param output_unit ...
1660!> \param para_env ...
1661!> \date 19.06.2003
1662!> \author MK
1663!> \version 1.0
1664! **************************************************************************************************
1665 SUBROUTINE write_neighbor_distribution(ab, qs_kind_set, output_unit, para_env)
1666 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1667 POINTER :: ab
1668 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
1669 INTEGER, INTENT(in) :: output_unit
1670 TYPE(mp_para_env_type), POINTER :: para_env
1671
1672 CHARACTER(len=*), PARAMETER :: routinen = 'write_neighbor_distribution'
1673 LOGICAL, PARAMETER :: full_output = .false.
1674
1675 INTEGER :: handle, ikind, inode, ipe, jkind, n, &
1676 nkind, nnode
1677 INTEGER(int_8) :: nblock_max, nblock_sum, nelement_max, &
1678 nelement_sum, tmp(2)
1679 INTEGER, ALLOCATABLE, DIMENSION(:) :: nblock, nelement, nnsgf
1680 TYPE(gto_basis_set_type), POINTER :: orb_basis_set
1682 DIMENSION(:), POINTER :: nl_iterator
1683
1684 CALL timeset(routinen, handle)
1685 associate(mype => para_env%mepos + 1, npe => para_env%num_pe)
1686
1687 ! Allocate work storage
1688 ALLOCATE (nblock(npe), nelement(npe))
1689 nblock(:) = 0
1690 nelement(:) = 0
1691 nkind = SIZE(qs_kind_set)
1692 ALLOCATE (nnsgf(nkind))
1693 nnsgf = 1
1694 DO ikind = 1, nkind
1695 CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set)
1696 IF (ASSOCIATED(orb_basis_set)) THEN
1697 CALL get_gto_basis_set(gto_basis_set=orb_basis_set, nsgf=nnsgf(ikind))
1698 END IF
1699 END DO
1700
1701 CALL neighbor_list_iterator_create(nl_iterator, ab)
1702 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1703 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, inode=inode, nnode=nnode)
1704 IF (inode == 1) THEN
1705 n = nnsgf(ikind)*nnsgf(jkind)
1706 nblock(mype) = nblock(mype) + nnode
1707 nelement(mype) = nelement(mype) + n*nnode
1708 END IF
1709 END DO
1710 CALL neighbor_list_iterator_release(nl_iterator)
1711
1712 IF (full_output) THEN
1713 ! XXXXXXXX should gather/scatter this on ionode
1714 CALL para_env%sum(nblock)
1715 CALL para_env%sum(nelement)
1716
1717 nblock_sum = sum(int(nblock, kind=int_8))
1718 nelement_sum = sum(int(nelement, kind=int_8))
1719 ELSE
1720 nblock_sum = nblock(mype)
1721 nblock_max = nblock(mype)
1722 nelement_sum = nelement(mype)
1723 nelement_max = nelement(mype)
1724 tmp = (/nblock_sum, nelement_sum/)
1725 CALL para_env%sum(tmp)
1726 nblock_sum = tmp(1); nelement_sum = tmp(2)
1727 tmp = (/nblock_max, nelement_max/)
1728 CALL para_env%max(tmp)
1729 nblock_max = tmp(1); nelement_max = tmp(2)
1730 END IF
1731
1732 IF (output_unit > 0) THEN
1733 IF (full_output) THEN
1734 WRITE (unit=output_unit, &
1735 fmt="(/,/,T2,A,/,/,T3,A,/,/,(T4,I6,T27,I10,T55,I10))") &
1736 "DISTRIBUTION OF THE NEIGHBOR LISTS", &
1737 "Process Number of particle pairs Number of matrix elements", &
1738 (ipe - 1, nblock(ipe), nelement(ipe), ipe=1, npe)
1739 WRITE (unit=output_unit, fmt="(/,T7,A3,T27,I10,T55,I10)") &
1740 "Sum", sum(nblock), sum(nelement)
1741 ELSE
1742 WRITE (unit=output_unit, fmt="(/,T2,A)") "DISTRIBUTION OF THE NEIGHBOR LISTS"
1743 WRITE (unit=output_unit, fmt="(T15,A,T68,I13)") "Total number of particle pairs:", nblock_sum
1744 WRITE (unit=output_unit, fmt="(T15,A,T68,I13)") "Total number of matrix elements:", nelement_sum
1745 WRITE (unit=output_unit, fmt="(T15,A,T68,I13)") "Average number of particle pairs:", (nblock_sum + npe - 1)/npe
1746 WRITE (unit=output_unit, fmt="(T15,A,T68,I13)") "Maximum number of particle pairs:", nblock_max
1747 WRITE (unit=output_unit, fmt="(T15,A,T68,I13)") "Average number of matrix element:", (nelement_sum + npe - 1)/npe
1748 WRITE (unit=output_unit, fmt="(T15,A,T68,I13)") "Maximum number of matrix elements:", nelement_max
1749 END IF
1750 END IF
1751 END associate
1752
1753 ! Release work storage
1754
1755 DEALLOCATE (nblock, nelement, nnsgf)
1756
1757 CALL timestop(handle)
1758
1759 END SUBROUTINE write_neighbor_distribution
1760
1761! **************************************************************************************************
1762!> \brief Write a set of neighbor lists to the output unit.
1763!> \param ab ...
1764!> \param particle_set ...
1765!> \param cell ...
1766!> \param para_env ...
1767!> \param neighbor_list_section ...
1768!> \param nl_type ...
1769!> \param middle_name ...
1770!> \param nlname ...
1771!> \date 04.03.2002
1772!> \par History
1773!> - Adapted to the new parallelized neighbor list version
1774!> (26.06.2003,MK)
1775!> \author MK
1776!> \version 1.0
1777! **************************************************************************************************
1778 SUBROUTINE write_neighbor_lists(ab, particle_set, cell, para_env, neighbor_list_section, &
1779 nl_type, middle_name, nlname)
1780
1781 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1782 POINTER :: ab
1783 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1784 TYPE(cell_type), POINTER :: cell
1785 TYPE(mp_para_env_type), POINTER :: para_env
1786 TYPE(section_vals_type), POINTER :: neighbor_list_section
1787 CHARACTER(LEN=*), INTENT(IN) :: nl_type, middle_name, nlname
1788
1789 CHARACTER(LEN=default_string_length) :: string, unit_str
1790 INTEGER :: iatom, inode, iw, jatom, nneighbor, nnode
1791 INTEGER, DIMENSION(3) :: cell_b
1792 REAL(dp) :: dab, unit_conv
1793 REAL(dp), DIMENSION(3) :: ra, rab, rb
1794 TYPE(cp_logger_type), POINTER :: logger
1795 TYPE(neighbor_list_iterator_p_type), &
1796 DIMENSION(:), POINTER :: nl_iterator
1797
1798 NULLIFY (logger)
1799 logger => cp_get_default_logger()
1800 IF (btest(cp_print_key_should_output(logger%iter_info, neighbor_list_section, &
1801 trim(nl_type)), &
1802 cp_p_file)) THEN
1803 iw = cp_print_key_unit_nr(logger=logger, &
1804 basis_section=neighbor_list_section, &
1805 print_key_path=trim(nl_type), &
1806 extension=".out", &
1807 middle_name=trim(middle_name), &
1808 local=.true., &
1809 log_filename=.false., &
1810 file_position="REWIND")
1811 associate(mype => para_env%mepos)
1812 CALL section_vals_val_get(neighbor_list_section, "UNIT", c_val=unit_str)
1813 unit_conv = cp_unit_from_cp2k(1.0_dp, trim(unit_str))
1814
1815 ! Print headline
1816 string = ""
1817 WRITE (unit=string, fmt="(A,I5,A)") &
1818 trim(nlname)//" IN "//trim(unit_str)//" (PROCESS", mype, ")"
1819 CALL compress(string)
1820 IF (iw > 0) WRITE (unit=iw, fmt="(/,/,T2,A)") trim(string)
1821
1822 nneighbor = 0
1823
1824 CALL neighbor_list_iterator_create(nl_iterator, ab)
1825 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1826 CALL get_iterator_info(nl_iterator, inode=inode, nnode=nnode, &
1827 iatom=iatom, jatom=jatom, cell=cell_b, r=rab)
1828 nneighbor = nneighbor + 1
1829 ra(:) = pbc(particle_set(iatom)%r, cell)
1830 rb(:) = ra(:) + rab(:)
1831 dab = sqrt(rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3))
1832 IF (iw > 0) THEN
1833 IF (inode == 1) THEN
1834 WRITE (unit=iw, fmt="(/,T2,I5,3X,I6,3X,3F12.6)") &
1835 iatom, nnode, ra(1:3)*unit_conv
1836 END IF
1837 WRITE (unit=iw, fmt="(T10,I6,3X,3I4,3F12.6,2X,F12.6)") &
1838 jatom, cell_b(1:3), rb(1:3)*unit_conv, dab*unit_conv
1839 END IF
1840 END DO
1841 CALL neighbor_list_iterator_release(nl_iterator)
1842
1843 string = ""
1844 WRITE (unit=string, fmt="(A,I12,A,I12)") &
1845 "Total number of neighbor interactions for process", mype, ":", &
1846 nneighbor
1847 CALL compress(string)
1848 IF (iw > 0) WRITE (unit=iw, fmt="(/,T2,A)") trim(string)
1849 CALL cp_print_key_finished_output(unit_nr=iw, &
1850 logger=logger, &
1851 basis_section=neighbor_list_section, &
1852 print_key_path=trim(nl_type), &
1853 local=.true.)
1854 END associate
1855 END IF
1856
1857 END SUBROUTINE write_neighbor_lists
1858
1859END 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, tb_tblite)
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 ...