dcdfort
Modern Fortran library for analyzing DCD trajectories
All Classes Namespaces Files Functions Variables Pages
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
152 
153  implicit none
154  class(dcdwriter), intent(inout) :: this
155 
156  close(this%u)
157 

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

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

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