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