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