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