LCOV - code coverage report
Current view: top level - src - fist_nonbond_env_types.F (source / functions) Coverage Total Hit
Test: CP2K Regtests (git:5064cfc) Lines: 96.3 % 218 210
Test Date: 2026-03-04 06:45:10 Functions: 41.7 % 12 5

            Line data    Source code
       1              : !--------------------------------------------------------------------------------------------------!
       2              : !   CP2K: A general program to perform molecular dynamics simulations                              !
       3              : !   Copyright 2000-2026 CP2K developers group <https://cp2k.org>                                   !
       4              : !                                                                                                  !
       5              : !   SPDX-License-Identifier: GPL-2.0-or-later                                                      !
       6              : !--------------------------------------------------------------------------------------------------!
       7              : 
       8              : ! **************************************************************************************************
       9              : !> \par History
      10              : !>      none
      11              : !> \author HAF
      12              : ! **************************************************************************************************
      13              : MODULE fist_nonbond_env_types
      14              :    USE ace_wrapper,                     ONLY: ace_model_release,&
      15              :                                               ace_model_type
      16              :    USE atomic_kind_types,               ONLY: atomic_kind_type
      17              :    USE cell_types,                      ONLY: cell_release,&
      18              :                                               cell_type
      19              :    USE deepmd_wrapper,                  ONLY: deepmd_model_release,&
      20              :                                               deepmd_model_type
      21              :    USE fist_neighbor_list_types,        ONLY: fist_neighbor_deallocate,&
      22              :                                               fist_neighbor_type
      23              :    USE kinds,                           ONLY: default_string_length,&
      24              :                                               dp
      25              :    USE pair_potential_types,            ONLY: &
      26              :         ace_type, allegro_type, gal21_type, gal_type, nequip_type, pair_potential_pp_release, &
      27              :         pair_potential_pp_type, siepmann_type, tersoff_type
      28              :    USE torch_api,                       ONLY: torch_model_release,&
      29              :                                               torch_model_type
      30              : #include "./base/base_uses.f90"
      31              : 
      32              :    IMPLICIT NONE
      33              :    PRIVATE
      34              : 
      35              :    CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'fist_nonbond_env_types'
      36              :    PUBLIC :: fist_nonbond_env_type, fist_nonbond_env_set, &
      37              :              fist_nonbond_env_get, fist_nonbond_env_create, &
      38              :              fist_nonbond_env_release, pos_type, eam_type, &
      39              :              nequip_data_type, &
      40              :              deepmd_data_type, ace_data_type
      41              : 
      42              : ! **************************************************************************************************
      43              :    TYPE pos_type
      44              :       REAL(KIND=dp) :: r(3) = 0.0_dp
      45              :    END TYPE pos_type
      46              : 
      47              :    TYPE eam_type
      48              :       REAL(KIND=dp) :: f_embed = 0.0_dp
      49              :       REAL(KIND=dp) :: rho = 0.0_dp
      50              :    END TYPE eam_type
      51              : 
      52              :    TYPE nequip_data_type
      53              :       INTEGER, POINTER        :: use_indices(:) => NULL()
      54              :       REAL(KIND=dp), POINTER  :: force(:, :) => NULL()
      55              :       REAL(KIND=dp)           :: virial(3, 3) = 0.0_dp
      56              :       TYPE(torch_model_type)  :: model
      57              :    END TYPE nequip_data_type
      58              : 
      59              :    TYPE deepmd_data_type
      60              :       INTEGER, POINTER        :: use_indices(:) => NULL()
      61              :       REAL(KIND=dp), POINTER  :: force(:, :) => NULL()
      62              :       REAL(KIND=dp)           :: virial(3, 3) = 0.0_dp
      63              :       TYPE(deepmd_model_type) :: model
      64              :    END TYPE deepmd_data_type
      65              : 
      66              :    TYPE ace_data_type
      67              :       INTEGER, ALLOCATABLE    :: use_indices(:)
      68              :       INTEGER, ALLOCATABLE    :: inverse_index_map(:)
      69              :       INTEGER                 :: natom = 0
      70              :       INTEGER                 :: nghost = 0
      71              :       INTEGER                 :: refupdate = 0
      72              :       INTEGER                 :: nei = 0
      73              :       INTEGER, ALLOCATABLE    :: uctype(:)
      74              :       INTEGER, ALLOCATABLE    :: attype(:)
      75              :       INTEGER, ALLOCATABLE    :: origin(:)
      76              :       INTEGER, ALLOCATABLE    :: shift(:, :)
      77              :       INTEGER, ALLOCATABLE    :: neiat(:)
      78              :       INTEGER, ALLOCATABLE    :: nlist(:)
      79              :       REAL(KIND=dp), ALLOCATABLE  :: force(:, :)
      80              :       REAL(KIND=dp), ALLOCATABLE  :: atpos(:, :)
      81              :       REAL(KIND=dp)           :: virial(3, 3) = 0.0_dp
      82              :       TYPE(ace_model_type)    :: model
      83              :    END TYPE ace_data_type
      84              : 
      85              : ! **************************************************************************************************
      86              :    TYPE fist_nonbond_env_type
      87              :       INTEGER                                    :: natom_types = -1
      88              :       INTEGER                                    :: counter = -1
      89              :       INTEGER                                    :: last_update = -1
      90              :       INTEGER                                    :: num_update = -1
      91              :       LOGICAL                                    :: do_nonbonded = .FALSE.
      92              :       LOGICAL                                    :: do_electrostatics = .FALSE.
      93              :       LOGICAL                                    :: shift_cutoff = .FALSE.
      94              :       CHARACTER(len=default_string_length)       :: unit_type = ""
      95              :       REAL(KIND=dp)                              :: lup = 0.0_dp
      96              :       REAL(KIND=dp)                              :: aup = 0.0_dp
      97              :       REAL(KIND=dp)                              :: ei_scale14 = 0.0_dp
      98              :       REAL(KIND=dp)                              :: vdw_scale14 = 0.0_dp
      99              :       REAL(KIND=dp)                              :: long_range_correction = 0.0_dp
     100              :       REAL(KIND=dp), DIMENSION(:, :), POINTER    :: rlist_cut => NULL()
     101              :       REAL(KIND=dp), DIMENSION(:, :), POINTER    :: rlist_lowsq => NULL()
     102              :       REAL(KIND=dp), DIMENSION(:, :), POINTER    :: ij_kind_full_fac => NULL()
     103              :       REAL(KIND=dp), DIMENSION(:), POINTER       :: charges => NULL()
     104              :       TYPE(fist_neighbor_type), POINTER          :: nonbonded => NULL()
     105              :       TYPE(pair_potential_pp_type), POINTER      :: potparm14 => NULL()
     106              :       TYPE(pair_potential_pp_type), POINTER      :: potparm => NULL()
     107              :       TYPE(cell_type), POINTER                   :: cell_last_update => NULL()
     108              :       TYPE(pos_type), DIMENSION(:), POINTER      :: r_last_update => NULL()
     109              :       TYPE(pos_type), DIMENSION(:), POINTER      :: r_last_update_pbc => NULL()
     110              :       TYPE(pos_type), DIMENSION(:), POINTER      :: rshell_last_update_pbc => NULL()
     111              :       TYPE(pos_type), DIMENSION(:), POINTER      :: rcore_last_update_pbc => NULL()
     112              :       TYPE(eam_type), DIMENSION(:), POINTER      :: eam_data => NULL()
     113              :       TYPE(deepmd_data_type), POINTER            :: deepmd_data => NULL()
     114              :       TYPE(ace_data_type), POINTER               :: ace_data => NULL()
     115              :       TYPE(nequip_data_type), POINTER    :: nequip_data => NULL()
     116              :    END TYPE fist_nonbond_env_type
     117              : 
     118              : CONTAINS
     119              : 
     120              : ! **************************************************************************************************
     121              : !> \brief sets a fist_nonbond_env
     122              : !> \param fist_nonbond_env the object to create
     123              : !> \param potparm14 ...
     124              : !> \param potparm ...
     125              : !> \param nonbonded ...
     126              : !> \param rlist_cut ...
     127              : !> \param rlist_lowsq ...
     128              : !> \param aup ...
     129              : !> \param lup ...
     130              : !> \param ei_scale14 ...
     131              : !> \param vdw_scale14 ...
     132              : !> \param shift_cutoff ...
     133              : !> \param do_electrostatics ...
     134              : !> \param r_last_update ...
     135              : !> \param r_last_update_pbc ...
     136              : !> \param rshell_last_update_pbc ...
     137              : !> \param rcore_last_update_pbc ...
     138              : !> \param cell_last_update ...
     139              : !> \param num_update ...
     140              : !> \param last_update ...
     141              : !> \param counter ...
     142              : !> \param natom_types ...
     143              : !> \param long_range_correction ...
     144              : !> \param ij_kind_full_fac ...
     145              : !> \param eam_data ...
     146              : !> \param nequip_data ...
     147              : !> \param deepmd_data ...
     148              : !> \param ace_data ...
     149              : !> \param charges ...
     150              : !> \par History
     151              : !>      12.2002 created [fawzi]
     152              : !> \author Fawzi Mohamed
     153              : ! **************************************************************************************************
     154       551147 :    SUBROUTINE fist_nonbond_env_get(fist_nonbond_env, potparm14, potparm, &
     155              :                                    nonbonded, rlist_cut, rlist_lowsq, aup, lup, ei_scale14, vdw_scale14, &
     156              :                                    shift_cutoff, do_electrostatics, r_last_update, r_last_update_pbc, rshell_last_update_pbc, &
     157              :                                    rcore_last_update_pbc, cell_last_update, num_update, last_update, &
     158              :                                    counter, natom_types, long_range_correction, ij_kind_full_fac, eam_data, &
     159              :                                    nequip_data, deepmd_data, ace_data, charges)
     160              : 
     161              :       TYPE(fist_nonbond_env_type), INTENT(IN)            :: fist_nonbond_env
     162              :       TYPE(pair_potential_pp_type), OPTIONAL, POINTER    :: potparm14, potparm
     163              :       TYPE(fist_neighbor_type), OPTIONAL, POINTER        :: nonbonded
     164              :       REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER  :: rlist_cut, rlist_lowsq
     165              :       REAL(KIND=dp), OPTIONAL                            :: aup, lup, ei_scale14, vdw_scale14
     166              :       LOGICAL, INTENT(OUT), OPTIONAL                     :: shift_cutoff, do_electrostatics
     167              :       TYPE(pos_type), DIMENSION(:), OPTIONAL, POINTER    :: r_last_update, r_last_update_pbc, &
     168              :                                                             rshell_last_update_pbc, &
     169              :                                                             rcore_last_update_pbc
     170              :       TYPE(cell_type), OPTIONAL, POINTER                 :: cell_last_update
     171              :       INTEGER, OPTIONAL                                  :: num_update, last_update, counter, &
     172              :                                                             natom_types
     173              :       REAL(KIND=dp), OPTIONAL                            :: long_range_correction
     174              :       REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER  :: ij_kind_full_fac
     175              :       TYPE(eam_type), DIMENSION(:), OPTIONAL, POINTER    :: eam_data
     176              :       TYPE(nequip_data_type), OPTIONAL, POINTER          :: nequip_data
     177              :       TYPE(deepmd_data_type), OPTIONAL, POINTER          :: deepmd_data
     178              :       TYPE(ace_data_type), OPTIONAL, POINTER             :: ace_data
     179              :       REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER     :: charges
     180              : 
     181       551147 :       IF (PRESENT(charges)) charges => fist_nonbond_env%charges
     182       551147 :       IF (PRESENT(potparm14)) potparm14 => fist_nonbond_env%potparm14
     183       551147 :       IF (PRESENT(eam_data)) eam_data => fist_nonbond_env%eam_data
     184       551147 :       IF (PRESENT(nequip_data)) nequip_data => fist_nonbond_env%nequip_data
     185       551147 :       IF (PRESENT(deepmd_data)) deepmd_data => fist_nonbond_env%deepmd_data
     186       551147 :       IF (PRESENT(ace_data)) ace_data => fist_nonbond_env%ace_data
     187       551147 :       IF (PRESENT(potparm)) potparm => fist_nonbond_env%potparm
     188       551147 :       IF (PRESENT(rlist_cut)) rlist_cut => fist_nonbond_env%rlist_cut
     189       551147 :       IF (PRESENT(rlist_lowsq)) rlist_lowsq => fist_nonbond_env%rlist_lowsq
     190       551147 :       IF (PRESENT(ij_kind_full_fac)) ij_kind_full_fac => fist_nonbond_env%ij_kind_full_fac
     191       551147 :       IF (PRESENT(nonbonded)) nonbonded => fist_nonbond_env%nonbonded
     192       551147 :       IF (PRESENT(r_last_update)) &
     193       251396 :          r_last_update => fist_nonbond_env%r_last_update
     194       551147 :       IF (PRESENT(r_last_update_pbc)) &
     195       406166 :          r_last_update_pbc => fist_nonbond_env%r_last_update_pbc
     196       551147 :       IF (PRESENT(rshell_last_update_pbc)) &
     197       165736 :          rshell_last_update_pbc => fist_nonbond_env%rshell_last_update_pbc
     198       551147 :       IF (PRESENT(rcore_last_update_pbc)) &
     199       165736 :          rcore_last_update_pbc => fist_nonbond_env%rcore_last_update_pbc
     200       551147 :       IF (PRESENT(cell_last_update)) &
     201        83764 :          cell_last_update => fist_nonbond_env%cell_last_update
     202       551147 :       IF (PRESENT(lup)) lup = fist_nonbond_env%lup
     203       551147 :       IF (PRESENT(aup)) aup = fist_nonbond_env%aup
     204       551147 :       IF (PRESENT(ei_scale14)) ei_scale14 = fist_nonbond_env%ei_scale14
     205       551147 :       IF (PRESENT(vdw_scale14)) vdw_scale14 = fist_nonbond_env%vdw_scale14
     206       551147 :       IF (PRESENT(shift_cutoff)) &
     207            0 :          shift_cutoff = fist_nonbond_env%shift_cutoff
     208       551147 :       IF (PRESENT(do_electrostatics)) do_electrostatics = fist_nonbond_env%do_electrostatics
     209       551147 :       IF (PRESENT(natom_types)) natom_types = fist_nonbond_env%natom_types
     210       551147 :       IF (PRESENT(counter)) counter = fist_nonbond_env%counter
     211       551147 :       IF (PRESENT(last_update)) last_update = fist_nonbond_env%last_update
     212       551147 :       IF (PRESENT(num_update)) num_update = fist_nonbond_env%num_update
     213       551147 :       IF (PRESENT(long_range_correction)) &
     214            0 :          long_range_correction = fist_nonbond_env%long_range_correction
     215       551147 :    END SUBROUTINE fist_nonbond_env_get
     216              : 
     217              : ! **************************************************************************************************
     218              : !> \brief sets a fist_nonbond_env
     219              : !> \param fist_nonbond_env the object to create
     220              : !> \param potparm14 ...
     221              : !> \param potparm ...
     222              : !> \param rlist_cut ...
     223              : !> \param rlist_lowsq ...
     224              : !> \param nonbonded ...
     225              : !> \param aup ...
     226              : !> \param lup ...
     227              : !> \param ei_scale14 ...
     228              : !> \param vdw_scale14 ...
     229              : !> \param shift_cutoff ...
     230              : !> \param do_electrostatics ...
     231              : !> \param r_last_update ...
     232              : !> \param r_last_update_pbc ...
     233              : !> \param rshell_last_update_pbc ...
     234              : !> \param rcore_last_update_pbc ...
     235              : !> \param cell_last_update ...
     236              : !> \param num_update ...
     237              : !> \param last_update ...
     238              : !> \param counter ...
     239              : !> \param natom_types ...
     240              : !> \param long_range_correction ...
     241              : !> \param eam_data ...
     242              : !> \param nequip_data ...
     243              : !> \param deepmd_data ...
     244              : !> \param ace_data ...
     245              : !> \param charges ...
     246              : !> \par History
     247              : !>      12.2002 created [fawzi]
     248              : !> \author Fawzi Mohamed
     249              : ! **************************************************************************************************
     250        97567 :    SUBROUTINE fist_nonbond_env_set(fist_nonbond_env, potparm14, potparm, &
     251              :                                    rlist_cut, rlist_lowsq, nonbonded, aup, lup, ei_scale14, vdw_scale14, &
     252              :                                    shift_cutoff, do_electrostatics, r_last_update, r_last_update_pbc, rshell_last_update_pbc, &
     253              :                                    rcore_last_update_pbc, cell_last_update, num_update, last_update, &
     254              :                                    counter, natom_types, long_range_correction, eam_data, &
     255              :                                    nequip_data, deepmd_data, ace_data, charges)
     256              : 
     257              :       TYPE(fist_nonbond_env_type), INTENT(INOUT)         :: fist_nonbond_env
     258              :       TYPE(pair_potential_pp_type), OPTIONAL, POINTER    :: potparm14, potparm
     259              :       REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER  :: rlist_cut, rlist_lowsq
     260              :       TYPE(fist_neighbor_type), OPTIONAL, POINTER        :: nonbonded
     261              :       REAL(KIND=dp), OPTIONAL                            :: aup, lup, ei_scale14, vdw_scale14
     262              :       LOGICAL, INTENT(IN), OPTIONAL                      :: shift_cutoff, do_electrostatics
     263              :       TYPE(pos_type), DIMENSION(:), OPTIONAL, POINTER    :: r_last_update, r_last_update_pbc, &
     264              :                                                             rshell_last_update_pbc, &
     265              :                                                             rcore_last_update_pbc
     266              :       TYPE(cell_type), OPTIONAL, POINTER                 :: cell_last_update
     267              :       INTEGER, OPTIONAL                                  :: num_update, last_update, counter, &
     268              :                                                             natom_types
     269              :       REAL(KIND=dp), OPTIONAL                            :: long_range_correction
     270              :       TYPE(eam_type), DIMENSION(:), OPTIONAL, POINTER    :: eam_data
     271              :       TYPE(nequip_data_type), OPTIONAL, POINTER          :: nequip_data
     272              :       TYPE(deepmd_data_type), OPTIONAL, POINTER          :: deepmd_data
     273              :       TYPE(ace_data_type), OPTIONAL, POINTER             :: ace_data
     274              :       REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER     :: charges
     275              : 
     276        97567 :       IF (PRESENT(potparm14)) fist_nonbond_env%potparm14 => potparm14
     277        97567 :       IF (PRESENT(eam_data)) fist_nonbond_env%eam_data => eam_data
     278        97567 :       IF (PRESENT(nequip_data)) fist_nonbond_env%nequip_data => nequip_data
     279        97567 :       IF (PRESENT(deepmd_data)) fist_nonbond_env%deepmd_data => deepmd_data
     280        97567 :       IF (PRESENT(ace_data)) fist_nonbond_env%ace_data => ace_data
     281        97567 :       IF (PRESENT(potparm)) fist_nonbond_env%potparm => potparm
     282        97567 :       IF (PRESENT(rlist_cut)) fist_nonbond_env%rlist_cut => rlist_cut
     283        97567 :       IF (PRESENT(charges)) fist_nonbond_env%charges => charges
     284        97567 :       IF (PRESENT(rlist_lowsq)) fist_nonbond_env%rlist_lowsq => rlist_lowsq
     285        97567 :       IF (PRESENT(nonbonded)) fist_nonbond_env%nonbonded => nonbonded
     286        97567 :       IF (PRESENT(r_last_update)) &
     287        11138 :          fist_nonbond_env%r_last_update => r_last_update
     288        97567 :       IF (PRESENT(r_last_update_pbc)) &
     289        11138 :          fist_nonbond_env%r_last_update_pbc => r_last_update_pbc
     290        97567 :       IF (PRESENT(rshell_last_update_pbc)) &
     291        11138 :          fist_nonbond_env%rshell_last_update_pbc => rshell_last_update_pbc
     292        97567 :       IF (PRESENT(rcore_last_update_pbc)) &
     293        11138 :          fist_nonbond_env%rcore_last_update_pbc => rcore_last_update_pbc
     294        97567 :       IF (PRESENT(cell_last_update)) &
     295        11138 :          fist_nonbond_env%cell_last_update => cell_last_update
     296        97567 :       IF (PRESENT(lup)) fist_nonbond_env%lup = lup
     297        97567 :       IF (PRESENT(aup)) fist_nonbond_env%aup = aup
     298        97567 :       IF (PRESENT(ei_scale14)) fist_nonbond_env%ei_scale14 = ei_scale14
     299        97567 :       IF (PRESENT(vdw_scale14)) fist_nonbond_env%vdw_scale14 = vdw_scale14
     300        97567 :       IF (PRESENT(shift_cutoff)) &
     301            0 :          fist_nonbond_env%shift_cutoff = shift_cutoff
     302        97567 :       IF (PRESENT(do_electrostatics)) fist_nonbond_env%do_electrostatics = do_electrostatics
     303        97567 :       IF (PRESENT(natom_types)) fist_nonbond_env%natom_types = natom_types
     304        97567 :       IF (PRESENT(counter)) fist_nonbond_env%counter = counter
     305        97567 :       IF (PRESENT(last_update)) fist_nonbond_env%last_update = last_update
     306        97567 :       IF (PRESENT(num_update)) fist_nonbond_env%num_update = num_update
     307        97567 :       IF (PRESENT(long_range_correction)) &
     308            0 :          fist_nonbond_env%long_range_correction = long_range_correction
     309        97567 :    END SUBROUTINE fist_nonbond_env_set
     310              : 
     311              : ! **************************************************************************************************
     312              : !> \brief allocates and intitializes a fist_nonbond_env
     313              : !> \param fist_nonbond_env the object to create
     314              : !> \param atomic_kind_set ...
     315              : !> \param potparm14 ...
     316              : !> \param potparm ...
     317              : !> \param do_nonbonded ...
     318              : !> \param do_electrostatics ...
     319              : !> \param verlet_skin ...
     320              : !> \param ewald_rcut ...
     321              : !> \param ei_scale14 ...
     322              : !> \param vdw_scale14 ...
     323              : !> \param shift_cutoff ...
     324              : !> \par History
     325              : !>      12.2002 created [fawzi]
     326              : !> \author Fawzi Mohamed
     327              : ! **************************************************************************************************
     328         2673 :    SUBROUTINE fist_nonbond_env_create(fist_nonbond_env, atomic_kind_set, &
     329              :                                       potparm14, potparm, do_nonbonded, do_electrostatics, verlet_skin, ewald_rcut, &
     330              :                                       ei_scale14, vdw_scale14, shift_cutoff)
     331              :       TYPE(fist_nonbond_env_type), INTENT(OUT)           :: fist_nonbond_env
     332              :       TYPE(atomic_kind_type), DIMENSION(:), POINTER      :: atomic_kind_set
     333              :       TYPE(pair_potential_pp_type), OPTIONAL, POINTER    :: potparm14, potparm
     334              :       LOGICAL, INTENT(IN)                                :: do_nonbonded, do_electrostatics
     335              :       REAL(KIND=dp), INTENT(IN)                          :: verlet_skin, ewald_rcut, ei_scale14, &
     336              :                                                             vdw_scale14
     337              :       LOGICAL, INTENT(IN)                                :: shift_cutoff
     338              : 
     339              :       NULLIFY (fist_nonbond_env%potparm14)
     340              :       NULLIFY (fist_nonbond_env%potparm)
     341              :       NULLIFY (fist_nonbond_env%rlist_cut)
     342              :       NULLIFY (fist_nonbond_env%rlist_lowsq)
     343              :       NULLIFY (fist_nonbond_env%ij_kind_full_fac)
     344              :       NULLIFY (fist_nonbond_env%nonbonded)
     345              :       NULLIFY (fist_nonbond_env%cell_last_update)
     346              :       NULLIFY (fist_nonbond_env%r_last_update)
     347              :       NULLIFY (fist_nonbond_env%r_last_update_pbc)
     348              :       NULLIFY (fist_nonbond_env%rshell_last_update_pbc)
     349              :       NULLIFY (fist_nonbond_env%rcore_last_update_pbc)
     350              :       NULLIFY (fist_nonbond_env%eam_data)
     351              :       NULLIFY (fist_nonbond_env%nequip_data)
     352              :       NULLIFY (fist_nonbond_env%deepmd_data)
     353              :       NULLIFY (fist_nonbond_env%ace_data)
     354              :       NULLIFY (fist_nonbond_env%charges)
     355              :       CALL init_fist_nonbond_env(fist_nonbond_env, atomic_kind_set, potparm14, &
     356              :                                  potparm, do_nonbonded, do_electrostatics, verlet_skin, ewald_rcut, ei_scale14, &
     357         2673 :                                  vdw_scale14, shift_cutoff)
     358         2673 :    END SUBROUTINE fist_nonbond_env_create
     359              : 
     360              : ! **************************************************************************************************
     361              : !> \brief Purpose: Initialise the FIST nonbond environment.
     362              : !> \param fist_nonbond_env the object to create
     363              : !> \param atomic_kind_set ...
     364              : !> \param potparm14 ...
     365              : !> \param potparm ...
     366              : !> \param do_nonbonded ...
     367              : !> \param do_electrostatics ...
     368              : !> \param verlet_skin ...
     369              : !> \param ewald_rcut ...
     370              : !> \param ei_scale14 ...
     371              : !> \param vdw_scale14 ...
     372              : !> \param shift_cutoff ...
     373              : ! **************************************************************************************************
     374         2673 :    SUBROUTINE init_fist_nonbond_env(fist_nonbond_env, atomic_kind_set, &
     375              :                                     potparm14, potparm, do_nonbonded, do_electrostatics, verlet_skin, ewald_rcut, ei_scale14, &
     376              :                                     vdw_scale14, shift_cutoff)
     377              : 
     378              :       TYPE(fist_nonbond_env_type), INTENT(INOUT)         :: fist_nonbond_env
     379              :       TYPE(atomic_kind_type), DIMENSION(:), POINTER      :: atomic_kind_set
     380              :       TYPE(pair_potential_pp_type), OPTIONAL, POINTER    :: potparm14, potparm
     381              :       LOGICAL, INTENT(IN)                                :: do_nonbonded, do_electrostatics
     382              :       REAL(KIND=dp), INTENT(IN)                          :: verlet_skin, ewald_rcut, ei_scale14, &
     383              :                                                             vdw_scale14
     384              :       LOGICAL, INTENT(IN)                                :: shift_cutoff
     385              : 
     386              :       INTEGER                                            :: idim, jdim, natom_types
     387              :       LOGICAL                                            :: check, use_potparm, use_potparm14
     388              :       REAL(KIND=dp)                                      :: fac, rcut, rlow
     389              : 
     390         2673 :       use_potparm14 = PRESENT(potparm14)
     391         2673 :       IF (use_potparm14) use_potparm14 = use_potparm14 .OR. ASSOCIATED(potparm14)
     392         2673 :       use_potparm = PRESENT(potparm)
     393         2673 :       IF (use_potparm) use_potparm = use_potparm .OR. ASSOCIATED(potparm)
     394         2673 :       NULLIFY (fist_nonbond_env%nonbonded)
     395         2673 :       NULLIFY (fist_nonbond_env%r_last_update)
     396         2673 :       NULLIFY (fist_nonbond_env%r_last_update_pbc)
     397         2673 :       NULLIFY (fist_nonbond_env%rshell_last_update_pbc)
     398         2673 :       NULLIFY (fist_nonbond_env%rcore_last_update_pbc)
     399         2673 :       NULLIFY (fist_nonbond_env%cell_last_update)
     400         2673 :       NULLIFY (fist_nonbond_env%rlist_cut)
     401         2673 :       NULLIFY (fist_nonbond_env%rlist_lowsq)
     402         2673 :       NULLIFY (fist_nonbond_env%ij_kind_full_fac)
     403         2673 :       fist_nonbond_env%unit_type = "ANGSTROM"
     404         2673 :       fist_nonbond_env%do_nonbonded = do_nonbonded
     405         2673 :       fist_nonbond_env%do_electrostatics = do_electrostatics
     406         2673 :       fist_nonbond_env%lup = 0
     407         2673 :       fist_nonbond_env%aup = 0
     408         2673 :       fist_nonbond_env%ei_scale14 = ei_scale14
     409         2673 :       fist_nonbond_env%vdw_scale14 = vdw_scale14
     410         2673 :       fist_nonbond_env%shift_cutoff = shift_cutoff
     411         2673 :       fist_nonbond_env%counter = 0
     412         2673 :       fist_nonbond_env%last_update = 0
     413         2673 :       fist_nonbond_env%num_update = 0
     414         2673 :       fist_nonbond_env%long_range_correction = 0
     415         2673 :       IF (do_nonbonded) THEN
     416         2657 :          natom_types = 1
     417              :          ! Determine size of kind arrays
     418         2657 :          natom_types = SIZE(atomic_kind_set)
     419         2657 :          IF (use_potparm14) THEN
     420         2625 :             check = (SIZE(potparm14%pot, 1) == natom_types)
     421         2625 :             CPASSERT(check)
     422              :          END IF
     423         2657 :          IF (use_potparm) THEN
     424         2625 :             check = (SIZE(potparm%pot, 1) == natom_types)
     425         2625 :             CPASSERT(check)
     426              :          END IF
     427        10628 :          ALLOCATE (fist_nonbond_env%rlist_cut(natom_types, natom_types))
     428         7971 :          ALLOCATE (fist_nonbond_env%rlist_lowsq(natom_types, natom_types))
     429         7971 :          ALLOCATE (fist_nonbond_env%ij_kind_full_fac(natom_types, natom_types))
     430       518241 :          fist_nonbond_env%ij_kind_full_fac = 1.0_dp
     431        13916 :          DO idim = 1, natom_types
     432       271708 :             DO jdim = idim, natom_types
     433       269051 :                IF ((use_potparm) .OR. (use_potparm14)) THEN
     434       257712 :                   IF (use_potparm) THEN
     435       257712 :                      rcut = SQRT(potparm%pot(idim, jdim)%pot%rcutsq)
     436       257712 :                      fac = potparm%pot(idim, jdim)%pot%spl_f%rscale(1)
     437       257712 :                      rlow = fac/(potparm%pot(idim, jdim)%pot%pair_spline_data(1)%spline_data%xn)
     438              :                   ELSE
     439            0 :                      rcut = SQRT(potparm14%pot(idim, jdim)%pot%rcutsq)
     440            0 :                      fac = potparm14%pot(idim, jdim)%pot%spl_f%rscale(1)
     441            0 :                      rlow = fac/(potparm14%pot(idim, jdim)%pot%pair_spline_data(1)%spline_data%xn)
     442              :                   END IF
     443              :                   ! Warning: rlist_rcut should only be used by the neighbor list
     444              :                   ! algorithm. It is not the cutoff for the evaluation of the
     445              :                   ! interactions because rlist_rcut includes the Verlet skin.
     446       257712 :                   rcut = MAX(rcut, ewald_rcut) + verlet_skin
     447       257712 :                   fist_nonbond_env%rlist_cut(idim, jdim) = rcut
     448       257712 :                   fist_nonbond_env%rlist_cut(jdim, idim) = rcut
     449       257712 :                   rlow = rlow*(1.06_dp)**2 ! 1.06_dp in order to have 1/2 Emax_spline
     450       257712 :                   fist_nonbond_env%rlist_lowsq(idim, jdim) = rlow
     451       257712 :                   fist_nonbond_env%rlist_lowsq(jdim, idim) = rlow
     452              :                   ! In case of manybody potential the neighbor list will be full.
     453              :                   ! This means that for each atom pair (a,b) of the current types,
     454              :                   ! atom a is in the neighbor list of b and b is in the neighbor
     455              :                   ! list of a. ij_kind_full_fac is used to correct for the double
     456              :                   ! counting in the conventional pair potentials cause by this
     457              :                   ! situation.
     458       515314 :                   IF (ANY(potparm%pot(idim, jdim)%pot%type == tersoff_type)) THEN
     459              :                      ! TODO: what if 14 is not of tersoff type while the normal
     460              :                      ! nonbond is? (or the reverse). We'd better impose
     461              :                      ! consistency.
     462          118 :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     463              :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     464              :                   END IF
     465       515427 :                   IF (ANY(potparm%pot(idim, jdim)%pot%type == siepmann_type)) THEN
     466              :                      ! TODO:see tersoff_type
     467            5 :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     468              :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     469              :                   END IF
     470       515414 :                   IF (ANY(potparm%pot(idim, jdim)%pot%type == ace_type)) THEN
     471           18 :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     472           18 :                      fist_nonbond_env%ij_kind_full_fac(jdim, idim) = 0.5_dp
     473              :                   END IF
     474       515431 :                   IF (ANY(potparm%pot(idim, jdim)%pot%type == gal_type)) THEN
     475            1 :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     476              :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     477              :                   END IF
     478       515431 :                   IF (ANY(potparm%pot(idim, jdim)%pot%type == gal21_type)) THEN
     479            1 :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     480              :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     481              :                   END IF
     482       515426 :                   IF (ANY(potparm%pot(idim, jdim)%pot%type == nequip_type)) THEN
     483            6 :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     484              :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     485              :                   END IF
     486       515426 :                   IF (ANY(potparm%pot(idim, jdim)%pot%type == allegro_type)) THEN
     487            6 :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     488              :                      fist_nonbond_env%ij_kind_full_fac(idim, jdim) = 0.5_dp
     489              :                   END IF
     490              :                ELSE
     491              :                   ! In case we don't use potparm for initialization let's account
     492              :                   ! only for the real-space part of the Ewald sum.
     493           80 :                   fist_nonbond_env%rlist_cut(idim, jdim) = ewald_rcut
     494           80 :                   fist_nonbond_env%rlist_cut(jdim, idim) = ewald_rcut
     495           80 :                   fist_nonbond_env%rlist_lowsq(idim, jdim) = 0.0_dp
     496           80 :                   fist_nonbond_env%rlist_lowsq(jdim, idim) = 0.0_dp
     497              :                END IF
     498              :             END DO
     499              :          END DO
     500         2657 :          IF (use_potparm14) fist_nonbond_env%potparm14 => potparm14
     501         2657 :          IF (use_potparm) fist_nonbond_env%potparm => potparm
     502         2657 :          fist_nonbond_env%natom_types = natom_types
     503              :       ELSE
     504           16 :          NULLIFY (fist_nonbond_env%potparm)
     505           16 :          NULLIFY (fist_nonbond_env%potparm14)
     506              :       END IF
     507         2673 :    END SUBROUTINE init_fist_nonbond_env
     508              : 
     509              : ! **************************************************************************************************
     510              : !> \brief releases the given fist_nonbond_env (see doc/ReferenceCounting.html)
     511              : !> \param fist_nonbond_env the object to release
     512              : !> \par History
     513              : !>      12.2002 created [fawzi]
     514              : !> \author Fawzi Mohamed
     515              : ! **************************************************************************************************
     516         2673 :    SUBROUTINE fist_nonbond_env_release(fist_nonbond_env)
     517              :       TYPE(fist_nonbond_env_type), INTENT(INOUT)         :: fist_nonbond_env
     518              : 
     519         2673 :       IF (ASSOCIATED(fist_nonbond_env%nonbonded)) THEN
     520         2503 :          CALL fist_neighbor_deallocate(fist_nonbond_env%nonbonded)
     521              :       END IF
     522              :       ! Release potparm
     523         2673 :       CALL pair_potential_pp_release(fist_nonbond_env%potparm)
     524              :       ! Release potparm14
     525         2673 :       CALL pair_potential_pp_release(fist_nonbond_env%potparm14)
     526         2673 :       IF (ASSOCIATED(fist_nonbond_env%r_last_update)) THEN
     527         2503 :          DEALLOCATE (fist_nonbond_env%r_last_update)
     528              :       END IF
     529         2673 :       IF (ASSOCIATED(fist_nonbond_env%r_last_update_pbc)) THEN
     530         2503 :          DEALLOCATE (fist_nonbond_env%r_last_update_pbc)
     531              :       END IF
     532         2673 :       IF (ASSOCIATED(fist_nonbond_env%charges)) THEN
     533            8 :          DEALLOCATE (fist_nonbond_env%charges)
     534              :       END IF
     535         2673 :       IF (ASSOCIATED(fist_nonbond_env%eam_data)) THEN
     536           12 :          DEALLOCATE (fist_nonbond_env%eam_data)
     537              :       END IF
     538         2673 :       IF (ASSOCIATED(fist_nonbond_env%nequip_data)) THEN
     539            4 :          IF (ASSOCIATED(fist_nonbond_env%nequip_data%force)) THEN
     540            4 :             DEALLOCATE (fist_nonbond_env%nequip_data%force)
     541              :          END IF
     542            4 :          IF (ASSOCIATED(fist_nonbond_env%nequip_data%use_indices)) THEN
     543            4 :             DEALLOCATE (fist_nonbond_env%nequip_data%use_indices)
     544              :          END IF
     545            4 :          CALL torch_model_release(fist_nonbond_env%nequip_data%model)
     546            4 :          DEALLOCATE (fist_nonbond_env%nequip_data)
     547              :       END IF
     548         2673 :       IF (ASSOCIATED(fist_nonbond_env%deepmd_data)) THEN
     549            2 :          IF (ASSOCIATED(fist_nonbond_env%deepmd_data%force)) THEN
     550            2 :             DEALLOCATE (fist_nonbond_env%deepmd_data%force)
     551              :          END IF
     552            2 :          IF (ASSOCIATED(fist_nonbond_env%deepmd_data%use_indices)) THEN
     553            2 :             DEALLOCATE (fist_nonbond_env%deepmd_data%use_indices)
     554              :          END IF
     555            2 :          CALL deepmd_model_release(fist_nonbond_env%deepmd_data%model)
     556            2 :          DEALLOCATE (fist_nonbond_env%deepmd_data)
     557              :       END IF
     558         2673 :       IF (ASSOCIATED(fist_nonbond_env%ace_data)) THEN
     559            6 :          IF (ALLOCATED(fist_nonbond_env%ace_data%use_indices)) THEN
     560            6 :             DEALLOCATE (fist_nonbond_env%ace_data%use_indices)
     561              :          END IF
     562            6 :          IF (ALLOCATED(fist_nonbond_env%ace_data%inverse_index_map)) THEN
     563            6 :             DEALLOCATE (fist_nonbond_env%ace_data%inverse_index_map)
     564              :          END IF
     565            6 :          IF (ALLOCATED(fist_nonbond_env%ace_data%force)) THEN
     566            6 :             DEALLOCATE (fist_nonbond_env%ace_data%force)
     567              :          END IF
     568            6 :          IF (ALLOCATED(fist_nonbond_env%ace_data%atpos)) THEN
     569            6 :             DEALLOCATE (fist_nonbond_env%ace_data%atpos)
     570              :          END IF
     571            6 :          IF (ALLOCATED(fist_nonbond_env%ace_data%uctype)) THEN
     572            6 :             DEALLOCATE (fist_nonbond_env%ace_data%uctype)
     573              :          END IF
     574            6 :          IF (ALLOCATED(fist_nonbond_env%ace_data%attype)) THEN
     575            6 :             DEALLOCATE (fist_nonbond_env%ace_data%attype)
     576              :          END IF
     577            6 :          IF (ALLOCATED(fist_nonbond_env%ace_data%origin)) THEN
     578            6 :             DEALLOCATE (fist_nonbond_env%ace_data%origin)
     579              :          END IF
     580            6 :          IF (ALLOCATED(fist_nonbond_env%ace_data%shift)) THEN
     581            6 :             DEALLOCATE (fist_nonbond_env%ace_data%shift)
     582              :          END IF
     583            6 :          IF (ALLOCATED(fist_nonbond_env%ace_data%neiat)) THEN
     584            6 :             DEALLOCATE (fist_nonbond_env%ace_data%neiat)
     585              :          END IF
     586            6 :          IF (ALLOCATED(fist_nonbond_env%ace_data%nlist)) THEN
     587            6 :             DEALLOCATE (fist_nonbond_env%ace_data%nlist)
     588              :          END IF
     589            6 :          CALL ace_model_release(fist_nonbond_env%ace_data%model)
     590            6 :          DEALLOCATE (fist_nonbond_env%ace_data)
     591              :       END IF
     592         2673 :       IF (ASSOCIATED(fist_nonbond_env%rshell_last_update_pbc)) THEN
     593          242 :          DEALLOCATE (fist_nonbond_env%rshell_last_update_pbc)
     594              :       END IF
     595         2673 :       IF (ASSOCIATED(fist_nonbond_env%rcore_last_update_pbc)) THEN
     596          242 :          DEALLOCATE (fist_nonbond_env%rcore_last_update_pbc)
     597              :       END IF
     598         2673 :       IF (ASSOCIATED(fist_nonbond_env%cell_last_update)) THEN
     599         2503 :          CALL cell_release(fist_nonbond_env%cell_last_update)
     600              :       END IF
     601         2673 :       IF (ASSOCIATED(fist_nonbond_env%ij_kind_full_fac)) THEN
     602         2657 :          DEALLOCATE (fist_nonbond_env%ij_kind_full_fac)
     603              :       END IF
     604         2673 :       IF (ASSOCIATED(fist_nonbond_env%rlist_cut)) THEN
     605         2657 :          DEALLOCATE (fist_nonbond_env%rlist_cut)
     606              :       END IF
     607         2673 :       IF (ASSOCIATED(fist_nonbond_env%rlist_lowsq)) THEN
     608         2657 :          DEALLOCATE (fist_nonbond_env%rlist_lowsq)
     609              :       END IF
     610         2673 :    END SUBROUTINE fist_nonbond_env_release
     611              : 
     612            0 : END MODULE fist_nonbond_env_types
        

Generated by: LCOV version 2.0-1