Add function to read single xtc frmames.

This commit is contained in:
Niels Müller
2017-08-18 13:22:58 +02:00
parent edc1bd759a
commit d178df5c7f
3 changed files with 57 additions and 5 deletions

View File

@ -8,7 +8,7 @@ Trajectories in trr format may be read with `.gromacs.reader.TRRReader`, which i
import os import os
from .tpxio import TPXReader from .tpxio import TPXReader
from .xtcio import XTCReader from .xtcio import XTCReader, read_xtcframe
from .enxio import EDRFile from .enxio import EDRFile
from .errors import FileTypeError from .errors import FileTypeError
from .gromacs.reader import index_filename_for_xtc from .gromacs.reader import index_filename_for_xtc

View File

@ -1,4 +1,4 @@
# from libc.stdio cimport FILE from libc.stdio cimport FILE
from utility cimport * from utility cimport *
@ -6,11 +6,36 @@ from utility cimport *
# ctypedef struct XDR: # ctypedef struct XDR:
# pass # pass
cdef extern from "gromacs/fileio/filetypes.h":
int fn2ftp(const char *fn)
#cdef extern from "gromacs/fileio/gmx_internal_xdr.h":
# ctypedef struct XDR:
# pass
#
# ctypedef enum xdr_op:
# XDR_ENCODE = 0
# XDR_DECODE = 1
# XDR_FREE = 2
#
# void xdrstdio_create (XDR *__xdrs, FILE *__file, xdr_op __xop)
cdef extern from "gromacs/fileio/gmxfio.h": cdef extern from "gromacs/fileio/gmxfio.h":
ctypedef struct t_fileio: ctypedef struct t_fileio:
pass pass
# FILE *fp # the file pointer */
# gmx_bool bRead # the file is open for reading */
# gmx_bool bDouble # write doubles instead of floats */
# gmx_bool bReadWrite # the file is open for reading and writing */
# char *fn # the file name */
# # XDR *xdr # the xdr data pointer */
# int xdrmode # the xdr mode */
# int iFTP # the file type identifier */
# t_fileio *next, *prev # next and previous file pointers in the linked list */
# # tMPI_Lock_t mtx; # content locking mutex. This is a fast lock
void gmx_fio_rewind(t_fileio *fio) void gmx_fio_rewind(t_fileio *fio)

View File

@ -16,9 +16,9 @@ from .errors import InvalidIndexException, InvalidMagicException, XTCError
cdef extern from "gromacs/fileio/xtcio.h": cdef extern from "gromacs/fileio/xtcio.h":
t_fileio *open_xtc(const char *filename, const char *mode) t_fileio *open_xtc(const char *filename, const char *mode) nogil
void close_xtc(t_fileio *fio) void close_xtc(t_fileio *fio) nogil
int read_first_xtc(t_fileio *fio, int read_first_xtc(t_fileio *fio,
int *natoms, gmx_int64_t *step, real *time, int *natoms, gmx_int64_t *step, real *time,
@ -26,7 +26,7 @@ cdef extern from "gromacs/fileio/xtcio.h":
int read_next_xtc(t_fileio *fio, int read_next_xtc(t_fileio *fio,
int natoms, gmx_int64_t *step, real *time, int natoms, gmx_int64_t *step, real *time,
matrix box, rvec *x, real *prec, gmx_bool *_bOK) matrix box, rvec *x, real *prec, gmx_bool *_bOK) nogil
if sizeof(real) == 4: if sizeof(real) == 4:
@ -69,6 +69,33 @@ cdef array get_xtc_index(t_fileio *fio):
return cache[:-1] return cache[:-1]
def read_xtcframe(fname, pos, natoms):
cdef t_fileio *fio = open_xtc(fname.encode(), 'r')
gmx_fio_seek(fio, pos)
cdef:
matrix box
gmx_bool _bOK
real time, prec
gmx_int64_t cur_step
np.ndarray[real, ndim=2] coords = np.empty((natoms, 3), dtype=np_real)
int nratms = natoms
with nogil:
read_next_xtc(fio, nratms, &cur_step, &time, box, <rvec *>coords.data, &prec, &_bOK)
close_xtc(fio)
if _bOK:
frame = XTCFrame()
frame._coordinates = coords
frame.index = cur_step
frame.time = time
frame.box = box
return frame
else:
raise
cdef class XTCReader: cdef class XTCReader:
cdef: cdef:
t_fileio *fio t_fileio *fio