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 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 struct 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 t_trxstatus* open_trx(const char* outfile, const char* filemode); 160 /* Open a TRX file and return an allocated status pointer */ 161 162 struct t_fileio* trx_get_fileio(t_trxstatus* status); 163 /* get a fileio from a trxstatus */ 164 165 float trx_get_time_of_final_frame(t_trxstatus* status); 166 /* get time of final frame. Only supported for TNG and XTC */ 167 168 gmx_bool bRmod_fd(double a, double b, double c, gmx_bool bDouble); 169 /* Returns TRUE when (a - b) MOD c = 0, using a margin which is slightly 170 * larger than the float/double precision. 171 */ 172 173 #if GMX_DOUBLE 174 # define bRmod(a, b, c) bRmod_fd(a, b, c, TRUE) 175 #else 176 # define bRmod(a, b, c) bRmod_fd(a, b, c, FALSE) 177 #endif 178 179 int check_times2(real t, real t0, gmx_bool bDouble); 180 /* This routine checkes if the read-in time is correct or not; 181 * returns -1 if t<tbegin or t MOD dt = t0, 182 * 0 if tbegin <= t <=tend+margin, 183 * 1 if t>tend 184 * where margin is 0.1*min(t-tp,tp-tpp), if this positive, 0 otherwise. 185 * tp and tpp should be the time of the previous frame and the one before. 186 * The mod is done with single or double precision accuracy depending 187 * on the value of bDouble. 188 */ 189 190 int check_times(real t); 191 /* This routine checkes if the read-in time is correct or not; 192 * returns -1 if t<tbegin, 193 * 0 if tbegin <= t <=tend, 194 * 1 if t>tend 195 */ 196 197 198 /* For trxframe.flags, used in trxframe read routines. 199 * When a READ flag is set, the field will be read when present, 200 * but a frame might be returned which does not contain the field. 201 * When a NEED flag is set, frames not containing the field will be skipped. 202 */ 203 #define TRX_READ_X (1u << 0u) 204 #define TRX_NEED_X (1u << 1u) 205 #define TRX_READ_V (1u << 2u) 206 #define TRX_NEED_V (1u << 3u) 207 #define TRX_READ_F (1u << 4u) 208 #define TRX_NEED_F (1u << 5u) 209 /* Useful for reading natoms from a trajectory without skipping */ 210 #define TRX_DONT_SKIP (1u << 6u) 211 212 /* For trxframe.not_ok */ 213 #define HEADER_NOT_OK (1u << 0u) 214 #define DATA_NOT_OK (1u << 1u) 215 #define FRAME_NOT_OK (HEADER_NOT_OK | DATA_NOT_OK) 216 217 bool read_first_frame(const gmx_output_env_t* oenv, 218 t_trxstatus** status, 219 const char* fn, 220 struct t_trxframe* fr, 221 int flags); 222 /* Read the first frame which is in accordance with flags, which are 223 * defined further up in this file. 224 * Memory will be allocated for flagged entries. 225 * The flags are copied to fr for subsequent calls to read_next_frame. 226 * Returns true when succeeded, false otherwise. 227 */ 228 229 bool read_next_frame(const gmx_output_env_t* oenv, t_trxstatus* status, struct t_trxframe* fr); 230 /* Reads the next frame which is in accordance with fr->flags. 231 * Returns true when succeeded, false otherwise. 232 */ 233 234 int read_first_x(const gmx_output_env_t* oenv, t_trxstatus** status, const char* fn, real* t, rvec** x, matrix box); 235 /* These routines read first coordinates and box, and allocates 236 * memory for the coordinates, for a trajectory file. 237 * The routine returns the number of atoms, or 0 when something is wrong. 238 * The integer in status should be passed to calls of read_next_x 239 */ 240 241 gmx_bool read_next_x(const gmx_output_env_t* oenv, t_trxstatus* status, real* t, rvec x[], matrix box); 242 /* Read coordinates and box from a trajectory file. Return TRUE when all well, 243 * or FALSE when end of file (or last frame requested by user). 244 * status is the integer set in read_first_x. 245 */ 246 247 void rewind_trj(t_trxstatus* status); 248 /* Rewind trajectory file as opened with read_first_x */ 249 250 struct t_topology* read_top(const char* fn, PbcType* pbcType); 251 /* Extract a topology data structure from a topology file. 252 * If pbcType!=NULL *pbcType gives the pbc type. 253 */ 254 255 #endif 256