Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Allocate set%ffnb dynamically. #1105

Merged
merged 2 commits into from
Sep 19, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
8 changes: 5 additions & 3 deletions src/gfnff/gfnff_ini.f90
Original file line number Diff line number Diff line change
Expand Up @@ -2450,10 +2450,12 @@ subroutine gfnff_topo_changes(env, neigh)
integer :: int_tmp(40)
integer :: i,j,idx,iTr,d1,d2,numnb,nnb

! check if hardcoded size of ffnb is still up to date
if (size(set%ffnb, dim=1).ne.neigh%numnb) call env%error('The array set%ffnb has not been adjusted to changes in neigh%numnb.', source)
! only do something if there are changes stored in set%ffnb
if(set%ffnb(1,1).ne.-1) then
if(allocated(set%ffnb)) then
! check if hardcoded size of ffnb is still up to date
if (size(set%ffnb, dim=1).ne.neigh%numnb) call env%error('The array set%ffnb has not been adjusted to changes in neigh%numnb.', source)
! there should not be any "-1" in set%ffnb anymore, if it was set up correctly
if (any(set%ffnb.eq.-1)) call env%error('GFN-FF neighbor list could not be adjusted!', source)
Thomas3R marked this conversation as resolved.
Show resolved Hide resolved
d2=size(set%ffnb, dim=2)
do i=1, d2
if (set%ffnb(1,i).eq.-1) exit
Expand Down
46 changes: 45 additions & 1 deletion src/set_module.f90
Original file line number Diff line number Diff line change
Expand Up @@ -781,7 +781,10 @@ subroutine rdcontrol(fname,env,copy_file)
case('cube' ); call rdblock(env,set_cube, line,id,copy,err,ncount)
case('write' ); call rdblock(env,set_write, line,id,copy,err,ncount)
case('gfn' ); call rdblock(env,set_gfn, line,id,copy,err,ncount)
case('ffnb' ); call rdblock(env,set_ffnb, line,id,copy,err,ncount)
case('ffnb' )
! dynamic allocation of ffnb array requires reading fname before calling rdblock
call alloc_ffnb(env, fname)
call rdblock(env,set_ffnb, line,id,copy,err,ncount)
case('scc' ); call rdblock(env,set_scc, line,id,copy,err,ncount)
case('oniom' ); call rdblock(env,set_oniom, line,id,copy,err,ncount)
case('opt' ); call rdblock(env,set_opt, line,id,copy,err,ncount)
Expand Down Expand Up @@ -1500,6 +1503,47 @@ subroutine set_gfn(env,key,val)
end select
end subroutine set_gfn

! determine number of GFN-FF neighbor list changes in control file
! and allocate set%ffnb accordingly
subroutine alloc_ffnb(env, fname)
type(TEnvironment), intent(inout) :: env
character(len=*),intent(in) :: fname
character(len=*), parameter :: source = 'alloc_ffnb'
character(len=:),allocatable :: line

integer :: id, ie, err
logical :: is_ffnb_block = .false.
integer :: copy = -1
integer :: n_changes = 0 ! number of atoms that neigh%nb should be adjusted for

call open_file(id,fname,'r')
if (id.eq.-1) then
call env%warning("could not find '"//fname//"'", source)
return
endif
! read first line
call mirror_line(id,copy,line,err)
! read control file and check
count_n:do
! check if the $ffnb block has been reached
if (line(1:5).eq."$ffnb".or.is_ffnb_block) then
is_ffnb_block = .true.
call mirror_line(id,copy,line,err)
if (is_iostat_end(err)) exit count_n ! check if EOF !
if (index(line,flag).ne.0) exit count_n ! check if new flag !
ie = index(line,equal) ! find the equal sign !
if (ie.eq.0) cycle ! cycle if there is no equal sign
n_changes = n_changes + 1
else ! otherwise read the next line
call mirror_line(id,copy,line,err)
if (is_iostat_end(err)) exit count_n ! check if EOF !
end if
end do count_n
if (.not.allocated(set%ffnb)) allocate(set%ffnb(42,n_changes), source=-1)

end subroutine alloc_ffnb


subroutine set_ffnb(env,key,val)
implicit none
character(len=*), parameter :: source = 'set_ffnb'
Expand Down
3 changes: 1 addition & 2 deletions src/setparam.f90
Original file line number Diff line number Diff line change
Expand Up @@ -517,9 +517,8 @@ module xtb_setparam
!> PTB settings
type(TPTBSetup) :: ptbsetup
!> GFN-FF manual setup of nb list via xcontrol
! allows a maximum of 164 atoms neighbors to be changed
! ffnb(42,i) stores the number of neighbors of atom i
integer :: ffnb(42,164) = -1
integer, allocatable :: ffnb(:,:)
end type TSet

type(TSet) :: set
Expand Down