(git:374b731)
Loading...
Searching...
No Matches
qs_neighbor_list_types.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Define the neighbor list data types and the corresponding functionality
10!> \par History
11!> - cleaned (23.07.2003,MK)
12!> - full refactoring, list iterators (20.10.2010, JGH)
13!> - add get_neighbor_list_set_p, return info for a set of neighborlists
14!> (07.2014,JGH)
15!> \author Matthias Krack (21.06.2000)
16! **************************************************************************************************
18
19 USE kinds, ONLY: dp
20 USE util, ONLY: locate,&
21 sort
22#include "./base/base_uses.f90"
23
24 IMPLICIT NONE
25
26 PRIVATE
27
28! *** Global parameters (in this module) ***
29
30 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_neighbor_list_types'
31
32! *** Definition of the data types for a linked list of neighbors ***
33
34! **************************************************************************************************
35 TYPE neighbor_node_type
36 PRIVATE
37 TYPE(neighbor_node_type), POINTER :: next_neighbor_node
38 REAL(dp), DIMENSION(3) :: r
39 INTEGER, DIMENSION(3) :: cell
40 INTEGER :: neighbor
41 END TYPE neighbor_node_type
42
43! **************************************************************************************************
44 TYPE neighbor_list_type
45 PRIVATE
46 TYPE(neighbor_list_type), POINTER :: next_neighbor_list
47 TYPE(neighbor_node_type), POINTER :: first_neighbor_node, &
48 last_neighbor_node
49 INTEGER :: atom, nnode
50 END TYPE neighbor_list_type
51
52! **************************************************************************************************
54 PRIVATE
55 TYPE(neighbor_list_type), POINTER :: first_neighbor_list, &
56 last_neighbor_list
57 INTEGER :: nlist
58 LOGICAL :: symmetric
60
61! **************************************************************************************************
63 TYPE(neighbor_list_type), POINTER :: neighbor_list
65
66! **************************************************************************************************
68 TYPE(neighbor_list_set_type), POINTER :: neighbor_list_set
69 INTEGER :: nl_size
70 INTEGER :: nl_start
71 INTEGER :: nl_end
72 TYPE(neighbor_list_task_type), DIMENSION(:), POINTER :: nlist_task
74
75! **************************************************************************************************
76 TYPE list_search_type
77 PRIVATE
78 INTEGER :: nlist
79 INTEGER, DIMENSION(:), POINTER :: atom_list
80 INTEGER, DIMENSION(:), POINTER :: atom_index
82 DIMENSION(:), POINTER :: neighbor_list
83 END TYPE list_search_type
84
85! **************************************************************************************************
87 INTEGER :: iatom, jatom, &
88 ikind, jkind, nkind, &
89 ilist, nlist, inode, nnode
90 REAL(kind=dp), DIMENSION(3) :: r
91 INTEGER, DIMENSION(3) :: cell
93 POINTER :: next ! Pointer for forming a linked list of tasks
95
97 MODULE PROCEDURE nl_sub_iterate
98 MODULE PROCEDURE nl_sub_iterate_ref
99 END INTERFACE
100
101! **************************************************************************************************
102! Neighbor List Iterator
103! **************************************************************************************************
104 TYPE neighbor_list_iterator_type
105 PRIVATE
106 INTEGER :: ikind, jkind, ilist, inode
107 INTEGER :: nkind, nlist, nnode
108 INTEGER :: iatom, jatom
109 TYPE(neighbor_list_set_p_type), &
110 DIMENSION(:), POINTER :: nl
111 TYPE(neighbor_list_type), POINTER :: neighbor_list
112 TYPE(neighbor_node_type), POINTER :: neighbor_node
113 TYPE(list_search_type), &
114 DIMENSION(:), POINTER :: list_search
115 END TYPE neighbor_list_iterator_type
116
118 PRIVATE
119 TYPE(neighbor_list_iterator_type), POINTER :: neighbor_list_iterator
120 INTEGER :: last
122! **************************************************************************************************
123
124! *** Public data types ***
125
126 PUBLIC :: neighbor_list_p_type, &
130
131! *** Public subroutines ***
132
133 PUBLIC :: add_neighbor_list, &
141
142! *** Iterator functions and types ***
143
151
152CONTAINS
153
154! **************************************************************************************************
155!> \brief Neighbor list iterator functions
156!> \param iterator_set ...
157!> \param nl ...
158!> \param search ...
159!> \param nthread ...
160!> \date 28.07.2010
161!> \author jhu
162!> \version 1.0
163! **************************************************************************************************
164 SUBROUTINE neighbor_list_iterator_create(iterator_set, nl, search, nthread)
166 DIMENSION(:), POINTER :: iterator_set
167 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
168 POINTER :: nl
169 LOGICAL, INTENT(IN), OPTIONAL :: search
170 INTEGER, INTENT(IN), OPTIONAL :: nthread
171
172 INTEGER :: iatom, il, ilist, mthread, nlist
173 TYPE(list_search_type), DIMENSION(:), POINTER :: list_search
174 TYPE(neighbor_list_iterator_type), POINTER :: iterator
175 TYPE(neighbor_list_type), POINTER :: neighbor_list
176
177 mthread = 1
178 IF (PRESENT(nthread)) mthread = nthread
179
180 ALLOCATE (iterator_set(0:mthread - 1))
181
182 DO il = 0, mthread - 1
183 ALLOCATE (iterator_set(il)%neighbor_list_iterator)
184
185 iterator => iterator_set(il)%neighbor_list_iterator
186
187 iterator%nl => nl
188
189 iterator%ikind = 0
190 iterator%jkind = 0
191 iterator%nkind = nint(sqrt(real(SIZE(nl), dp)))
192
193 iterator%ilist = 0
194 iterator%nlist = 0
195 iterator%inode = 0
196 iterator%nnode = 0
197
198 iterator%iatom = 0
199 iterator%jatom = 0
200
201 NULLIFY (iterator%neighbor_list)
202 NULLIFY (iterator%neighbor_node)
203 NULLIFY (iterator%list_search)
204 END DO
205
206 iterator_set(:)%last = 0
207
208 IF (PRESENT(search)) THEN
209 IF (search) THEN
210 ALLOCATE (list_search(SIZE(nl)))
211 DO il = 1, SIZE(nl)
212 IF (ASSOCIATED(nl(il)%neighbor_list_set)) THEN
213 CALL get_neighbor_list_set(neighbor_list_set=nl(il)%neighbor_list_set, nlist=nlist)
214 list_search(il)%nlist = nlist
215 ALLOCATE (list_search(il)%atom_list(nlist))
216 ALLOCATE (list_search(il)%atom_index(nlist))
217 ALLOCATE (list_search(il)%neighbor_list(nlist))
218
219 NULLIFY (neighbor_list)
220 DO ilist = 1, nlist
221 IF (.NOT. ASSOCIATED(neighbor_list)) THEN
222 neighbor_list => first_list(nl(il)%neighbor_list_set)
223 ELSE
224 neighbor_list => neighbor_list%next_neighbor_list
225 END IF
226 CALL get_neighbor_list(neighbor_list=neighbor_list, atom=iatom)
227 list_search(il)%atom_list(ilist) = iatom
228 list_search(il)%neighbor_list(ilist)%neighbor_list => neighbor_list
229 END DO
230 CALL sort(list_search(il)%atom_list, nlist, list_search(il)%atom_index)
231
232 ELSE
233 list_search(il)%nlist = -1
234 NULLIFY (list_search(il)%atom_list, list_search(il)%atom_index, list_search(il)%neighbor_list)
235 END IF
236 END DO
237 DO il = 0, mthread - 1
238 iterator => iterator_set(il)%neighbor_list_iterator
239 iterator%list_search => list_search
240 END DO
241 END IF
242 END IF
243
244 END SUBROUTINE neighbor_list_iterator_create
245
246! **************************************************************************************************
247!> \brief ...
248!> \param iterator_set ...
249! **************************************************************************************************
250 SUBROUTINE neighbor_list_iterator_release(iterator_set)
252 DIMENSION(:), POINTER :: iterator_set
253
254 INTEGER :: il, mthread
255 TYPE(neighbor_list_iterator_type), POINTER :: iterator
256
257!all threads have the same search list
258
259 iterator => iterator_set(0)%neighbor_list_iterator
260 IF (ASSOCIATED(iterator%list_search)) THEN
261 DO il = 1, SIZE(iterator%list_search)
262 IF (iterator%list_search(il)%nlist >= 0) THEN
263 DEALLOCATE (iterator%list_search(il)%atom_list)
264 DEALLOCATE (iterator%list_search(il)%atom_index)
265 DEALLOCATE (iterator%list_search(il)%neighbor_list)
266 END IF
267 END DO
268 DEALLOCATE (iterator%list_search)
269 END IF
270
271 mthread = SIZE(iterator_set)
272 DO il = 0, mthread - 1
273 DEALLOCATE (iterator_set(il)%neighbor_list_iterator)
274 END DO
275 DEALLOCATE (iterator_set)
276
277 END SUBROUTINE neighbor_list_iterator_release
278
279! **************************************************************************************************
280!> \brief ...
281!> \param iterator_set ...
282!> \param ikind ...
283!> \param jkind ...
284!> \param iatom ...
285!> \param mepos ...
286! **************************************************************************************************
287 SUBROUTINE nl_set_sub_iterator(iterator_set, ikind, jkind, iatom, mepos)
289 DIMENSION(:), POINTER :: iterator_set
290 INTEGER, INTENT(IN) :: ikind, jkind, iatom
291 INTEGER, INTENT(IN), OPTIONAL :: mepos
292
293 INTEGER :: i, ij, ilist, me, nlist, nnode
294 TYPE(list_search_type), POINTER :: list_search
295 TYPE(neighbor_list_iterator_type), POINTER :: iterator
296 TYPE(neighbor_list_type), POINTER :: neighbor_list
297
298 IF (PRESENT(mepos)) THEN
299 me = mepos
300 ELSE
301 me = 0
302 END IF
303
304 ! Set up my thread-local iterator for the list of iatom / jkind nodes
305
306 iterator => iterator_set(me)%neighbor_list_iterator
307 ij = ikind + iterator%nkind*(jkind - 1)
308 IF (ASSOCIATED(iterator%list_search)) THEN
309 list_search => iterator%list_search(ij)
310 nlist = list_search%nlist
311 ilist = 0
312 NULLIFY (neighbor_list)
313 IF (nlist > 0) THEN
314 i = locate(list_search%atom_list, iatom)
315 i = list_search%atom_index(i)
316 IF (i > 0) neighbor_list => list_search%neighbor_list(i)%neighbor_list
317 ilist = i
318 END IF
319 IF (ASSOCIATED(neighbor_list)) THEN
320 CALL get_neighbor_list(neighbor_list=neighbor_list, nnode=nnode)
321 ELSE
322 nnode = 0
323 END IF
324 ELSE
325 cpabort("")
326 END IF
327
328 iterator%ikind = ikind
329 iterator%jkind = jkind
330
331 iterator%ilist = ilist
332 iterator%nlist = nlist
333 iterator%inode = 0
334 iterator%nnode = nnode
335
336 iterator%iatom = iatom
337 iterator%jatom = 0
338
339 iterator%neighbor_list => neighbor_list
340 NULLIFY (iterator%neighbor_node)
341
342 END SUBROUTINE nl_set_sub_iterator
343
344! **************************************************************************************************
345!> \brief ...
346!> \param iterator_set ...
347!> \param mepos ...
348!> \return ...
349! **************************************************************************************************
350 FUNCTION neighbor_list_iterate(iterator_set, mepos) RESULT(istat)
352 DIMENSION(:), POINTER :: iterator_set
353 INTEGER, OPTIONAL :: mepos
354 INTEGER :: istat
355
356 INTEGER :: iab, last, me
357 TYPE(neighbor_list_iterator_type), POINTER :: iterator
358 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
359 POINTER :: nl
360
361 IF (SIZE(iterator_set) .NE. 1 .AND. .NOT. PRESENT(mepos)) &
362 cpabort("Parallel iterator calls must include 'mepos'")
363
364 IF (PRESENT(mepos)) THEN
365 me = mepos
366 ELSE
367 me = 0
368 END IF
369
370 istat = 0
371
372!$OMP CRITICAL(neighbour_list_iterate_critical)
373 last = iterator_set(0)%last
374 IF (last /= me) THEN
375 iterator_set(me)%neighbor_list_iterator = iterator_set(last)%neighbor_list_iterator
376 END IF
377 iterator => iterator_set(me)%neighbor_list_iterator
378 nl => iterator%nl
379
380 IF (iterator%inode < iterator%nnode) THEN
381 ! we can be sure that there is another node in this list
382 iterator%inode = iterator%inode + 1
383 iterator%neighbor_node => iterator%neighbor_node%next_neighbor_node
384 ELSE
385 iab = max(iterator%ikind + iterator%nkind*(iterator%jkind - 1), 0)
386 kindloop: DO ! look for the next list with nnode /= 0
387 listloop: DO
388 IF (iterator%ilist >= iterator%nlist) EXIT listloop
389 iterator%ilist = iterator%ilist + 1
390 IF (ASSOCIATED(iterator%neighbor_list)) THEN
391 iterator%neighbor_list => iterator%neighbor_list%next_neighbor_list
392 ELSE
393 iterator%neighbor_list => first_list(nl(iab)%neighbor_list_set)
394 END IF
395 CALL get_neighbor_list(neighbor_list=iterator%neighbor_list, atom=iterator%iatom, &
396 nnode=iterator%nnode)
397 IF (iterator%nnode > 0) EXIT kindloop
398 END DO listloop
399 IF (iab >= iterator%nkind**2) THEN
400 istat = 1
401 EXIT kindloop
402 ELSE
403 iab = iab + 1
404 iterator%jkind = (iab - 1)/iterator%nkind + 1
405 iterator%ikind = iab - iterator%nkind*(iterator%jkind - 1)
406 iterator%ilist = 0
407 IF (.NOT. ASSOCIATED(nl(iab)%neighbor_list_set)) THEN
408 iterator%ilist = 0
409 iterator%nlist = 0
410 ELSE
411 CALL get_neighbor_list_set(neighbor_list_set= &
412 nl(iab)%neighbor_list_set, nlist=iterator%nlist)
413 iterator%ilist = 0
414 END IF
415 NULLIFY (iterator%neighbor_list)
416 END IF
417 END DO kindloop
418 IF (istat == 0) THEN
419 iterator%inode = 1
420 iterator%neighbor_node => first_node(iterator%neighbor_list)
421 END IF
422 END IF
423 IF (istat == 0) THEN
424 CALL get_neighbor_node(neighbor_node=iterator%neighbor_node, neighbor=iterator%jatom)
425 END IF
426
427 ! mark the last iterator updated
428 iterator_set(:)%last = me
429!$OMP END CRITICAL(neighbour_list_iterate_critical)
430
431 END FUNCTION neighbor_list_iterate
432
433! **************************************************************************************************
434!> \brief ...
435!> \param iterator_set ...
436!> \param mepos ...
437!> \return ...
438! **************************************************************************************************
439 FUNCTION nl_sub_iterate(iterator_set, mepos) RESULT(istat)
441 DIMENSION(:), POINTER :: iterator_set
442 INTEGER, INTENT(IN), OPTIONAL :: mepos
443 INTEGER :: istat
444
445 INTEGER :: me
446 TYPE(neighbor_list_iterator_type), POINTER :: iterator
447
448 ! Each thread's sub-iterator are independent, no need to synchronise with other threads
449
450 IF (PRESENT(mepos)) THEN
451 me = mepos
452 ELSE
453 me = 0
454 END IF
455
456 istat = 0
457
458 iterator => iterator_set(me)%neighbor_list_iterator
459
460 IF (ASSOCIATED(iterator%neighbor_list)) THEN
461 IF (iterator%inode >= iterator%nnode) THEN
462 ! end of loop
463 istat = 1
464 ELSEIF (iterator%inode == 0) THEN
465 iterator%inode = 1
466 iterator%neighbor_node => first_node(iterator%neighbor_list)
467 ELSEIF (iterator%inode > 0) THEN
468 ! we can be sure that there is another node in this list
469 iterator%inode = iterator%inode + 1
470 iterator%neighbor_node => iterator%neighbor_node%next_neighbor_node
471 ELSE
472 cpabort("wrong")
473 END IF
474 ELSE
475 ! no list available
476 istat = 1
477 END IF
478 IF (istat == 0) THEN
479 CALL get_neighbor_node(neighbor_node=iterator%neighbor_node, neighbor=iterator%jatom)
480 END IF
481
482 END FUNCTION nl_sub_iterate
483
484! **************************************************************************************************
485!> \brief wrap nl_sub_iterate s.t. external loop over kinds and calls to nl_set_sub_iterator
486!> are no longer needed. This fixes first atom of iter_sub to second atom of iter_ref.
487!> \param iter_sub ...
488!> \param iter_ref ...
489!> \param mepos ...
490!> \return ...
491! **************************************************************************************************
492 RECURSIVE FUNCTION nl_sub_iterate_ref(iter_sub, iter_ref, mepos) RESULT(iter_stat)
494 DIMENSION(:), POINTER :: iter_sub, iter_ref
495 INTEGER, INTENT(IN), OPTIONAL :: mepos
496 INTEGER :: iter_stat
497
498 INTEGER :: atom_ref, kind_ref, kind_sub, me, nkind
499 TYPE(neighbor_list_iterator_type), POINTER :: iterator
500
501 IF (PRESENT(mepos)) THEN
502 me = mepos
503 ELSE
504 me = 0
505 END IF
506
507 iterator => iter_sub(me)%neighbor_list_iterator
508 kind_sub = iterator%jkind
509
510 CALL get_iterator_info(iter_ref, jatom=atom_ref, jkind=kind_ref)
511
512 IF (iterator%inode == 0) THEN
513 CALL nl_set_sub_iterator(iter_sub, kind_ref, max(kind_sub, 1), atom_ref)
514 END IF
515 iter_stat = nl_sub_iterate(iter_sub)
516 IF (iter_stat == 0) RETURN
517
518 nkind = iterator%nkind
519
520 IF (kind_sub == nkind) THEN
521 CALL nl_set_sub_iterator(iter_sub, kind_ref, 1, atom_ref)
522 RETURN
523 ELSE
524 kind_sub = kind_sub + 1
525 CALL nl_set_sub_iterator(iter_sub, kind_ref, kind_sub, atom_ref)
526 iter_stat = nl_sub_iterate_ref(iter_sub, iter_ref)
527 END IF
528
529 END FUNCTION
530
531! **************************************************************************************************
532!> \brief ...
533!> \param iterator_set ...
534!> \param mepos ...
535!> \param ikind ...
536!> \param jkind ...
537!> \param nkind ...
538!> \param ilist ...
539!> \param nlist ...
540!> \param inode ...
541!> \param nnode ...
542!> \param iatom ...
543!> \param jatom ...
544!> \param r ...
545!> \param cell ...
546! **************************************************************************************************
547 SUBROUTINE get_iterator_info(iterator_set, mepos, &
548 ikind, jkind, nkind, ilist, nlist, inode, nnode, iatom, jatom, r, cell)
550 DIMENSION(:), POINTER :: iterator_set
551 INTEGER, OPTIONAL :: mepos, ikind, jkind, nkind, ilist, &
552 nlist, inode, nnode, iatom, jatom
553 REAL(dp), DIMENSION(3), OPTIONAL :: r
554 INTEGER, DIMENSION(3), OPTIONAL :: cell
555
556 INTEGER :: me
557 TYPE(neighbor_list_iterator_type), POINTER :: iterator
558
559 IF (SIZE(iterator_set) .NE. 1 .AND. .NOT. PRESENT(mepos)) &
560 cpabort("Parallel iterator calls must include 'mepos'")
561
562 IF (PRESENT(mepos)) THEN
563 me = mepos
564 ELSE
565 me = 0
566 END IF
567 iterator => iterator_set(me)%neighbor_list_iterator
568
569 IF (PRESENT(ikind)) ikind = iterator%ikind
570 IF (PRESENT(jkind)) jkind = iterator%jkind
571 IF (PRESENT(nkind)) nkind = iterator%nkind
572 IF (PRESENT(ilist)) ilist = iterator%ilist
573 IF (PRESENT(nlist)) nlist = iterator%nlist
574 IF (PRESENT(inode)) inode = iterator%inode
575 IF (PRESENT(nnode)) nnode = iterator%nnode
576 IF (PRESENT(iatom)) iatom = iterator%iatom
577 IF (PRESENT(jatom)) jatom = iterator%jatom
578 IF (PRESENT(r)) THEN
579 CALL get_neighbor_node(neighbor_node=iterator%neighbor_node, r=r)
580 END IF
581 IF (PRESENT(cell)) THEN
582 CALL get_neighbor_node(neighbor_node=iterator%neighbor_node, cell=cell)
583 END IF
584
585 END SUBROUTINE get_iterator_info
586
587! **************************************************************************************************
588!> \brief Captures the current state of the iterator in a neighbor_list_task_type
589!> \param iterator_set the iterator / array of iterators (for multiple threads)
590!> \param task the task structure which is returned
591!> \param mepos OpenMP thread index
592! **************************************************************************************************
593 SUBROUTINE get_iterator_task(iterator_set, task, mepos)
595 DIMENSION(:), POINTER :: iterator_set
596 TYPE(neighbor_list_task_type), INTENT(OUT) :: task
597 INTEGER, OPTIONAL :: mepos
598
599 IF (PRESENT(mepos)) THEN
600 CALL get_iterator_info(iterator_set, mepos=mepos, ikind=task%ikind, jkind=task%jkind, &
601 nkind=task%nkind, &
602 ilist=task%ilist, nlist=task%nlist, &
603 inode=task%inode, nnode=task%nnode, &
604 iatom=task%iatom, jatom=task%jatom, &
605 r=task%r, cell=task%cell)
606 ELSE
607 CALL get_iterator_info(iterator_set, ikind=task%ikind, jkind=task%jkind, &
608 nkind=task%nkind, &
609 ilist=task%ilist, nlist=task%nlist, &
610 inode=task%inode, nnode=task%nnode, &
611 iatom=task%iatom, jatom=task%jatom, &
612 r=task%r, cell=task%cell)
613 END IF
614
615 NULLIFY (task%next)
616
617 END SUBROUTINE get_iterator_task
618
619! **************************************************************************************************
620!> \brief Add a new neighbor list to a neighbor list set.
621!> \param neighbor_list_set ...
622!> \param atom ...
623!> \param neighbor_list ...
624!> \date 13.09.2000
625!> \author MK
626!> \version 1.0
627! **************************************************************************************************
628 SUBROUTINE add_neighbor_list(neighbor_list_set, atom, neighbor_list)
629
630 TYPE(neighbor_list_set_type), POINTER :: neighbor_list_set
631 INTEGER, INTENT(IN) :: atom
632 TYPE(neighbor_list_type), POINTER :: neighbor_list
633
634 TYPE(neighbor_list_type), POINTER :: new_neighbor_list
635
636 IF (ASSOCIATED(neighbor_list_set)) THEN
637
638 IF (ASSOCIATED(neighbor_list_set%last_neighbor_list)) THEN
639
640 new_neighbor_list => &
641 neighbor_list_set%last_neighbor_list%next_neighbor_list
642
643 IF (.NOT. ASSOCIATED(new_neighbor_list)) THEN
644
645! *** Allocate a new neighbor list ***
646
647 ALLOCATE (new_neighbor_list)
648
649 NULLIFY (new_neighbor_list%next_neighbor_list)
650 NULLIFY (new_neighbor_list%first_neighbor_node)
651
652! *** Link the new neighbor list to the neighbor list set ***
653
654 neighbor_list_set%last_neighbor_list%next_neighbor_list => new_neighbor_list
655
656 END IF
657
658 ELSE
659
660 new_neighbor_list => neighbor_list_set%first_neighbor_list
661
662 IF (.NOT. ASSOCIATED(new_neighbor_list)) THEN
663
664! *** Allocate a new first neighbor list ***
665
666 ALLOCATE (new_neighbor_list)
667
668 NULLIFY (new_neighbor_list%next_neighbor_list)
669 NULLIFY (new_neighbor_list%first_neighbor_node)
670
671! *** Link the new first neighbor list to the neighbor list set ***
672
673 neighbor_list_set%first_neighbor_list => new_neighbor_list
674
675 END IF
676
677 END IF
678
679! *** Store the data set of the new neighbor list ***
680
681 NULLIFY (new_neighbor_list%last_neighbor_node)
682 new_neighbor_list%atom = atom
683 new_neighbor_list%nnode = 0
684
685! *** Update the pointer to the last neighbor ***
686! *** list of the neighbor list set ***
687
688 neighbor_list_set%last_neighbor_list => new_neighbor_list
689
690! *** Increment the neighbor list counter ***
691
692 neighbor_list_set%nlist = neighbor_list_set%nlist + 1
693
694! *** Return a pointer to the new neighbor list ***
695
696 neighbor_list => new_neighbor_list
697
698 ELSE
699
700 cpabort("The requested neighbor list set is not associated")
701
702 END IF
703
704 END SUBROUTINE add_neighbor_list
705
706! **************************************************************************************************
707!> \brief Add a new neighbor list node to a neighbor list.
708!> \param neighbor_list ...
709!> \param neighbor ...
710!> \param cell ...
711!> \param r ...
712!> \param exclusion_list ...
713!> \param nkind ...
714!> \date 23.06.2000
715!> \author MK
716!> \version 1.0
717! **************************************************************************************************
718 SUBROUTINE add_neighbor_node(neighbor_list, neighbor, cell, r, exclusion_list, nkind)
719
720 TYPE(neighbor_list_type), POINTER :: neighbor_list
721 INTEGER, INTENT(IN) :: neighbor
722 INTEGER, DIMENSION(3), INTENT(IN) :: cell
723 REAL(dp), DIMENSION(3), INTENT(IN) :: r
724 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: exclusion_list
725 INTEGER, INTENT(IN), OPTIONAL :: nkind
726
727 INTEGER :: iatom, my_nkind
728 TYPE(neighbor_node_type), POINTER :: new_neighbor_node
729
730 IF (ASSOCIATED(neighbor_list)) THEN
731
732! *** Check for exclusions ***
733
734 IF (PRESENT(exclusion_list)) THEN
735 IF (ASSOCIATED(exclusion_list)) THEN
736 DO iatom = 1, SIZE(exclusion_list)
737 IF (exclusion_list(iatom) == 0) EXIT
738 IF (exclusion_list(iatom) == neighbor) RETURN
739 END DO
740 END IF
741 END IF
742
743 my_nkind = 0
744 IF (PRESENT(nkind)) my_nkind = nkind
745
746 IF (ASSOCIATED(neighbor_list%last_neighbor_node)) THEN
747
748 new_neighbor_node => neighbor_list%last_neighbor_node%next_neighbor_node
749
750 IF (.NOT. ASSOCIATED(new_neighbor_node)) THEN
751
752! *** Allocate a new neighbor node ***
753
754 ALLOCATE (new_neighbor_node)
755
756 NULLIFY (new_neighbor_node%next_neighbor_node)
757
758! *** Link the new neighbor node to the neighbor list ***
759
760 neighbor_list%last_neighbor_node%next_neighbor_node => new_neighbor_node
761
762 END IF
763
764 ELSE
765
766 new_neighbor_node => neighbor_list%first_neighbor_node
767
768 IF (.NOT. ASSOCIATED(new_neighbor_node)) THEN
769
770! *** Allocate a new first neighbor node ***
771
772 ALLOCATE (new_neighbor_node)
773
774 NULLIFY (new_neighbor_node%next_neighbor_node)
775
776! *** Link the new first neighbor node to the neighbor list ***
777
778 neighbor_list%first_neighbor_node => new_neighbor_node
779
780 END IF
781
782 END IF
783
784! *** Store the data set of the new neighbor ***
785
786 new_neighbor_node%neighbor = neighbor
787 new_neighbor_node%cell(:) = cell(:)
788 new_neighbor_node%r(:) = r(:)
789
790! *** Update the pointer to the last neighbor node of the neighbor list ***
791
792 neighbor_list%last_neighbor_node => new_neighbor_node
793
794! *** Increment the neighbor node counter ***
795
796 neighbor_list%nnode = neighbor_list%nnode + 1
797
798 ELSE
799
800 cpabort("The requested neighbor list is not associated")
801
802 END IF
803
804 END SUBROUTINE add_neighbor_node
805
806! **************************************************************************************************
807!> \brief Allocate and initialize a set of neighbor lists.
808!> \param neighbor_list_set ...
809!> \param symmetric ...
810!> \date 23.06.2000
811!> \author MK
812!> \version 1.0
813! **************************************************************************************************
814 SUBROUTINE allocate_neighbor_list_set(neighbor_list_set, symmetric)
815
816 TYPE(neighbor_list_set_type), POINTER :: neighbor_list_set
817 LOGICAL, INTENT(IN) :: symmetric
818
819! *** Deallocate the old neighbor list set ***
820
821 IF (ASSOCIATED(neighbor_list_set)) THEN
822 CALL deallocate_neighbor_list_set(neighbor_list_set)
823 END IF
824
825! *** Allocate a set of neighbor lists ***
826
827 ALLOCATE (neighbor_list_set)
828
829 NULLIFY (neighbor_list_set%first_neighbor_list)
830
831! *** Initialize the pointers to the first neighbor list ***
832
833 CALL init_neighbor_list_set(neighbor_list_set, symmetric)
834
835 END SUBROUTINE allocate_neighbor_list_set
836
837! **************************************************************************************************
838!> \brief Deallocate a neighbor list.
839!> \param neighbor_list ...
840!> \date 20.09.2002
841!> \author MK
842!> \version 1.0
843! **************************************************************************************************
844 SUBROUTINE deallocate_neighbor_list(neighbor_list)
845
846 TYPE(neighbor_list_type), POINTER :: neighbor_list
847
848 TYPE(neighbor_node_type), POINTER :: neighbor_node, next_neighbor_node
849
850 IF (ASSOCIATED(neighbor_list)) THEN
851
852 neighbor_node => neighbor_list%first_neighbor_node
853
854 DO WHILE (ASSOCIATED(neighbor_node))
855 next_neighbor_node => neighbor_node%next_neighbor_node
856 DEALLOCATE (neighbor_node)
857 neighbor_node => next_neighbor_node
858 END DO
859
860 DEALLOCATE (neighbor_list)
861
862 END IF
863
864 END SUBROUTINE deallocate_neighbor_list
865
866! **************************************************************************************************
867!> \brief Deallocate a neighbor list set.
868!> \param neighbor_list_set ...
869!> \date 03.11.2000
870!> \author MK
871!> \version 1.0
872! **************************************************************************************************
873 SUBROUTINE deallocate_neighbor_list_set(neighbor_list_set)
874 TYPE(neighbor_list_set_type), POINTER :: neighbor_list_set
875
876 TYPE(neighbor_list_type), POINTER :: neighbor_list, next_neighbor_list
877
878 IF (ASSOCIATED(neighbor_list_set)) THEN
879
880 neighbor_list => neighbor_list_set%first_neighbor_list
881
882 DO WHILE (ASSOCIATED(neighbor_list))
883 next_neighbor_list => neighbor_list%next_neighbor_list
884 CALL deallocate_neighbor_list(neighbor_list)
885 neighbor_list => next_neighbor_list
886 END DO
887
888 DEALLOCATE (neighbor_list_set)
889
890 END IF
891
892 END SUBROUTINE deallocate_neighbor_list_set
893
894! **************************************************************************************************
895!> \brief Return a pointer to the first neighbor list of a neighbor list set.
896!> \param neighbor_list_set ...
897!> \return ...
898!> \date 13.09.2000
899!> \author MK
900!> \version 1.0
901! **************************************************************************************************
902 FUNCTION first_list(neighbor_list_set) RESULT(first_neighbor_list)
903
904 TYPE(neighbor_list_set_type), POINTER :: neighbor_list_set
905 TYPE(neighbor_list_type), POINTER :: first_neighbor_list
906
907 first_neighbor_list => neighbor_list_set%first_neighbor_list
908
909 END FUNCTION first_list
910
911! **************************************************************************************************
912!> \brief Return a pointer to the first neighbor node of a neighbor list.
913!> \param neighbor_list ...
914!> \return ...
915!> \date 23.06.2000,
916!> \author MK
917!> \version 1.0
918! **************************************************************************************************
919 FUNCTION first_node(neighbor_list) RESULT(first_neighbor_node)
920
921 TYPE(neighbor_list_type), POINTER :: neighbor_list
922 TYPE(neighbor_node_type), POINTER :: first_neighbor_node
923
924 first_neighbor_node => neighbor_list%first_neighbor_node
925
926 END FUNCTION first_node
927
928! **************************************************************************************************
929!> \brief Return the requested data of a neighbor list.
930!> \param neighbor_list ...
931!> \param atom ...
932!> \param nnode ...
933!> \date 13.09.2000
934!> \author MK
935!> \version 1.0
936! **************************************************************************************************
937 SUBROUTINE get_neighbor_list(neighbor_list, atom, nnode)
938
939 TYPE(neighbor_list_type), POINTER :: neighbor_list
940 INTEGER, INTENT(OUT), OPTIONAL :: atom, nnode
941
942 IF (ASSOCIATED(neighbor_list)) THEN
943
944 IF (PRESENT(atom)) atom = neighbor_list%atom
945 IF (PRESENT(nnode)) nnode = neighbor_list%nnode
946
947 ELSE
948
949 cpabort("The requested neighbor list is not associated")
950
951 END IF
952
953 END SUBROUTINE get_neighbor_list
954
955! **************************************************************************************************
956!> \brief Return the components of a neighbor list set.
957!> \param neighbor_list_set ...
958!> \param nlist ...
959!> \param symmetric ...
960!> \date 10.11.2000
961!> \author MK
962!> \version 1.0
963! **************************************************************************************************
964 SUBROUTINE get_neighbor_list_set(neighbor_list_set, nlist, symmetric)
965
966 TYPE(neighbor_list_set_type), POINTER :: neighbor_list_set
967 INTEGER, INTENT(OUT), OPTIONAL :: nlist
968 LOGICAL, INTENT(OUT), OPTIONAL :: symmetric
969
970 IF (ASSOCIATED(neighbor_list_set)) THEN
971
972 IF (PRESENT(nlist)) nlist = neighbor_list_set%nlist
973 IF (PRESENT(symmetric)) symmetric = neighbor_list_set%symmetric
974
975 ELSE
976
977 cpabort("The requested neighbor list set is not associated")
978
979 END IF
980
981 END SUBROUTINE get_neighbor_list_set
982
983! **************************************************************************************************
984!> \brief Return the components of the first neighbor list set.
985!> \param neighbor_list_sets ...
986!> \param nlist ...
987!> \param symmetric ...
988!> \date 07.2014
989!> \author JGH
990!> \version 1.0
991! **************************************************************************************************
992 SUBROUTINE get_neighbor_list_set_p(neighbor_list_sets, nlist, symmetric)
993
994 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
995 POINTER :: neighbor_list_sets
996 INTEGER, INTENT(OUT), OPTIONAL :: nlist
997 LOGICAL, INTENT(OUT), OPTIONAL :: symmetric
998
999 INTEGER :: i
1000 TYPE(neighbor_list_set_type), POINTER :: neighbor_list_set
1001
1002 IF (ASSOCIATED(neighbor_list_sets)) THEN
1003
1004 NULLIFY (neighbor_list_set)
1005 DO i = 1, SIZE(neighbor_list_sets)
1006 neighbor_list_set => neighbor_list_sets(i)%neighbor_list_set
1007 IF (ASSOCIATED(neighbor_list_set)) EXIT
1008 END DO
1009
1010 IF (ASSOCIATED(neighbor_list_set)) THEN
1011 IF (PRESENT(nlist)) nlist = neighbor_list_set%nlist
1012 IF (PRESENT(symmetric)) symmetric = neighbor_list_set%symmetric
1013 ELSE
1014 CALL cp_abort(__location__, "No neighbor list set is associated. "// &
1015 "Did you specify *all* required basis-sets, eg. for ADMM?")
1016 END IF
1017
1018 ELSE
1019
1020 cpabort("The requested neighbor list sets are not associated")
1021
1022 END IF
1023
1024 END SUBROUTINE get_neighbor_list_set_p
1025
1026! **************************************************************************************************
1027!> \brief Return the requested data of a neighbor node.
1028!> \param neighbor_node ...
1029!> \param neighbor ...
1030!> \param cell ...
1031!> \param r ...
1032!> \date 23.06.2000
1033!> \author MK
1034!> \version 1.0
1035! **************************************************************************************************
1036 SUBROUTINE get_neighbor_node(neighbor_node, neighbor, cell, r)
1037
1038 TYPE(neighbor_node_type), POINTER :: neighbor_node
1039 INTEGER, INTENT(OUT), OPTIONAL :: neighbor
1040 INTEGER, DIMENSION(3), INTENT(OUT), OPTIONAL :: cell
1041 REAL(dp), DIMENSION(3), INTENT(OUT), OPTIONAL :: r
1042
1043 IF (ASSOCIATED(neighbor_node)) THEN
1044
1045 IF (PRESENT(neighbor)) neighbor = neighbor_node%neighbor
1046 IF (PRESENT(r)) r(:) = neighbor_node%r(:)
1047 IF (PRESENT(cell)) cell(:) = neighbor_node%cell(:)
1048
1049 ELSE
1050
1051 cpabort("The requested neighbor node is not associated")
1052
1053 END IF
1054
1055 END SUBROUTINE get_neighbor_node
1056
1057! **************************************************************************************************
1058!> \brief Initialize a neighbor list set. Nothing is (de)allocated here.
1059!> This routine is also used to prepare a neighbor list set for
1060!> overwriting.
1061!> \param neighbor_list_set ...
1062!> \param symmetric ...
1063!> \date 20.09.2002
1064!> \author MK
1065!> \version 1.0
1066! **************************************************************************************************
1067 SUBROUTINE init_neighbor_list_set(neighbor_list_set, symmetric)
1068
1069 TYPE(neighbor_list_set_type), POINTER :: neighbor_list_set
1070 LOGICAL, INTENT(IN) :: symmetric
1071
1072 IF (ASSOCIATED(neighbor_list_set)) THEN
1073
1074 ! *** Initialize the pointers to the last neighbor list ***
1075 NULLIFY (neighbor_list_set%last_neighbor_list)
1076
1077 ! *** Initialize the neighbor list counter ***
1078 neighbor_list_set%nlist = 0
1079
1080 ! *** Initialize the neighbor list build properties
1081 neighbor_list_set%symmetric = symmetric
1082
1083 ELSE
1084
1085 cpabort("The requested neighbor list set is not associated")
1086
1087 END IF
1088
1089 END SUBROUTINE init_neighbor_list_set
1090
1091! **************************************************************************************************
1092!> \brief releases an array of neighbor_list_sets
1093!> \param nlists ...
1094!> \author Ole Schuett
1095! **************************************************************************************************
1097 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1098 POINTER :: nlists
1099
1100 INTEGER :: i
1101
1102 IF (ASSOCIATED(nlists)) THEN
1103 DO i = 1, SIZE(nlists)
1104 CALL deallocate_neighbor_list_set(nlists(i)%neighbor_list_set)
1105 END DO
1106 IF (ASSOCIATED(nlists(1)%nlist_task)) THEN
1107 DEALLOCATE (nlists(1)%nlist_task)
1108 END IF
1109 DEALLOCATE (nlists)
1110 END IF
1111 END SUBROUTINE release_neighbor_list_sets
1112
1113END MODULE qs_neighbor_list_types
Definition atom.F:9
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
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 deallocate_neighbor_list_set(neighbor_list_set)
Deallocate a neighbor list set.
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 nl_set_sub_iterator(iterator_set, ikind, jkind, iatom, mepos)
...
subroutine, public get_neighbor_list_set(neighbor_list_set, nlist, symmetric)
Return the components of a neighbor list set.
subroutine, public neighbor_list_iterator_release(iterator_set)
...
subroutine, public get_neighbor_list_set_p(neighbor_list_sets, nlist, symmetric)
Return the components of the first neighbor list 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)
...
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