2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016,2017, 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
46 #include "gromacs/math/utilities.h"
47 #include "gromacs/utility/arraysize.h"
48 #include "gromacs/utility/basedefinitions.h"
49 #include "gromacs/utility/exceptions.h"
50 #include "gromacs/utility/smalloc.h"
52 #include "selmethod.h"
54 /** Defines the comparison operator for comparison expressions. */
57 CMP_INVALID, /**< Indicates an error */
62 CMP_EQUAL, /**< '==' */
66 /** The operand has a single value. */
67 #define CMP_SINGLEVAL 1
68 /** The operand value is dynamic. */
69 #define CMP_DYNAMICVAL 2
70 /** The value is real. */
72 /** The integer array is allocated. */
73 #define CMP_ALLOCINT 16
74 /** The real array is allocated. */
75 #define CMP_ALLOCREAL 32
78 * Data structure for comparison expression operand values.
82 /** Flags that describe the type of the operand. */
84 /** (Array of) integer value(s). */
86 /** (Array of) real value(s). */
91 * Data structure for comparison expression evaluation.
95 /** Comparison operator as a string. */
97 /** Comparison operator type. */
100 t_compare_value left;
102 t_compare_value right;
103 } t_methoddata_compare;
106 * Allocates data for comparison expression evaluation.
108 * \param[in] npar Not used (should be 5).
109 * \param[in,out] param Method parameters (should point to a copy of
110 * \ref smparams_compare).
111 * \returns Pointer to the allocated data (\c t_methoddata_compare).
113 * Allocates memory for a \c t_methoddata_compare structure.
116 init_data_compare(int npar, gmx_ana_selparam_t *param);
118 * Initializes data for comparison expression evaluation.
120 * \param[in] top Not used.
121 * \param[in] npar Not used (should be 5).
122 * \param[in] param Method parameters (should point to \ref smparams_compare).
123 * \param[in] data Should point to a \c t_methoddata_compare.
124 * \returns 0 if the input data is valid, -1 on error.
127 init_compare(const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data);
128 /** Frees the memory allocated for comparison expression evaluation. */
130 free_data_compare(void *data);
132 * Evaluates comparison expressions.
134 * \param[in] context Not used.
135 * \param[in] g Evaluation index group.
136 * \param[out] out Output data structure (\p out->u.g is used).
137 * \param[in] data Should point to a \c t_methoddata_compare.
140 evaluate_compare(const gmx::SelMethodEvalContext &context,
141 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
143 /** Parameters for comparison expression evaluation. */
144 static gmx_ana_selparam_t smparams_compare[] = {
145 {"int1", {INT_VALUE, -1, {nullptr}}, nullptr,
146 SPAR_OPTIONAL | SPAR_DYNAMIC | SPAR_ATOMVAL},
147 {"real1", {REAL_VALUE, -1, {nullptr}}, nullptr,
148 SPAR_OPTIONAL | SPAR_DYNAMIC | SPAR_ATOMVAL},
149 {"op", {STR_VALUE, 1, {nullptr}}, nullptr, 0},
150 {"int2", {INT_VALUE, -1, {nullptr}}, nullptr,
151 SPAR_OPTIONAL | SPAR_DYNAMIC | SPAR_ATOMVAL},
152 {"real2", {REAL_VALUE, -1, {nullptr}}, nullptr,
153 SPAR_OPTIONAL | SPAR_DYNAMIC | SPAR_ATOMVAL},
156 /** \internal Selection method data for comparison expression evaluation. */
157 gmx_ana_selmethod_t sm_compare = {
158 "cmp", GROUP_VALUE, SMETH_SINGLEVAL,
159 asize(smparams_compare), smparams_compare,
168 {nullptr, nullptr, 0, nullptr},
172 * Returns a \c e_comparison_t value corresponding to an operator.
174 * \param[in] str String to process.
175 * \returns The comparison type corresponding to the first one or two
176 * characters of \p str.
178 * \p str can contain any number of characters; only the first two
180 * If the beginning of \p str does not match any of the recognized types,
181 * \ref CMP_INVALID is returned.
183 static e_comparison_t
184 comparison_type(char *str)
188 case '<': return (str[1] == '=') ? CMP_LEQ : CMP_LESS;
189 case '>': return (str[1] == '=') ? CMP_GEQ : CMP_GTR;
190 case '=': return (str[1] == '=') ? CMP_EQUAL : CMP_INVALID;
191 case '!': return (str[1] == '=') ? CMP_NEQ : CMP_INVALID;
197 * Returns a string corresponding to a \c e_comparison_t value.
199 * \param[in] cmpt Comparison type to convert.
200 * \returns Pointer to a string that corresponds to \p cmpt.
202 * The return value points to a string constant and should not be \p free'd.
204 * The function returns NULL if \p cmpt is not one of the valid values.
207 comparison_type_str(e_comparison_t cmpt)
209 const char *p = nullptr;
212 case CMP_INVALID: p = "INVALID"; break;
213 case CMP_LESS: p = "<"; break;
214 case CMP_LEQ: p = "<="; break;
215 case CMP_GTR: p = ">"; break;
216 case CMP_GEQ: p = ">="; break;
217 case CMP_EQUAL: p = "=="; break;
218 case CMP_NEQ: p = "!="; break;
219 // No default clause so we intentionally get compiler errors
220 // if new selection choices are added later.
226 * \param[in] fp File to receive the output.
227 * \param[in] data Should point to a \c t_methoddata_compare.
230 _gmx_selelem_print_compare_info(FILE *fp, void *data)
232 t_methoddata_compare *d = (t_methoddata_compare *)data;
235 /* Print the left value */
236 if ((d->left.flags & CMP_SINGLEVAL) && !(d->left.flags & CMP_DYNAMICVAL))
238 if (d->left.flags & CMP_REALVAL)
240 fprintf(fp, "%f ", d->left.r[0]);
244 fprintf(fp, "%d ", d->left.i[0]);
247 /* Print the operator */
248 if (d->cmpt != CMP_INVALID)
250 fprintf(fp, "%s", comparison_type_str(d->cmpt));
254 fprintf(fp, "%s", d->cmpop);
256 /* Print the right value */
257 if ((d->right.flags & CMP_SINGLEVAL) && !(d->right.flags & CMP_DYNAMICVAL))
259 if (d->right.flags & CMP_REALVAL)
261 fprintf(fp, " %f", d->right.r[0]);
265 fprintf(fp, " %d", d->right.i[0]);
272 init_data_compare(int /* npar */, gmx_ana_selparam_t *param)
274 t_methoddata_compare *data;
277 param[2].val.u.s = &data->cmpop;
282 * Reverses a comparison operator.
284 * \param[in] type Comparison operator to reverse.
285 * \returns The correct comparison operator that equals \p type when the
286 * left and right sides are interchanged.
288 static e_comparison_t
289 reverse_comparison_type(e_comparison_t type)
293 case CMP_LESS: return CMP_GTR;
294 case CMP_LEQ: return CMP_GEQ;
295 case CMP_GTR: return CMP_LESS;
296 case CMP_GEQ: return CMP_LEQ;
303 * Initializes the value storage for comparison expression.
305 * \param[out] val Value structure to initialize.
306 * \param[in] param Parameters to use for initialization.
307 * \returns The number of values provided for the value, 0 on error.
310 init_comparison_value(t_compare_value *val, gmx_ana_selparam_t param[2])
315 if (param[0].flags & SPAR_SET)
317 val->flags |= (param[0].flags & SPAR_DYNAMIC) ? CMP_DYNAMICVAL : 0;
318 val->flags |= !(param[0].flags & SPAR_ATOMVAL) ? CMP_SINGLEVAL : 0;
320 val->i = param[0].val.u.i;
322 else if (param[1].flags & SPAR_SET)
324 val->flags |= (param[1].flags & SPAR_DYNAMIC) ? CMP_DYNAMICVAL : 0;
325 val->flags |= !(param[1].flags & SPAR_ATOMVAL) ? CMP_SINGLEVAL : 0;
326 val->flags |= CMP_REALVAL;
328 val->r = param[1].val.u.r;
340 * Converts an integer value to floating point.
342 * \param[in] n Number of values in the \p val->u array.
343 * \param[in,out] val Value to convert.
346 convert_int_real(int n, t_compare_value *val)
352 for (i = 0; i < n; ++i)
354 rv[i] = (real)val->i[i];
356 /* Free the previous value if one is present. */
359 val->flags |= CMP_REALVAL | CMP_ALLOCREAL;
363 * Converts a floating point value to integer.
365 * \param[in] n Number of values in the \p val->u array.
366 * \param[in,out] val Value to convert.
367 * \param[in] cmpt Comparison operator type.
368 * \param[in] bRight true if \p val appears on the right hand size of
370 * \returns 0 on success, EINVAL on error.
372 * The values are rounded such that the same comparison operator can be used.
375 convert_real_int(int n, t_compare_value *val, e_comparison_t cmpt, bool bRight)
382 cmpt = reverse_comparison_type(cmpt);
385 /* Round according to the comparison type */
386 for (i = 0; i < n; ++i)
392 iv[i] = static_cast<int>(std::ceil(val->r[i]));
396 iv[i] = static_cast<int>(std::floor(val->r[i]));
401 /* TODO: Implement, although it isn't typically very useful.
402 * Implementation is only a matter of proper initialization,
403 * the evaluation function can already handle this case with
404 * proper preparations. */
405 GMX_THROW(gmx::NotImplementedError("Equality comparison between dynamic integer and static real expressions not implemented"));
406 case CMP_INVALID: /* Should not be reached */
408 GMX_THROW(gmx::InternalError("Invalid comparison type"));
411 /* Free the previous value if one is present. */
414 val->flags &= ~CMP_REALVAL;
415 val->flags |= CMP_ALLOCINT;
419 init_compare(const gmx_mtop_t * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
421 t_methoddata_compare *d = (t_methoddata_compare *)data;
424 /* Store the values */
425 n1 = init_comparison_value(&d->left, ¶m[0]);
426 n2 = init_comparison_value(&d->right, ¶m[3]);
427 /* Store the comparison type */
428 d->cmpt = comparison_type(d->cmpop);
429 if (d->cmpt == CMP_INVALID)
431 GMX_THROW(gmx::InternalError("Invalid comparison type"));
433 /* Convert the values to the same type */
434 /* TODO: Currently, there are no dynamic integer-valued selection methods,
435 * which means that only the branches with convert_int_real() will ever be
436 * taken. It should be considered whether it is necessary to support these
437 * other cases at all.
439 if ((d->left.flags & CMP_REALVAL) && !(d->right.flags & CMP_REALVAL))
441 if (d->left.flags & d->right.flags & CMP_DYNAMICVAL)
443 /* Nothing can be done */
445 else if (!(d->right.flags & CMP_DYNAMICVAL))
447 convert_int_real(n2, &d->right);
449 else /* d->left is static */
451 convert_real_int(n1, &d->left, d->cmpt, false);
454 else if (!(d->left.flags & CMP_REALVAL) && (d->right.flags & CMP_REALVAL))
456 if (d->left.flags & d->right.flags & CMP_DYNAMICVAL)
458 /* Reverse the sides to place the integer on the right */
460 d->left.r = d->right.r;
461 d->right.r = nullptr;
462 d->right.i = d->left.i;
464 flags = d->left.flags;
465 d->left.flags = d->right.flags;
466 d->right.flags = flags;
467 d->cmpt = reverse_comparison_type(d->cmpt);
469 else if (!(d->left.flags & CMP_DYNAMICVAL))
471 convert_int_real(n1, &d->left);
473 else /* d->right is static */
475 convert_real_int(n2, &d->right, d->cmpt, true);
481 * \param data Data to free (should point to a \c t_methoddata_compare).
483 * Frees the memory allocated for \c t_methoddata_compare.
486 free_data_compare(void *data)
488 t_methoddata_compare *d = (t_methoddata_compare *)data;
491 if (d->left.flags & CMP_ALLOCINT)
495 if (d->left.flags & CMP_ALLOCREAL)
499 if (d->right.flags & CMP_ALLOCINT)
503 if (d->right.flags & CMP_ALLOCREAL)
511 * Implementation for evaluate_compare() for integer values.
513 * \param[in] g Evaluation index group.
514 * \param[out] out Output data structure (\p out->u.g is used).
515 * \param[in] data Should point to a \c t_methoddata_compare.
518 evaluate_compare_int(gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
520 t_methoddata_compare *d = (t_methoddata_compare *)data;
525 for (i = i1 = i2 = ig = 0; i < g->isize; ++i)
532 case CMP_INVALID: break;
533 case CMP_LESS: bAccept = a < b; break;
534 case CMP_LEQ: bAccept = a <= b; break;
535 case CMP_GTR: bAccept = a > b; break;
536 case CMP_GEQ: bAccept = a >= b; break;
537 case CMP_EQUAL: bAccept = a == b; break;
538 case CMP_NEQ: bAccept = a != b; break;
542 out->u.g->index[ig++] = g->index[i];
544 if (!(d->left.flags & CMP_SINGLEVAL))
548 if (!(d->right.flags & CMP_SINGLEVAL))
553 out->u.g->isize = ig;
557 * Implementation for evaluate_compare() if either value is non-integer.
559 * \param[in] g Evaluation index group.
560 * \param[out] out Output data structure (\p out->u.g is used).
561 * \param[in] data Should point to a \c t_methoddata_compare.
563 * Left value is assumed to be real-valued; right value can be either.
564 * This is ensured by the initialization method.
567 evaluate_compare_real(gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
569 t_methoddata_compare *d = (t_methoddata_compare *)data;
574 for (i = i1 = i2 = ig = 0; i < g->isize; ++i)
577 b = (d->right.flags & CMP_REALVAL) ? d->right.r[i2] : d->right.i[i2];
581 case CMP_INVALID: break;
582 case CMP_LESS: bAccept = a < b; break;
583 case CMP_LEQ: bAccept = a <= b; break;
584 case CMP_GTR: bAccept = a > b; break;
585 case CMP_GEQ: bAccept = a >= b; break;
586 case CMP_EQUAL: bAccept = gmx_within_tol(a, b, GMX_REAL_EPS); break;
587 case CMP_NEQ: bAccept = !gmx_within_tol(a, b, GMX_REAL_EPS); break;
591 out->u.g->index[ig++] = g->index[i];
593 if (!(d->left.flags & CMP_SINGLEVAL))
597 if (!(d->right.flags & CMP_SINGLEVAL))
602 out->u.g->isize = ig;
606 evaluate_compare(const gmx::SelMethodEvalContext & /*context*/,
607 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
609 t_methoddata_compare *d = (t_methoddata_compare *)data;
611 if (!((d->left.flags | d->right.flags) & CMP_REALVAL))
613 evaluate_compare_int(g, out, data);
617 evaluate_compare_real(g, out, data);