More topology handling.

This commit is contained in:
Niels Müller
2016-11-09 12:57:49 +01:00
parent 108200cea8
commit b0bc90fa6c
3 changed files with 47 additions and 10 deletions

View File

@ -1,4 +1,6 @@
from libc.stdio cimport FILE
from utility cimport *
from math cimport *
from mdtypes cimport *
@ -60,11 +62,32 @@ cdef extern from "gromacs/topology/block.h":
pass
cdef extern from "gromacs/topology/idef.h":
ctypedef struct t_ilist:
pass
ctypedef struct t_ilist:
pass
ctypedef struct gmx_ffparams_t:
pass
cdef enum t_ft_enum:
F_LJ
ctypedef union t_iparams:
pass
ctypedef int t_functype
ctypedef struct gmx_cmap_t:
pass
ctypedef struct gmx_ffparams_t:
int ntypes;
int atnr;
t_functype *functype;
t_iparams *iparams;
double reppow; # The repulsion power for VdW: C12*r^-reppow */
real fudgeQQ; # The scaling factor for Coulomb 1-4: f*q1*q2 */
gmx_cmap_t cmap_grid; # The dihedral correction maps */
void pr_ffparams(FILE *fp, int indent, const char *title,
const gmx_ffparams_t *ffparams, gmx_bool bShowNumbers)
cdef extern from "gromacs/topology/topology.h":
ctypedef struct gmx_moltype_t:

View File

@ -102,6 +102,17 @@ cdef open_tpx(const char* filename, t_inputrec *ir, matrix box, int *natoms, gmx
stdio.setbuf(stdio.stderr, NULL)
return return_code
cdef read_ffparams(gmx_ffparams_t *ffparams, gmx_bool bShowNumbers):
cdef char buffer[stdio.BUFSIZ]
stdio.setbuf(stdio.stdout, buffer)
pr_ffparams(stdio.stdout, 0, '', ffparams, bShowNumbers)
stdio.fflush(stdio.stderr)
stdio.fseek(stdio.stderr, 0, stdio.SEEK_END)
stdio.setbuf(stdio.stderr, NULL)
return buffer
cdef class TPXReader:
cdef:
t_tpxheader header
@ -208,6 +219,9 @@ cdef class TPXReader:
return self.input_record.nstenergy
def read_ff(self):
return read_ffparams(&self.topology.ffparams, True)
def __cinit__(self, filename):
filename = cstr(filename)

View File

@ -22,11 +22,11 @@ if 'gromacs' in os.environ.get('LD_LIBRARY_PATH', ''):
if 'gromacs' in p:
library_dirs.append(p)
lib = p
gmx_root = lib.split('lib')[0]
include = os.path.join(gmx_root, 'include')
if os.path.exists(include):
include_dirs.append(include)
check_header_version(include)
gmx_root = lib.split('lib')[0]
include = os.path.join(gmx_root, 'include')
if os.path.exists(include):
include_dirs.append(include)
check_header_version(include)
extensions = [
Extension(