Include selection keyword details in user guide
[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 *const help_same[] = {
174     "::",
175     "",
176     "  same KEYWORD as ATOM_EXPR",
177     "",
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",
197      "Extending selections", asize(help_same), help_same},
198 };
199
200 /*! \brief
201  * Selection method data for the \p same method.
202  *
203  * This selection method is used for matching string keywords. The parser
204  * never sees this method; _gmx_selelem_custom_init_same() replaces sm_same
205  * with this method in cases where it is required.
206  */
207 static gmx_ana_selmethod_t sm_same_str = {
208     "same", GROUP_VALUE, SMETH_SINGLEVAL,
209     asize(smparams_same_str), smparams_same_str,
210     &init_data_same,
211     NULL,
212     &init_same,
213     NULL,
214     &free_data_same,
215     &init_frame_same_str,
216     &evaluate_same_str,
217     NULL,
218     {NULL, NULL, 0, NULL},
219 };
220
221 static void *
222 init_data_same(int /* npar */, gmx_ana_selparam_t *param)
223 {
224     t_methoddata_same *data;
225
226     snew(data, 1);
227     data->as_s_sorted = NULL;
228     param[1].nvalptr  = &data->nas;
229     return data;
230 }
231
232 /*!
233  * \param[in,out] method  The method to initialize.
234  * \param[in,out] params  Pointer to the first parameter.
235  * \param[in]     scanner Scanner data structure.
236  *
237  * If \p *method is not a \c same method, this function returns
238  * immediately.
239  */
240 void
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;
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_THROW(gmx::InvalidInputError("'same' should be followed by a single keyword"));
257     }
258     gmx_ana_selmethod_t *kwmethod = kwvalues.front().expr->u.expr.method;
259     if (kwmethod->type == STR_VALUE)
260     {
261         *method = &sm_same_str;
262     }
263
264     /* We do custom processing for the "as" parameter. */
265     gmx::SelectionParserParameterList::iterator asparam = ++params->begin();
266     if (asparam != params->end() && asparam->name() == sm_same.param[1].name)
267     {
268         const gmx::SelectionParserValueList &asvalues = asparam->values();
269         if (asvalues.size() != 1 || !asvalues.front().hasExpressionValue())
270         {
271             // TODO: Think about providing more informative context.
272             GMX_THROW(gmx::InvalidInputError("'same ... as' should be followed by a single expression"));
273         }
274         const gmx::SelectionTreeElementPointer &child = asvalues.front().expr;
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, child, scanner);
280         /* Replace the second parameter with one with a value from \p kwelem. */
281         std::string pname = asparam->name();
282         *asparam = gmx::SelectionParserParameter::createFromExpression(pname, kwelem);
283     }
284 }
285
286 static void
287 init_same(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
288 {
289     t_methoddata_same *d = (t_methoddata_same *)data;
290
291     d->val.ptr = param[0].val.u.ptr;
292     d->as.ptr  = param[1].val.u.ptr;
293     if (param[1].val.type == STR_VALUE)
294     {
295         snew(d->as_s_sorted, d->nas);
296     }
297     if (!(param[0].flags & SPAR_ATOMVAL))
298     {
299         GMX_THROW(gmx::InvalidInputError(
300                           "The 'same' selection keyword combined with a "
301                           "non-keyword does not make sense"));
302     }
303 }
304
305 /*!
306  * \param data Data to free (should point to a \ref t_methoddata_same).
307  */
308 static void
309 free_data_same(void *data)
310 {
311     t_methoddata_same *d = (t_methoddata_same *)data;
312
313     sfree(d->as_s_sorted);
314     sfree(d);
315 }
316
317 /*! \brief
318  * Helper function for comparison of two integers.
319  */
320 static int
321 cmp_int(const void *a, const void *b)
322 {
323     if (*(int *)a < *(int *)b)
324     {
325         return -1;
326     }
327     if (*(int *)a > *(int *)b)
328     {
329         return 1;
330     }
331     return 0;
332 }
333
334 static void
335 init_frame_same_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */, void *data)
336 {
337     t_methoddata_same *d = (t_methoddata_same *)data;
338     int                i, j;
339
340     /* Collapse adjacent values, and check whether the array is sorted. */
341     d->bSorted = true;
342     if (d->nas == 0)
343     {
344         return;
345     }
346     for (i = 1, j = 0; i < d->nas; ++i)
347     {
348         if (d->as.i[i] != d->as.i[j])
349         {
350             if (d->as.i[i] < d->as.i[j])
351             {
352                 d->bSorted = false;
353             }
354             ++j;
355             d->as.i[j] = d->as.i[i];
356         }
357     }
358     d->nas = j + 1;
359
360     if (!d->bSorted)
361     {
362         qsort(d->as.i, d->nas, sizeof(d->as.i[0]), &cmp_int);
363         /* More identical values may become adjacent after sorting. */
364         for (i = 1, j = 0; i < d->nas; ++i)
365         {
366             if (d->as.i[i] != d->as.i[j])
367             {
368                 ++j;
369                 d->as.i[j] = d->as.i[i];
370             }
371         }
372         d->nas = j + 1;
373     }
374 }
375
376 /*!
377  * See sel_updatefunc() for description of the parameters.
378  * \p data should point to a \c t_methoddata_same.
379  *
380  * Calculates which values in \c data->val.i can be found in \c data->as.i
381  * (assumed sorted), and writes the corresponding atoms to output.
382  * If \c data->val is sorted, uses a linear scan of both arrays, otherwise a
383  * binary search of \c data->as is performed for each block of values in
384  * \c data->val.
385  */
386 static void
387 evaluate_same_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
388                   gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
389 {
390     t_methoddata_same     *d = (t_methoddata_same *)data;
391     int                    i, j;
392
393     out->u.g->isize = 0;
394     i               = j = 0;
395     while (j < g->isize)
396     {
397         if (d->bSorted)
398         {
399             /* If we are sorted, we can do a simple linear scan. */
400             while (i < d->nas && d->as.i[i] < d->val.i[j])
401             {
402                 ++i;
403             }
404         }
405         else
406         {
407             /* If not, we must do a binary search of all the values. */
408             int i1, i2;
409
410             i1 = 0;
411             i2 = d->nas;
412             while (i2 - i1 > 1)
413             {
414                 int itry = (i1 + i2) / 2;
415                 if (d->as.i[itry] <= d->val.i[j])
416                 {
417                     i1 = itry;
418                 }
419                 else
420                 {
421                     i2 = itry;
422                 }
423             }
424             i = (d->as.i[i1] == d->val.i[j] ? i1 : d->nas);
425         }
426         /* Check whether the value was found in the as list. */
427         if (i == d->nas || d->as.i[i] != d->val.i[j])
428         {
429             /* If not, skip all atoms with the same value. */
430             int tmpval = d->val.i[j];
431             ++j;
432             while (j < g->isize && d->val.i[j] == tmpval)
433             {
434                 ++j;
435             }
436         }
437         else
438         {
439             /* Copy all the atoms with this value to the output. */
440             while (j < g->isize && d->val.i[j] == d->as.i[i])
441             {
442                 out->u.g->index[out->u.g->isize++] = g->index[j];
443                 ++j;
444             }
445         }
446         if (j < g->isize && d->val.i[j] < d->val.i[j - 1])
447         {
448             d->bSorted = false;
449         }
450     }
451 }
452
453 /*! \brief
454  * Helper function for comparison of two strings.
455  */
456 static int
457 cmp_str(const void *a, const void *b)
458 {
459     return strcmp(*(char **)a, *(char **)b);
460 }
461
462 static void
463 init_frame_same_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */, void *data)
464 {
465     t_methoddata_same *d = (t_methoddata_same *)data;
466     int                i, j;
467
468     /* Collapse adjacent values.
469      * For strings, it's unlikely that the values would be sorted originally,
470      * so set bSorted always to false. */
471     d->bSorted        = false;
472     if (d->nas == 0)
473     {
474         return;
475     }
476     d->as_s_sorted[0] = d->as.s[0];
477     for (i = 1, j = 0; i < d->nas; ++i)
478     {
479         if (strcmp(d->as.s[i], d->as_s_sorted[j]) != 0)
480         {
481             ++j;
482             d->as_s_sorted[j] = d->as.s[i];
483         }
484     }
485     d->nas = j + 1;
486
487     qsort(d->as_s_sorted, d->nas, sizeof(d->as_s_sorted[0]), &cmp_str);
488     /* More identical values may become adjacent after sorting. */
489     for (i = 1, j = 0; i < d->nas; ++i)
490     {
491         if (strcmp(d->as_s_sorted[i], d->as_s_sorted[j]) != 0)
492         {
493             ++j;
494             d->as_s_sorted[j] = d->as_s_sorted[i];
495         }
496     }
497     d->nas = j + 1;
498 }
499
500 /*!
501  * See sel_updatefunc() for description of the parameters.
502  * \p data should point to a \c t_methoddata_same.
503  *
504  * Calculates which strings in \c data->val.s can be found in \c data->as.s
505  * (assumed sorted), and writes the corresponding atoms to output.
506  * A binary search of \c data->as is performed for each block of values in
507  * \c data->val.
508  */
509 static void
510 evaluate_same_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
511                   gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
512 {
513     t_methoddata_same     *d = (t_methoddata_same *)data;
514     int                    j;
515
516     out->u.g->isize = 0;
517     j               = 0;
518     while (j < g->isize)
519     {
520         /* Do a binary search of the strings. */
521         void *ptr;
522         ptr = bsearch(&d->val.s[j], d->as_s_sorted, d->nas,
523                       sizeof(d->as_s_sorted[0]), &cmp_str);
524         /* Check whether the value was found in the as list. */
525         if (ptr == NULL)
526         {
527             /* If not, skip all atoms with the same value. */
528             const char *tmpval = d->val.s[j];
529             ++j;
530             while (j < g->isize && strcmp(d->val.s[j], tmpval) == 0)
531             {
532                 ++j;
533             }
534         }
535         else
536         {
537             const char *tmpval = d->val.s[j];
538             /* Copy all the atoms with this value to the output. */
539             while (j < g->isize && strcmp(d->val.s[j], tmpval) == 0)
540             {
541                 out->u.g->index[out->u.g->isize++] = g->index[j];
542                 ++j;
543             }
544         }
545     }
546 }