2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018,2019,2020, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "trajectoryframe.h"
43 #include "gromacs/math/veccompare.h"
44 #include "gromacs/topology/atoms.h"
45 #include "gromacs/utility/compare.h"
46 #include "gromacs/utility/exceptions.h"
47 #include "gromacs/utility/smalloc.h"
48 #include "gromacs/utility/stringutil.h"
50 void comp_frame(FILE* fp, t_trxframe* fr1, t_trxframe* fr2, gmx_bool bRMSD, real ftol, real abstol)
53 cmp_int(fp, "not_ok", -1, fr1->not_ok, fr2->not_ok);
54 cmp_int(fp, "natoms", -1, fr1->natoms, fr2->natoms);
55 if (cmp_bool(fp, "bStep", -1, fr1->bStep, fr2->bStep))
57 cmp_int(fp, "step", -1, fr1->step, fr2->step);
59 cmp_int(fp, "step", -1, fr1->step, fr2->step);
60 if (cmp_bool(fp, "bTime", -1, fr1->bTime, fr2->bTime))
62 cmp_real(fp, "time", -1, fr1->time, fr2->time, ftol, abstol);
64 if (cmp_bool(fp, "bLambda", -1, fr1->bLambda, fr2->bLambda))
66 cmp_real(fp, "lambda", -1, fr1->lambda, fr2->lambda, ftol, abstol);
68 if (cmp_bool(fp, "bAtoms", -1, fr1->bAtoms, fr2->bAtoms))
70 compareAtoms(fp, fr1->atoms, fr2->atoms, ftol, abstol);
72 if (cmp_bool(fp, "bPrec", -1, fr1->bPrec, fr2->bPrec))
74 cmp_real(fp, "prec", -1, fr1->prec, fr2->prec, ftol, abstol);
76 if (cmp_bool(fp, "bX", -1, fr1->bX, fr2->bX))
78 cmp_rvecs(fp, "x", std::min(fr1->natoms, fr2->natoms), fr1->x, fr2->x, bRMSD, ftol, abstol);
80 if (cmp_bool(fp, "bV", -1, fr1->bV, fr2->bV))
82 cmp_rvecs(fp, "v", std::min(fr1->natoms, fr2->natoms), fr1->v, fr2->v, bRMSD, ftol, abstol);
84 if (cmp_bool(fp, "bF", -1, fr1->bF, fr2->bF))
86 cmp_rvecs(fp, "f", std::min(fr1->natoms, fr2->natoms), fr1->f, fr2->f, bRMSD, ftol, abstol);
88 if (cmp_bool(fp, "bBox", -1, fr1->bBox, fr2->bBox))
90 cmp_rvecs(fp, "box", 3, fr1->box, fr2->box, FALSE, ftol, abstol);
94 void done_frame(t_trxframe* frame)
98 done_atom(frame->atoms);
109 TrajectoryFrame::TrajectoryFrame(const t_trxframe& frame) : frame_(frame)
111 // This would be nicer as an initializer, but once uncrustify is
112 // happy, Doxygen can't parse it.
113 box_ = { { { { 0 } } } };
117 GMX_THROW(APIError("Cannot handle trajectory frame that lacks a step number"));
121 GMX_THROW(APIError("Cannot handle trajectory frame that lacks a time"));
125 for (int d = 0; d < DIM; ++d)
127 for (int dd = 0; dd < DIM; ++dd)
129 box_[d][dd] = frame.box[d][dd];
135 std::string TrajectoryFrame::frameName() const
137 return formatString("Time %f Step %" PRId64, frame_.time, frame_.step);
140 std::int64_t TrajectoryFrame::step() const
145 double TrajectoryFrame::time() const
150 PbcType TrajectoryFrame::pbc() const
152 return frame_.pbcType;
155 ArrayRef<const RVec> TrajectoryFrame::x() const
159 return arrayRefFromArray(reinterpret_cast<RVec*>(frame_.x), frame_.natoms);
167 ArrayRef<const RVec> TrajectoryFrame::v() const
171 return arrayRefFromArray(reinterpret_cast<RVec*>(frame_.v), frame_.natoms);
179 ArrayRef<const RVec> TrajectoryFrame::f() const
183 return arrayRefFromArray(reinterpret_cast<RVec*>(frame_.f), frame_.natoms);
191 bool TrajectoryFrame::hasBox() const
196 const BoxMatrix& TrajectoryFrame::box() const