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

Module that contains dcdwriter class. More...

Data Types

type  dcdwriter
 dcdwriter class More...
 

Functions/Subroutines

subroutine dcdwriter_open (this, filename)
 Opens new file to write to. More...
 
subroutine dcdwriter_write_header (this, istart, nevery, timestep, natoms)
 Writes header to new DCD file. More...
 
subroutine dcdwriter_close (this)
 Closes DCD file. More...
 
subroutine dcdwriter_write_next (this, xyz, box_in)
 Writes snapshot to an open DCD file. More...
 

Detailed Description

Module that contains dcdwriter class.

Function/Subroutine Documentation

◆ dcdwriter_close()

subroutine dcdfort_writer::dcdwriter_close ( class(dcdwriter), intent(inout)  this)

Closes DCD file.

Parameters
[in,out]thisdcdwriter object
149 
150  implicit none
151  class(dcdwriter), intent(inout) :: this
152 
153  close(this%u)
154 

◆ dcdwriter_open()

subroutine dcdfort_writer::dcdwriter_open ( class(dcdwriter), intent(inout)  this,
character (len=*), intent(in)  filename 
)

Opens new file to write to.

Parameters
[in,out]thisdcdwriter object
[in]filenamename of new DCD file to write to
59 
60  implicit none
61  character (len=*), intent(in) :: filename
62  class(dcdwriter), intent(inout) :: this
63  integer :: filesize
64  logical :: ex
65 
66  open(newunit=this%u, file=trim(filename), form="unformatted", access="stream", status="new")
67 

◆ dcdwriter_write_header()

subroutine dcdfort_writer::dcdwriter_write_header ( class(dcdwriter), intent(inout)  this,
integer, intent(in)  istart,
integer, intent(in)  nevery,
real, intent(in)  timestep,
integer, intent(in)  natoms 
)

Writes header to new DCD file.

Parameters
[in,out]thisdcdwriter object
[in]istartfirst timestep in file
[in]neveryhow often snapshots written (in timesteps)
[in]timestepsimulation timestep
[in]natomsnumber of atoms in each snapshot
77 
78  implicit none
79 
80  integer :: dummy, i, n
81  integer, intent(in) :: istart, nevery, natoms
82  real, intent(in) :: timestep
83  class(dcdwriter), intent(inout) :: this
84  character (len=79) :: remarks1, remarks2
85  character (len=8) :: date
86  character (len=10) :: time
87 
88  call date_and_time(date=date,time=time)
89  remarks1 = "Created by libdcdfort"
90  remarks2 = "REMARK Created on "//date//" "//time
91 
92  this%nframes = 0
93  this%iend = istart
94  this%nevery = nevery
95 
96  write(this%u) 84
97  write(this%u) "CORD"
98 
99  inquire(unit=this%u, pos=this%nframes_pos)
100  ! Number of snapshots in file
101  write(this%u) this%nframes
102 
103  ! Timestep of first snapshot
104  write(this%u) istart
105 
106  ! Save snapshots every this many steps
107  write(this%u) nevery
108 
109  inquire(unit=this%u, pos=this%iend_pos)
110  ! Timestep of last snapshot
111  write(this%u) this%iend
112 
113  do i = 1, 4
114  write(this%u) 0
115  end do
116 
117  ! Has unit cell
118  write(this%u) 1
119 
120  ! Simulation timestep
121  write(this%u) timestep
122 
123  do i = 1, 9
124  write(this%u) 0
125  end do
126 
127  ! Pretend to be CHARMM version 24
128  write(this%u) 24
129  write(this%u) 84
130  write(this%u) 164
131 
132  write(this%u) 2
133  write(this%u) remarks1//c_null_char
134  write(this%u) remarks2//c_null_char
135 
136  write(this%u) 164
137  write(this%u) 4
138 
139  ! Number of atoms in each snapshot
140  write(this%u) natoms
141 
142  write(this%u) 4
143 

◆ dcdwriter_write_next()

subroutine dcdfort_writer::dcdwriter_write_next ( class(dcdwriter), intent(inout)  this,
real, dimension(:,:), intent(in)  xyz,
real(8), dimension(6), intent(in)  box_in 
)

Writes snapshot to an open DCD file.

Writes a new snapshot to a DCD file. Header should have already been written.

Parameters
[in,out]thisdcdwriter object
[in]xyzcoordinates of all atoms in this snapshot
[in]box_inbox dimensions for this snapshot
163 
164  implicit none
165  real, intent(in) :: xyz(:,:)
166  real(8), intent(in) :: box_in(6)
167  real(8) :: box(6)
168  integer :: dummy
169  class(dcdwriter), intent(inout) :: this
170  integer :: coord_size
171 
172  coord_size = size(xyz,2)*4
173  box = box_in
174 
175  ! Should be 48 (6 double precision floats)
176  write(this%u) 48
177 
178  box(4) = (90.0 - box(4)) * dsin(pi/180.0)
179  box(5) = (90.0 - box(5)) * dsin(pi/180.0)
180  box(6) = (90.0 - box(6)) * dsin(pi/180.0)
181  write(this%u) box(1) ! A
182  write(this%u) box(6) ! gamma
183  write(this%u) box(2) ! B
184  write(this%u) box(5) ! beta
185  write(this%u) box(4) ! alpha
186  write(this%u) box(3) ! C
187 
188  write(this%u) 48
189  write(this%u) coord_size
190 
191  write(this%u) xyz(1,:)
192 
193  write(this%u) 48
194  write(this%u) coord_size
195 
196  write(this%u) xyz(2,:)
197 
198  write(this%u) 48
199  write(this%u) coord_size
200 
201  write(this%u) xyz(3,:)
202 
203  write(this%u) 48
204 
205  inquire(unit=this%u, pos=this%curr_pos)
206 
207  this%nframes = this%nframes+1
208  this%iend = this%iend + this%nevery
209 
210  ! Go back and update header
211  write(this%u, pos=this%nframes_pos) this%nframes
212  write(this%u, pos=this%iend_pos) this%iend
213  write(this%u, pos=this%curr_pos)
214