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