Skip to content

Commit 91dbeea

Browse files
Enable hyphen separated range of atoms to modify mass (#1179)
* enable hyphen separated range of atoms to modify mass Signed-off-by: Marcel Müller <[email protected]> * Update src/constrain_param.f90 Co-authored-by: Marvin Friede <[email protected]> * cycle instead of return for single syntax error in comma-separated configuration Signed-off-by: Marcel Müller <[email protected]> * rename ati to iat Signed-off-by: Marcel Müller <[email protected]> * forgot to change from return to cycle on both positions Signed-off-by: Marcel Müller <[email protected]> --------- Signed-off-by: Marcel Müller <[email protected]> Co-authored-by: Marvin Friede <[email protected]>
1 parent 5f7a2e2 commit 91dbeea

File tree

1 file changed

+41
-20
lines changed

1 file changed

+41
-20
lines changed

src/constrain_param.f90

Lines changed: 41 additions & 20 deletions
Original file line numberDiff line numberDiff line change
@@ -1227,6 +1227,7 @@ subroutine set_split(env,key,val,nat,at,idMap,xyz)
12271227
end subroutine set_split
12281228

12291229
subroutine set_hess(env,key,val,nat,at,idMap,xyz)
1230+
use xtb_type_atomlist, only : TAtomList
12301231
use xtb_splitparam
12311232
implicit none
12321233
character(len=*), parameter :: source = 'userdata_hess'
@@ -1238,14 +1239,12 @@ subroutine set_hess(env,key,val,nat,at,idMap,xyz)
12381239
type(TIdentityMap), intent(in) :: idMap
12391240
real(wp),intent(in) :: xyz(3,nat)
12401241

1241-
integer :: idum
1242-
real(wp) :: ddum
1243-
logical :: ldum
1244-
integer :: i,j
1242+
type(TAtomList) :: atl
12451243
integer, allocatable :: list(:)
1246-
1247-
integer :: narg
1244+
real(wp) :: ddum
1245+
integer :: i,j,idum,iat,narg
12481246
character(len=p_str_length),dimension(p_arg_length) :: argv
1247+
character(len=256) :: warningstring
12491248

12501249
call parse(val,comma,argv,narg)
12511250
if (set%verbose) then
@@ -1279,15 +1278,26 @@ subroutine set_hess(env,key,val,nat,at,idMap,xyz)
12791278
endif
12801279
do i = 1, narg, 2
12811280
j = i+1
1282-
if (getValue(env,trim(argv(i)),idum).and.&
1283-
getValue(env,trim(argv(j)),ddum)) then
1284-
if (idum.gt.nat) then
1285-
call env%warning('Attempted setting atom mass not present in system.',source)
1281+
if (getValue(env,trim(argv(j)),ddum)) then
1282+
call atl%new(argv(i))
1283+
if (atl%get_error()) then
1284+
call env%warning('something is wrong in the mass list',source)
12861285
cycle
12871286
endif
1288-
atmass(idum) = ddum
1289-
write(env%unit,'(a,1x,i0,1x,a,1x,g0)') &
1290-
'mass of atom ',idum,' changed to',atmass(idum)
1287+
call atl%to_list(list)
1288+
do idum = 1, size(list)
1289+
iat = list(idum)
1290+
if (iat.gt.nat) then
1291+
write(warningstring, '(a, i0, a)') 'Attempted setting atom mass for atom ', &
1292+
& iat, ' that is not present in system.'
1293+
call env%warning(trim(warningstring), source)
1294+
cycle
1295+
endif
1296+
atmass(iat) = ddum
1297+
write(env%unit,'(a,1x,i0,1x,a,1x,g0)') &
1298+
& 'mass of atom ',iat,' changed to',atmass(iat)
1299+
enddo
1300+
call atl%destroy()
12911301
endif
12921302
enddo
12931303
case('scale mass')
@@ -1296,15 +1306,26 @@ subroutine set_hess(env,key,val,nat,at,idMap,xyz)
12961306
endif
12971307
do i = 1, narg, 2
12981308
j = i+1
1299-
if (getValue(env,trim(argv(i)),idum).and.&
1300-
getValue(env,trim(argv(j)),ddum)) then
1301-
if (idum.gt.nat) then
1302-
call env%warning('Attempted scaling atom not present in system.',source)
1309+
if (getValue(env,trim(argv(j)),ddum)) then
1310+
call atl%new(argv(i))
1311+
if (atl%get_error()) then
1312+
call env%warning('something is wrong in the mass list',source)
13031313
cycle
13041314
endif
1305-
atmass(idum) = atmass(idum)*ddum
1306-
write(env%unit,'(a,1x,i0,1x,a,1x,g0)') &
1307-
'mass of atom ',idum,' changed to',atmass(idum)
1315+
call atl%to_list(list)
1316+
do idum = 1, size(list)
1317+
iat = list(idum)
1318+
if (iat.gt.nat) then
1319+
write(warningstring, '(a, i0, a)') 'Attempted setting atom mass for atom ', &
1320+
& iat, ' that is not present in system.'
1321+
call env%warning(trim(warningstring), source)
1322+
cycle
1323+
endif
1324+
atmass(iat) = atmass(iat)*ddum
1325+
write(env%unit,'(a,1x,i0,1x,a,1x,g0)') &
1326+
'mass of atom ',iat,' changed to',atmass(iat)
1327+
enddo
1328+
call atl%destroy()
13081329
endif
13091330
enddo
13101331
end select

0 commit comments

Comments
 (0)