2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014, 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 "gromacs/legacyheaders/typedefs.h"
40 #include "tng/tng_io_fwd.h"
49 /*! \brief Open a TNG trajectory file
51 * \param filename Name of file to open
52 * \param mode Can be set to 'r', 'w' or 'a' for reading, writing or appending respectively.
53 * \param tng_data_p Pointer to an allocated tng_trajectory_t into which a handle to a TNG trajectory will be stored.
55 * Handles all I/O errors internally via fatal error
57 void gmx_tng_open(const char *filename,
59 tng_trajectory_t *tng_data_p);
61 /*! \brief Finish writing a TNG trajectory file */
62 void gmx_tng_close(tng_trajectory_t *tng);
64 /*!\brief Add molecular topology information to TNG output (if
67 * \param tng Valid handle to a TNG trajectory
68 * \param mtop Pointer to a topology (can be NULL)
70 void gmx_tng_add_mtop(tng_trajectory_t tng,
71 const gmx_mtop_t *mtop);
73 /*! \brief Do all TNG preparation for full-precision whole-system
74 * trajectory writing during MD simulations.
76 * \param tng Valid handle to a TNG trajectory
77 * \param mtop Global topology
78 * \param ir Input settings (for writing frequencies)
80 void gmx_tng_prepare_md_writing(tng_trajectory_t tng,
81 const gmx_mtop_t *mtop,
82 const t_inputrec *ir);
84 /*! \brief Set the default compression precision for TNG writing
86 * \param tng Valid handle to a TNG trajectory
87 * \param prec GROMACS-style precision setting (i.e. 1000 for 3 digits of precision) */
88 void gmx_tng_set_compression_precision(tng_trajectory_t tng,
91 /*! \brief Do all TNG preparation for low-precision selection-based
92 * trajectory writing during MD simulations.
94 * \param tng Valid handle to a TNG trajectory
95 * \param mtop Global topology
96 * \param ir Input settings (for writing frequencies)
98 void gmx_tng_prepare_low_prec_writing(tng_trajectory_t tng,
99 const gmx_mtop_t *mtop,
100 const t_inputrec *ir);
102 /*! \brief Write a frame to a TNG file
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
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(tng_trajectory_t tng,
118 const gmx_bool bUseLossyCompression,
120 real elapsedPicoSeconds,
128 /*! \brief Write the current frame set to disk. Perform compression
131 * \param tng Valid handle to a TNG trajectory
133 void fflush_tng(tng_trajectory_t tng);
135 /*! \brief Get the time (in picoseconds) of the final frame in the
138 * \param tng Valid handle to a TNG trajectory
140 float gmx_tng_get_time_of_final_frame(tng_trajectory_t tng);
146 #endif /* GMX_FILEIO_TNGIO_H */