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