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)
private

Closes DCD file.

Parameters
[in,out]thisdcdwriter object
154 
155  implicit none
156  class(dcdwriter), intent(inout) :: this
157 
158  close(this%u)
159 

◆ 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
60 
61  implicit none
62  character (len=*), intent(in) :: filename
63  class(dcdwriter), intent(inout) :: this
64 
65  open(newunit=this%u, file=trim(filename), form="unformatted", access="stream", status="replace")
66 

◆ dcdwriter_write_header()

subroutine dcdfort_writer::dcdwriter_write_header ( class(dcdwriter), intent(inout)  this,
integer(kind=int32), intent(in)  istart,
integer(kind=int32), intent(in)  nevery,
real(kind=real32), intent(in)  timestep,
integer(kind=int32), intent(in)  natoms 
)
private

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

◆ dcdwriter_write_next()

subroutine dcdfort_writer::dcdwriter_write_next ( class(dcdwriter), intent(inout)  this,
real(kind=real32), dimension(:,:), intent(in)  xyz,
real(kind=real64), dimension(6), intent(in)  box_in 
)
private

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
168 
169  implicit none
170  real(kind=real32), intent(in) :: xyz(:,:)
171  real(kind=real64), intent(in) :: box_in(6)
172  real(kind=real64) :: box(6)
173  class(dcdwriter), intent(inout) :: this
174  integer(kind=int32) :: coord_size
175 
176  coord_size = size(xyz,1)*4
177  box = box_in
178 
179  ! Should be 48 (6 double precision floats)
180  write(this%u) 48
181 
182  write(this%u) box(1) ! A
183  write(this%u) box(6) ! gamma
184  write(this%u) box(2) ! B
185  write(this%u) box(5) ! beta
186  write(this%u) box(4) ! alpha
187  write(this%u) box(3) ! C
188 
189  write(this%u) 48
190  write(this%u) coord_size
191 
192  write(this%u) xyz(:,1)
193 
194  write(this%u) coord_size
195  write(this%u) coord_size
196 
197  write(this%u) xyz(:,2)
198 
199  write(this%u) coord_size
200  write(this%u) coord_size
201 
202  write(this%u) xyz(:,3)
203 
204  write(this%u) coord_size
205 
206  inquire(unit=this%u, pos=this%curr_pos)
207 
208  this%nframes = this%nframes+1
209  this%iend = this%iend + this%nevery
210 
211  ! Go back and update header
212  write(this%u, pos=this%nframes_pos) this%nframes
213  write(this%u, pos=this%iend_pos) this%iend
214  write(this%u, pos=this%curr_pos)
215 
216  flush(this%u)
217