Update energy and trajectory frame types
[gromacs.git] / src / gromacs / trajectory / trajectoryframe.h
blob04fbf601ca4ef2a8ba83c5915e5322c9c5c320d8
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,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.
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 <array>
49 #include "gromacs/math/vectypes.h"
50 #include "gromacs/utility/arrayref.h"
51 #include "gromacs/utility/basedefinitions.h"
52 #include "gromacs/utility/real.h"
54 struct t_atoms;
56 typedef struct t_trxframe
58 int not_ok; /* integrity flags */
59 gmx_bool bDouble; /* Double precision? */
60 int natoms; /* number of atoms (atoms, x, v, f, index) */
61 gmx_bool bStep;
62 gmx_int64_t step; /* MD step number */
63 gmx_bool bTime;
64 real time; /* time of the frame */
65 gmx_bool bLambda;
66 gmx_bool bFepState; /* does it contain fep_state? */
67 real lambda; /* free energy perturbation lambda */
68 int fep_state; /* which fep state are we in? */
69 gmx_bool bAtoms;
70 t_atoms *atoms; /* atoms struct (natoms) */
71 gmx_bool bPrec;
72 real prec; /* precision of x, fraction of 1 nm */
73 gmx_bool bX;
74 rvec *x; /* coordinates (natoms) */
75 gmx_bool bV;
76 rvec *v; /* velocities (natoms) */
77 gmx_bool bF;
78 rvec *f; /* forces (natoms) */
79 gmx_bool bBox;
80 matrix box; /* the 3 box vectors */
81 gmx_bool bPBC;
82 int ePBC; /* the type of pbc */
83 gmx_bool bIndex;
84 int *index; /* atom indices of contained coordinates */
85 } t_trxframe;
87 void comp_frame(FILE *fp, t_trxframe *fr1, t_trxframe *fr2,
88 gmx_bool bRMSD, real ftol, real abstol);
90 void done_frame(t_trxframe *frame);
92 namespace gmx
95 /*!\brief A 3x3 matrix data type useful for simulation boxes
97 * \todo Implement a full replacement for C-style real[DIM][DIM] */
98 using BoxMatrix = std::array <std::array<real, DIM>, DIM>;
100 /*! \internal
101 * \brief Contains a valid trajectory frame.
103 * Valid frames have a step and time, but need not have any particular
104 * other fields.
106 * \todo Eventually t_trxframe should be replaced by a class such as
107 * this. Currently we need to introduce BoxMatrix so that we can have
108 * a normal C++ getter that returns the contents of a box matrix,
109 * since you cannot use a real[DIM][DIM] as a function return type.
111 * \todo Consider a std::optional work-alike type for expressing that
112 * a field may or may not have content. */
113 class TrajectoryFrame
115 public:
116 /*! \brief Constructor
118 * \throws APIError If \c frame lacks either step or time.
120 explicit TrajectoryFrame(const t_trxframe &frame);
121 /*! \brief Return a string that helps users identify this frame, containing time and step number.
123 * \throws std::bad_alloc when out of memory */
124 std::string frameName() const;
125 //! Step number read from the trajectory file frame.
126 std::int64_t step() const;
127 //! Time read from the trajectory file frame.
128 double time() const;
129 //! The PBC characteristics of the box.
130 int pbc() const;
131 //! Get a view of position coordinates of the frame (which could be empty).
132 ArrayRef<const RVec> x() const;
133 //! Get a view of velocity coordinates of the frame (which could be empty).
134 ArrayRef<const RVec> v() const;
135 //! Get a view of force coordinates of the frame (which could be empty).
136 ArrayRef<const RVec> f() const;
137 //! Return whether the frame has a box.
138 bool hasBox() const;
139 //! Return a handle to the frame's box, which is all zero if the frame has no box.
140 const BoxMatrix &box() const;
141 // TODO make this private when updating trajectory comparison code
142 //! Handle to trajectory data
143 const t_trxframe &frame_;
144 private:
145 //! Box matrix data from the frame_.
146 BoxMatrix box_;
149 } // namespace gmx
151 #endif