4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
17 * GROMACS: A message-passing parallel molecular dynamics implementation
18 * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19 * Comp. Phys. Comm. 91, 43-56 (1995)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
27 * Good ROcking Metal Altar for Chronical Sinners
33 static char *SRCID_confio_h
= "$Id$";
36 #ident "@(#) confio.h 1.16 11/23/92"
37 #endif /* HAVE_IDENT */
41 /* For reading coordinate files it is assumed that enough memory
42 * has been allocated beforehand.
59 extern void init_t_atoms(t_atoms
*atoms
, int natoms
, bool bPdbinfo
);
60 /* allocate memory for the arrays, set nr to natoms and nres to 0
61 * set pdbinfo to NULL or allocate memory for it */
63 extern void free_t_atoms(t_atoms
*atoms
);
64 /* free all the arrays and set the nr and nres to 0 */
66 void clear_g96info(t_g96info
*info
);
67 /* set all bools in the info struct to FALSE */
69 int read_g96_conf(FILE *fp
,char *infile
,int nwanted
,t_g96info
*info
,
70 char *title
,t_atoms
*atoms
,rvec
*x
, rvec
*v
,matrix box
);
71 /* read a Gromos96 coordinate or trajectory file, *
72 * returns the number of atoms *
73 * sets what's in the frame in info *
74 * read from fp, infile is only needed for error messages *
75 * nwanted is the number of wanted coordinates, *
76 * set this to -1 if you want to know the number of atoms in the file *
77 * title, atoms, x, v can all be NULL, in which case they won't be read */
79 void write_g96_conf(FILE *out
,char *title
,t_atoms
*atoms
,
80 rvec
*x
,rvec
*v
,matrix box
,
81 int nindex
,atom_id
*index
);
82 /* write a Gromos96 coordinate file *
83 * x, v and index can be NULL */
85 extern bool gro_next_x(FILE *status
,real
*t
,int natoms
,rvec x
[],matrix box
);
86 extern int gro_first_x(FILE *status
, real
*t
, rvec
**x
, matrix box
);
87 /* read first/next x frame from gro file */
89 extern bool gro_next_x_or_v(FILE *status
,real
*t
,int natoms
,
90 rvec x
[],rvec
*v
,matrix box
);
91 extern int gro_first_x_or_v(FILE *status
, real
*t
,
92 rvec
**x
, rvec
**v
, matrix box
);
93 /* read first/next x and/or v frame from gro file */
95 extern bool gro_next_v(FILE *status
,real
*t
,int natoms
,rvec
*v
,matrix box
);
96 extern int gro_first_v(FILE *status
, real
*t
, rvec
**v
, matrix box
);
97 /* read first/next v frame from gro file */
99 extern bool gro_next_x_v(FILE *status
,real
*t
,int natoms
,
100 rvec x
[],rvec
*v
,matrix box
);
101 extern int gro_first_x_v(FILE *status
, real
*t
,
102 rvec
**x
, rvec
**v
, matrix box
);
103 /* read first/next x and v frame from gro file */
105 extern void write_hconf(FILE *out
,char *title
,
106 t_atoms
*atoms
,rvec
*x
,
109 extern void write_hconf_indexed(FILE *out
,char *title
,t_atoms
*atoms
,
110 int nx
,atom_id index
[],
111 rvec
*x
,rvec
*v
,matrix box
);
113 extern void write_hconf_p(FILE *out
,char *title
,t_atoms
*atoms
, int pr
,
114 rvec
*x
,rvec
*v
,matrix box
);
115 /* Write a Gromos file with precision pr: number of decimal places in x,
116 * v has one place more. */
118 extern void write_xdr_conf(char *outfile
,char *title
,t_atoms
*atoms
,
119 rvec x
[],rvec
*v
,matrix box
);
121 extern void read_xdr_coordnum(char *infile
,int *natoms
);
123 extern void read_xdr_conf(char *infile
,char *title
,t_atoms
*atoms
,
124 rvec x
[],rvec
*v
,matrix box
);
126 void write_sto_conf_indexed(char *outfile
,char *title
,t_atoms
*atoms
,
127 rvec x
[],rvec
*v
,matrix box
,
128 atom_id nindex
,atom_id index
[]);
129 /* like write_sto_conf, but indexed */
131 extern void write_sto_conf(char *outfile
, char *title
,t_atoms
*atoms
,
132 rvec x
[],rvec
*v
, matrix box
);
133 /* write atoms, x, v (if .gro and not NULL) and box (if not NULL)
134 * to an STO (.gro or .pdb) file */
136 extern void get_stx_coordnum (char *infile
,int *natoms
);
137 /* read the number of atoms from an STX file */
139 extern void read_stx_conf(char *infile
, char *title
,t_atoms
*atoms
,
140 rvec x
[],rvec
*v
, matrix box
);
141 /* read atoms, x, v and box from an STX file */
147 #endif /* _confio_h */