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