2 * This file is part of the GROMACS molecular simulation package.
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.
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.
37 * Implements internal selection method for comparison expressions.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
44 #include "gromacs/legacyheaders/macros.h"
46 #include "gromacs/math/utilities.h"
47 #include "gromacs/selection/selmethod.h"
48 #include "gromacs/utility/common.h"
49 #include "gromacs/utility/exceptions.h"
50 #include "gromacs/utility/smalloc.h"
52 /** Defines the comparison operator for comparison expressions. */
55 CMP_INVALID, /**< Indicates an error */
60 CMP_EQUAL, /**< '==' */
64 /** The operand has a single value. */
65 #define CMP_SINGLEVAL 1
66 /** The operand value is dynamic. */
67 #define CMP_DYNAMICVAL 2
68 /** The value is real. */
70 /** The integer array is allocated. */
71 #define CMP_ALLOCINT 16
72 /** The real array is allocated. */
73 #define CMP_ALLOCREAL 32
76 * Data structure for comparison expression operand values.
80 /** Flags that describe the type of the operand. */
82 /** (Array of) integer value(s). */
84 /** (Array of) real value(s). */
89 * Data structure for comparison expression evaluation.
93 /** Comparison operator as a string. */
95 /** Comparison operator type. */
100 t_compare_value right;
101 } t_methoddata_compare;
104 * Allocates data for comparison expression evaluation.
106 * \param[in] npar Not used (should be 5).
107 * \param[in,out] param Method parameters (should point to a copy of
108 * \ref smparams_compare).
109 * \returns Pointer to the allocated data (\c t_methoddata_compare).
111 * Allocates memory for a \c t_methoddata_compare structure.
114 init_data_compare(int npar, gmx_ana_selparam_t *param);
116 * Initializes data for comparison expression evaluation.
118 * \param[in] top Not used.
119 * \param[in] npar Not used (should be 5).
120 * \param[in] param Method parameters (should point to \ref smparams_compare).
121 * \param[in] data Should point to a \c t_methoddata_compare.
122 * \returns 0 if the input data is valid, -1 on error.
125 init_compare(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
126 /** Frees the memory allocated for comparison expression evaluation. */
128 free_data_compare(void *data);
129 /** Evaluates comparison expressions. */
131 evaluate_compare(t_topology *top, t_trxframe *fr, t_pbc *pbc,
132 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
134 /** Parameters for comparison expression evaluation. */
135 static gmx_ana_selparam_t smparams_compare[] = {
136 {"int1", {INT_VALUE, -1, {NULL}}, NULL,
137 SPAR_OPTIONAL | SPAR_DYNAMIC | SPAR_ATOMVAL},
138 {"real1", {REAL_VALUE, -1, {NULL}}, NULL,
139 SPAR_OPTIONAL | SPAR_DYNAMIC | SPAR_ATOMVAL},
140 {"op", {STR_VALUE, 1, {NULL}}, NULL, 0},
141 {"int2", {INT_VALUE, -1, {NULL}}, NULL,
142 SPAR_OPTIONAL | SPAR_DYNAMIC | SPAR_ATOMVAL},
143 {"real2", {REAL_VALUE, -1, {NULL}}, NULL,
144 SPAR_OPTIONAL | SPAR_DYNAMIC | SPAR_ATOMVAL},
147 /** \internal Selection method data for comparison expression evaluation. */
148 gmx_ana_selmethod_t sm_compare = {
149 "cmp", GROUP_VALUE, SMETH_SINGLEVAL,
150 asize(smparams_compare), smparams_compare,
163 * Returns a \c e_comparison_t value corresponding to an operator.
165 * \param[in] str String to process.
166 * \returns The comparison type corresponding to the first one or two
167 * characters of \p str.
169 * \p str can contain any number of characters; only the first two
171 * If the beginning of \p str does not match any of the recognized types,
172 * \ref CMP_INVALID is returned.
174 static e_comparison_t
175 comparison_type(char *str)
179 case '<': return (str[1] == '=') ? CMP_LEQ : CMP_LESS;
180 case '>': return (str[1] == '=') ? CMP_GEQ : CMP_GTR;
181 case '=': return (str[1] == '=') ? CMP_EQUAL : CMP_INVALID;
182 case '!': return (str[1] == '=') ? CMP_NEQ : CMP_INVALID;
188 * Returns a string corresponding to a \c e_comparison_t value.
190 * \param[in] cmpt Comparison type to convert.
191 * \returns Pointer to a string that corresponds to \p cmpt.
193 * The return value points to a string constant and should not be \p free'd.
195 * The function returns NULL if \p cmpt is not one of the valid values.
198 comparison_type_str(e_comparison_t cmpt)
200 const char *p = NULL;
203 case CMP_INVALID: p = "INVALID"; break;
204 case CMP_LESS: p = "<"; break;
205 case CMP_LEQ: p = "<="; break;
206 case CMP_GTR: p = ">"; break;
207 case CMP_GEQ: p = ">="; break;
208 case CMP_EQUAL: p = "=="; break;
209 case CMP_NEQ: p = "!="; break;
210 // No default clause so we intentionally get compiler errors
211 // if new selection choices are added later.
217 * \param[in] fp File to receive the output.
218 * \param[in] data Should point to a \c t_methoddata_compare.
221 _gmx_selelem_print_compare_info(FILE *fp, void *data)
223 t_methoddata_compare *d = (t_methoddata_compare *)data;
226 /* Print the left value */
227 if ((d->left.flags & CMP_SINGLEVAL) && !(d->left.flags & CMP_DYNAMICVAL))
229 if (d->left.flags & CMP_REALVAL)
231 fprintf(fp, "%f ", d->left.r[0]);
235 fprintf(fp, "%d ", d->left.i[0]);
238 /* Print the operator */
239 if (d->cmpt != CMP_INVALID)
241 fprintf(fp, "%s", comparison_type_str(d->cmpt));
245 fprintf(fp, "%s", d->cmpop);
247 /* Print the right value */
248 if ((d->right.flags & CMP_SINGLEVAL) && !(d->right.flags & CMP_DYNAMICVAL))
250 if (d->right.flags & CMP_REALVAL)
252 fprintf(fp, " %f", d->right.r[0]);
256 fprintf(fp, " %d", d->right.i[0]);
263 init_data_compare(int /* npar */, gmx_ana_selparam_t *param)
265 t_methoddata_compare *data;
268 param[2].val.u.s = &data->cmpop;
273 * Reverses a comparison operator.
275 * \param[in] type Comparison operator to reverse.
276 * \returns The correct comparison operator that equals \p type when the
277 * left and right sides are interchanged.
279 static e_comparison_t
280 reverse_comparison_type(e_comparison_t type)
284 case CMP_LESS: return CMP_GTR;
285 case CMP_LEQ: return CMP_GEQ;
286 case CMP_GTR: return CMP_LESS;
287 case CMP_GEQ: return CMP_LEQ;
294 * Initializes the value storage for comparison expression.
296 * \param[out] val Value structure to initialize.
297 * \param[in] param Parameters to use for initialization.
298 * \returns The number of values provided for the value, 0 on error.
301 init_comparison_value(t_compare_value *val, gmx_ana_selparam_t param[2])
306 if (param[0].flags & SPAR_SET)
308 val->flags |= (param[0].flags & SPAR_DYNAMIC) ? CMP_DYNAMICVAL : 0;
309 val->flags |= !(param[0].flags & SPAR_ATOMVAL) ? CMP_SINGLEVAL : 0;
311 val->i = param[0].val.u.i;
313 else if (param[1].flags & SPAR_SET)
315 val->flags |= (param[1].flags & SPAR_DYNAMIC) ? CMP_DYNAMICVAL : 0;
316 val->flags |= !(param[1].flags & SPAR_ATOMVAL) ? CMP_SINGLEVAL : 0;
317 val->flags |= CMP_REALVAL;
319 val->r = param[1].val.u.r;
331 * Converts an integer value to floating point.
333 * \param[in] n Number of values in the \p val->u array.
334 * \param[in,out] val Value to convert.
337 convert_int_real(int n, t_compare_value *val)
343 for (i = 0; i < n; ++i)
345 rv[i] = (real)val->i[i];
347 /* Free the previous value if one is present. */
350 val->flags |= CMP_REALVAL | CMP_ALLOCREAL;
354 * Converts a floating point value to integer.
356 * \param[in] n Number of values in the \p val->u array.
357 * \param[in,out] val Value to convert.
358 * \param[in] cmpt Comparison operator type.
359 * \param[in] bRight true if \p val appears on the right hand size of
361 * \returns 0 on success, EINVAL on error.
363 * The values are rounded such that the same comparison operator can be used.
366 convert_real_int(int n, t_compare_value *val, e_comparison_t cmpt, bool bRight)
373 cmpt = reverse_comparison_type(cmpt);
376 /* Round according to the comparison type */
377 for (i = 0; i < n; ++i)
383 iv[i] = static_cast<int>(std::ceil(val->r[i]));
387 iv[i] = static_cast<int>(std::floor(val->r[i]));
392 /* TODO: Implement, although it isn't typically very useful.
393 * Implementation is only a matter of proper initialization,
394 * the evaluation function can already handle this case with
395 * proper preparations. */
396 GMX_THROW(gmx::NotImplementedError("Equality comparison between dynamic integer and static real expressions not implemented"));
397 case CMP_INVALID: /* Should not be reached */
399 GMX_THROW(gmx::InternalError("Invalid comparison type"));
402 /* Free the previous value if one is present. */
405 val->flags &= ~CMP_REALVAL;
406 val->flags |= CMP_ALLOCINT;
410 init_compare(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
412 t_methoddata_compare *d = (t_methoddata_compare *)data;
415 /* Store the values */
416 n1 = init_comparison_value(&d->left, ¶m[0]);
417 n2 = init_comparison_value(&d->right, ¶m[3]);
418 /* Store the comparison type */
419 d->cmpt = comparison_type(d->cmpop);
420 if (d->cmpt == CMP_INVALID)
422 GMX_THROW(gmx::InternalError("Invalid comparison type"));
424 /* Convert the values to the same type */
425 /* TODO: Currently, there are no dynamic integer-valued selection methods,
426 * which means that only the branches with convert_int_real() will ever be
427 * taken. It should be considered whether it is necessary to support these
428 * other cases at all.
430 if ((d->left.flags & CMP_REALVAL) && !(d->right.flags & CMP_REALVAL))
432 if (d->left.flags & d->right.flags & CMP_DYNAMICVAL)
434 /* Nothing can be done */
436 else if (!(d->right.flags & CMP_DYNAMICVAL))
438 convert_int_real(n2, &d->right);
440 else /* d->left is static */
442 convert_real_int(n1, &d->left, d->cmpt, false);
445 else if (!(d->left.flags & CMP_REALVAL) && (d->right.flags & CMP_REALVAL))
447 if (d->left.flags & d->right.flags & CMP_DYNAMICVAL)
449 /* Reverse the sides to place the integer on the right */
451 d->left.r = d->right.r;
453 d->right.i = d->left.i;
455 flags = d->left.flags;
456 d->left.flags = d->right.flags;
457 d->right.flags = flags;
458 d->cmpt = reverse_comparison_type(d->cmpt);
460 else if (!(d->left.flags & CMP_DYNAMICVAL))
462 convert_int_real(n1, &d->left);
464 else /* d->right is static */
466 convert_real_int(n2, &d->right, d->cmpt, true);
472 * \param data Data to free (should point to a \c t_methoddata_compare).
474 * Frees the memory allocated for \c t_methoddata_compare.
477 free_data_compare(void *data)
479 t_methoddata_compare *d = (t_methoddata_compare *)data;
482 if (d->left.flags & CMP_ALLOCINT)
486 if (d->left.flags & CMP_ALLOCREAL)
490 if (d->right.flags & CMP_ALLOCINT)
494 if (d->right.flags & CMP_ALLOCREAL)
502 * Implementation for evaluate_compare() for integer values.
504 * \param[in] top Not used.
505 * \param[in] fr Not used.
506 * \param[in] pbc Not used.
507 * \param[in] g Evaluation index group.
508 * \param[out] out Output data structure (\p out->u.g is used).
509 * \param[in] data Should point to a \c t_methoddata_compare.
512 evaluate_compare_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
513 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
515 t_methoddata_compare *d = (t_methoddata_compare *)data;
520 GMX_UNUSED_VALUE(top);
521 GMX_UNUSED_VALUE(fr);
522 GMX_UNUSED_VALUE(pbc);
523 for (i = i1 = i2 = ig = 0; i < g->isize; ++i)
530 case CMP_INVALID: break;
531 case CMP_LESS: bAccept = a < b; break;
532 case CMP_LEQ: bAccept = a <= b; break;
533 case CMP_GTR: bAccept = a > b; break;
534 case CMP_GEQ: bAccept = a >= b; break;
535 case CMP_EQUAL: bAccept = a == b; break;
536 case CMP_NEQ: bAccept = a != b; break;
540 out->u.g->index[ig++] = g->index[i];
542 if (!(d->left.flags & CMP_SINGLEVAL))
546 if (!(d->right.flags & CMP_SINGLEVAL))
551 out->u.g->isize = ig;
555 * Implementation for evaluate_compare() if either value is non-integer.
557 * \param[in] top Not used.
558 * \param[in] fr Not used.
559 * \param[in] pbc Not used.
560 * \param[in] g Evaluation index group.
561 * \param[out] out Output data structure (\p out->u.g is used).
562 * \param[in] data Should point to a \c t_methoddata_compare.
564 * Left value is assumed to be real-valued; right value can be either.
565 * This is ensured by the initialization method.
568 evaluate_compare_real(t_topology *top, t_trxframe *fr, t_pbc *pbc,
569 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
571 t_methoddata_compare *d = (t_methoddata_compare *)data;
576 GMX_UNUSED_VALUE(top);
577 GMX_UNUSED_VALUE(fr);
578 GMX_UNUSED_VALUE(pbc);
579 for (i = i1 = i2 = ig = 0; i < g->isize; ++i)
582 b = (d->right.flags & CMP_REALVAL) ? d->right.r[i2] : d->right.i[i2];
586 case CMP_INVALID: break;
587 case CMP_LESS: bAccept = a < b; break;
588 case CMP_LEQ: bAccept = a <= b; break;
589 case CMP_GTR: bAccept = a > b; break;
590 case CMP_GEQ: bAccept = a >= b; break;
591 case CMP_EQUAL: bAccept = gmx_within_tol(a, b, GMX_REAL_EPS); break;
592 case CMP_NEQ: bAccept = !gmx_within_tol(a, b, GMX_REAL_EPS); break;
596 out->u.g->index[ig++] = g->index[i];
598 if (!(d->left.flags & CMP_SINGLEVAL))
602 if (!(d->right.flags & CMP_SINGLEVAL))
607 out->u.g->isize = ig;
611 * \param[in] top Not used.
612 * \param[in] fr Not used.
613 * \param[in] pbc Not used.
614 * \param[in] g Evaluation index group.
615 * \param[out] out Output data structure (\p out->u.g is used).
616 * \param[in] data Should point to a \c t_methoddata_compare.
619 evaluate_compare(t_topology *top, t_trxframe *fr, t_pbc *pbc,
620 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
622 t_methoddata_compare *d = (t_methoddata_compare *)data;
624 if (!((d->left.flags | d->right.flags) & CMP_REALVAL))
626 evaluate_compare_int(top, fr, pbc, g, out, data);
630 evaluate_compare_real(top, fr, pbc, g, out, data);