Rename all source files from - to _ for consistency.
[gromacs.git] / src / gromacs / selection / sm_same.cpp
blob26d7c5d16c63d95eb73a70a58731c8d213c6253e
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
35 /*! \internal \file
36 * \brief
37 * Implements the \p same selection method.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
42 #include "gmxpre.h"
44 #include <cstdlib>
45 #include <cstring>
47 #include "gromacs/utility/arraysize.h"
48 #include "gromacs/utility/exceptions.h"
49 #include "gromacs/utility/smalloc.h"
51 #include "keywords.h"
52 #include "parsetree.h"
53 #include "selelem.h"
54 #include "selmethod.h"
55 #include "selmethod_impl.h"
57 /*! \internal
58 * \brief
59 * Data structure for the \p same selection method.
61 * To avoid duplicate initialization code, the same data structure is used
62 * for matching both integer and string keywords; hence the unions.
64 * \ingroup module_selection
66 typedef struct
68 /** Value for each atom to match. */
69 union
71 int *i;
72 char **s;
73 void *ptr;
74 } val;
75 /*! \brief
76 * Number of values in the \p as array.
78 * For string values, this is actually the number of values in the
79 * \p as_s_sorted array.
81 int nas;
82 /** Values to match against. */
83 union
85 int *i;
86 char **s;
87 void *ptr;
88 } as;
89 /*! \brief
90 * Separate array for sorted \p as.s array.
92 * The array of strings returned as the output value of a parameter should
93 * not be messed with to avoid memory corruption (the pointers in the array
94 * may be reused for several evaluations), so we keep our own copy for
95 * modifications.
97 char **as_s_sorted;
98 /** Whether simple matching can be used. */
99 bool bSorted;
100 } t_methoddata_same;
102 /*! \brief
103 * Allocates data for the \p same selection method.
105 * \param[in] npar Not used (should be 2).
106 * \param[in,out] param Method parameters (should point to a copy of
107 * ::smparams_same_int or ::smparams_same_str).
108 * \returns Pointer to the allocated data (\ref t_methoddata_same).
110 static void *
111 init_data_same(int npar, gmx_ana_selparam_t *param);
112 /*! \brief
113 * Initializes the \p same selection method.
115 * \param top Not used.
116 * \param npar Not used (should be 2).
117 * \param param Initialized method parameters (should point to a copy of
118 * ::smparams_same_int or ::smparams_same_str).
119 * \param data Pointer to \ref t_methoddata_same to initialize.
120 * \returns 0 on success, -1 on failure.
122 static void
123 init_same(const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data);
124 /** Frees the data allocated for the \p same selection method. */
125 static void
126 free_data_same(void *data);
127 /*! \brief
128 * Initializes the evaluation of the \p same selection method for a frame.
130 * \param[in] context Not used.
131 * \param data Should point to a \ref t_methoddata_same.
133 * Sorts the \c data->as.i array and removes identical values for faster and
134 * simpler lookup.
136 static void
137 init_frame_same_int(const gmx::SelMethodEvalContext &context, void *data);
138 /** Evaluates the \p same selection method. */
139 static void
140 evaluate_same_int(const gmx::SelMethodEvalContext &context,
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.
145 * \param[in] context Not used.
146 * \param data Should point to a \ref t_methoddata_same.
148 * Sorts the \c data->as.s array and removes identical values for faster and
149 * simpler lookup.
151 static void
152 init_frame_same_str(const gmx::SelMethodEvalContext &context, void *data);
153 /** Evaluates the \p same selection method. */
154 static void
155 evaluate_same_str(const gmx::SelMethodEvalContext &context,
156 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
158 /** Parameters for the \p same selection method. */
159 static gmx_ana_selparam_t smparams_same_int[] = {
160 {nullptr, {INT_VALUE, -1, {nullptr}}, nullptr, SPAR_DYNAMIC | SPAR_ATOMVAL},
161 {"as", {INT_VALUE, -1, {nullptr}}, nullptr, SPAR_DYNAMIC | SPAR_VARNUM},
164 /** Parameters for the \p same selection method. */
165 static gmx_ana_selparam_t smparams_same_str[] = {
166 {nullptr, {STR_VALUE, -1, {nullptr}}, nullptr, SPAR_DYNAMIC | SPAR_ATOMVAL},
167 {"as", {STR_VALUE, -1, {nullptr}}, nullptr, SPAR_DYNAMIC | SPAR_VARNUM},
170 /** Help text for the \p same selection method. */
171 static const char *const help_same[] = {
172 "::",
174 " same KEYWORD as ATOM_EXPR",
177 "The keyword [TT]same[tt] can be used to select all atoms for which",
178 "the given [TT]KEYWORD[tt] matches any of the atoms in [TT]ATOM_EXPR[tt].",
179 "Keywords that evaluate to integer or string values are supported.",
182 /*! \internal \brief Selection method data for the \p same method. */
183 gmx_ana_selmethod_t sm_same = {
184 "same", GROUP_VALUE, 0,
185 asize(smparams_same_int), smparams_same_int,
186 &init_data_same,
187 nullptr,
188 &init_same,
189 nullptr,
190 &free_data_same,
191 &init_frame_same_int,
192 &evaluate_same_int,
193 nullptr,
194 {"same KEYWORD as ATOM_EXPR",
195 "Extending selections", asize(help_same), help_same},
198 /*! \brief
199 * Selection method data for the \p same method.
201 * This selection method is used for matching string keywords. The parser
202 * never sees this method; _gmx_selelem_custom_init_same() replaces sm_same
203 * with this method in cases where it is required.
205 static gmx_ana_selmethod_t sm_same_str = {
206 "same", GROUP_VALUE, SMETH_SINGLEVAL,
207 asize(smparams_same_str), 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},
219 static void *
220 init_data_same(int /* npar */, gmx_ana_selparam_t *param)
222 t_methoddata_same *data;
224 snew(data, 1);
225 data->as_s_sorted = nullptr;
226 param[1].nvalptr = &data->nas;
227 return data;
231 * \param[in,out] method The method to initialize.
232 * \param[in,out] params Pointer to the first parameter.
233 * \param[in] scanner Scanner data structure.
235 * If \p *method is not a \c same method, this function returns
236 * immediately.
238 void
239 _gmx_selelem_custom_init_same(gmx_ana_selmethod_t **method,
240 const gmx::SelectionParserParameterListPointer &params,
241 void *scanner)
244 /* Do nothing if this is not a same method. */
245 if (!*method || (*method)->name != sm_same.name || params->empty())
247 return;
250 const gmx::SelectionParserValueList &kwvalues = params->front().values();
251 if (kwvalues.size() != 1 || !kwvalues.front().hasExpressionValue()
252 || kwvalues.front().expr->type != SEL_EXPRESSION)
254 GMX_THROW(gmx::InvalidInputError("'same' should be followed by a single keyword"));
256 gmx_ana_selmethod_t *kwmethod = kwvalues.front().expr->u.expr.method;
257 if (kwmethod->type == STR_VALUE)
259 *method = &sm_same_str;
262 /* We do custom processing for the "as" parameter. */
263 gmx::SelectionParserParameterList::iterator asparam = ++params->begin();
264 if (asparam != params->end() && asparam->name() == sm_same.param[1].name)
266 const gmx::SelectionParserValueList &asvalues = asparam->values();
267 if (asvalues.size() != 1 || !asvalues.front().hasExpressionValue())
269 // TODO: Think about providing more informative context.
270 GMX_THROW(gmx::InvalidInputError("'same ... as' should be followed by a single expression"));
272 const gmx::SelectionTreeElementPointer &child = asvalues.front().expr;
273 /* Create a second keyword evaluation element for the keyword given as
274 * the first parameter, evaluating the keyword in the group given by the
275 * second parameter. */
276 gmx::SelectionTreeElementPointer kwelem
277 = _gmx_sel_init_keyword_evaluator(kwmethod, child, scanner);
278 /* Replace the second parameter with one with a value from \p kwelem. */
279 std::string pname = asparam->name();
280 *asparam = gmx::SelectionParserParameter::createFromExpression(pname, kwelem);
284 static void
285 init_same(const gmx_mtop_t * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
287 t_methoddata_same *d = static_cast<t_methoddata_same *>(data);
289 d->val.ptr = param[0].val.u.ptr;
290 d->as.ptr = param[1].val.u.ptr;
291 if (param[1].val.type == STR_VALUE)
293 snew(d->as_s_sorted, d->nas);
295 if (!(param[0].flags & SPAR_ATOMVAL))
297 GMX_THROW(gmx::InvalidInputError(
298 "The 'same' selection keyword combined with a "
299 "non-keyword does not make sense"));
304 * \param data Data to free (should point to a \ref t_methoddata_same).
306 static void
307 free_data_same(void *data)
309 t_methoddata_same *d = static_cast<t_methoddata_same *>(data);
311 sfree(d->as_s_sorted);
312 sfree(d);
315 /*! \brief
316 * Helper function for comparison of two integers.
318 static int
319 cmp_int(const void *a, const void *b)
321 if (*reinterpret_cast<const int*>(a) < *reinterpret_cast<const int*>(b))
323 return -1;
325 if (*reinterpret_cast<const int*>(a) > *reinterpret_cast<const int*>(b))
327 return 1;
329 return 0;
332 static void
333 init_frame_same_int(const gmx::SelMethodEvalContext & /*context*/, void *data)
335 t_methoddata_same *d = static_cast<t_methoddata_same *>(data);
336 int i, j;
338 /* Collapse adjacent values, and check whether the array is sorted. */
339 d->bSorted = true;
340 if (d->nas == 0)
342 return;
344 for (i = 1, j = 0; i < d->nas; ++i)
346 if (d->as.i[i] != d->as.i[j])
348 if (d->as.i[i] < d->as.i[j])
350 d->bSorted = false;
352 ++j;
353 d->as.i[j] = d->as.i[i];
356 d->nas = j + 1;
358 if (!d->bSorted)
360 qsort(d->as.i, d->nas, sizeof(d->as.i[0]), &cmp_int);
361 /* More identical values may become adjacent after sorting. */
362 for (i = 1, j = 0; i < d->nas; ++i)
364 if (d->as.i[i] != d->as.i[j])
366 ++j;
367 d->as.i[j] = d->as.i[i];
370 d->nas = j + 1;
375 * See sel_updatefunc() for description of the parameters.
376 * \p data should point to a \c t_methoddata_same.
378 * Calculates which values in \c data->val.i can be found in \c data->as.i
379 * (assumed sorted), and writes the corresponding atoms to output.
380 * If \c data->val is sorted, uses a linear scan of both arrays, otherwise a
381 * binary search of \c data->as is performed for each block of values in
382 * \c data->val.
384 static void
385 evaluate_same_int(const gmx::SelMethodEvalContext & /*context*/,
386 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
388 t_methoddata_same *d = static_cast<t_methoddata_same *>(data);
389 int i, j;
391 out->u.g->isize = 0;
392 i = j = 0;
393 while (j < g->isize)
395 if (d->bSorted)
397 /* If we are sorted, we can do a simple linear scan. */
398 while (i < d->nas && d->as.i[i] < d->val.i[j])
400 ++i;
403 else
405 /* If not, we must do a binary search of all the values. */
406 int i1, i2;
408 i1 = 0;
409 i2 = d->nas;
410 while (i2 - i1 > 1)
412 int itry = (i1 + i2) / 2;
413 if (d->as.i[itry] <= d->val.i[j])
415 i1 = itry;
417 else
419 i2 = itry;
422 i = (d->as.i[i1] == d->val.i[j] ? i1 : d->nas);
424 /* Check whether the value was found in the as list. */
425 if (i == d->nas || d->as.i[i] != d->val.i[j])
427 /* If not, skip all atoms with the same value. */
428 int tmpval = d->val.i[j];
429 ++j;
430 while (j < g->isize && d->val.i[j] == tmpval)
432 ++j;
435 else
437 /* Copy all the atoms with this value to the output. */
438 while (j < g->isize && d->val.i[j] == d->as.i[i])
440 out->u.g->index[out->u.g->isize++] = g->index[j];
441 ++j;
444 if (j < g->isize && d->val.i[j] < d->val.i[j - 1])
446 d->bSorted = false;
451 /*! \brief
452 * Helper function for comparison of two strings.
454 static int
455 cmp_str(const void *a, const void *b)
457 return strcmp(*static_cast<char *const*>(a), *static_cast<char *const*>(b));
460 static void
461 init_frame_same_str(const gmx::SelMethodEvalContext & /*context*/, void *data)
463 t_methoddata_same *d = static_cast<t_methoddata_same *>(data);
464 int i, j;
466 /* Collapse adjacent values.
467 * For strings, it's unlikely that the values would be sorted originally,
468 * so set bSorted always to false. */
469 d->bSorted = false;
470 if (d->nas == 0)
472 return;
474 d->as_s_sorted[0] = d->as.s[0];
475 for (i = 1, j = 0; i < d->nas; ++i)
477 if (strcmp(d->as.s[i], d->as_s_sorted[j]) != 0)
479 ++j;
480 d->as_s_sorted[j] = d->as.s[i];
483 d->nas = j + 1;
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)
489 if (strcmp(d->as_s_sorted[i], d->as_s_sorted[j]) != 0)
491 ++j;
492 d->as_s_sorted[j] = d->as_s_sorted[i];
495 d->nas = j + 1;
499 * See sel_updatefunc() for description of the parameters.
500 * \p data should point to a \c t_methoddata_same.
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.
507 static void
508 evaluate_same_str(const gmx::SelMethodEvalContext & /*context*/,
509 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
511 t_methoddata_same *d = static_cast<t_methoddata_same *>(data);
512 int j;
514 out->u.g->isize = 0;
515 j = 0;
516 while (j < g->isize)
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 == nullptr)
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)
530 ++j;
533 else
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)
539 out->u.g->index[out->u.g->isize++] = g->index[j];
540 ++j;