f0e77d69ea11a021b2c876c930a7461dd6ba9223
[alexxy/gromacs.git] / src / gromacs / trajectory / trajectoryframe.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 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 /* The gmx_bools indicate whether a field was read from the trajectory.
40  * Do not try to use a pointer when its gmx_bool is FALSE, as memory might
41  * not be allocated.
42  */
43 #ifndef GMX_TRAJECTORY_TRX_H
44 #define GMX_TRAJECTORY_TRX_H
45
46 #include <cstdio>
47
48 #include <array>
49
50 #include "gromacs/math/vectypes.h"
51 #include "gromacs/utility/arrayref.h"
52 #include "gromacs/utility/basedefinitions.h"
53 #include "gromacs/utility/real.h"
54
55 struct t_atoms;
56
57 typedef struct t_trxframe // NOLINT (clang-analyzer-optin.performance.Padding)
58 {
59     int      not_ok;  /* integrity flags                  */
60     gmx_bool bDouble; /* Double precision?                */
61     int      natoms;  /* number of atoms (atoms, x, v, f, index) */
62     gmx_bool bStep;
63     int64_t  step; /* MD step number                   */
64     gmx_bool bTime;
65     real     time; /* time of the frame                */
66     gmx_bool bLambda;
67     gmx_bool bFepState; /* does it contain fep_state?       */
68     real     lambda;    /* free energy perturbation lambda  */
69     int      fep_state; /* which fep state are we in? */
70     gmx_bool bAtoms;
71     t_atoms* atoms; /* atoms struct (natoms)            */
72     gmx_bool bPrec;
73     real     prec; /* precision of x, fraction of 1 nm */
74     gmx_bool bX;
75     rvec*    x; /* coordinates (natoms)             */
76     gmx_bool bV;
77     rvec*    v; /* velocities (natoms)              */
78     gmx_bool bF;
79     rvec*    f; /* forces (natoms)                  */
80     gmx_bool bBox;
81     matrix   box; /* the 3 box vectors                */
82     gmx_bool bPBC;
83     int      ePBC; /* the type of pbc                  */
84     gmx_bool bIndex;
85     int*     index; /* atom indices of contained coordinates */
86 } t_trxframe;
87
88 void comp_frame(FILE* fp, t_trxframe* fr1, t_trxframe* fr2, gmx_bool bRMSD, real ftol, real abstol);
89
90 void done_frame(t_trxframe* frame);
91
92 namespace gmx
93 {
94
95 /*!\brief A 3x3 matrix data type useful for simulation boxes
96  *
97  * \todo Implement a full replacement for C-style real[DIM][DIM] */
98 using BoxMatrix = std::array<std::array<real, DIM>, DIM>;
99
100 /*! \internal
101  * \brief Contains a valid trajectory frame.
102  *
103  * Valid frames have a step and time, but need not have any particular
104  * other fields.
105  *
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.
110  *
111  * \todo Consider a std::optional work-alike type for expressing that
112  * a field may or may not have content. */
113 class TrajectoryFrame
114 {
115 public:
116     /*! \brief Constructor
117      *
118      * \throws APIError If \c frame lacks either step or time.
119      */
120     explicit TrajectoryFrame(const t_trxframe& frame);
121     /*! \brief Return a string that helps users identify this frame, containing time and step number.
122      *
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
142 private:
143     //! Handle to trajectory data
144     const t_trxframe& frame_;
145     //! Box matrix data from the frame_.
146     BoxMatrix box_;
147 };
148
149 } // namespace gmx
150
151 #endif