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