dcdfort
Modern Fortran library for analyzing DCD trajectories
Data Types | Functions/Subroutines
dcdfort_index Module Reference

Module that contains IndexFile type. More...

Data Types

type  indexfile
 IndexFile class. More...
 
type  ndxgroups
 ndxgroups class More...
 

Functions/Subroutines

subroutine indexfile_read (this, filename, N)
 Open, read in, and process the GROMACS-style index file, storing information in memory. More...
 
integer function indexfile_get (this, group_name, I)
 Gets the number of atoms in a group. More...
 

Detailed Description

Module that contains IndexFile type.

This module contains the IndexFile class which can be used to read in and process GROMACS-style index files. It should not be necessary to create an object from this class independently from a Trajectory object, since it is included in the Trajectory class.

Function/Subroutine Documentation

◆ indexfile_get()

integer function dcdfort_index::indexfile_get ( class(indexfile), intent(inout)  this,
character (len=*), intent(in)  group_name,
integer, intent(in), optional  I 
)
private

Gets the number of atoms in a group.

Parameters
[in,out]thisIndexFile class
[in]group_nameindex group title or name (in brackets in the index file)
[in]Ithe location in the group
Returns
If an atom is specified, integer returns the overall index for that atom; otherwise, returns number of atoms in group
222 
223  implicit none
224  integer :: indexfile_get
225  class(indexfile), intent(inout) :: this
226  character (len=*), intent(in) :: group_name
227  integer, intent(in), optional :: i
228  integer :: j
229 
230  do j = 1, size(this%group)
231 
232  if (trim(this%group(j)%title) .eq. trim(group_name)) then
233 
234  indexfile_get = merge(this%group(j)%LOC(i), this%group(j)%NUMATOMS, present(i))
235  return
236 
237  end if
238 
239  end do
240 
241  ! If user asked to get the number of atoms in an index group, and that index group is not
242  ! in the index file, just return 0. If the user specified an atom number, then throw an error,
243  ! since the overall index cannot be returned in that case
244  if (.not. present(i)) then
245  if (this%group_warning) then
246  write(error_unit, '(a)') "LIBDCDFORT WARNING: No atoms found in index group '"//trim(group_name)//"'."
247  write(error_unit, '(a)') "This warning will not appear again for any other index groups."
248  this%group_warning = .false.
249  end if
250  indexfile_get = 0
251  else
252  indexfile_get = -1
253  write(error_unit, '(a)') "LIBDCDFORT ERROR: "//trim(group_name)//" is not in index file. The groups available are:"
254  do j = 1, size(this%group)
255  write(error_unit,'(a10,a,i0,a)') this%group(j)%title, " (", this%group(j)%NUMATOMS, ")"
256  end do
257  stop 1
258  end if
259 

◆ indexfile_read()

subroutine dcdfort_index::indexfile_read ( class(indexfile), intent(inout)  this,
character (len=*), intent(in)  filename,
integer, intent(in), optional  N 
)

Open, read in, and process the GROMACS-style index file, storing information in memory.

Parameters
[in,out]thisIndexFile class
[in]filenamename of GROMACS-style index file
[in]Nnumber of atoms in the entire system; used as a check
66 
67  implicit none
68  class(indexfile), intent(inout) :: this
69  character (len=*), intent(in) :: filename
70  character (len=2048) :: line, ncols_string, fmt_string
71  integer :: index_file_unit, io_status, ngrps, i, j, ncols
72  integer, allocatable :: indices_tmp(:), title_loc(:), num_array(:)
73  logical :: ex
74  integer, intent(in), optional :: n
75 
76  ! Does the file exist?
77  inquire(file=trim(filename), exist=ex)
78  if (ex .eqv. .false.) call error_stop_program(trim(filename)//" does not exist.")
79 
80  ! Is in index file?
81  open(newunit=index_file_unit, file=trim(filename), status="old")
82  read(index_file_unit, '(a)', iostat=io_status) line
83  if (index(line, "[") .eq. 0) call error_stop_program(trim(filename)//" is not a valid index file.")
84 
85  ! How many groups are in it?
86  rewind index_file_unit
87  io_status = 0
88  ngrps = 0
89  do while (io_status .eq. 0)
90  read(index_file_unit, '(a)', iostat=io_status) line
91  if (io_status .ne. 0) goto 100
92  if (index(line, "[") .ne. 0) ngrps = ngrps + 1
93  end do
94 
95 100 continue
96 
97  if (allocated(this%group)) deallocate(this%group)
98  allocate(this%group(ngrps), title_loc(ngrps+1)) ! Add one to include end of file
99 
100  ! Now find the title locations and save their names
101  rewind index_file_unit
102  i = 1
103  j = 1
104  io_status = 0
105  do while (io_status .eq. 0)
106 
107  read(index_file_unit, '(a)', iostat=io_status) line
108  if (io_status .ne. 0) goto 200
109  if (index(line, "[") .ne. 0) then
110  this%group(i)%title = trim(line(index(line, "[")+2:index(line, "]")-2))
111  title_loc(i) = j
112  i = i + 1
113  end if
114  j = j + 1
115 
116  end do
117 
118 200 continue
119 
120  if (this%group(1)%title .ne. "System") call error_stop_program("Index file does not have 'System' group as first group.")
121 
122  ! Index files can have a varying number of columns. This attempts to
123  ! detect the correct number by reading in the second line of the file,
124  ! which should be a list of indices for the "System" group.
125  ncols = 100
126  allocate(num_array(ncols))
127  io_status = 5000
128  do while (io_status .ne. 0)
129  ncols = ncols - 1
130  write(ncols_string, '(i0)') ncols
131  write(fmt_string, '(a)') '('//trim(ncols_string)//'i0)'
132  rewind index_file_unit
133  read(index_file_unit, '(a)', iostat=io_status) line
134  read(index_file_unit, '(a)', iostat=io_status) line
135  read(line, *, iostat=io_status) num_array(1:ncols)
136  end do
137 
138  title_loc(i) = j ! End of file location
139 
140  ! Now finally get all of the indices for each group
141  ! Allocate for total number of atoms in system, since that is the maximum
142  allocate(indices_tmp(n))
143  do i = 1, ngrps
144 
145  ! Initial guess only how many items are in the group
146  ! Add 1, bc loop subtracts 1 at the beginning
147  this%group(i)%NUMATOMS = (title_loc(i+1)-title_loc(i)-1)*ncols + 1
148 
149  if (n < this%group(i)%NUMATOMS) this%group(i)%NUMATOMS = n + 1
150  io_status = 5000
151 
152  do while (io_status .ne. 0)
153 
154  ! Our guess was too large if we made it back here, go to the beginning and reduce our guess by 1, try again
155  rewind index_file_unit
156  this%group(i)%NUMATOMS = this%group(i)%NUMATOMS - 1
157  if (this%group(i)%NUMATOMS .le. 0) then
158  this%group(i)%NUMATOMS = 0
159  goto 300
160  end if
161 
162  ! Read all the way to the group
163  do j = 1, title_loc(i); read(index_file_unit, '(a)', iostat=io_status) line; end do
164 
165  ! Attempt to read into array
166  read(index_file_unit, *, iostat=io_status) indices_tmp(1:this%group(i)%NUMATOMS)
167 
168  end do
169 
170  ! Specifying array bounds for array to be allocated is not required for F2008 but is required for F2003
171  allocate(this%group(i)%LOC(1:this%group(i)%NUMATOMS), source=indices_tmp(1:this%group(i)%NUMATOMS))
172 
173 300 cycle
174 
175  end do
176  deallocate(indices_tmp)
177 
178  do i = 1, ngrps-1
179  do j = i+1, ngrps
180  if (i .ne. j) then
181  if (this%group(i)%title .eq. this%group(j)%title) then
182  write(error_unit,*)
183  write(error_unit,'(a, a, a)') "LIBDCDFORT WARNING: Index group ", this%group(i)%title, &
184  " was specified more than once in index file."
185  write(error_unit,*)
186  end if
187  end if
188  end do
189  end do
190 
191  ! Error checking to see if index file goes with the trajectory file
192 
193  ! If the number of atoms is not the same in the System group (first group) and dcd file
194  if (this%group(1)%numatoms .ne. n .or. this%group(1)%loc(this%group(1)%numatoms) .ne. n) then
195  call error_stop_program("Index file does not match dcd file.")
196  end if
197 
198  if (present(n)) then
199  do i = 1, ngrps
200 
201  ! If number of atoms in index group is larger than number of atoms in dcd file
202  if (this%group(i)%numatoms .gt. n) call error_stop_program("Index file does not match dcd file.")
203 
204  ! If a location number is greater than number of atoms in dcd file
205  do j = 1, this%group(i)%numatoms
206  if (this%group(i)%loc(j) .gt. n) call error_stop_program("Index file does not match dcd file.")
207  end do
208 
209  end do
210  end if
211 
212  close(index_file_unit)
213