Avoid crash with empty reference for 'same as' selection
[alexxy/gromacs.git] / src / gromacs / selection / sm_same.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012,2013,2014,2015, 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 /*! \internal \file
36  * \brief
37  * Implements the \p same selection method.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_selection
41  */
42 #include <stdlib.h>
43
44 #include "gromacs/legacyheaders/macros.h"
45
46 #include "gromacs/selection/selmethod.h"
47 #include "gromacs/utility/exceptions.h"
48 #include "gromacs/utility/smalloc.h"
49
50 #include "keywords.h"
51 #include "parsetree.h"
52 #include "selelem.h"
53
54 /*! \internal
55  * \brief
56  * Data structure for the \p same selection method.
57  *
58  * To avoid duplicate initialization code, the same data structure is used
59  * for matching both integer and string keywords; hence the unions.
60  *
61  * \ingroup module_selection
62  */
63 typedef struct
64 {
65     /** Value for each atom to match. */
66     union
67     {
68         int                 *i;
69         char               **s;
70         void                *ptr;
71     }                        val;
72     /*! \brief
73      * Number of values in the \p as array.
74      *
75      * For string values, this is actually the number of values in the
76      * \p as_s_sorted array.
77      */
78     int                      nas;
79     /** Values to match against. */
80     union
81     {
82         int                 *i;
83         char               **s;
84         void                *ptr;
85     }                        as;
86     /*! \brief
87      * Separate array for sorted \p as.s array.
88      *
89      * The array of strings returned as the output value of a parameter should
90      * not be messed with to avoid memory corruption (the pointers in the array
91      * may be reused for several evaluations), so we keep our own copy for
92      * modifications.
93      */
94     char                   **as_s_sorted;
95     /** Whether simple matching can be used. */
96     bool                     bSorted;
97 } t_methoddata_same;
98
99 /*! \brief
100  * Allocates data for the \p same selection method.
101  *
102  * \param[in]     npar  Not used (should be 2).
103  * \param[in,out] param Method parameters (should point to a copy of
104  *      ::smparams_same_int or ::smparams_same_str).
105  * \returns Pointer to the allocated data (\ref t_methoddata_same).
106  */
107 static void *
108 init_data_same(int npar, gmx_ana_selparam_t *param);
109 /*! \brief
110  * Initializes the \p same selection method.
111  *
112  * \param   top   Not used.
113  * \param   npar  Not used (should be 2).
114  * \param   param Initialized method parameters (should point to a copy of
115  *      ::smparams_same_int or ::smparams_same_str).
116  * \param   data  Pointer to \ref t_methoddata_same to initialize.
117  * \returns 0 on success, -1 on failure.
118  */
119 static void
120 init_same(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
121 /** Frees the data allocated for the \p same selection method. */
122 static void
123 free_data_same(void *data);
124 /*! \brief
125  * Initializes the evaluation of the \p same selection method for a frame.
126  *
127  * \param[in]  top  Not used.
128  * \param[in]  fr   Current frame.
129  * \param[in]  pbc  PBC structure.
130  * \param      data Should point to a \ref t_methoddata_same.
131  *
132  * Sorts the \c data->as.i array and removes identical values for faster and
133  * simpler lookup.
134  */
135 static void
136 init_frame_same_int(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
137 /** Evaluates the \p same selection method. */
138 static void
139 evaluate_same_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
140                   gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
141 /*! \brief
142  * Initializes the evaluation of the \p same selection method for a frame.
143  *
144  * \param[in]  top  Not used.
145  * \param[in]  fr   Current frame.
146  * \param[in]  pbc  PBC structure.
147  * \param      data Should point to a \ref t_methoddata_same.
148  *
149  * Sorts the \c data->as.s array and removes identical values for faster and
150  * simpler lookup.
151  */
152 static void
153 init_frame_same_str(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
154 /** Evaluates the \p same selection method. */
155 static void
156 evaluate_same_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
157                   gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
158
159 /** Parameters for the \p same selection method. */
160 static gmx_ana_selparam_t smparams_same_int[] = {
161     {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_ATOMVAL},
162     {"as", {INT_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
163 };
164
165 /** Parameters for the \p same selection method. */
166 static gmx_ana_selparam_t smparams_same_str[] = {
167     {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_ATOMVAL},
168     {"as", {STR_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
169 };
170
171 /** Help text for the \p same selection method. */
172 static const char *help_same[] = {
173     "EXTENDING SELECTIONS[PAR]",
174
175     "[TT]same KEYWORD as ATOM_EXPR[tt][PAR]",
176
177     "The keyword [TT]same[tt] can be used to select all atoms for which",
178     "the given [TT]KEYWORD[tt] matches any of the atoms in [TT]ATOM_EXPR[tt].",
179     "Keywords that evaluate to integer or string values are supported.",
180 };
181
182 /*! \internal \brief Selection method data for the \p same method. */
183 gmx_ana_selmethod_t sm_same = {
184     "same", GROUP_VALUE, 0,
185     asize(smparams_same_int), smparams_same_int,
186     &init_data_same,
187     NULL,
188     &init_same,
189     NULL,
190     &free_data_same,
191     &init_frame_same_int,
192     &evaluate_same_int,
193     NULL,
194     {"same KEYWORD as ATOM_EXPR", asize(help_same), help_same},
195 };
196
197 /*! \brief
198  * Selection method data for the \p same method.
199  *
200  * This selection method is used for matching string keywords. The parser
201  * never sees this method; _gmx_selelem_custom_init_same() replaces sm_same
202  * with this method in cases where it is required.
203  */
204 static gmx_ana_selmethod_t sm_same_str = {
205     "same", GROUP_VALUE, SMETH_SINGLEVAL,
206     asize(smparams_same_str), smparams_same_str,
207     &init_data_same,
208     NULL,
209     &init_same,
210     NULL,
211     &free_data_same,
212     &init_frame_same_str,
213     &evaluate_same_str,
214     NULL,
215     {"same KEYWORD as ATOM_EXPR", asize(help_same), help_same},
216 };
217
218 static void *
219 init_data_same(int /* npar */, gmx_ana_selparam_t *param)
220 {
221     t_methoddata_same *data;
222
223     snew(data, 1);
224     data->as_s_sorted = NULL;
225     param[1].nvalptr  = &data->nas;
226     return data;
227 }
228
229 /*!
230  * \param[in,out] method  The method to initialize.
231  * \param[in,out] params  Pointer to the first parameter.
232  * \param[in]     scanner Scanner data structure.
233  * \returns       0 on success, a non-zero error code on error.
234  *
235  * If \p *method is not a \c same method, this function returns zero
236  * immediately.
237  */
238 int
239 _gmx_selelem_custom_init_same(gmx_ana_selmethod_t                           **method,
240                               const gmx::SelectionParserParameterListPointer &params,
241                               void                                           *scanner)
242 {
243
244     /* Do nothing if this is not a same method. */
245     if (!*method || (*method)->name != sm_same.name || params->empty())
246     {
247         return 0;
248     }
249
250     const gmx::SelectionParserValueList &kwvalues = params->front().values();
251     if (kwvalues.size() != 1 || !kwvalues.front().hasExpressionValue()
252         || kwvalues.front().expr->type != SEL_EXPRESSION)
253     {
254         _gmx_selparser_error(scanner, "'same' should be followed by a single keyword");
255         return -1;
256     }
257     gmx_ana_selmethod_t *kwmethod = kwvalues.front().expr->u.expr.method;
258     if (kwmethod->type == STR_VALUE)
259     {
260         *method = &sm_same_str;
261     }
262
263     /* We do custom processing for the "as" parameter. */
264     gmx::SelectionParserParameterList::iterator asparam = ++params->begin();
265     if (asparam != params->end() && asparam->name() == sm_same.param[1].name)
266     {
267         gmx::SelectionParserParameterList    kwparams;
268         gmx::SelectionParserValueListPointer values(
269                 new gmx::SelectionParserValueList(asparam->values()));
270         kwparams.push_back(
271                 gmx::SelectionParserParameter::create(NULL, move(values)));
272
273         /* Create a second keyword evaluation element for the keyword given as
274          * the first parameter, evaluating the keyword in the group given by the
275          * second parameter. */
276         gmx::SelectionTreeElementPointer kwelem
277             = _gmx_sel_init_keyword_evaluator(kwmethod, kwparams, scanner);
278         // FIXME: Use exceptions.
279         if (!kwelem)
280         {
281             return -1;
282         }
283         /* Replace the second parameter with one with a value from \p kwelem. */
284         std::string pname = asparam->name();
285         *asparam = gmx::SelectionParserParameter::createFromExpression(pname, kwelem);
286     }
287     return 0;
288 }
289
290 static void
291 init_same(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
292 {
293     t_methoddata_same *d = (t_methoddata_same *)data;
294
295     d->val.ptr = param[0].val.u.ptr;
296     d->as.ptr  = param[1].val.u.ptr;
297     if (param[1].val.type == STR_VALUE)
298     {
299         snew(d->as_s_sorted, d->nas);
300     }
301     if (!(param[0].flags & SPAR_ATOMVAL))
302     {
303         GMX_THROW(gmx::InvalidInputError(
304                           "The 'same' selection keyword combined with a "
305                           "non-keyword does not make sense"));
306     }
307 }
308
309 /*!
310  * \param data Data to free (should point to a \ref t_methoddata_same).
311  */
312 static void
313 free_data_same(void *data)
314 {
315     t_methoddata_same *d = (t_methoddata_same *)data;
316
317     sfree(d->as_s_sorted);
318     sfree(d);
319 }
320
321 /*! \brief
322  * Helper function for comparison of two integers.
323  */
324 static int
325 cmp_int(const void *a, const void *b)
326 {
327     if (*(int *)a < *(int *)b)
328     {
329         return -1;
330     }
331     if (*(int *)a > *(int *)b)
332     {
333         return 1;
334     }
335     return 0;
336 }
337
338 static void
339 init_frame_same_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */, void *data)
340 {
341     t_methoddata_same *d = (t_methoddata_same *)data;
342     int                i, j;
343
344     /* Collapse adjacent values, and check whether the array is sorted. */
345     d->bSorted = true;
346     if (d->nas == 0)
347     {
348         return;
349     }
350     for (i = 1, j = 0; i < d->nas; ++i)
351     {
352         if (d->as.i[i] != d->as.i[j])
353         {
354             if (d->as.i[i] < d->as.i[j])
355             {
356                 d->bSorted = false;
357             }
358             ++j;
359             d->as.i[j] = d->as.i[i];
360         }
361     }
362     d->nas = j + 1;
363
364     if (!d->bSorted)
365     {
366         qsort(d->as.i, d->nas, sizeof(d->as.i[0]), &cmp_int);
367         /* More identical values may become adjacent after sorting. */
368         for (i = 1, j = 0; i < d->nas; ++i)
369         {
370             if (d->as.i[i] != d->as.i[j])
371             {
372                 ++j;
373                 d->as.i[j] = d->as.i[i];
374             }
375         }
376         d->nas = j + 1;
377     }
378 }
379
380 /*!
381  * See sel_updatefunc() for description of the parameters.
382  * \p data should point to a \c t_methoddata_same.
383  *
384  * Calculates which values in \c data->val.i can be found in \c data->as.i
385  * (assumed sorted), and writes the corresponding atoms to output.
386  * If \c data->val is sorted, uses a linear scan of both arrays, otherwise a
387  * binary search of \c data->as is performed for each block of values in
388  * \c data->val.
389  */
390 static void
391 evaluate_same_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
392                   gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
393 {
394     t_methoddata_same     *d = (t_methoddata_same *)data;
395     int                    i, j;
396
397     out->u.g->isize = 0;
398     i               = j = 0;
399     while (j < g->isize)
400     {
401         if (d->bSorted)
402         {
403             /* If we are sorted, we can do a simple linear scan. */
404             while (i < d->nas && d->as.i[i] < d->val.i[j])
405             {
406                 ++i;
407             }
408         }
409         else
410         {
411             /* If not, we must do a binary search of all the values. */
412             int i1, i2;
413
414             i1 = 0;
415             i2 = d->nas;
416             while (i2 - i1 > 1)
417             {
418                 int itry = (i1 + i2) / 2;
419                 if (d->as.i[itry] <= d->val.i[j])
420                 {
421                     i1 = itry;
422                 }
423                 else
424                 {
425                     i2 = itry;
426                 }
427             }
428             i = (d->as.i[i1] == d->val.i[j] ? i1 : d->nas);
429         }
430         /* Check whether the value was found in the as list. */
431         if (i == d->nas || d->as.i[i] != d->val.i[j])
432         {
433             /* If not, skip all atoms with the same value. */
434             int tmpval = d->val.i[j];
435             ++j;
436             while (j < g->isize && d->val.i[j] == tmpval)
437             {
438                 ++j;
439             }
440         }
441         else
442         {
443             /* Copy all the atoms with this value to the output. */
444             while (j < g->isize && d->val.i[j] == d->as.i[i])
445             {
446                 out->u.g->index[out->u.g->isize++] = g->index[j];
447                 ++j;
448             }
449         }
450         if (j < g->isize && d->val.i[j] < d->val.i[j - 1])
451         {
452             d->bSorted = false;
453         }
454     }
455 }
456
457 /*! \brief
458  * Helper function for comparison of two strings.
459  */
460 static int
461 cmp_str(const void *a, const void *b)
462 {
463     return strcmp(*(char **)a, *(char **)b);
464 }
465
466 static void
467 init_frame_same_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */, void *data)
468 {
469     t_methoddata_same *d = (t_methoddata_same *)data;
470     int                i, j;
471
472     /* Collapse adjacent values.
473      * For strings, it's unlikely that the values would be sorted originally,
474      * so set bSorted always to false. */
475     d->bSorted        = false;
476     if (d->nas == 0)
477     {
478         return;
479     }
480     d->as_s_sorted[0] = d->as.s[0];
481     for (i = 1, j = 0; i < d->nas; ++i)
482     {
483         if (strcmp(d->as.s[i], d->as_s_sorted[j]) != 0)
484         {
485             ++j;
486             d->as_s_sorted[j] = d->as.s[i];
487         }
488     }
489     d->nas = j + 1;
490
491     qsort(d->as_s_sorted, d->nas, sizeof(d->as_s_sorted[0]), &cmp_str);
492     /* More identical values may become adjacent after sorting. */
493     for (i = 1, j = 0; i < d->nas; ++i)
494     {
495         if (strcmp(d->as_s_sorted[i], d->as_s_sorted[j]) != 0)
496         {
497             ++j;
498             d->as_s_sorted[j] = d->as_s_sorted[i];
499         }
500     }
501     d->nas = j + 1;
502 }
503
504 /*!
505  * See sel_updatefunc() for description of the parameters.
506  * \p data should point to a \c t_methoddata_same.
507  *
508  * Calculates which strings in \c data->val.s can be found in \c data->as.s
509  * (assumed sorted), and writes the corresponding atoms to output.
510  * A binary search of \c data->as is performed for each block of values in
511  * \c data->val.
512  */
513 static void
514 evaluate_same_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
515                   gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
516 {
517     t_methoddata_same     *d = (t_methoddata_same *)data;
518     int                    j;
519
520     out->u.g->isize = 0;
521     j               = 0;
522     while (j < g->isize)
523     {
524         /* Do a binary search of the strings. */
525         void *ptr;
526         ptr = bsearch(&d->val.s[j], d->as_s_sorted, d->nas,
527                       sizeof(d->as_s_sorted[0]), &cmp_str);
528         /* Check whether the value was found in the as list. */
529         if (ptr == NULL)
530         {
531             /* If not, skip all atoms with the same value. */
532             const char *tmpval = d->val.s[j];
533             ++j;
534             while (j < g->isize && strcmp(d->val.s[j], tmpval) == 0)
535             {
536                 ++j;
537             }
538         }
539         else
540         {
541             const char *tmpval = d->val.s[j];
542             /* Copy all the atoms with this value to the output. */
543             while (j < g->isize && strcmp(d->val.s[j], tmpval) == 0)
544             {
545                 out->u.g->index[out->u.g->isize++] = g->index[j];
546                 ++j;
547             }
548         }
549     }
550 }