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