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