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, 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, every)
 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, every, skip)
 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
379 
380  implicit none
381  class(trajectory), intent(inout) :: this
382 
383  call this%dcd%close()
384 

◆ 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
box dimensions of specified snapshot. Array containing 6 elements, ordered as A, B, C, alpha, beta, gamma. A = length of unit cell vector along x-axis; B = length of unit cell vector in xy-plane; C = length of unit cell vector in yz-plane; alpha = cosine of angle between B and C; beta = cosine of angle between A and C; gamma = cosine of angle between A and B;
347 
348  implicit none
349  real(8) :: trajectory_get_box(6)
350  class(trajectory), intent(in) :: this
351  integer, intent(in) :: frame
352  integer :: i
353 
354  call trajectory_check_frame(this, frame)
355  do i = 1, 6
356  trajectory_get_box(i) = this%frameArray(frame)%box(i)
357  end do
358 

◆ 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
121 
122  implicit none
123  integer :: trajectory_get_natoms
124  class(trajectory), intent(inout) :: this
125  character (len=*), intent(in), optional :: group
126 
127  if (this%read_only_index_group .and. present(group)) then
128  call error_stop_program("Do not specify an index group in natoms() when already specifying an &
129  &index group with read() or read_next().")
130  end if
131 
132  trajectory_get_natoms = merge(this%ndx%get_natoms(group), this%NUMATOMS, present(group))
133 

◆ 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
267 
268  implicit none
269  real :: trajectory_get_xyz(3)
270  integer, intent(in) :: frame, atom
271  integer :: atom_tmp, natoms
272  class(trajectory), intent(inout) :: this
273  character (len=1024) :: message
274  character (len=*), intent(in), optional :: group
275 
276  call trajectory_check_frame(this, frame)
277 
278  if (this%read_only_index_group .and. present(group)) then
279  call error_stop_program("Do not specify an index group in x() when already specifying an &
280  &index group with read() or read_next().")
281  end if
282 
283  atom_tmp = merge(this%ndx%get(group, atom), atom, present(group))
284  natoms = merge(this%natoms(group), this%natoms(), present(group))
285 
286  if (atom > natoms .or. atom < 1) then
287  write(message, "(a,i0,a,i0,a)") "Tried to access atom number ", atom_tmp, " when there are ", &
288  natoms, ". Note that Fortran uses one-based indexing."
289  call error_stop_program(trim(message))
290  end if
291 
292  ! Faster to list individually
293  trajectory_get_xyz(1) = this%frameArray(frame)%xyz(atom_tmp,1)
294  trajectory_get_xyz(2) = this%frameArray(frame)%xyz(atom_tmp,2)
295  trajectory_get_xyz(3) = this%frameArray(frame)%xyz(atom_tmp,3)
296 

◆ trajectory_open()

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

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

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

◆ trajectory_read()

subroutine dcdfort_trajectory::trajectory_read ( class(trajectory), intent(inout)  this,
character (len=*)  dcdfile,
character (len=*), optional  ndxfile,
character (len=*), optional  ndxgrp,
integer, intent(in), optional  every,
integer, intent(in), optional  skip 
)
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; otherwise, all particles read into memory.
[in]everyRead in every this many frames; default is 1
[in]skipSkip this many frames at the beginning of the trajectory file; default is 0
308 
309  implicit none
310  class(trajectory), intent(inout) :: this
311  character (len=*), optional :: ndxfile, ndxgrp
312  character (len=*) :: dcdfile
313  integer :: n
314  integer, intent(in), optional :: every, skip
315 
316  call this%open(dcdfile, ndxfile)
317 
318  if (present(every)) then
319  write(error_unit,'(a,i0,a)') prompt//"Saving ", every, " snapshots into memory."
320  end if
321  if (present(skip)) then
322  write(error_unit,'(a,i0,a)') prompt//"Skipping first ", skip, " snapshots."
323  end if
324 
325  if (present(skip)) then
326  n = this%skip_next(skip)
327  n = this%read_next(this%NFRAMES-n, ndxgrp, every)
328  else
329  n = this%read_next(this%NFRAMES, ndxgrp, every)
330  end if
331 
332  call this%close()
333 

◆ trajectory_read_next()

integer function dcdfort_trajectory::trajectory_read_next ( class(trajectory), intent(inout)  this,
integer, intent(in), optional  F,
character (len=*), optional  ndxgrp,
integer, intent(in), optional  every 
)
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
[in]everyRead in every this many frames; default is to read in every frame
Returns
number of frames read in
168 
169  implicit none
170  integer :: trajectory_read_next
171  class(trajectory), intent(inout) :: this
172  integer, intent(in), optional :: f, every
173  character (len=*), optional :: ndxgrp
174  integer :: i, j, n, stat, s, k
175  real, allocatable :: xyz(:,:)
176 
177  ! If the user specified how many frames to read and it is greater than one, use it
178  n = merge(f, 1, present(f))
179  s = merge(every, 1, present(every)) - 1
180 
181  ! Are we near the end of the file?
182  n = min(this%FRAMES_REMAINING, n)
183  this%FRAMES_REMAINING = this%FRAMES_REMAINING - n
184 
185  if (present(every)) n = n / every
186 
187  if (allocated(this%frameArray)) deallocate(this%frameArray)
188  allocate(this%frameArray(n))
189 
190  write(error_unit,*)
191 
192  this%read_only_index_group = .false.
193 
194  if (present(ndxgrp)) then
195 
196  allocate(xyz(this%N,3))
197  this%NUMATOMS = this%natoms(trim(ndxgrp))
198  do i = 1, n
199 
200  if (modulo(i*(s+1), 1000) .eq. 0) call print_frames_saved(i, ndxgrp)
201 
202  allocate(this%frameArray(i)%xyz(this%NUMATOMS,3))
203  call this%dcd%skip_next(s)
204  call this%dcd%read_next(xyz, this%frameArray(i)%box)
205 
206  do j = 1, size(this%ndx%group)
207  if (trim(this%ndx%group(j)%title) .eq. trim(ndxgrp)) then
208  ! Tends to be faster to list this out individually
209  do k = 1, this%NUMATOMS
210  this%frameArray(i)%xyz(k,1) = xyz(this%ndx%group(j)%LOC(k),1)
211  this%frameArray(i)%xyz(k,2) = xyz(this%ndx%group(j)%LOC(k),2)
212  this%frameArray(i)%xyz(k,3) = xyz(this%ndx%group(j)%LOC(k),3)
213  end do
214  exit
215  end if
216  end do
217 
218  end do
219  deallocate(xyz)
220 
221  this%read_only_index_group = .true.
222 
223  else
224 
225  do i = 1, n
226 
227  if (modulo(i*(s+1), 1000) .eq. 0) call print_frames_saved(i)
228 
229  allocate(this%frameArray(i)%xyz(this%NUMATOMS,3))
230  call this%dcd%skip_next(s)
231  call this%dcd%read_next(this%frameArray(i)%xyz, this%frameArray(i)%box)
232 
233  end do
234 
235  end if
236 
237  this%frames_read = n
238  trajectory_read_next = n
239  call print_frames_saved(n, ndxgrp)
240 

◆ 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
142 
143  class(trajectory), intent(inout) :: this
144  integer, intent(in), optional :: f
145  integer :: trajectory_skip_next, i, stat, n
146 
147  ! If the user specified how many frames to read and it is greater than one, use it
148  n = merge(f, 1, present(f))
149 
150  ! Are we near the end of the file?
151  n = min(this%FRAMES_REMAINING, n)
152  this%FRAMES_REMAINING = this%FRAMES_REMAINING - n
153 
154  call this%dcd%skip_next(n)
155  write(error_unit,'(a,i0,a)') prompt//"Skipped ", n, " frames."
156  trajectory_skip_next = n
157