Skip to content

Commit fc1f226

Browse files
committed
correct removal of inertia during MD for default group
1 parent 536ed22 commit fc1f226

File tree

2 files changed

+28
-31
lines changed

2 files changed

+28
-31
lines changed

source/cluster.f

Lines changed: 22 additions & 25 deletions
Original file line numberDiff line numberDiff line change
@@ -154,35 +154,32 @@ subroutine cluster
154154
c
155155
c pack atoms of each group into a contiguous indexed list
156156
c
157-
if (use_group) then
158-
do i = 1, n
159-
list(i) = grplist(i)
160-
end do
161-
call sort3 (n,list,kgrp)
157+
do i = 1, n
158+
list(i) = grplist(i)
159+
end do
160+
call sort3 (n,list,kgrp)
162161
c
163162
c find the first and last atom in each of the groups
164163
c
165-
k = list(1)
166-
igrp(1,k) = 1
167-
do i = 1, n
168-
j = list(i)
169-
if (j .ne. k) then
170-
igrp(2,k) = i - 1
171-
igrp(1,j) = i
172-
k = j
173-
end if
174-
ngrp = max(j,ngrp)
175-
end do
176-
igrp(2,j) = n
164+
k = list(1)
165+
igrp(1,k) = 1
166+
do i = 1, n
167+
j = list(i)
168+
if (j .ne. k) then
169+
igrp(2,k) = i - 1
170+
igrp(1,j) = i
171+
k = j
172+
end if
173+
ngrp = max(j,ngrp)
174+
end do
175+
igrp(2,j) = n
177176
c
178177
c sort the list of atoms in each group by atom number
179178
c
180-
do i = 0, ngrp
181-
size = igrp(2,i) - igrp(1,i) + 1
182-
if (igrp(1,i) .ne. 0)
183-
& call sort (size,kgrp(igrp(1,i)))
184-
end do
185-
end if
179+
do i = 0, ngrp
180+
size = igrp(2,i) - igrp(1,i) + 1
181+
if (igrp(1,i) .ne. 0) call sort (size,kgrp(igrp(1,i)))
182+
end do
186183
c
187184
c perform deallocation of some local arrays
188185
c
@@ -229,7 +226,7 @@ subroutine cluster
229226
c
230227
c compute the total mass of all atoms in each group
231228
c
232-
do i = 1, ngrp
229+
do i = 0, ngrp
233230
grpmass(i) = 0.0d0
234231
do j = igrp(1,i), igrp(2,i)
235232
grpmass(i) = grpmass(i) + mass(kgrp(j))
@@ -239,7 +236,7 @@ subroutine cluster
239236
c output the final list of atoms in each group
240237
c
241238
if (use_group .and. debug) then
242-
do i = 1, ngrp
239+
do i = 0, ngrp
243240
size = igrp(2,i) - igrp(1,i) + 1
244241
if (size .ne. 0) then
245242
write (iout,50) i

source/image.f

Lines changed: 6 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -12,23 +12,23 @@
1212
c ################################################################
1313
c
1414
c
15-
c "image" takes the components of pairwise distance between
16-
c two points in a periodic box and converts to the components
17-
c of the minimum image distance
15+
c "image" takes the components of pairwise distance between two
16+
c points in a periodic box and converts to the components of the
17+
c minimum image distance
1818
c
1919
c literature reference:
2020
c
2121
c U. K. Deiters, "Efficient Coding of the Minimum Image Convention",
2222
c Zeitschrift fur Physikalische Chemie, 227, 345-352 (2013)
2323
c
24-
c note the "do while" clause below can be written using the "nint"
25-
c intrinsic, and the two forms give equivalent values:
24+
c note the "do while" clauses below can be written using the "nint"
25+
c intrinsic, and both forms give the same value; for example:
2626
c
2727
c do while (abs(xr) .gt. xbox2)
2828
c xr = xr - sign(xbox,xr) vs. xr = xr - xbox*nint(xr/xbox)
2929
c end do
3030
c
31-
c which one is faster depends upon specific machine and compiler
31+
c which form is faster depends upon specific machine and compiler
3232
c combinations, and other implementations are also possible
3333
c
3434
c

0 commit comments

Comments
 (0)