2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
38 #ifndef GMX_FILEIO_TRXIO_H
39 #define GMX_FILEIO_TRXIO_H
41 #include "../legacyheaders/readinp.h"
42 #include "../legacyheaders/oenv.h"
51 #if 0 /* avoid screwing up indentation */
60 /* a dedicated status type contains fp, etc. */
61 typedef struct t_trxstatus t_trxstatus;
63 /* I/O function types */
65 /************************************************
66 * Trajectory functions
67 ************************************************/
69 int prec2ndec(real prec);
70 /* Convert precision in 1/(nm) to number of decimal places */
72 /*! \brief Convert number of decimal places to trajectory precision in
74 real ndec2prec(int ndec);
76 void clear_trxframe(struct t_trxframe *fr, gmx_bool bFirst);
77 /* Set all content gmx_booleans to FALSE.
78 * When bFirst = TRUE, set natoms=-1, all pointers to NULL
79 * and all data to zero.
82 void set_trxframe_ePBC(struct t_trxframe *fr, int ePBC);
83 /* Set the type of periodic boundary conditions, ePBC=-1 is not set */
85 int nframes_read(t_trxstatus *status);
86 /* Returns the number of frames read from the trajectory */
88 int write_trxframe_indexed(t_trxstatus *status, struct t_trxframe *fr, int nind,
89 const atom_id *ind, gmx_conect gc);
90 /* Write an indexed frame to a TRX file, see write_trxframe. gc may be NULL */
92 int write_trxframe(t_trxstatus *status, struct t_trxframe *fr, gmx_conect gc);
93 /* Write a frame to a TRX file.
94 * Only entries for which the gmx_boolean is TRUE will be written,
95 * except for step, time, lambda and/or box, which may not be
96 * omitted for certain trajectory formats.
97 * The precision for .xtc and .gro is fr->prec, when fr->bPrec=FALSE,
98 * the precision is set to 1000.
99 * gc is important for pdb file writing only and may be NULL.
102 int write_trx(t_trxstatus *status, int nind, const atom_id *ind, struct t_atoms *atoms,
103 int step, real time, matrix box, rvec x[], rvec *v,
105 /* Write an indexed frame to a TRX file.
107 * atoms can be NULL for file types which don't need atom names.
110 void trjtools_gmx_prepare_tng_writing(const char *filename,
116 const struct gmx_mtop_t *mtop,
117 const atom_id *index,
118 const char *index_group_name);
119 /* Sets up *out for writing TNG. If *in != NULL and contains a TNG trajectory
120 * some data, e.g. molecule system, will be copied over from *in to *out.
121 * If *in == NULL a file name (infile) of a TNG file can be provided instead
122 * and used for copying data to *out.
123 * If there is no TNG input natoms is used to create "implicit atoms" (no atom
124 * or molecular data present). If natoms == -1 the number of atoms are
125 * not known (or there is already a TNG molecule system to copy, in which case
126 * natoms is not required anyhow). If an group of indexed atoms are written
127 * natoms must be the length of index. index_group_name is the name of the
131 /*! \brief Write a trxframe to the TNG file in status.
133 * This function is needed because both t_trxstatus and
134 * tng_trajectory_t are encapsulated, so client trajectory-writing
135 * code with a t_trxstatus can't just call the TNG writing
137 void write_tng_frame(t_trxstatus *status,
138 struct t_trxframe *fr);
140 void close_trx(t_trxstatus *status);
141 /* Close trj file as opened with read_first_x, read_frist_frame
142 * or open_trx. Identical to close_trj.
145 t_trxstatus *open_trx(const char *outfile, const char *filemode);
146 /* Open a TRX file and return an allocated status pointer */
148 t_fileio *trx_get_fileio(t_trxstatus *status);
149 /* get a fileio from a trxstatus */
151 float trx_get_time_of_final_frame(t_trxstatus *status);
152 /* get time of final frame. Only supported for TNG and XTC */
154 gmx_bool bRmod_fd(double a, double b, double c, gmx_bool bDouble);
155 /* Returns TRUE when (a - b) MOD c = 0, using a margin which is slightly
156 * larger than the float/double precision.
160 #define bRmod(a, b, c) bRmod_fd(a, b, c, TRUE)
162 #define bRmod(a, b, c) bRmod_fd(a, b, c, FALSE)
165 int check_times2(real t, real t0, gmx_bool bDouble);
166 /* This routine checkes if the read-in time is correct or not;
167 * returns -1 if t<tbegin or t MOD dt = t0,
168 * 0 if tbegin <= t <=tend+margin,
170 * where margin is 0.1*min(t-tp,tp-tpp), if this positive, 0 otherwise.
171 * tp and tpp should be the time of the previous frame and the one before.
172 * The mod is done with single or double precision accuracy depending
173 * on the value of bDouble.
176 int check_times(real t);
177 /* This routine checkes if the read-in time is correct or not;
178 * returns -1 if t<tbegin,
179 * 0 if tbegin <= t <=tend,
187 /* For trxframe.flags, used in trxframe read routines.
188 * When a READ flag is set, the field will be read when present,
189 * but a frame might be returned which does not contain the field.
190 * When a NEED flag is set, frames not containing the field will be skipped.
192 #define TRX_READ_X (1<<0)
193 #define TRX_NEED_X (1<<1)
194 #define TRX_READ_V (1<<2)
195 #define TRX_NEED_V (1<<3)
196 #define TRX_READ_F (1<<4)
197 #define TRX_NEED_F (1<<5)
198 /* Useful for reading natoms from a trajectory without skipping */
199 #define TRX_DONT_SKIP (1<<6)
201 /* For trxframe.not_ok */
202 #define HEADER_NOT_OK (1<<0)
203 #define DATA_NOT_OK (1<<1)
204 #define FRAME_NOT_OK (HEADER_NOT_OK | DATA_NOT_OK)
206 int read_first_frame(const output_env_t oenv, t_trxstatus **status,
207 const char *fn, struct t_trxframe *fr, int flags);
208 /* Read the first frame which is in accordance with flags, which are
209 * defined further up in this file.
210 * Returns natoms when succeeded, 0 otherwise.
211 * Memory will be allocated for flagged entries.
212 * The flags are copied to fr for subsequent calls to read_next_frame.
213 * Returns TRUE when succeeded, FALSE otherwise.
216 gmx_bool read_next_frame(const output_env_t oenv, t_trxstatus *status,
217 struct t_trxframe *fr);
218 /* Reads the next frame which is in accordance with fr->flags.
219 * Returns TRUE when succeeded, FALSE otherwise.
222 int read_first_x(const output_env_t oenv, t_trxstatus **status,
223 const char *fn, real *t, rvec **x, matrix box);
224 /* These routines read first coordinates and box, and allocates
225 * memory for the coordinates, for a trajectory file.
226 * The routine returns the number of atoms, or 0 when something is wrong.
227 * The integer in status should be passed to calls of read_next_x
230 gmx_bool read_next_x(const output_env_t oenv, t_trxstatus *status, real *t, rvec x[], matrix box);
231 /* Read coordinates and box from a trajectory file. Return TRUE when all well,
232 * or FALSE when end of file (or last frame requested by user).
233 * status is the integer set in read_first_x.
236 void close_trj(t_trxstatus *status);
237 /* Close trj file as opened with read_first_x, read_frist_frame
238 * or open_trx. Identical to close_trx.
241 void rewind_trj(t_trxstatus *status);
242 /* Rewind trj file as opened with read_first_x */
244 struct t_topology *read_top(const char *fn, int *ePBC);
245 /* Extract a topology data structure from a topology file.
246 * If ePBC!=NULL *ePBC gives the pbc type.