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