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

Module that contains Trajectory type. More...

Data Types

type  trajectory
 Trajectory class. More...
 

Functions/Subroutines

subroutine trajectory_open (this, filename_in, ndxfile)
 Trajectory class method which opens DCD file and optionally index file. More...
 
integer function trajectory_get_natoms (this, group)
 Trajectory class method which gets the number of atoms in the system or an index group. More...
 
integer function trajectory_skip_next (this, F)
 Trajectory class method which skips a specified number of frames. More...
 
integer function trajectory_read_next (this, F, ndxgrp)
 Trajectory class method which reads a specified number of frames into memory after using the open() method. More...
 
real function, dimension(3) trajectory_get_xyz (this, frame, atom, group)
 Trajectory class method which returns the coordinates of a particle. More...
 
subroutine trajectory_read (this, dcdfile, ndxfile, ndxgrp)
 Trajectory class method which opens, reads, and closes a trajectory file. More...
 
real(8) function, dimension(6) trajectory_get_box (this, frame)
 Trajectory class method which returns the box from a simulation frame. More...
 
subroutine trajectory_close (this)
 Trajectory class method which closes a DCD file which was opened with open() More...
 

Detailed Description

Module that contains Trajectory type.

Function/Subroutine Documentation

◆ trajectory_close()

subroutine dcdfort_trajectory::trajectory_close ( class(trajectory), intent(inout)  this)
private

Trajectory class method which closes a DCD file which was opened with open()

Parameters
[in,out]thisTrajectory class
342 
343  implicit none
344  class(trajectory), intent(inout) :: this
345 
346  call this%dcd%close()
347 

◆ trajectory_get_box()

real(8) function, dimension(6) dcdfort_trajectory::trajectory_get_box ( class(trajectory), intent(in)  this,
integer, intent(in)  frame 
)
private

Trajectory class method which returns the box from a simulation frame.

Parameters
[in,out]thisTrajectory class
[in]frameWhich frame to return the box dimensions
Returns
A real(8) array with 6 elements containing the x, y, and z dimensions of the box as well as the alpha, beta, and gamma angles
313 
314  implicit none
315  real(8) :: trajectory_get_box(6)
316  class(trajectory), intent(in) :: this
317  integer, intent(in) :: frame
318 
319  call trajectory_check_frame(this, frame)
320  trajectory_get_box = this%frameArray(frame)%box
321 

◆ trajectory_get_natoms()

integer function dcdfort_trajectory::trajectory_get_natoms ( class(trajectory), intent(inout)  this,
character (len=*), intent(in), optional  group 
)
private

Trajectory class method which gets the number of atoms in the system or an index group.

Parameters
[in,out]thisthe Trajectory object
[in]groupname of index group
Returns
number of atoms in index group (if specified) or in system
122 
123  implicit none
124  integer :: trajectory_get_natoms
125  class(trajectory), intent(inout) :: this
126  character (len=*), intent(in), optional :: group
127 
128  if (this%read_only_index_group .and. present(group)) then
129  call error_stop_program("Do not specify an index group in natoms() when already specifying an &
130  &index group with read() or read_next().")
131  end if
132 
133  trajectory_get_natoms = merge(this%ndx%get_natoms(group), this%NUMATOMS, present(group))
134 

◆ trajectory_get_xyz()

real function, dimension(3) dcdfort_trajectory::trajectory_get_xyz ( class(trajectory), intent(inout)  this,
integer, intent(in)  frame,
integer, intent(in)  atom,
character (len=*), intent(in), optional  group 
)
private

Trajectory class method which returns the coordinates of a particle.

Gets the coordinates of a particle from the read in trajectory. Returns a real array with 3 elements (x, y, and z). An optional index group can be specified; if so, the atomid is in relationship to the group.

Parameters
[in,out]thisTrajectory class
[in]frame
[in]atomatomid of particle to get
[in]groupoptional index group
Returns
coordinate of the particle
257 
258  implicit none
259  real :: trajectory_get_xyz(3)
260  integer, intent(in) :: frame, atom
261  integer :: atom_tmp, natoms
262  class(trajectory), intent(inout) :: this
263  character (len=1024) :: message
264  character (len=*), intent(in), optional :: group
265 
266  call trajectory_check_frame(this, frame)
267 
268  if (this%read_only_index_group .and. present(group)) then
269  call error_stop_program("Do not specify an index group in x() when already specifying an &
270  &index group with read() or read_next().")
271  end if
272 
273  atom_tmp = merge(this%ndx%get(group, atom), atom, present(group))
274  natoms = merge(this%natoms(group), this%natoms(), present(group))
275 
276  if (atom > natoms .or. atom < 1) then
277  write(message, "(a,i0,a,i0,a)") "Tried to access atom number ", atom_tmp, " when there are ", &
278  natoms, ". Note that Fortran uses one-based indexing."
279  call error_stop_program(trim(message))
280  end if
281 
282  trajectory_get_xyz = this%frameArray(frame)%xyz(:,atom_tmp)
283 

◆ trajectory_open()

subroutine dcdfort_trajectory::trajectory_open ( class(trajectory), intent(inout)  this,
character (len=*), intent(in)  filename_in,
character (len=*), intent(in), optional  ndxfile 
)

Trajectory class method which opens DCD file and optionally index file.

Parameters
[in,out]thisthe Trajectory object
[in]filename_inname of DCD file
[in]ndxfilename of GROMACS style index file
88 
89  implicit none
90  class(trajectory), intent(inout) :: this
91  character (len=*), intent(in) :: filename_in
92  character (len=*), intent(in), optional :: ndxfile
93  character (len=206) :: filetype
94  logical :: ex
95  integer :: i, j
96 
97  call this%dcd%open(trim(filename_in))
98 
99  write(error_unit,'(a)') "Opened "//trim(filename_in)//" for reading."
100  call this%dcd%read_header(this%nframes, this%istart, this%nevery, this%iend, this%timestep, this%NUMATOMS)
101 
102  write(error_unit,'(i0,a)') this%NUMATOMS, " atoms present in system."
103  write(error_unit,'(i0,a)') this%NFRAMES, " frames present in trajectory file."
104  write(error_unit,'(a,i0)') "First timestep saved: ", this%ISTART
105  write(error_unit,'(a,i0)') "Last timestep saved: ", this%IEND
106  write(error_unit,'(a,f12.6)') "Simulation timestep: ", this%timestep
107  write(error_unit,'(a,i0,a)') "Trajectory written every ", this%nevery, " timesteps."
108 
109  this%FRAMES_REMAINING = this%NFRAMES
110 
111  this%N = this%NUMATOMS ! Save for use when user selects just one group
112  if (present(ndxfile)) call this%ndx%read(ndxfile, this%NUMATOMS)
113 

◆ trajectory_read()

subroutine dcdfort_trajectory::trajectory_read ( class(trajectory), intent(inout)  this,
character (len=*)  dcdfile,
character (len=*), optional  ndxfile,
character (len=*), optional  ndxgrp 
)
private

Trajectory class method which opens, reads, and closes a trajectory file.

Parameters
[in,out]thisTrajectory class
[in]dcdfileName of DCD trajectory file
[in]ndxfileName of optional index file
[in]ndxgrpName of optional group. If specified, only that group will be read into memory.
292 
293  implicit none
294  class(trajectory), intent(inout) :: this
295  character (len=*), optional :: ndxfile, ndxgrp
296  character (len=*) :: dcdfile
297  integer :: n
298 
299  call this%open(dcdfile, ndxfile)
300 
301  n = this%read_next(this%NFRAMES, ndxgrp)
302 
303  call this%close()
304 

◆ trajectory_read_next()

integer function dcdfort_trajectory::trajectory_read_next ( class(trajectory), intent(inout)  this,
integer, intent(in), optional  F,
character (len=*), optional  ndxgrp 
)
private

Trajectory class method which reads a specified number of frames into memory after using the open() method.

Parameters
[in,out]thisTrajectory class
[in]Fnumber of frames to read in; if not specified, 1 frame is read
[in]ndxgrpread only this index group into memory
Returns
number of frames read in
175 
176  implicit none
177  integer :: trajectory_read_next
178  class(trajectory), intent(inout) :: this
179  integer, intent(in), optional :: f
180  character (len=*), optional :: ndxgrp
181  integer :: i, j, n, stat
182  real, allocatable :: xyz(:,:)
183 
184  ! If the user specified how many frames to read and it is greater than one, use it
185  n = merge(f, 1, present(f))
186 
187  ! Are we near the end of the file?
188  n = min(this%FRAMES_REMAINING, n)
189  this%FRAMES_REMAINING = this%FRAMES_REMAINING - n
190 
191  if (allocated(this%frameArray)) deallocate(this%frameArray)
192  allocate(this%frameArray(n))
193 
194  write(error_unit,*)
195 
196  this%read_only_index_group = .false.
197 
198  if (present(ndxgrp)) then
199 
200  allocate(xyz(3,this%N))
201  this%NUMATOMS = this%natoms(trim(ndxgrp))
202  do i = 1, n
203 
204  if (modulo(i, 1000) .eq. 0) call print_frames_saved(i)
205 
206  allocate(this%frameArray(i)%xyz(3,this%NUMATOMS))
207  call this%dcd%read_next(xyz, this%frameArray(i)%box)
208 
209  do j = 1, size(this%ndx%group)
210  if (trim(this%ndx%group(j)%title) .eq. trim(ndxgrp)) then
211  this%frameArray(i)%xyz = xyz(:,this%ndx%group(j)%LOC)
212  exit
213  end if
214  end do
215 
216  end do
217  deallocate(xyz)
218 
219  this%read_only_index_group = .true.
220 
221  else
222 
223  do i = 1, n
224 
225  if (modulo(i, 1000) .eq. 0) call print_frames_saved(i)
226 
227  allocate(this%frameArray(i)%xyz(3,this%NUMATOMS))
228  call this%dcd%read_next(this%frameArray(i)%xyz, this%frameArray(i)%box)
229 
230  end do
231 
232  end if
233 
234  this%frames_read = n
235  trajectory_read_next = n
236  call print_frames_saved(n)
237 

◆ trajectory_skip_next()

integer function dcdfort_trajectory::trajectory_skip_next ( class(trajectory), intent(inout)  this,
integer, intent(in), optional  F 
)
private

Trajectory class method which skips a specified number of frames.

Parameters
[in,out]thisTrajectory class
[in]Fnumber of frames to skip; if not specified, 1 frame is skipped
Returns
number of frames skipped
143 
144  class(trajectory), intent(inout) :: this
145  integer, intent(in), optional :: f
146  integer :: trajectory_skip_next, i, stat, n
147 
148  ! If the user specified how many frames to read and it is greater than one, use it
149  n = merge(f, 1, present(f))
150 
151  ! Are we near the end of the file?
152  n = min(this%FRAMES_REMAINING, n)
153  this%FRAMES_REMAINING = this%FRAMES_REMAINING - n
154 
155  do i = 1, n
156  call this%dcd%skip_next()
157  if (stat .ne. 0) then
158  write(error_unit,'(a,i0,a)') "Skipped ", i-1, " frames."
159  trajectory_skip_next = i-1
160  exit
161  end if
162  end do
163  write(error_unit,'(a,i0,a)') "Skipped ", n, " frames."
164  trajectory_skip_next = n
165