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