Merge release-4-6 into release-5-0
[alexxy/gromacs.git] / src / gromacs / fileio / trxio.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 #ifndef GMX_FILEIO_TRXIO_H
39 #define GMX_FILEIO_TRXIO_H
40
41 #include "../legacyheaders/typedefs.h"
42 #include "filenm.h"
43 #include "../legacyheaders/readinp.h"
44 #include "pdbio.h"
45 #include "../legacyheaders/oenv.h"
46 #include "gmxfio.h"
47
48 #ifdef __cplusplus
49 extern "C" {
50 #endif
51 #if 0 /* avoid screwing up indentation */
52 }
53 #endif
54
55 /* a dedicated status type contains fp, etc. */
56 typedef struct t_trxstatus t_trxstatus;
57
58 /* I/O function types */
59
60 /************************************************
61  *             Trajectory functions
62  ************************************************/
63
64 int prec2ndec(real prec);
65 /* Convert precision in 1/(nm) to number of decimal places */
66
67 /*! \brief Convert number of decimal places to trajectory precision in
68  * 1/(nm) */
69 real ndec2prec(int ndec);
70
71 void clear_trxframe(t_trxframe *fr, gmx_bool bFirst);
72 /* Set all content gmx_booleans to FALSE.
73  * When bFirst = TRUE, set natoms=-1, all pointers to NULL
74  *                     and all data to zero.
75  */
76
77 void set_trxframe_ePBC(t_trxframe *fr, int ePBC);
78 /* Set the type of periodic boundary conditions, ePBC=-1 is not set */
79
80 int nframes_read(t_trxstatus *status);
81 /* Returns the number of frames read from the trajectory */
82
83 int write_trxframe_indexed(t_trxstatus *status, t_trxframe *fr, int nind,
84                            const atom_id *ind, gmx_conect gc);
85 /* Write an indexed frame to a TRX file, see write_trxframe. gc may be NULL */
86
87 int write_trxframe(t_trxstatus *status, t_trxframe *fr, gmx_conect gc);
88 /* Write a frame to a TRX file.
89  * Only entries for which the gmx_boolean is TRUE will be written,
90  * except for step, time, lambda and/or box, which may not be
91  * omitted for certain trajectory formats.
92  * The precision for .xtc and .gro is fr->prec, when fr->bPrec=FALSE,
93  * the precision is set to 1000.
94  * gc is important for pdb file writing only and may be NULL.
95  */
96
97 int write_trx(t_trxstatus *status, int nind, const atom_id *ind, t_atoms *atoms,
98               int step, real time, matrix box, rvec x[], rvec *v,
99               gmx_conect gc);
100 /* Write an indexed frame to a TRX file.
101  * v can be NULL.
102  * atoms can be NULL for file types which don't need atom names.
103  */
104
105 void trjtools_gmx_prepare_tng_writing(const char       *filename,
106                                       char              filemode,
107                                       t_trxstatus      *in,
108                                       t_trxstatus     **out,
109                                       const char       *infile,
110                                       const int         natoms,
111                                       const gmx_mtop_t *mtop,
112                                       const atom_id    *index,
113                                       const char       *index_group_name);
114 /* Sets up *out for writing TNG. If *in != NULL and contains a TNG trajectory
115  * some data, e.g. molecule system, will be copied over from *in to *out.
116  * If *in == NULL a file name (infile) of a TNG file can be provided instead
117  * and used for copying data to *out.
118  * If there is no TNG input natoms is used to create "implicit atoms" (no atom
119  * or molecular data present). If natoms == -1 the number of atoms are
120  * not known (or there is already a TNG molecule system to copy, in which case
121  * natoms is not required anyhow). If an group of indexed atoms are written
122  * natoms must be the length of index. index_group_name is the name of the
123  * index group.
124  */
125
126 /*! \brief Write a trxframe to the TNG file in status.
127  *
128  * This function is needed because both t_trxstatus and
129  * tng_trajectory_t are encapsulated, so client trajectory-writing
130  * code with a t_trxstatus can't just call the TNG writing
131  * function. */
132 void write_tng_frame(t_trxstatus *status,
133                      t_trxframe  *fr);
134
135 void close_trx(t_trxstatus *status);
136 /* Close trj file as opened with read_first_x, read_frist_frame
137  * or open_trx. Identical to close_trj.
138  */
139
140 t_trxstatus *open_trx(const char *outfile, const char *filemode);
141 /* Open a TRX file and return an allocated status pointer */
142
143 t_fileio *trx_get_fileio(t_trxstatus *status);
144 /* get a fileio from a trxstatus */
145
146 float trx_get_time_of_final_frame(t_trxstatus *status);
147 /* get time of final frame. Only supported for TNG and XTC */
148
149 gmx_bool bRmod_fd(double a, double b, double c, gmx_bool bDouble);
150 /* Returns TRUE when (a - b) MOD c = 0, using a margin which is slightly
151  * larger than the float/double precision.
152  */
153
154 #ifdef GMX_DOUBLE
155 #define bRmod(a, b, c) bRmod_fd(a, b, c, TRUE)
156 #else
157 #define bRmod(a, b, c) bRmod_fd(a, b, c, FALSE)
158 #endif
159
160 int check_times2(real t, real t0, gmx_bool bDouble);
161 /* This routine checkes if the read-in time is correct or not;
162  * returns -1 if t<tbegin or t MOD dt = t0,
163  *          0 if tbegin <= t <=tend+margin,
164  *          1 if t>tend
165  * where margin is 0.1*min(t-tp,tp-tpp), if this positive, 0 otherwise.
166  * tp and tpp should be the time of the previous frame and the one before.
167  * The mod is done with single or double precision accuracy depending
168  * on the value of bDouble.
169  */
170
171 int check_times(real t);
172 /* This routine checkes if the read-in time is correct or not;
173  * returns -1 if t<tbegin,
174  *          0 if tbegin <= t <=tend,
175  *          1 if t>tend
176  */
177
178
179
180
181
182 /* For trxframe.flags, used in trxframe read routines.
183  * When a READ flag is set, the field will be read when present,
184  * but a frame might be returned which does not contain the field.
185  * When a NEED flag is set, frames not containing the field will be skipped.
186  */
187 #define TRX_READ_X    (1<<0)
188 #define TRX_NEED_X    (1<<1)
189 #define TRX_READ_V    (1<<2)
190 #define TRX_NEED_V    (1<<3)
191 #define TRX_READ_F    (1<<4)
192 #define TRX_NEED_F    (1<<5)
193 /* Useful for reading natoms from a trajectory without skipping */
194 #define TRX_DONT_SKIP (1<<6)
195
196 /* For trxframe.not_ok */
197 #define HEADER_NOT_OK (1<<0)
198 #define DATA_NOT_OK   (1<<1)
199 #define FRAME_NOT_OK  (HEADER_NOT_OK | DATA_NOT_OK)
200
201 int read_first_frame(const output_env_t oenv, t_trxstatus **status,
202                      const char *fn, t_trxframe *fr, int flags);
203 /* Read the first frame which is in accordance with flags, which are
204  * defined further up in this file.
205  * Returns natoms when succeeded, 0 otherwise.
206  * Memory will be allocated for flagged entries.
207  * The flags are copied to fr for subsequent calls to read_next_frame.
208  * Returns TRUE when succeeded, FALSE otherwise.
209  */
210
211 gmx_bool read_next_frame(const output_env_t oenv, t_trxstatus *status,
212                          t_trxframe *fr);
213 /* Reads the next frame which is in accordance with fr->flags.
214  * Returns TRUE when succeeded, FALSE otherwise.
215  */
216
217 int read_first_x(const output_env_t oenv, t_trxstatus **status,
218                  const char *fn, real *t, rvec **x, matrix box);
219 /* These routines read first coordinates and box, and allocates
220  * memory for the coordinates, for a trajectory file.
221  * The routine returns the number of atoms, or 0 when something is wrong.
222  * The integer in status should be passed to calls of read_next_x
223  */
224
225 gmx_bool read_next_x(const output_env_t oenv, t_trxstatus *status, real *t, rvec x[], matrix box);
226 /* Read coordinates and box from a trajectory file. Return TRUE when all well,
227  * or FALSE when end of file (or last frame requested by user).
228  * status is the integer set in read_first_x.
229  */
230
231 void close_trj(t_trxstatus *status);
232 /* Close trj file as opened with read_first_x, read_frist_frame
233  * or open_trx. Identical to close_trx.
234  */
235
236 void rewind_trj(t_trxstatus *status);
237 /* Rewind trj file as opened with read_first_x */
238
239 t_topology *read_top(const char *fn, int *ePBC);
240 /* Extract a topology data structure from a topology file.
241  * If ePBC!=NULL *ePBC gives the pbc type.
242  */
243
244 #ifdef __cplusplus
245 }
246 #endif
247
248 #endif