dcdfort
Modern Fortran library for analyzing DCD trajectories
All Classes Namespaces Files Functions Variables Pages
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...
 
pure integer(kind=int32) 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()

pure integer(kind=int32) function dcdfort_index::indexfile_get ( class(indexfile), intent(in)  this,
character (len=*), intent(in)  group_name,
integer(kind=int32), 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
228 
229  implicit none
230  integer(kind=int32) :: indexfile_get
231  class(IndexFile), intent(in) :: this
232  character (len=*), intent(in) :: group_name
233  integer(kind=int32), intent(in), optional :: I
234  integer(kind=int32) :: J
235 
236  do j = 1, size(this%group)
237 
238  if (trim(this%group(j)%title) .eq. trim(group_name)) then
239 
240  if (present(i)) then
241  indexfile_get = this%group(j)%LOC(i)
242  else
243  indexfile_get = this%group(j)%NUMATOMS
244  end if
245  return
246 
247  end if
248 
249  end do
250 
251  ! If user asked to get the number of atoms in an index group, and that index group is not
252  ! in the index file, just return 0. If the user specified an atom number, then throw an error,
253  ! since the overall index cannot be returned in that case
254  if (.not. present(i)) then
255  indexfile_get = 0
256  else
257  indexfile_get = -1
258  call error_stop_program(trim(group_name)//" is not in the index file.")
259  end if
260 

◆ indexfile_read()

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