SYCL: Avoid using no_init read accessor in rocFFT
[alexxy/gromacs.git] / src / gromacs / fileio / tngio.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
5  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 #ifndef GMX_FILEIO_TNGIO_H
38 #define GMX_FILEIO_TNGIO_H
39
40 #include <cstdio>
41
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/utility/basedefinitions.h"
44 #include "gromacs/utility/real.h"
45
46 struct gmx_mtop_t;
47 struct t_inputrec;
48 struct gmx_tng_trajectory;
49 typedef struct gmx_tng_trajectory* gmx_tng_trajectory_t;
50 struct t_trxframe;
51
52 namespace gmx
53 {
54 template<typename>
55 class ArrayRef;
56 }
57 /*! \brief Open a TNG trajectory file
58  *
59  * \param filename   Name of file to open
60  * \param mode       Can be set to 'r', 'w' or 'a' for reading, writing or appending respectively.
61  * \param tng_data_p Pointer to an allocated gmx_tng_trajectory_t into which a handle to a TNG trajectory will be stored.
62  *
63  * Handles all I/O errors internally via fatal error
64  */
65 void gmx_tng_open(const char* filename, char mode, gmx_tng_trajectory_t* tng_data_p);
66
67 /*! \brief Finish writing a TNG trajectory file */
68 void gmx_tng_close(gmx_tng_trajectory_t* tng);
69
70 /*!\brief Add molecular topology information to TNG output (if
71  * available)
72  *
73  * \param tng   Valid handle to a TNG trajectory
74  * \param mtop  Pointer to a topology (can be NULL)
75  */
76 void gmx_tng_add_mtop(gmx_tng_trajectory_t tng, const gmx_mtop_t* mtop);
77
78 /*! \brief Do all TNG preparation for full-precision whole-system
79  * trajectory writing during MD simulations.
80  *
81  * \param tng   Valid handle to a TNG trajectory
82  * \param mtop  Global topology
83  * \param ir    Input settings (for writing frequencies)
84  */
85 void gmx_tng_prepare_md_writing(gmx_tng_trajectory_t tng, const gmx_mtop_t* mtop, const t_inputrec* ir);
86
87 /*! \brief Set the default compression precision for TNG writing
88  *
89  * \param tng   Valid handle to a TNG trajectory
90  * \param prec  GROMACS-style precision setting (i.e. 1000 for 3 digits of precision) */
91 void gmx_tng_set_compression_precision(gmx_tng_trajectory_t tng, real prec);
92
93 /*! \brief Do all TNG preparation for low-precision selection-based
94  * trajectory writing during MD simulations.
95  *
96  * \param tng   Valid handle to a TNG trajectory
97  * \param mtop  Global topology
98  * \param ir    Input settings (for writing frequencies)
99  */
100 void gmx_tng_prepare_low_prec_writing(gmx_tng_trajectory_t tng, const gmx_mtop_t* mtop, const t_inputrec* ir);
101
102 /*! \brief Write a frame to a TNG file
103  *
104  * \param tng                  Valid handle to a TNG trajectory
105  * \param bUseLossyCompression Whether to use lossy compression
106  * \param step                 MD step number
107  * \param elapsedPicoSeconds   Elapsed MD time
108  * \param lambda               Free-energy lambda value
109  * \param box                  Simulation box
110  * \param nAtoms               Number of atoms (i.e. vector lengths)
111  * \param x                    Vector of position coordinates
112  * \param v                    Vector of elocities
113  * \param f                    Vector of forces
114  *
115  * The pointers tng, x, v, f may be NULL, which triggers not writing
116  * (that component). box can only be NULL if x is also NULL. */
117 void gmx_fwrite_tng(gmx_tng_trajectory_t tng,
118                     gmx_bool             bUseLossyCompression,
119                     int64_t              step,
120                     real                 elapsedPicoSeconds,
121                     real                 lambda,
122                     const rvec*          box,
123                     int                  nAtoms,
124                     const rvec*          x,
125                     const rvec*          v,
126                     const rvec*          f);
127
128 /*! \brief Write the current frame set to disk. Perform compression
129  * etc.
130  *
131  * \param tng Valid handle to a TNG trajectory
132  */
133 void fflush_tng(gmx_tng_trajectory_t tng);
134
135 /*! \brief Get the time (in picoseconds) of the final frame in the
136  * trajectory.
137  *
138  * \param tng Valid handle to a TNG trajectory
139  */
140 float gmx_tng_get_time_of_final_frame(gmx_tng_trajectory_t tng);
141
142 /*! \brief Prepare to write TNG output from trajectory conversion tools */
143 void gmx_prepare_tng_writing(const char*              filename,
144                              char                     mode,
145                              gmx_tng_trajectory_t*    in,
146                              gmx_tng_trajectory_t*    out,
147                              int                      nAtoms,
148                              const struct gmx_mtop_t* mtop,
149                              gmx::ArrayRef<const int> index,
150                              const char*              indexGroupName);
151
152 /*! \brief Write a trxframe to a TNG file
153  *
154  * \param output Trajectory to write to
155  * \param frame  Frame data to write
156  * \param natoms Number of atoms to actually write
157  *
158  * The natoms field in frame is the number of atoms in the system. The
159  * parameter natoms supports writing an index-group subset of the
160  * atoms.
161  */
162 void gmx_write_tng_from_trxframe(gmx_tng_trajectory_t output, const t_trxframe* frame, int natoms);
163
164 /*! \brief Creates a molecule containing only the indexed atoms and sets
165  * the number of all other molecules to 0. Works similar to a
166  * selection group. */
167 void gmx_tng_setup_atom_subgroup(gmx_tng_trajectory_t tng, gmx::ArrayRef<const int> ind, const char* name);
168
169 /*! \brief Read the first/next TNG frame. */
170 gmx_bool gmx_read_next_tng_frame(gmx_tng_trajectory_t input,
171                                  struct t_trxframe*   fr,
172                                  int64_t*             requestedIds,
173                                  int                  numRequestedIds);
174
175 /*! \brief Print the molecule system to stream */
176 void gmx_print_tng_molecule_system(gmx_tng_trajectory_t input, FILE* stream);
177
178 /*! \brief Get a list of block IDs present in the next frame with data. */
179 gmx_bool gmx_get_tng_data_block_types_of_next_frame(gmx_tng_trajectory_t input,
180                                                     int                  frame,
181                                                     int                  nRequestedIds,
182                                                     int64_t*             requestedIds,
183                                                     int64_t*             nextFrame,
184                                                     int64_t*             nBlocks,
185                                                     int64_t**            blockIds);
186
187 /*! \brief Get data of the next frame with data from the data block
188  * with the specified block ID. */
189 gmx_bool gmx_get_tng_data_next_frame_of_block_type(gmx_tng_trajectory_t input,
190                                                    int64_t              blockId,
191                                                    real**               values,
192                                                    int64_t*             frameNumber,
193                                                    double*              frameTime,
194                                                    int64_t*             nValuesPerFrame,
195                                                    int64_t*             nAtoms,
196                                                    real*                prec,
197                                                    char*                name,
198                                                    int                  maxLen,
199                                                    gmx_bool*            bOK);
200
201 /*! \brief Get the output interval of box size.
202  *
203  * \return The box output interval, or -1 when TNG support is not available. */
204 int gmx_tng_get_box_output_interval(gmx_tng_trajectory_t gmx_tng);
205
206 /*! \brief Get the output interval of lambda.
207  *
208  * \return The box output interval, or -1 when TNG support is not available. */
209 int gmx_tng_get_lambda_output_interval(gmx_tng_trajectory_t gmx_tng);
210
211 #endif /* GMX_FILEIO_TNGIO_H */