68500a628b39fc64e50e8719643f786f13e138c8
[alexxy/gromacs.git] / src / gromacs / trajectory / trajectoryframe.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include "trajectoryframe.h"
38
39 #include <cstdio>
40
41 #include <algorithm>
42
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"
49
50 void comp_frame(FILE *fp, t_trxframe *fr1, t_trxframe *fr2,
51                 gmx_bool bRMSD, real ftol, real abstol)
52 {
53     fprintf(fp, "\n");
54     cmp_int(fp, "not_ok", -1, fr1->not_ok, fr2->not_ok);
55     cmp_int(fp, "natoms", -1, fr1->natoms, fr2->natoms);
56     if (cmp_bool(fp, "bStep", -1, fr1->bStep, fr2->bStep))
57     {
58         cmp_int(fp, "step", -1, fr1->step, fr2->step);
59     }
60     cmp_int(fp, "step", -1, fr1->step, fr2->step);
61     if (cmp_bool(fp, "bTime", -1, fr1->bTime, fr2->bTime))
62     {
63         cmp_real(fp, "time", -1, fr1->time, fr2->time, ftol, abstol);
64     }
65     if (cmp_bool(fp, "bLambda", -1, fr1->bLambda, fr2->bLambda))
66     {
67         cmp_real(fp, "lambda", -1, fr1->lambda, fr2->lambda, ftol, abstol);
68     }
69     if (cmp_bool(fp, "bAtoms", -1, fr1->bAtoms, fr2->bAtoms))
70     {
71         cmp_atoms(fp, fr1->atoms, fr2->atoms, ftol, abstol);
72     }
73     if (cmp_bool(fp, "bPrec", -1, fr1->bPrec, fr2->bPrec))
74     {
75         cmp_real(fp, "prec", -1, fr1->prec, fr2->prec, ftol, abstol);
76     }
77     if (cmp_bool(fp, "bX", -1, fr1->bX, fr2->bX))
78     {
79         cmp_rvecs(fp, "x", std::min(fr1->natoms, fr2->natoms), fr1->x, fr2->x, bRMSD, ftol, abstol);
80     }
81     if (cmp_bool(fp, "bV", -1, fr1->bV, fr2->bV))
82     {
83         cmp_rvecs(fp, "v", std::min(fr1->natoms, fr2->natoms), fr1->v, fr2->v, bRMSD, ftol, abstol);
84     }
85     if (cmp_bool(fp, "bF", -1, fr1->bF, fr2->bF))
86     {
87         cmp_rvecs(fp, "f", std::min(fr1->natoms, fr2->natoms), fr1->f, fr2->f, bRMSD, ftol, abstol);
88     }
89     if (cmp_bool(fp, "bBox", -1, fr1->bBox, fr2->bBox))
90     {
91         cmp_rvecs(fp, "box", 3, fr1->box, fr2->box, FALSE, ftol, abstol);
92     }
93 }
94
95 void done_frame(t_trxframe *frame)
96 {
97     if (frame->atoms)
98     {
99         done_atom(frame->atoms);
100         sfree(frame->atoms);
101     }
102     sfree(frame->x);
103     sfree(frame->v);
104     sfree(frame->f);
105 }
106
107 namespace gmx
108 {
109
110 TrajectoryFrame::TrajectoryFrame(const t_trxframe &frame)
111     : frame_(frame)
112 {
113     // This would be nicer as an initializer, but once uncrustify is
114     // happy, Doxygen can't parse it.
115     box_ = {{{{0}}}};
116
117     if (!frame.bStep)
118     {
119         GMX_THROW(APIError("Cannot handle trajectory frame that lacks a step number"));
120     }
121     if (!frame.bTime)
122     {
123         GMX_THROW(APIError("Cannot handle trajectory frame that lacks a time"));
124     }
125     if (frame.bBox)
126     {
127         for (int d = 0; d < DIM; ++d)
128         {
129             for (int dd = 0; dd < DIM; ++dd)
130             {
131                 box_[d][dd] = frame.box[d][dd];
132             }
133         }
134     }
135 }
136
137 std::string TrajectoryFrame::frameName() const
138 {
139     return formatString("Time %f Step %" GMX_PRId64, frame_.time, frame_.step);
140 }
141
142 std::int64_t TrajectoryFrame::step() const
143 {
144     return frame_.step;
145 }
146
147 double TrajectoryFrame::time() const
148 {
149     return frame_.time;
150 }
151
152 int TrajectoryFrame::pbc() const
153 {
154     return frame_.ePBC;
155 }
156
157 ArrayRef<const RVec> TrajectoryFrame::x() const
158 {
159     return arrayRefFromArray(reinterpret_cast<RVec *>(frame_.x),
160                              frame_.natoms);
161 }
162
163 ArrayRef<const RVec> TrajectoryFrame::v() const
164 {
165     if (frame_.bV)
166     {
167         return arrayRefFromArray(reinterpret_cast<RVec *>(frame_.v),
168                                  frame_.natoms);
169     }
170     else
171     {
172         return EmptyArrayRef();
173     }
174 }
175
176 ArrayRef<const RVec> TrajectoryFrame::f() const
177 {
178     if (frame_.bF)
179     {
180         return arrayRefFromArray(reinterpret_cast<RVec *>(frame_.f),
181                                  frame_.natoms);
182     }
183     else
184     {
185         return EmptyArrayRef();
186     }
187 }
188
189 bool TrajectoryFrame::hasBox() const
190 {
191     return frame_.bBox;
192 }
193
194 const BoxMatrix &TrajectoryFrame::box() const
195 {
196     return box_;
197 }
198
199 } // namespace gmx