Uncrustify template.cpp and selection.cpp.
[alexxy/gromacs.git] / src / gromacs / selection / selection.cpp
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11  * Copyright (c) 2001-2009, The GROMACS development team,
12  * check out http://www.gromacs.org for more information.
13
14  * This program is free software; you can redistribute it and/or
15  * modify it under the terms of the GNU General Public License
16  * as published by the Free Software Foundation; either version 2
17  * of the License, or (at your option) any later version.
18  *
19  * If you want to redistribute modifications, please consider that
20  * scientific software is very special. Version control is crucial -
21  * bugs must be traceable. We will be happy to consider code for
22  * inclusion in the official distribution, but derived work must not
23  * be called official GROMACS. Details are found in the README & COPYING
24  * files - if they are missing, get the official version at www.gromacs.org.
25  *
26  * To help us fund GROMACS development, we humbly ask that you cite
27  * the papers on the package - you can find them in the top README file.
28  *
29  * For more info, check our website at http://www.gromacs.org
30  */
31 /*! \internal \file
32  * \brief
33  * Implements classes in selection.h.
34  *
35  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36  * \ingroup module_selection
37  */
38 #include "selection.h"
39
40 #include "position.h"
41 #include "selelem.h"
42 #include "selvalue.h"
43
44 namespace gmx
45 {
46
47 namespace internal
48 {
49
50 SelectionData::SelectionData(SelectionTreeElement *elem,
51                              const char           *selstr)
52     : name_(elem->name()), selectionText_(selstr),
53       rootElement_(*elem), coveredFractionType_(CFRAC_NONE),
54       coveredFraction_(1.0), averageCoveredFraction_(1.0),
55       bDynamic_(false), bDynamicCoveredFraction_(false)
56 {
57     gmx_ana_pos_clear(&rawPositions_);
58
59     if (elem->child->type == SEL_CONST)
60     {
61         // TODO: This is not exception-safe if any called function throws.
62         gmx_ana_pos_copy(&rawPositions_, elem->child->v.u.p, true);
63     }
64     else
65     {
66         SelectionTreeElementPointer child = elem->child;
67         child->flags     &= ~SEL_ALLOCVAL;
68         _gmx_selvalue_setstore(&child->v, &rawPositions_);
69         /* We should also skip any modifiers to determine the dynamic
70          * status. */
71         while (child->type == SEL_MODIFIER)
72         {
73             child = child->child;
74             if (child->type == SEL_SUBEXPRREF)
75             {
76                 child = child->child;
77                 /* Because most subexpression elements are created
78                  * during compilation, we need to check for them
79                  * explicitly here.
80                  */
81                 if (child->type == SEL_SUBEXPR)
82                 {
83                     child = child->child;
84                 }
85             }
86         }
87         /* For variable references, we should skip the
88          * SEL_SUBEXPRREF and SEL_SUBEXPR elements. */
89         if (child->type == SEL_SUBEXPRREF)
90         {
91             child = child->child->child;
92         }
93         bDynamic_ = (child->child->flags & SEL_DYNAMIC);
94     }
95     initCoveredFraction(CFRAC_NONE);
96 }
97
98
99 SelectionData::~SelectionData()
100 {
101     gmx_ana_pos_deinit(&rawPositions_);
102 }
103
104
105 bool
106 SelectionData::initCoveredFraction(e_coverfrac_t type)
107 {
108     coveredFractionType_ = type;
109     if (type == CFRAC_NONE)
110     {
111         bDynamicCoveredFraction_ = false;
112     }
113     else if (!_gmx_selelem_can_estimate_cover(rootElement()))
114     {
115         coveredFractionType_     = CFRAC_NONE;
116         bDynamicCoveredFraction_ = false;
117     }
118     else
119     {
120         bDynamicCoveredFraction_ = true;
121     }
122     coveredFraction_        = bDynamicCoveredFraction_ ? 0.0 : 1.0;
123     averageCoveredFraction_ = coveredFraction_;
124     return type == CFRAC_NONE || coveredFractionType_ != CFRAC_NONE;
125 }
126
127 namespace
128 {
129
130 /*! \brief
131  * Helper function to compute total masses and charges for positions.
132  *
133  * \param[in]  top     Topology to take atom masses from.
134  * \param[in]  pos     Positions to compute masses and charges for.
135  * \param[out] masses  Output masses.
136  * \param[out] charges Output charges.
137  *
138  * Does not throw if enough space has been reserved for the output vectors.
139  */
140 void computeMassesAndCharges(const t_topology *top, const gmx_ana_pos_t &pos,
141                              std::vector<real> *masses,
142                              std::vector<real> *charges)
143 {
144     GMX_ASSERT(top != NULL, "Should not have been called with NULL topology");
145     masses->clear();
146     charges->clear();
147     for (int b = 0; b < pos.nr; ++b)
148     {
149         real mass   = 0.0;
150         real charge = 0.0;
151         for (int i = pos.m.mapb.index[b]; i < pos.m.mapb.index[b+1]; ++i)
152         {
153             int index = pos.g->index[i];
154             mass   += top->atoms.atom[index].m;
155             charge += top->atoms.atom[index].q;
156         }
157         masses->push_back(mass);
158         charges->push_back(charge);
159     }
160 }
161
162 }       // namespace
163
164 void
165 SelectionData::initializeMassesAndCharges(const t_topology *top)
166 {
167     GMX_ASSERT(posMass_.empty() && posCharge_.empty(),
168                "Should not be called more than once");
169     posMass_.reserve(posCount());
170     posCharge_.reserve(posCount());
171     if (top == NULL)
172     {
173         posMass_.resize(posCount(), 1.0);
174         posCharge_.resize(posCount(), 0.0);
175     }
176     else
177     {
178         computeMassesAndCharges(top, rawPositions_, &posMass_, &posCharge_);
179     }
180 }
181
182
183 void
184 SelectionData::refreshMassesAndCharges(const t_topology *top)
185 {
186     if (top != NULL && isDynamic() && !hasFlag(efSelection_DynamicMask))
187     {
188         computeMassesAndCharges(top, rawPositions_, &posMass_, &posCharge_);
189     }
190 }
191
192
193 void
194 SelectionData::updateCoveredFractionForFrame()
195 {
196     if (isCoveredFractionDynamic())
197     {
198         real cfrac = _gmx_selelem_estimate_coverfrac(rootElement());
199         coveredFraction_         = cfrac;
200         averageCoveredFraction_ += cfrac;
201     }
202 }
203
204
205 void
206 SelectionData::computeAverageCoveredFraction(int nframes)
207 {
208     if (isCoveredFractionDynamic() && nframes > 0)
209     {
210         averageCoveredFraction_ /= nframes;
211     }
212 }
213
214
215 void
216 SelectionData::restoreOriginalPositions(const t_topology *top)
217 {
218     if (isDynamic())
219     {
220         gmx_ana_pos_t &p = rawPositions_;
221         gmx_ana_index_copy(p.g, rootElement().v.u.g, false);
222         p.g->name = NULL;
223         gmx_ana_indexmap_update(&p.m, p.g, hasFlag(gmx::efSelection_DynamicMask));
224         p.nr = p.m.nr;
225         refreshMassesAndCharges(top);
226     }
227 }
228
229 }   // namespace internal
230
231
232 void
233 Selection::printInfo(FILE *fp) const
234 {
235     fprintf(fp, "\"%s\" (%d position%s, %d atom%s%s)", name(),
236             posCount(),  posCount()  == 1 ? "" : "s",
237             atomCount(), atomCount() == 1 ? "" : "s",
238             isDynamic() ? ", dynamic" : "");
239     fprintf(fp, "\n");
240 }
241
242
243 void
244 Selection::printDebugInfo(FILE *fp, int nmaxind) const
245 {
246     const gmx_ana_pos_t &p = data().rawPositions_;
247
248     fprintf(fp, "  ");
249     printInfo(fp);
250     fprintf(fp, "    ");
251     gmx_ana_index_dump(fp, p.g, -1, nmaxind);
252
253     fprintf(fp, "    Block (size=%d):", p.m.mapb.nr);
254     if (!p.m.mapb.index)
255     {
256         fprintf(fp, " (null)");
257     }
258     else
259     {
260         int n = p.m.mapb.nr;
261         if (nmaxind >= 0 && n > nmaxind)
262         {
263             n = nmaxind;
264         }
265         for (int i = 0; i <= n; ++i)
266         {
267             fprintf(fp, " %d", p.m.mapb.index[i]);
268         }
269         if (n < p.m.mapb.nr)
270         {
271             fprintf(fp, " ...");
272         }
273     }
274     fprintf(fp, "\n");
275
276     int n = posCount();
277     if (nmaxind >= 0 && n > nmaxind)
278     {
279         n = nmaxind;
280     }
281     fprintf(fp, "    RefId:");
282     if (!p.m.refid)
283     {
284         fprintf(fp, " (null)");
285     }
286     else
287     {
288         for (int i = 0; i < n; ++i)
289         {
290             fprintf(fp, " %d", p.m.refid[i]);
291         }
292         if (n < posCount())
293         {
294             fprintf(fp, " ...");
295         }
296     }
297     fprintf(fp, "\n");
298
299     fprintf(fp, "    MapId:");
300     if (!p.m.mapid)
301     {
302         fprintf(fp, " (null)");
303     }
304     else
305     {
306         for (int i = 0; i < n; ++i)
307         {
308             fprintf(fp, " %d", p.m.mapid[i]);
309         }
310         if (n < posCount())
311         {
312             fprintf(fp, " ...");
313         }
314     }
315     fprintf(fp, "\n");
316 }
317
318 } // namespace gmx