ee035a0fd8a95c1624be08e598f4235d9d1bbc8e
[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,2015,2016,2017,2018, 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 "gromacs/fileio/pdbio.h"
42
43 #ifdef __cplusplus
44 extern "C" {
45 #endif
46
47 struct gmx_mtop_t;
48 struct gmx_output_env_t;
49 struct t_atoms;
50 struct t_fileio;
51 struct t_topology;
52 struct t_trxframe;
53
54 /* a dedicated status type contains fp, etc. */
55 typedef struct t_trxstatus t_trxstatus;
56
57 /* I/O function types */
58
59 /************************************************
60  *             Trajectory functions
61  ************************************************/
62
63 int prec2ndec(real prec);
64 /* Convert precision in 1/(nm) to number of decimal places */
65
66 /*! \brief Convert number of decimal places to trajectory precision in
67  * 1/(nm) */
68 real ndec2prec(int ndec);
69
70 void clear_trxframe(struct t_trxframe *fr, gmx_bool bFirst);
71 /* Set all content gmx_booleans to FALSE.
72  * When bFirst = TRUE, set natoms=-1, all pointers to NULL
73  *                     and all data to zero.
74  */
75
76 void set_trxframe_ePBC(struct t_trxframe *fr, int ePBC);
77 /* Set the type of periodic boundary conditions, ePBC=-1 is not set */
78
79 int nframes_read(t_trxstatus *status);
80 /* Returns the number of frames read from the trajectory */
81
82 int write_trxframe_indexed(t_trxstatus *status, const t_trxframe *fr, int nind,
83                            const int *ind, gmx_conect gc);
84 /* Write an indexed frame to a TRX file, see write_trxframe. gc may be NULL */
85
86 int write_trxframe(t_trxstatus *status, struct t_trxframe *fr, gmx_conect gc);
87 /* Write a frame to a TRX file.
88  * Only entries for which the gmx_boolean is TRUE will be written,
89  * except for step, time, lambda and/or box, which may not be
90  * omitted for certain trajectory formats.
91  * The precision for .xtc and .gro is fr->prec, when fr->bPrec=FALSE,
92  * the precision is set to 1000.
93  * gc is important for pdb file writing only and may be NULL.
94  */
95
96 int write_trx(t_trxstatus *status, int nind, const int *ind, const t_atoms *atoms,
97               int step, real time, matrix box, rvec x[], rvec *v,
98               gmx_conect gc);
99 /* Write an indexed frame to a TRX file.
100  * v can be NULL.
101  * atoms can be NULL for file types which don't need atom names.
102  */
103 t_trxstatus *
104 trjtools_gmx_prepare_tng_writing(const char               *filename,
105                                  char                      filemode,
106                                  t_trxstatus              *in,
107                                  const char               *infile,
108                                  const int                 natoms,
109                                  const struct gmx_mtop_t  *mtop,
110                                  const int                *index,
111                                  const char               *index_group_name);
112 /* Sets up *out for writing TNG. If *in != NULL and contains a TNG trajectory
113  * some data, e.g. molecule system, will be copied over from *in to the return value.
114  * If *in == NULL a file name (infile) of a TNG file can be provided instead
115  * and used for copying data to the return value.
116  * If there is no TNG input natoms is used to create "implicit atoms" (no atom
117  * or molecular data present). If natoms == -1 the number of atoms are
118  * not known (or there is already a TNG molecule system to copy, in which case
119  * natoms is not required anyhow). If an group of indexed atoms are written
120  * natoms must be the length of index. index_group_name is the name of the
121  * index group.
122  */
123
124 /*! \brief Write a trxframe to the TNG file in status.
125  *
126  * This function is needed because both t_trxstatus and
127  * gmx_tng_trajectory_t are encapsulated, so client trajectory-writing
128  * code with a t_trxstatus can't just call the TNG writing
129  * function. */
130 void write_tng_frame(t_trxstatus        *status,
131                      struct t_trxframe  *fr);
132
133 void close_trx(t_trxstatus *status);
134 /* Close trajectory file as opened with read_first_x, read_first_frame
135  * or open_trx.
136  * Also frees memory in the structure.
137  */
138
139 t_trxstatus *open_trx(const char *outfile, const char *filemode);
140 /* Open a TRX file and return an allocated status pointer */
141
142 struct t_fileio *trx_get_fileio(t_trxstatus *status);
143 /* get a fileio from a trxstatus */
144
145 float trx_get_time_of_final_frame(t_trxstatus *status);
146 /* get time of final frame. Only supported for TNG and XTC */
147
148 gmx_bool bRmod_fd(double a, double b, double c, gmx_bool bDouble);
149 /* Returns TRUE when (a - b) MOD c = 0, using a margin which is slightly
150  * larger than the float/double precision.
151  */
152
153 #if GMX_DOUBLE
154 #define bRmod(a, b, c) bRmod_fd(a, b, c, TRUE)
155 #else
156 #define bRmod(a, b, c) bRmod_fd(a, b, c, FALSE)
157 #endif
158
159 int check_times2(real t, real t0, gmx_bool bDouble);
160 /* This routine checkes if the read-in time is correct or not;
161  * returns -1 if t<tbegin or t MOD dt = t0,
162  *          0 if tbegin <= t <=tend+margin,
163  *          1 if t>tend
164  * where margin is 0.1*min(t-tp,tp-tpp), if this positive, 0 otherwise.
165  * tp and tpp should be the time of the previous frame and the one before.
166  * The mod is done with single or double precision accuracy depending
167  * on the value of bDouble.
168  */
169
170 int check_times(real t);
171 /* This routine checkes if the read-in time is correct or not;
172  * returns -1 if t<tbegin,
173  *          0 if tbegin <= t <=tend,
174  *          1 if t>tend
175  */
176
177
178
179
180
181 /* For trxframe.flags, used in trxframe read routines.
182  * When a READ flag is set, the field will be read when present,
183  * but a frame might be returned which does not contain the field.
184  * When a NEED flag is set, frames not containing the field will be skipped.
185  */
186 #define TRX_READ_X    (1<<0)
187 #define TRX_NEED_X    (1<<1)
188 #define TRX_READ_V    (1<<2)
189 #define TRX_NEED_V    (1<<3)
190 #define TRX_READ_F    (1<<4)
191 #define TRX_NEED_F    (1<<5)
192 /* Useful for reading natoms from a trajectory without skipping */
193 #define TRX_DONT_SKIP (1<<6)
194
195 /* For trxframe.not_ok */
196 #define HEADER_NOT_OK (1<<0)
197 #define DATA_NOT_OK   (1<<1)
198 #define FRAME_NOT_OK  (HEADER_NOT_OK | DATA_NOT_OK)
199
200 bool read_first_frame(const gmx_output_env_t *oenv, t_trxstatus **status,
201                       const char *fn, struct t_trxframe *fr, int flags);
202 /* Read the first frame which is in accordance with flags, which are
203  * defined further up in this file.
204  * Memory will be allocated for flagged entries.
205  * The flags are copied to fr for subsequent calls to read_next_frame.
206  * Returns true when succeeded, false otherwise.
207  */
208
209 bool read_next_frame(const gmx_output_env_t *oenv, t_trxstatus *status,
210                      struct t_trxframe *fr);
211 /* Reads the next frame which is in accordance with fr->flags.
212  * Returns true when succeeded, false otherwise.
213  */
214
215 int read_first_x(const gmx_output_env_t *oenv, t_trxstatus **status,
216                  const char *fn, real *t, rvec **x, matrix box);
217 /* These routines read first coordinates and box, and allocates
218  * memory for the coordinates, for a trajectory file.
219  * The routine returns the number of atoms, or 0 when something is wrong.
220  * The integer in status should be passed to calls of read_next_x
221  */
222
223 gmx_bool read_next_x(const gmx_output_env_t *oenv, t_trxstatus *status, real *t, rvec x[], matrix box);
224 /* Read coordinates and box from a trajectory file. Return TRUE when all well,
225  * or FALSE when end of file (or last frame requested by user).
226  * status is the integer set in read_first_x.
227  */
228
229 void rewind_trj(t_trxstatus *status);
230 /* Rewind trajectory file as opened with read_first_x */
231
232 struct t_topology *read_top(const char *fn, int *ePBC);
233 /* Extract a topology data structure from a topology file.
234  * If ePBC!=NULL *ePBC gives the pbc type.
235  */
236
237 #ifdef __cplusplus
238 }
239 #endif
240
241 #endif