Use smart pointers to manage the selection tree.
[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
128 void
129 SelectionData::initializeMassesAndCharges(const t_topology *top)
130 {
131     posInfo_.reserve(posCount());
132     for (int b = 0; b < posCount(); ++b)
133     {
134         real mass   = 1.0;
135         real charge = 0.0;
136         if (top != NULL)
137         {
138             mass = 0.0;
139             for (int i = rawPositions_.m.mapb.index[b];
140                      i < rawPositions_.m.mapb.index[b+1];
141                      ++i)
142             {
143                 int index = rawPositions_.g->index[i];
144                 mass   += top->atoms.atom[index].m;
145                 charge += top->atoms.atom[index].q;
146             }
147         }
148         posInfo_.push_back(PositionInfo(mass, charge));
149     }
150     if (isDynamic() && !hasFlag(efSelection_DynamicMask))
151     {
152         originalPosInfo_ = posInfo_;
153     }
154 }
155
156
157 void
158 SelectionData::refreshMassesAndCharges()
159 {
160     if (!originalPosInfo_.empty())
161     {
162         posInfo_.clear();
163         for (int i = 0; i < posCount(); ++i)
164         {
165             int refid  = rawPositions_.m.refid[i];
166             posInfo_.push_back(originalPosInfo_[refid]);
167         }
168     }
169 }
170
171
172 void
173 SelectionData::updateCoveredFractionForFrame()
174 {
175     if (isCoveredFractionDynamic())
176     {
177         real cfrac = _gmx_selelem_estimate_coverfrac(rootElement());
178         coveredFraction_ = cfrac;
179         averageCoveredFraction_ += cfrac;
180     }
181 }
182
183
184 void
185 SelectionData::computeAverageCoveredFraction(int nframes)
186 {
187     if (isCoveredFractionDynamic() && nframes > 0)
188     {
189         averageCoveredFraction_ /= nframes;
190     }
191 }
192
193
194 void
195 SelectionData::restoreOriginalPositions()
196 {
197     if (isDynamic())
198     {
199         gmx_ana_pos_t &p = rawPositions_;
200         gmx_ana_index_copy(p.g, rootElement().v.u.g, false);
201         p.g->name = NULL;
202         gmx_ana_indexmap_update(&p.m, p.g, hasFlag(gmx::efSelection_DynamicMask));
203         p.nr = p.m.nr;
204         refreshMassesAndCharges();
205     }
206 }
207
208 } // namespace internal
209
210
211 void
212 Selection::printInfo(FILE *fp) const
213 {
214     fprintf(fp, "\"%s\" (%d position%s, %d atom%s%s)", name(),
215             posCount(),  posCount()  == 1 ? "" : "s",
216             atomCount(), atomCount() == 1 ? "" : "s",
217             isDynamic() ? ", dynamic" : "");
218     fprintf(fp, "\n");
219 }
220
221
222 void
223 Selection::printDebugInfo(FILE *fp, int nmaxind) const
224 {
225     const gmx_ana_pos_t &p = data().rawPositions_;
226
227     fprintf(fp, "  ");
228     printInfo(fp);
229     fprintf(fp, "    ");
230     gmx_ana_index_dump(fp, p.g, -1, nmaxind);
231
232     fprintf(fp, "    Block (size=%d):", p.m.mapb.nr);
233     if (!p.m.mapb.index)
234     {
235         fprintf(fp, " (null)");
236     }
237     else
238     {
239         int n = p.m.mapb.nr;
240         if (nmaxind >= 0 && n > nmaxind)
241             n = nmaxind;
242         for (int i = 0; i <= n; ++i)
243             fprintf(fp, " %d", p.m.mapb.index[i]);
244         if (n < p.m.mapb.nr)
245             fprintf(fp, " ...");
246     }
247     fprintf(fp, "\n");
248
249     int n = posCount();
250     if (nmaxind >= 0 && n > nmaxind)
251         n = nmaxind;
252     fprintf(fp, "    RefId:");
253     if (!p.m.refid)
254     {
255         fprintf(fp, " (null)");
256     }
257     else
258     {
259         for (int i = 0; i < n; ++i)
260             fprintf(fp, " %d", p.m.refid[i]);
261         if (n < posCount())
262             fprintf(fp, " ...");
263     }
264     fprintf(fp, "\n");
265
266     fprintf(fp, "    MapId:");
267     if (!p.m.mapid)
268     {
269         fprintf(fp, " (null)");
270     }
271     else
272     {
273         for (int i = 0; i < n; ++i)
274             fprintf(fp, " %d", p.m.mapid[i]);
275         if (n < posCount())
276             fprintf(fp, " ...");
277     }
278     fprintf(fp, "\n");
279 }
280
281 } // namespace gmx