CP2K 2.4 (Revision 12889)

qs_neighbor_list_types.f90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------------!
00002 !   CP2K: A general program to perform molecular dynamics simulations         !
00003 !   Copyright (C) 2000 - 2013  CP2K developers group                          !
00004 !-----------------------------------------------------------------------------!
00005 
00006 ! *****************************************************************************
00013 MODULE qs_neighbor_list_types
00014 
00015   USE kinds,                           ONLY: dp,&
00016                                              int_size
00017   USE termination,                     ONLY: stop_memory,&
00018                                              stop_program
00019   USE util,                            ONLY: locate,&
00020                                              sort
00021 #include "cp_common_uses.h"
00022 
00023   IMPLICIT NONE
00024 
00025   PRIVATE
00026 
00027 ! *** Global parameters (in this module) ***
00028 
00029   CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_neighbor_list_types'
00030 
00031 ! *** Definition of the data types for a linked list of neighbors ***
00032 
00033 ! *****************************************************************************
00034   TYPE neighbor_node_type
00035     PRIVATE
00036     TYPE(neighbor_node_type), POINTER :: next_neighbor_node
00037     REAL(dp), DIMENSION(3)            :: r
00038     INTEGER, DIMENSION(3)             :: cell
00039     INTEGER                           :: neighbor
00040   END TYPE neighbor_node_type
00041 
00042 ! *****************************************************************************
00043   TYPE neighbor_list_type
00044     PRIVATE
00045     TYPE(neighbor_list_type), POINTER :: next_neighbor_list
00046     TYPE(neighbor_node_type), POINTER :: first_neighbor_node,
00047                                          last_neighbor_node
00048     INTEGER                           :: atom,nnode
00049   END TYPE neighbor_list_type
00050 
00051 ! *****************************************************************************
00052   TYPE neighbor_list_set_type
00053     PRIVATE
00054     TYPE(neighbor_list_type), POINTER :: first_neighbor_list,
00055                                          last_neighbor_list
00056     REAL(dp)                          :: r_max
00057     INTEGER                           :: nlist
00058     LOGICAL                           :: symmetric
00059     LOGICAL                           :: mic
00060     LOGICAL                           :: molecular
00061   END TYPE neighbor_list_set_type
00062 
00063 ! *****************************************************************************
00064   TYPE neighbor_list_p_type
00065     TYPE(neighbor_list_type), POINTER :: neighbor_list
00066   END TYPE neighbor_list_p_type
00067 
00068 ! *****************************************************************************
00069   TYPE neighbor_list_set_p_type
00070     TYPE(neighbor_list_set_type), POINTER :: neighbor_list_set
00071   END TYPE neighbor_list_set_p_type
00072 
00073 ! *****************************************************************************
00074   TYPE list_search_type
00075     PRIVATE
00076     INTEGER                               :: nlist
00077     INTEGER, DIMENSION(:), POINTER        :: atom_list
00078     INTEGER, DIMENSION(:), POINTER        :: atom_index
00079     TYPE(neighbor_list_p_type), 
00080       DIMENSION(:), POINTER               :: neighbor_list
00081   END TYPE list_search_type
00082 
00083 ! *****************************************************************************
00084 ! Neighbor List Iterator
00085 ! *****************************************************************************
00086   TYPE neighbor_list_iterator_type
00087     PRIVATE
00088     INTEGER                               :: ikind, jkind, ilist, inode
00089     INTEGER                               :: nkind, nlist, nnode
00090     INTEGER                               :: iatom, jatom
00091     TYPE(neighbor_list_set_p_type), 
00092       DIMENSION(:), POINTER               :: nl
00093     TYPE(neighbor_list_type), POINTER     :: neighbor_list
00094     TYPE(neighbor_node_type), POINTER     :: neighbor_node
00095     TYPE(list_search_type), 
00096       DIMENSION(:), POINTER               :: list_search
00097   END TYPE neighbor_list_iterator_type
00098 
00099   TYPE neighbor_list_iterator_p_type
00100     PRIVATE
00101     TYPE(neighbor_list_iterator_type), POINTER :: neighbor_list_iterator
00102     INTEGER                                    :: last
00103   END TYPE neighbor_list_iterator_p_type
00104 ! *****************************************************************************
00105 
00106 ! *** Public data types ***
00107 
00108   PUBLIC :: neighbor_list_p_type,&
00109             neighbor_list_set_type,&
00110             neighbor_list_set_p_type,&
00111             neighbor_list_type
00112 
00113 ! *** Public subroutines ***
00114 
00115   PUBLIC :: add_neighbor_list,&
00116             add_neighbor_node,&
00117             allocate_neighbor_list_set,&
00118             deallocate_neighbor_list_set,&
00119             get_neighbor_list_set
00120 
00121 ! *** Iterator functions and types ***
00122 
00123   PUBLIC :: neighbor_list_iterator_p_type,&
00124             neighbor_list_iterator_type,&
00125             neighbor_list_iterator_create,&
00126             neighbor_list_iterator_release,&
00127             neighbor_list_iterate,&
00128             nl_set_sub_iterator,&
00129             nl_sub_iterate,&
00130             get_iterator_info
00131 
00132 CONTAINS
00133 
00134 ! *****************************************************************************
00140   SUBROUTINE neighbor_list_iterator_create(iterator_set,nl,search,nthread)
00141     TYPE(neighbor_list_iterator_p_type), 
00142       DIMENSION(:), POINTER                  :: iterator_set
00143     TYPE(neighbor_list_set_p_type), 
00144       DIMENSION(:), POINTER                  :: nl
00145     LOGICAL, INTENT(IN), OPTIONAL            :: search
00146     INTEGER, INTENT(IN), OPTIONAL            :: nthread
00147 
00148     CHARACTER(len=*), PARAMETER :: 
00149       routineN = 'neighbor_list_iterator_create', 
00150       routineP = moduleN//':'//routineN
00151 
00152     INTEGER                                  :: iatom, il, ilist, istat, 
00153                                                 mthread, nlist
00154     TYPE(list_search_type), DIMENSION(:), 
00155       POINTER                                :: list_search
00156     TYPE(neighbor_list_iterator_type), 
00157       POINTER                                :: iterator
00158     TYPE(neighbor_list_type), POINTER        :: neighbor_list
00159 
00160     mthread = 1
00161     IF(PRESENT(nthread)) mthread=nthread
00162 
00163     ALLOCATE(iterator_set(0:mthread-1),STAT=istat)
00164     IF (istat /= 0) CALL stop_memory(routineN,moduleN,__LINE__,&
00165                                      "iterator_set",int_size*mthread)
00166 
00167     DO il=0,mthread-1
00168        ALLOCATE (iterator_set(il)%neighbor_list_iterator,STAT=istat)
00169        IF (istat /= 0) CALL stop_memory(routineN,moduleN,__LINE__,&
00170                                         "neighbor_list_iterator",0)
00171 
00172        iterator => iterator_set(il)%neighbor_list_iterator
00173 
00174        iterator%nl => nl
00175 
00176        iterator%ikind = 0
00177        iterator%jkind = 0
00178        iterator%nkind = NINT(SQRT(REAL(SIZE(nl),dp)))
00179 
00180        iterator%ilist = 0
00181        iterator%nlist = 0
00182        iterator%inode = 0
00183        iterator%nnode = 0
00184 
00185        iterator%iatom = 0
00186        iterator%jatom = 0
00187 
00188        NULLIFY(iterator%neighbor_list)
00189        NULLIFY(iterator%neighbor_node)
00190        NULLIFY(iterator%list_search)
00191     END DO
00192 
00193     iterator_set(:)%last = 0
00194 
00195     IF (PRESENT(search)) THEN
00196       IF (search) THEN
00197          ALLOCATE(list_search(SIZE(nl)),STAT=istat)
00198          IF (istat /= 0) CALL stop_memory(routineN,moduleN,__LINE__,&
00199                                           "list_search",0)
00200          DO il=1,SIZE(nl)
00201            IF (ASSOCIATED(nl(il)%neighbor_list_set)) THEN
00202               CALL get_neighbor_list_set(neighbor_list_set=nl(il)%neighbor_list_set,nlist=nlist)
00203               list_search(il)%nlist = nlist
00204               ALLOCATE (list_search(il)%atom_list(nlist),STAT=istat)
00205               IF (istat /= 0) CALL stop_memory(routineN,moduleN,__LINE__,&
00206                                                "list_search(il)%atom_list",int_size*nlist)
00207               ALLOCATE (list_search(il)%atom_index(nlist),STAT=istat)
00208               IF (istat /= 0) CALL stop_memory(routineN,moduleN,__LINE__,&
00209                                                "list_search(il)%atom_index",int_size*nlist)
00210               ALLOCATE (list_search(il)%neighbor_list(nlist),STAT=istat)
00211               IF (istat /= 0) CALL stop_memory(routineN,moduleN,__LINE__,&
00212                                                "list_search(il)%neighbor_list",int_size*nlist)
00213 
00214               NULLIFY ( neighbor_list )
00215               DO ilist=1,nlist
00216                  IF ( .NOT. ASSOCIATED(neighbor_list) ) THEN
00217                     neighbor_list => first_list(nl(il)%neighbor_list_set)
00218                  ELSE
00219                     neighbor_list => neighbor_list%next_neighbor_list
00220                  END IF
00221                  CALL get_neighbor_list(neighbor_list=neighbor_list,atom=iatom)
00222                  list_search(il)%atom_list(ilist) = iatom
00223                  list_search(il)%neighbor_list(ilist)%neighbor_list => neighbor_list
00224               END DO
00225               CALL sort (list_search(il)%atom_list,nlist,list_search(il)%atom_index)
00226 
00227            ELSE
00228               list_search(il)%nlist = -1
00229               NULLIFY(list_search(il)%atom_list,list_search(il)%atom_index,list_search(il)%neighbor_list)
00230            END IF
00231          END DO
00232          DO il=0,mthread-1
00233             iterator => iterator_set(il)%neighbor_list_iterator
00234             iterator%list_search => list_search
00235          END DO
00236       END IF
00237     END IF
00238 
00239   END SUBROUTINE neighbor_list_iterator_create
00240 
00241   SUBROUTINE neighbor_list_iterator_release(iterator_set)
00242     TYPE(neighbor_list_iterator_p_type), 
00243       DIMENSION(:), POINTER                  :: iterator_set
00244 
00245     CHARACTER(len=*), PARAMETER :: 
00246       routineN = 'neighbor_list_iterator_release', 
00247       routineP = moduleN//':'//routineN
00248 
00249     INTEGER                                  :: il, istat, mthread
00250     TYPE(neighbor_list_iterator_type), 
00251       POINTER                                :: iterator
00252 
00253 !all threads have the same search list
00254 
00255     iterator => iterator_set(0)%neighbor_list_iterator
00256     IF (ASSOCIATED(iterator%list_search)) THEN
00257        DO il=1,SIZE(iterator%list_search)
00258           IF (iterator%list_search(il)%nlist >= 0) THEN
00259              DEALLOCATE (iterator%list_search(il)%atom_list,STAT=istat)
00260              IF (istat /= 0) CALL stop_memory(routineN,moduleN,__LINE__,&
00261                                               "iterator%list_search(il)%atom_list")
00262              DEALLOCATE (iterator%list_search(il)%atom_index,STAT=istat)
00263              IF (istat /= 0) CALL stop_memory(routineN,moduleN,__LINE__,&
00264                                               "iterator%list_search(il)%atom_index")
00265              DEALLOCATE (iterator%list_search(il)%neighbor_list,STAT=istat)
00266              IF (istat /= 0) CALL stop_memory(routineN,moduleN,__LINE__,&
00267                                               "iterator%list_search(il)%neighbor_list")
00268           END IF
00269        END DO
00270        DEALLOCATE (iterator%list_search,STAT=istat)
00271        IF (istat /= 0) CALL stop_memory(routineN,moduleN,__LINE__,"iterator%list_search")
00272     END IF
00273 
00274     mthread = SIZE(iterator_set)
00275     DO il=0,mthread-1
00276        DEALLOCATE (iterator_set(il)%neighbor_list_iterator,STAT=istat)
00277        IF (istat /= 0) CALL stop_memory(routineN,moduleN,__LINE__,"neighbor_list_iterator")
00278     END DO
00279     DEALLOCATE (iterator_set,STAT=istat)
00280     IF (istat /= 0) CALL stop_memory(routineN,moduleN,__LINE__,"iterator_set")
00281 
00282   END SUBROUTINE neighbor_list_iterator_release
00283 
00284   SUBROUTINE nl_set_sub_iterator(iterator_set,ikind,jkind,iatom)
00285     TYPE(neighbor_list_iterator_p_type), 
00286       DIMENSION(:), POINTER                  :: iterator_set
00287     INTEGER, INTENT(IN)                      :: ikind, jkind, iatom
00288 
00289     CHARACTER(len=*), PARAMETER :: routineN = 'nl_set_sub_iterator', 
00290       routineP = moduleN//':'//routineN
00291 
00292     INTEGER                                  :: i, ij, il, ilist, mthread, 
00293                                                 nlist, nnode
00294     TYPE(list_search_type), POINTER          :: list_search
00295     TYPE(neighbor_list_iterator_type), 
00296       POINTER                                :: iterator
00297     TYPE(neighbor_list_type), POINTER        :: neighbor_list
00298 
00299     iterator => iterator_set(0)%neighbor_list_iterator
00300     ij = ikind + iterator%nkind*(jkind - 1)
00301     IF (ASSOCIATED(iterator%list_search)) THEN
00302        list_search => iterator%list_search(ij)
00303        nlist = list_search%nlist
00304        ilist = 0
00305        NULLIFY(neighbor_list)
00306        IF (nlist > 0) THEN
00307          i = locate(list_search%atom_list,iatom)
00308          i = list_search%atom_index(i)
00309          IF (i > 0) neighbor_list => list_search%neighbor_list(i)%neighbor_list
00310          ilist = i
00311        END IF
00312        IF (ASSOCIATED(neighbor_list)) THEN
00313           CALL get_neighbor_list(neighbor_list=neighbor_list,nnode=nnode)
00314        ELSE
00315           nnode=0
00316        END IF
00317     ELSE
00318        STOP
00319     END IF
00320 
00321     mthread = SIZE(iterator_set)
00322     DO il=0,mthread-1
00323        iterator => iterator_set(il)%neighbor_list_iterator
00324 
00325        iterator%ikind = ikind
00326        iterator%jkind = jkind
00327 
00328        iterator%ilist = ilist
00329        iterator%nlist = nlist
00330        iterator%inode = 0
00331        iterator%nnode = nnode
00332 
00333        iterator%iatom = iatom
00334        iterator%jatom = 0
00335 
00336        iterator%neighbor_list => neighbor_list
00337        NULLIFY(iterator%neighbor_node)
00338     END DO
00339 
00340     iterator_set(:)%last = 0
00341 
00342   END SUBROUTINE nl_set_sub_iterator
00343 
00344   FUNCTION neighbor_list_iterate(iterator_set,mepos) RESULT(istat)
00345     TYPE(neighbor_list_iterator_p_type), 
00346       DIMENSION(:), POINTER                  :: iterator_set
00347     INTEGER, OPTIONAL                        :: mepos
00348     INTEGER                                  :: istat
00349 
00350     CHARACTER(len=*), PARAMETER :: routineN = 'neighbor_list_iterate', 
00351       routineP = moduleN//':'//routineN
00352 
00353     INTEGER                                  :: iab, last, me
00354     TYPE(neighbor_list_iterator_type), 
00355       POINTER                                :: iterator
00356     TYPE(neighbor_list_set_p_type), 
00357       DIMENSION(:), POINTER                  :: nl
00358 
00359     IF(PRESENT(mepos)) THEN
00360       me=mepos
00361     ELSE
00362       me=0
00363     END IF
00364 
00365     istat  = 0
00366 
00367 !$omp critical
00368     last = iterator_set(0)%last
00369     IF(last /= me) THEN
00370       iterator_set(me)%neighbor_list_iterator = iterator_set(last)%neighbor_list_iterator
00371     END IF
00372     iterator => iterator_set(me)%neighbor_list_iterator
00373     nl => iterator%nl
00374 
00375     IF(iterator%inode < iterator%nnode) THEN
00376       ! we can be sure that there is another node in this list
00377       iterator%inode = iterator%inode + 1
00378       iterator%neighbor_node => iterator%neighbor_node%next_neighbor_node
00379     ELSE
00380       iab = MAX(iterator%ikind + iterator%nkind*(iterator%jkind - 1),0)
00381       kindloop : DO ! look for the next list with nnode /= 0
00382          listloop : DO
00383             IF(iterator%ilist >= iterator%nlist) EXIT listloop
00384             iterator%ilist = iterator%ilist + 1
00385             IF(ASSOCIATED(iterator%neighbor_list)) THEN
00386                iterator%neighbor_list => iterator%neighbor_list%next_neighbor_list
00387             ELSE
00388                iterator%neighbor_list => first_list(nl(iab)%neighbor_list_set)
00389             END IF
00390             CALL get_neighbor_list(neighbor_list=iterator%neighbor_list,atom=iterator%iatom,&
00391                                  nnode=iterator%nnode)
00392             IF(iterator%nnode > 0) EXIT kindloop
00393          END DO listloop
00394          IF (iab >= iterator%nkind**2) THEN
00395            istat = 1
00396            EXIT kindloop
00397          ELSE
00398            iab = iab + 1
00399            iterator%jkind = (iab-1)/iterator%nkind + 1
00400            iterator%ikind = iab - iterator%nkind*(iterator%jkind - 1)
00401            iterator%ilist = 0
00402            IF (.NOT.ASSOCIATED(nl(iab)%neighbor_list_set)) THEN
00403              iterator%ilist = 0
00404              iterator%nlist = 0
00405            ELSE
00406              CALL get_neighbor_list_set(neighbor_list_set=&
00407                      nl(iab)%neighbor_list_set,nlist=iterator%nlist)
00408              iterator%ilist = 0
00409            END IF
00410            NULLIFY(iterator%neighbor_list)
00411          END IF
00412       END DO kindloop
00413       iterator%inode = 1
00414       IF (istat==0) iterator%neighbor_node => first_node(iterator%neighbor_list)
00415     END IF
00416     IF (istat==0) THEN
00417       CALL get_neighbor_node(neighbor_node=iterator%neighbor_node,neighbor=iterator%jatom)
00418     END IF
00419 
00420     ! mark the last iterator updated
00421     iterator_set(:)%last = me
00422 !$omp end critical
00423 
00424   END FUNCTION neighbor_list_iterate
00425 
00426   FUNCTION nl_sub_iterate(iterator_set,mepos) RESULT(istat)
00427     TYPE(neighbor_list_iterator_p_type), 
00428       DIMENSION(:), POINTER                  :: iterator_set
00429     INTEGER, OPTIONAL                        :: mepos
00430     INTEGER                                  :: istat
00431 
00432     CHARACTER(len=*), PARAMETER :: routineN = 'nl_sub_iterate', 
00433       routineP = moduleN//':'//routineN
00434 
00435     INTEGER                                  :: last, me
00436     TYPE(neighbor_list_iterator_type), 
00437       POINTER                                :: iterator
00438 
00439     IF(PRESENT(mepos)) THEN
00440       me=mepos
00441     ELSE
00442       me=0
00443     END IF
00444 
00445     istat  = 0
00446 
00447 !$omp critical
00448     last = iterator_set(0)%last
00449     IF(last /= me) THEN
00450       iterator_set(me)%neighbor_list_iterator = iterator_set(last)%neighbor_list_iterator
00451     END IF
00452     iterator => iterator_set(me)%neighbor_list_iterator
00453 
00454     IF (ASSOCIATED(iterator%neighbor_list)) THEN
00455        IF (iterator%inode >= iterator%nnode) THEN
00456           ! end of loop
00457           istat  = 1
00458        ELSEIF (iterator%inode == 0) THEN
00459           iterator%inode = 1
00460           iterator%neighbor_node => first_node(iterator%neighbor_list)
00461        ELSEIF (iterator%inode > 0) THEN
00462           ! we can be sure that there is another node in this list
00463           iterator%inode = iterator%inode + 1
00464           iterator%neighbor_node => iterator%neighbor_node%next_neighbor_node
00465        ELSE
00466           STOP "wrong"
00467        END IF
00468     ELSE
00469        ! no list available
00470        istat  = 1
00471     END IF
00472     IF (istat==0) THEN
00473       CALL get_neighbor_node(neighbor_node=iterator%neighbor_node,neighbor=iterator%jatom)
00474     END IF
00475 
00476     ! mark the last iterator updated
00477     iterator_set(:)%last = me
00478 !$omp end critical
00479 
00480   END FUNCTION nl_sub_iterate
00481 
00482   SUBROUTINE get_iterator_info(iterator_set,mepos,&
00483     ikind,jkind,nkind,ilist,nlist,inode,nnode,iatom,jatom,r,cell)
00484     TYPE(neighbor_list_iterator_p_type), 
00485       DIMENSION(:), POINTER                  :: iterator_set
00486     INTEGER, OPTIONAL                        :: mepos, ikind, jkind, nkind, 
00487                                                 ilist, nlist, inode, nnode, 
00488                                                 iatom, jatom
00489     REAL(dp), DIMENSION(3), OPTIONAL         :: r
00490     INTEGER, DIMENSION(3), OPTIONAL          :: cell
00491 
00492     CHARACTER(len=*), PARAMETER :: routineN = 'get_iterator_info', 
00493       routineP = moduleN//':'//routineN
00494 
00495     INTEGER                                  :: me
00496     TYPE(neighbor_list_iterator_type), 
00497       POINTER                                :: iterator
00498 
00499     IF(PRESENT(mepos)) THEN
00500       me=mepos
00501     ELSE
00502       me=0
00503     END IF
00504     iterator => iterator_set(me)%neighbor_list_iterator
00505 
00506     IF(PRESENT(ikind)) ikind=iterator%ikind
00507     IF(PRESENT(jkind)) jkind=iterator%jkind
00508     IF(PRESENT(nkind)) nkind=iterator%nkind
00509     IF(PRESENT(ilist)) ilist=iterator%ilist
00510     IF(PRESENT(nlist)) nlist=iterator%nlist
00511     IF(PRESENT(inode)) inode=iterator%inode
00512     IF(PRESENT(nnode)) nnode=iterator%nnode
00513     IF(PRESENT(iatom)) iatom=iterator%iatom
00514     IF(PRESENT(jatom)) jatom=iterator%jatom
00515     IF(PRESENT(r)) THEN
00516       CALL get_neighbor_node(neighbor_node=iterator%neighbor_node,r=r)
00517     END IF
00518     IF(PRESENT(cell)) THEN
00519       CALL get_neighbor_node(neighbor_node=iterator%neighbor_node,cell=cell)
00520     END IF
00521 
00522   END SUBROUTINE get_iterator_info
00523 
00524 ! *****************************************************************************
00530   SUBROUTINE add_neighbor_list(neighbor_list_set,atom,neighbor_list)
00531 
00532     TYPE(neighbor_list_set_type), POINTER    :: neighbor_list_set
00533     INTEGER, INTENT(IN)                      :: atom
00534     TYPE(neighbor_list_type), POINTER        :: neighbor_list
00535 
00536     CHARACTER(len=*), PARAMETER :: routineN = 'add_neighbor_list', 
00537       routineP = moduleN//':'//routineN
00538 
00539     INTEGER                                  :: istat
00540     TYPE(neighbor_list_type), POINTER        :: new_neighbor_list
00541 
00542     IF (ASSOCIATED(neighbor_list_set)) THEN
00543 
00544       IF (ASSOCIATED(neighbor_list_set%last_neighbor_list)) THEN
00545 
00546         new_neighbor_list =>&
00547           neighbor_list_set%last_neighbor_list%next_neighbor_list
00548 
00549         IF (.NOT.ASSOCIATED(new_neighbor_list)) THEN
00550 
00551 !         *** Allocate a new neighbor list ***
00552 
00553           ALLOCATE (new_neighbor_list,STAT=istat)
00554           IF (istat /= 0) THEN
00555             CALL stop_memory(routineN,moduleN,__LINE__,&
00556                              "new_neighbor_list",0)
00557           END IF
00558 
00559           NULLIFY (new_neighbor_list%next_neighbor_list)
00560           NULLIFY (new_neighbor_list%first_neighbor_node)
00561 
00562 !         *** Link the new neighbor list to the neighbor list set ***
00563 
00564           neighbor_list_set%last_neighbor_list%next_neighbor_list => new_neighbor_list
00565 
00566         END IF
00567 
00568       ELSE
00569 
00570         new_neighbor_list => neighbor_list_set%first_neighbor_list
00571 
00572         IF (.NOT.ASSOCIATED(new_neighbor_list)) THEN
00573 
00574 !         *** Allocate a new first neighbor list ***
00575 
00576           ALLOCATE (new_neighbor_list,STAT=istat)
00577           IF (istat /= 0) THEN
00578             CALL stop_memory(routineN,moduleN,__LINE__,&
00579                              "new_neighbor_list",0)
00580           END IF
00581 
00582           NULLIFY (new_neighbor_list%next_neighbor_list)
00583           NULLIFY (new_neighbor_list%first_neighbor_node)
00584 
00585 !         *** Link the new first neighbor list to the neighbor list set ***
00586 
00587           neighbor_list_set%first_neighbor_list => new_neighbor_list
00588 
00589         END IF
00590 
00591       END IF
00592 
00593 !     *** Store the data set of the new neighbor list ***
00594 
00595       NULLIFY (new_neighbor_list%last_neighbor_node)
00596       new_neighbor_list%atom = atom
00597       new_neighbor_list%nnode = 0
00598 
00599 !     *** Update the pointer to the last neighbor ***
00600 !     *** list of the neighbor list set           ***
00601 
00602       neighbor_list_set%last_neighbor_list => new_neighbor_list
00603 
00604 !     *** Increment the neighbor list counter ***
00605 
00606       neighbor_list_set%nlist = neighbor_list_set%nlist + 1
00607 
00608 !     *** Return a pointer to the new neighbor list ***
00609 
00610       neighbor_list => new_neighbor_list
00611 
00612     ELSE
00613 
00614       CALL stop_program(routineN,moduleN,__LINE__,&
00615                         "The requested neighbor list set is not associated")
00616 
00617     END IF
00618 
00619   END SUBROUTINE add_neighbor_list
00620 
00621 ! *****************************************************************************
00627   SUBROUTINE add_neighbor_node(neighbor_list,neighbor,cell,r,exclusion_list,nkind)
00628 
00629     TYPE(neighbor_list_type), POINTER        :: neighbor_list
00630     INTEGER, INTENT(IN)                      :: neighbor
00631     INTEGER, DIMENSION(3), INTENT(IN)        :: cell
00632     REAL(dp), DIMENSION(3), INTENT(IN)       :: r
00633     INTEGER, DIMENSION(:), OPTIONAL, POINTER :: exclusion_list
00634     INTEGER, INTENT(IN), OPTIONAL            :: nkind
00635 
00636     CHARACTER(len=*), PARAMETER :: routineN = 'add_neighbor_node', 
00637       routineP = moduleN//':'//routineN
00638 
00639     INTEGER                                  :: iatom, istat, my_nkind
00640     TYPE(neighbor_node_type), POINTER        :: new_neighbor_node
00641 
00642     IF (ASSOCIATED(neighbor_list)) THEN
00643 
00644 !     *** Check for exclusions ***
00645 
00646       IF (PRESENT(exclusion_list)) THEN
00647         IF ( ASSOCIATED ( exclusion_list ) ) THEN
00648           DO iatom=1,SIZE(exclusion_list)
00649             IF (exclusion_list(iatom) == 0) EXIT
00650             IF (exclusion_list(iatom) == neighbor) RETURN
00651           END DO
00652         END IF
00653       END IF
00654 
00655       my_nkind = 0
00656       IF(PRESENT(nkind)) my_nkind = nkind
00657 
00658       IF (ASSOCIATED(neighbor_list%last_neighbor_node)) THEN
00659 
00660         new_neighbor_node => neighbor_list%last_neighbor_node%next_neighbor_node
00661 
00662         IF (.NOT.ASSOCIATED(new_neighbor_node)) THEN
00663 
00664 !         *** Allocate a new neighbor node ***
00665 
00666           ALLOCATE (new_neighbor_node,STAT=istat)
00667           IF (istat /= 0) THEN
00668             CALL stop_memory(routineN,moduleN,__LINE__,&
00669                              "new_neighbor_node",0)
00670           END IF
00671 
00672           NULLIFY (new_neighbor_node%next_neighbor_node)
00673 
00674 !         *** Link the new neighbor node to the neighbor list ***
00675 
00676           neighbor_list%last_neighbor_node%next_neighbor_node => new_neighbor_node
00677 
00678         END IF
00679 
00680       ELSE
00681 
00682         new_neighbor_node => neighbor_list%first_neighbor_node
00683 
00684         IF (.NOT.ASSOCIATED(new_neighbor_node)) THEN
00685 
00686 !         *** Allocate a new first neighbor node ***
00687 
00688           ALLOCATE (new_neighbor_node,STAT=istat)
00689           IF (istat /= 0) THEN
00690             CALL stop_memory(routineN,moduleN,__LINE__,&
00691                              "new_neighbor_node",0)
00692           END IF
00693 
00694           NULLIFY (new_neighbor_node%next_neighbor_node)
00695 
00696 !         *** Link the new first neighbor node to the neighbor list ***
00697 
00698           neighbor_list%first_neighbor_node => new_neighbor_node
00699 
00700         END IF
00701 
00702       END IF
00703 
00704 !     *** Store the data set of the new neighbor ***
00705 
00706       new_neighbor_node%neighbor = neighbor
00707       new_neighbor_node%cell(:) = cell(:)
00708       new_neighbor_node%r(:) = r(:)
00709 
00710 !     *** Update the pointer to the last neighbor node of the neighbor list ***
00711 
00712       neighbor_list%last_neighbor_node => new_neighbor_node
00713 
00714 !     *** Increment the neighbor node counter ***
00715 
00716       neighbor_list%nnode = neighbor_list%nnode + 1
00717 
00718     ELSE
00719 
00720       CALL stop_program(routineN,moduleN,__LINE__,&
00721                         "The requested neighbor list is not associated")
00722 
00723     END IF
00724 
00725   END SUBROUTINE add_neighbor_node
00726 
00727 ! *****************************************************************************
00733   SUBROUTINE allocate_neighbor_list_set(neighbor_list_set,r_max,symmetric,mic,molecular)
00734 
00735     TYPE(neighbor_list_set_type), POINTER    :: neighbor_list_set
00736     REAL(dp), INTENT(IN)                     :: r_max
00737     LOGICAL, INTENT(IN)                      :: symmetric, mic, molecular
00738 
00739     CHARACTER(len=*), PARAMETER :: routineN = 'allocate_neighbor_list_set', 
00740       routineP = moduleN//':'//routineN
00741 
00742     INTEGER                                  :: istat
00743 
00744 !   *** Deallocate the old neighbor list set ***
00745 
00746     IF (ASSOCIATED(neighbor_list_set)) THEN
00747       CALL deallocate_neighbor_list_set(neighbor_list_set)
00748     END IF
00749 
00750 !   *** Allocate a set of neighbor lists ***
00751 
00752     ALLOCATE (neighbor_list_set,STAT=istat)
00753     IF (istat /= 0) THEN
00754       CALL stop_memory(routineN,moduleN,__LINE__,&
00755                        "neighbor_list_set",0)
00756     END IF
00757 
00758     NULLIFY (neighbor_list_set%first_neighbor_list)
00759 
00760 !   *** Initialize the pointers to the first neighbor list ***
00761 
00762     CALL init_neighbor_list_set(neighbor_list_set,r_max,symmetric,mic,molecular)
00763 
00764   END SUBROUTINE allocate_neighbor_list_set
00765 
00766 ! *****************************************************************************
00772   SUBROUTINE deallocate_neighbor_list(neighbor_list)
00773 
00774     TYPE(neighbor_list_type), POINTER        :: neighbor_list
00775 
00776     CHARACTER(len=*), PARAMETER :: routineN = 'deallocate_neighbor_list', 
00777       routineP = moduleN//':'//routineN
00778 
00779     INTEGER                                  :: istat
00780     TYPE(neighbor_node_type), POINTER        :: neighbor_node, 
00781                                                 next_neighbor_node
00782 
00783     IF (ASSOCIATED(neighbor_list)) THEN
00784 
00785       neighbor_node => neighbor_list%first_neighbor_node
00786 
00787       DO WHILE (ASSOCIATED(neighbor_node))
00788         next_neighbor_node => neighbor_node%next_neighbor_node
00789         DEALLOCATE (neighbor_node,STAT=istat)
00790         IF (istat /= 0) THEN
00791           CALL stop_memory(routineN,moduleN,__LINE__,"neighbor_node")
00792         END IF
00793         neighbor_node => next_neighbor_node
00794       END DO
00795 
00796       DEALLOCATE (neighbor_list,STAT=istat)
00797       IF (istat /= 0) THEN
00798         CALL stop_memory(routineN,moduleN,__LINE__,"neighbor_list")
00799       END IF
00800 
00801     END IF
00802 
00803   END SUBROUTINE deallocate_neighbor_list
00804 
00805 ! *****************************************************************************
00811   SUBROUTINE deallocate_neighbor_list_set(neighbor_list_set)
00812     TYPE(neighbor_list_set_type), POINTER    :: neighbor_list_set
00813 
00814     CHARACTER(len=*), PARAMETER :: routineN = 'deallocate_neighbor_list_set', 
00815       routineP = moduleN//':'//routineN
00816 
00817     INTEGER                                  :: istat
00818     TYPE(neighbor_list_type), POINTER        :: neighbor_list, 
00819                                                 next_neighbor_list
00820 
00821     IF (ASSOCIATED(neighbor_list_set)) THEN
00822 
00823       neighbor_list => neighbor_list_set%first_neighbor_list
00824 
00825       DO WHILE (ASSOCIATED(neighbor_list))
00826         next_neighbor_list => neighbor_list%next_neighbor_list
00827         CALL  deallocate_neighbor_list(neighbor_list)
00828         neighbor_list => next_neighbor_list
00829       END DO
00830 
00831       DEALLOCATE (neighbor_list_set,STAT=istat)
00832       IF (istat /= 0) THEN
00833         CALL stop_memory(routineN,moduleN,__LINE__,"neighbor_list_set")
00834       END IF
00835 
00836     END IF
00837 
00838   END SUBROUTINE deallocate_neighbor_list_set
00839 
00840 ! *****************************************************************************
00846   FUNCTION first_list(neighbor_list_set) RESULT(first_neighbor_list)
00847 
00848     TYPE(neighbor_list_set_type), POINTER    :: neighbor_list_set
00849     TYPE(neighbor_list_type), POINTER        :: first_neighbor_list
00850 
00851     first_neighbor_list => neighbor_list_set%first_neighbor_list
00852 
00853   END FUNCTION first_list
00854 
00855 ! *****************************************************************************
00861   FUNCTION first_node(neighbor_list) RESULT(first_neighbor_node)
00862 
00863     TYPE(neighbor_list_type), POINTER        :: neighbor_list
00864     TYPE(neighbor_node_type), POINTER        :: first_neighbor_node
00865 
00866     first_neighbor_node => neighbor_list%first_neighbor_node
00867 
00868   END FUNCTION first_node
00869 
00870 ! *****************************************************************************
00876   SUBROUTINE get_neighbor_list(neighbor_list,atom,nnode)
00877 
00878     TYPE(neighbor_list_type), POINTER        :: neighbor_list
00879     INTEGER, INTENT(OUT), OPTIONAL           :: atom, nnode
00880 
00881     CHARACTER(len=*), PARAMETER :: routineN = 'get_neighbor_list', 
00882       routineP = moduleN//':'//routineN
00883 
00884     IF (ASSOCIATED(neighbor_list)) THEN
00885 
00886       IF (PRESENT(atom)) atom = neighbor_list%atom
00887       IF (PRESENT(nnode)) nnode = neighbor_list%nnode
00888 
00889     ELSE
00890 
00891       CALL stop_program(routineN,moduleN,__LINE__,&
00892                         "The requested neighbor list is not associated")
00893 
00894     END IF
00895 
00896   END SUBROUTINE get_neighbor_list
00897 
00898 ! *****************************************************************************
00904   SUBROUTINE get_neighbor_list_set(neighbor_list_set,r_max,nlist,symmetric,mic,molecular)
00905 
00906     TYPE(neighbor_list_set_type), POINTER    :: neighbor_list_set
00907     REAL(dp), INTENT(OUT), OPTIONAL          :: r_max
00908     INTEGER, INTENT(OUT), OPTIONAL           :: nlist
00909     LOGICAL, INTENT(OUT), OPTIONAL           :: symmetric, mic, molecular
00910 
00911     CHARACTER(len=*), PARAMETER :: routineN = 'get_neighbor_list_set', 
00912       routineP = moduleN//':'//routineN
00913 
00914     IF (ASSOCIATED(neighbor_list_set)) THEN
00915 
00916       IF (PRESENT(r_max)) r_max = neighbor_list_set%r_max
00917       IF (PRESENT(nlist)) nlist = neighbor_list_set%nlist
00918       IF (PRESENT(symmetric)) symmetric = neighbor_list_set%symmetric
00919       IF (PRESENT(mic)) mic = neighbor_list_set%mic
00920       IF (PRESENT(molecular)) molecular = neighbor_list_set%molecular
00921 
00922     ELSE
00923 
00924       CALL stop_program(routineN,moduleN,__LINE__,&
00925                         "The requested neighbor list set is not associated")
00926 
00927     END IF
00928 
00929   END SUBROUTINE get_neighbor_list_set
00930 
00931 ! *****************************************************************************
00937   SUBROUTINE get_neighbor_node(neighbor_node,neighbor,cell,r)
00938 
00939     TYPE(neighbor_node_type), POINTER        :: neighbor_node
00940     INTEGER, INTENT(OUT), OPTIONAL           :: neighbor
00941     INTEGER, DIMENSION(3), INTENT(OUT), 
00942       OPTIONAL                               :: cell
00943     REAL(dp), DIMENSION(3), INTENT(OUT), 
00944       OPTIONAL                               :: r
00945 
00946     CHARACTER(len=*), PARAMETER :: routineN = 'get_neighbor_node', 
00947       routineP = moduleN//':'//routineN
00948 
00949     IF (ASSOCIATED(neighbor_node)) THEN
00950 
00951       IF (PRESENT(neighbor)) neighbor = neighbor_node%neighbor
00952       IF (PRESENT(r)) r(:) = neighbor_node%r(:)
00953       IF (PRESENT(cell)) cell(:) = neighbor_node%cell(:)
00954 
00955     ELSE
00956 
00957       CALL stop_program(routineN,moduleN,__LINE__,&
00958                         "The requested neighbor node is not associated")
00959 
00960     END IF
00961 
00962   END SUBROUTINE get_neighbor_node
00963 
00964 ! *****************************************************************************
00972   SUBROUTINE init_neighbor_list_set(neighbor_list_set,r_max,symmetric,mic,molecular)
00973 
00974     TYPE(neighbor_list_set_type), POINTER    :: neighbor_list_set
00975     REAL(dp), INTENT(IN)                     :: r_max
00976     LOGICAL, INTENT(IN)                      :: symmetric, mic, molecular
00977 
00978     CHARACTER(len=*), PARAMETER :: routineN = 'init_neighbor_list_set', 
00979       routineP = moduleN//':'//routineN
00980 
00981     IF (ASSOCIATED(neighbor_list_set)) THEN
00982 
00983       ! *** Initialize the pointers to the last neighbor list ***
00984       NULLIFY (neighbor_list_set%last_neighbor_list)
00985 
00986       ! *** Initialize the neighbor list counter ***
00987       neighbor_list_set%nlist = 0
00988 
00989       ! *** Initialize the maximum interaction radius (optional) ***
00990       neighbor_list_set%r_max = r_max
00991 
00992       ! *** Initialize the neighbor list build properties
00993       neighbor_list_set%symmetric = symmetric
00994       neighbor_list_set%mic       = mic
00995       neighbor_list_set%molecular = molecular
00996 
00997     ELSE
00998 
00999       CALL stop_program(routineN,moduleN,__LINE__,&
01000                         "The requested neighbor list set is not associated")
01001 
01002     END IF
01003 
01004   END SUBROUTINE init_neighbor_list_set
01005 
01006 END MODULE qs_neighbor_list_types