Import cmake Modules/FindCUDA.cmake
[gromacs.git] / src / gromacs / trajectory / trajectoryframe.h
blob3901eaaef71daabab692c0447ca2c539cf6ed9ea
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
38 /* The gmx_bools indicate whether a field was read from the trajectory.
39 * Do not try to use a pointer when its gmx_bool is FALSE, as memory might
40 * not be allocated.
42 #ifndef GMX_TRAJECTORY_TRX_H
43 #define GMX_TRAJECTORY_TRX_H
45 #include <cstdio>
47 #include "gromacs/math/vectypes.h"
48 #include "gromacs/utility/basedefinitions.h"
49 #include "gromacs/utility/real.h"
51 struct t_atoms;
53 typedef struct t_trxframe
55 int not_ok; /* integrity flags */
56 gmx_bool bDouble; /* Double precision? */
57 int natoms; /* number of atoms (atoms, x, v, f, index) */
58 gmx_bool bStep;
59 gmx_int64_t step; /* MD step number */
60 gmx_bool bTime;
61 real time; /* time of the frame */
62 gmx_bool bLambda;
63 gmx_bool bFepState; /* does it contain fep_state? */
64 real lambda; /* free energy perturbation lambda */
65 int fep_state; /* which fep state are we in? */
66 gmx_bool bAtoms;
67 t_atoms *atoms; /* atoms struct (natoms) */
68 gmx_bool bPrec;
69 real prec; /* precision of x, fraction of 1 nm */
70 gmx_bool bX;
71 rvec *x; /* coordinates (natoms) */
72 gmx_bool bV;
73 rvec *v; /* velocities (natoms) */
74 gmx_bool bF;
75 rvec *f; /* forces (natoms) */
76 gmx_bool bBox;
77 matrix box; /* the 3 box vectors */
78 gmx_bool bPBC;
79 int ePBC; /* the type of pbc */
80 gmx_bool bIndex;
81 int *index; /* atom indices of contained coordinates */
82 } t_trxframe;
84 void comp_frame(FILE *fp, t_trxframe *fr1, t_trxframe *fr2,
85 gmx_bool bRMSD, real ftol, real abstol);
87 void done_frame(t_trxframe *frame);
89 #endif