Compile nonbonded kernels as C++
[alexxy/gromacs.git] / src / 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,2018, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 #ifndef GMX_FILEIO_TRXIO_H
39 #define GMX_FILEIO_TRXIO_H
40
41 #include "gromacs/fileio/pdbio.h"
42
43 struct gmx_mtop_t;
44 struct gmx_output_env_t;
45 struct t_atoms;
46 struct t_fileio;
47 struct t_topology;
48 struct t_trxframe;
49
50 /* a dedicated status type contains fp, etc. */
51 typedef struct t_trxstatus t_trxstatus;
52
53 /* I/O function types */
54
55 /************************************************
56  *             Trajectory functions
57  ************************************************/
58
59 int prec2ndec(real prec);
60 /* Convert precision in 1/(nm) to number of decimal places */
61
62 /*! \brief Convert number of decimal places to trajectory precision in
63  * 1/(nm) */
64 real ndec2prec(int ndec);
65
66 void clear_trxframe(struct t_trxframe *fr, gmx_bool bFirst);
67 /* Set all content gmx_booleans to FALSE.
68  * When bFirst = TRUE, set natoms=-1, all pointers to NULL
69  *                     and all data to zero.
70  */
71
72 void set_trxframe_ePBC(struct t_trxframe *fr, int ePBC);
73 /* Set the type of periodic boundary conditions, ePBC=-1 is not set */
74
75 int nframes_read(t_trxstatus *status);
76 /* Returns the number of frames read from the trajectory */
77
78 int write_trxframe_indexed(t_trxstatus *status, const t_trxframe *fr, int nind,
79                            const int *ind, gmx_conect gc);
80 /* Write an indexed frame to a TRX file, see write_trxframe. gc may be NULL */
81
82 int write_trxframe(t_trxstatus *status, struct t_trxframe *fr, gmx_conect gc);
83 /* Write a frame to a TRX file.
84  * Only entries for which the gmx_boolean is TRUE will be written,
85  * except for step, time, lambda and/or box, which may not be
86  * omitted for certain trajectory formats.
87  * The precision for .xtc and .gro is fr->prec, when fr->bPrec=FALSE,
88  * the precision is set to 1000.
89  * gc is important for pdb file writing only and may be NULL.
90  */
91
92 int write_trx(t_trxstatus *status, int nind, const int *ind, const t_atoms *atoms,
93               int step, real time, matrix box, rvec x[], rvec *v,
94               gmx_conect gc);
95 /* Write an indexed frame to a TRX file.
96  * v can be NULL.
97  * atoms can be NULL for file types which don't need atom names.
98  */
99 t_trxstatus *
100 trjtools_gmx_prepare_tng_writing(const char               *filename,
101                                  char                      filemode,
102                                  t_trxstatus              *in,
103                                  const char               *infile,
104                                  int                       natoms,
105                                  const struct gmx_mtop_t  *mtop,
106                                  const int                *index,
107                                  const char               *index_group_name);
108 /* Sets up *out for writing TNG. If *in != NULL and contains a TNG trajectory
109  * some data, e.g. molecule system, will be copied over from *in to the return value.
110  * If *in == NULL a file name (infile) of a TNG file can be provided instead
111  * and used for copying data to the return value.
112  * If there is no TNG input natoms is used to create "implicit atoms" (no atom
113  * or molecular data present). If natoms == -1 the number of atoms are
114  * not known (or there is already a TNG molecule system to copy, in which case
115  * natoms is not required anyhow). If an group of indexed atoms are written
116  * natoms must be the length of index. index_group_name is the name of the
117  * index group.
118  */
119
120 /*! \brief Write a trxframe to the TNG file in status.
121  *
122  * This function is needed because both t_trxstatus and
123  * gmx_tng_trajectory_t are encapsulated, so client trajectory-writing
124  * code with a t_trxstatus can't just call the TNG writing
125  * function. */
126 void write_tng_frame(t_trxstatus        *status,
127                      struct t_trxframe  *fr);
128
129 void close_trx(t_trxstatus *status);
130 /* Close trajectory file as opened with read_first_x, read_first_frame
131  * or open_trx.
132  * Also frees memory in the structure.
133  */
134
135 t_trxstatus *open_trx(const char *outfile, const char *filemode);
136 /* Open a TRX file and return an allocated status pointer */
137
138 struct t_fileio *trx_get_fileio(t_trxstatus *status);
139 /* get a fileio from a trxstatus */
140
141 float trx_get_time_of_final_frame(t_trxstatus *status);
142 /* get time of final frame. Only supported for TNG and XTC */
143
144 gmx_bool bRmod_fd(double a, double b, double c, gmx_bool bDouble);
145 /* Returns TRUE when (a - b) MOD c = 0, using a margin which is slightly
146  * larger than the float/double precision.
147  */
148
149 #if GMX_DOUBLE
150 #define bRmod(a, b, c) bRmod_fd(a, b, c, TRUE)
151 #else
152 #define bRmod(a, b, c) bRmod_fd(a, b, c, FALSE)
153 #endif
154
155 int check_times2(real t, real t0, gmx_bool bDouble);
156 /* This routine checkes if the read-in time is correct or not;
157  * returns -1 if t<tbegin or t MOD dt = t0,
158  *          0 if tbegin <= t <=tend+margin,
159  *          1 if t>tend
160  * where margin is 0.1*min(t-tp,tp-tpp), if this positive, 0 otherwise.
161  * tp and tpp should be the time of the previous frame and the one before.
162  * The mod is done with single or double precision accuracy depending
163  * on the value of bDouble.
164  */
165
166 int check_times(real t);
167 /* This routine checkes if the read-in time is correct or not;
168  * returns -1 if t<tbegin,
169  *          0 if tbegin <= t <=tend,
170  *          1 if t>tend
171  */
172
173
174
175
176
177 /* For trxframe.flags, used in trxframe read routines.
178  * When a READ flag is set, the field will be read when present,
179  * but a frame might be returned which does not contain the field.
180  * When a NEED flag is set, frames not containing the field will be skipped.
181  */
182 #define TRX_READ_X    (1<<0)
183 #define TRX_NEED_X    (1<<1)
184 #define TRX_READ_V    (1<<2)
185 #define TRX_NEED_V    (1<<3)
186 #define TRX_READ_F    (1<<4)
187 #define TRX_NEED_F    (1<<5)
188 /* Useful for reading natoms from a trajectory without skipping */
189 #define TRX_DONT_SKIP (1<<6)
190
191 /* For trxframe.not_ok */
192 #define HEADER_NOT_OK (1<<0)
193 #define DATA_NOT_OK   (1<<1)
194 #define FRAME_NOT_OK  (HEADER_NOT_OK | DATA_NOT_OK)
195
196 bool read_first_frame(const gmx_output_env_t *oenv, t_trxstatus **status,
197                       const char *fn, struct t_trxframe *fr, int flags);
198 /* Read the first frame which is in accordance with flags, which are
199  * defined further up in this file.
200  * Memory will be allocated for flagged entries.
201  * The flags are copied to fr for subsequent calls to read_next_frame.
202  * Returns true when succeeded, false otherwise.
203  */
204
205 bool read_next_frame(const gmx_output_env_t *oenv, t_trxstatus *status,
206                      struct t_trxframe *fr);
207 /* Reads the next frame which is in accordance with fr->flags.
208  * Returns true when succeeded, false otherwise.
209  */
210
211 int read_first_x(const gmx_output_env_t *oenv, t_trxstatus **status,
212                  const char *fn, real *t, rvec **x, matrix box);
213 /* These routines read first coordinates and box, and allocates
214  * memory for the coordinates, for a trajectory file.
215  * The routine returns the number of atoms, or 0 when something is wrong.
216  * The integer in status should be passed to calls of read_next_x
217  */
218
219 gmx_bool read_next_x(const gmx_output_env_t *oenv, t_trxstatus *status, real *t, rvec x[], matrix box);
220 /* Read coordinates and box from a trajectory file. Return TRUE when all well,
221  * or FALSE when end of file (or last frame requested by user).
222  * status is the integer set in read_first_x.
223  */
224
225 void rewind_trj(t_trxstatus *status);
226 /* Rewind trajectory file as opened with read_first_x */
227
228 struct t_topology *read_top(const char *fn, int *ePBC);
229 /* Extract a topology data structure from a topology file.
230  * If ePBC!=NULL *ePBC gives the pbc type.
231  */
232
233 #endif