3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2009, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
32 * \brief Implementation of the \p same selection method.
44 #include <selmethod.h>
47 #include "parsetree.h"
51 * Data structure for the \p same selection method.
53 * To avoid duplicate initialization code, the same data structure is used
54 * for matching both integer and string keywords; hence the unions.
58 /** Value for each atom to match. */
66 * Number of values in the \p as array.
68 * For string values, this is actually the number of values in the
69 * \p as_s_sorted array.
72 /** Values to match against. */
80 * Separate array for sorted \p as.s array.
82 * The array of strings returned as the output value of a parameter should
83 * not be messed with to avoid memory corruption (the pointers in the array
84 * may be reused for several evaluations), so we keep our own copy for
88 /** Whether simple matching can be used. */
92 /** Allocates data for the \p same selection method. */
94 init_data_same(int npar
, gmx_ana_selparam_t
*param
);
95 /** Initializes the \p same selection method. */
97 init_same(t_topology
*top
, int npar
, gmx_ana_selparam_t
*param
, void *data
);
98 /** Frees the data allocated for the \p same selection method. */
100 free_data_same(void *data
);
101 /** Initializes the evaluation of the \p same selection method for a frame. */
103 init_frame_same_int(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
, void *data
);
104 /** Evaluates the \p same selection method. */
106 evaluate_same_int(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
,
107 gmx_ana_index_t
*g
, gmx_ana_selvalue_t
*out
, void *data
);
108 /** Initializes the evaluation of the \p same selection method for a frame. */
110 init_frame_same_str(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
, void *data
);
111 /** Evaluates the \p same selection method. */
113 evaluate_same_str(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
,
114 gmx_ana_index_t
*g
, gmx_ana_selvalue_t
*out
, void *data
);
116 /** Parameters for the \p same selection method. */
117 static gmx_ana_selparam_t smparams_same_int
[] = {
118 {NULL
, {INT_VALUE
, -1, {NULL
}}, NULL
, SPAR_DYNAMIC
| SPAR_ATOMVAL
},
119 {"as", {INT_VALUE
, -1, {NULL
}}, NULL
, SPAR_DYNAMIC
| SPAR_VARNUM
},
122 /** Parameters for the \p same selection method. */
123 static gmx_ana_selparam_t smparams_same_str
[] = {
124 {NULL
, {STR_VALUE
, -1, {NULL
}}, NULL
, SPAR_DYNAMIC
| SPAR_ATOMVAL
},
125 {"as", {STR_VALUE
, -1, {NULL
}}, NULL
, SPAR_DYNAMIC
| SPAR_VARNUM
},
128 /** Help text for the \p same selection method. */
129 static const char *help_same
[] = {
130 "EXTENDING SELECTIONS[PAR]",
132 "[TT]same KEYWORD as ATOM_EXPR[tt][PAR]",
134 "The keyword [TT]same[tt] can be used to select all atoms for which",
135 "the given [TT]KEYWORD[tt] matches any of the atoms in [TT]ATOM_EXPR[tt].",
136 "Keywords that evaluate to integer or string values are supported.",
139 /*! \internal \brief Selection method data for the \p same method. */
140 gmx_ana_selmethod_t sm_same
= {
141 "same", GROUP_VALUE
, 0,
142 asize(smparams_same_int
), smparams_same_int
,
148 &init_frame_same_int
,
151 {"same KEYWORD as ATOM_EXPR", asize(help_same
), help_same
},
155 * Selection method data for the \p same method.
157 * This selection method is used for matching string keywords. The parser
158 * never sees this method; _gmx_selelem_custom_init_same() replaces sm_same
159 * with this method in cases where it is required.
161 static gmx_ana_selmethod_t sm_same_str
= {
162 "same", GROUP_VALUE
, SMETH_SINGLEVAL
,
163 asize(smparams_same_str
), smparams_same_str
,
169 &init_frame_same_str
,
172 {"same KEYWORD as ATOM_EXPR", asize(help_same
), help_same
},
176 * \param[in] npar Not used (should be 2).
177 * \param[in,out] param Method parameters (should point to
178 * \ref smparams_same).
179 * \returns Pointer to the allocated data (\ref t_methoddata_same).
182 init_data_same(int npar
, gmx_ana_selparam_t
*param
)
184 t_methoddata_same
*data
;
187 data
->as_s_sorted
= NULL
;
188 param
[1].nvalptr
= &data
->nas
;
193 * \param[in,out] method The method to initialize.
194 * \param[in,out] params Pointer to the first parameter.
195 * \param[in] scanner Scanner data structure.
196 * \returns 0 on success, a non-zero error code on error.
198 * If \p *method is not a \c same method, this function returns zero
202 _gmx_selelem_custom_init_same(gmx_ana_selmethod_t
**method
,
203 t_selexpr_param
*params
,
206 gmx_ana_selmethod_t
*kwmethod
;
208 t_selexpr_param
*param
;
212 /* Do nothing if this is not a same method. */
213 if (!*method
|| (*method
)->name
!= sm_same
.name
)
218 if (params
->nval
!= 1 || !params
->value
->bExpr
219 || params
->value
->u
.expr
->type
!= SEL_EXPRESSION
)
221 _gmx_selparser_error("error: 'same' should be followed by a single keyword");
224 kwmethod
= params
->value
->u
.expr
->u
.expr
.method
;
226 if (kwmethod
->type
== STR_VALUE
)
228 *method
= &sm_same_str
;
231 /* We do custom processing with the second parameter, so remove it from
232 * the params list, but save the name for later. */
233 param
= params
->next
;
237 /* Create a second keyword evaluation element for the keyword given as
238 * the first parameter, evaluating the keyword in the group given by the
239 * second parameter. */
240 rc
= _gmx_sel_init_keyword_evaluator(&kwelem
, kwmethod
, param
, scanner
);
246 /* Replace the second parameter with one with a value from \p kwelem. */
247 param
= _gmx_selexpr_create_param(pname
);
249 param
->value
= _gmx_selexpr_create_value_expr(kwelem
);
250 params
->next
= param
;
255 * \param top Not used.
256 * \param npar Not used (should be 2).
257 * \param param Initialized method parameters (should point to a copy of
258 * \ref smparams_same).
259 * \param data Pointer to \ref t_methoddata_same to initialize.
260 * \returns 0 on success, -1 on failure.
263 init_same(t_topology
*top
, int npar
, gmx_ana_selparam_t
*param
, void *data
)
265 t_methoddata_same
*d
= (t_methoddata_same
*)data
;
267 d
->val
.ptr
= param
[0].val
.u
.ptr
;
268 d
->as
.ptr
= param
[1].val
.u
.ptr
;
269 if (param
[1].val
.type
== STR_VALUE
)
271 snew(d
->as_s_sorted
, d
->nas
);
273 if (!(param
[0].flags
& SPAR_ATOMVAL
))
275 fprintf(stderr
, "ERROR: the same selection keyword combined with a "
276 "non-keyword does not make sense\n");
283 * \param data Data to free (should point to a \ref t_methoddata_same).
286 free_data_same(void *data
)
288 t_methoddata_same
*d
= (t_methoddata_same
*)data
;
290 sfree(d
->as_s_sorted
);
294 * Helper function for comparison of two integers.
297 cmp_int(const void *a
, const void *b
)
299 if (*(int *)a
< *(int *)b
)
303 if (*(int *)a
> *(int *)b
)
311 * \param[in] top Not used.
312 * \param[in] fr Current frame.
313 * \param[in] pbc PBC structure.
314 * \param data Should point to a \ref t_methoddata_same.
315 * \returns 0 on success, a non-zero error code on error.
317 * Sorts the \c data->as.i array and removes identical values for faster and
321 init_frame_same_int(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
, void *data
)
323 t_methoddata_same
*d
= (t_methoddata_same
*)data
;
326 /* Collapse adjacent values, and check whether the array is sorted. */
328 for (i
= 1, j
= 0; i
< d
->nas
; ++i
)
330 if (d
->as
.i
[i
] != d
->as
.i
[j
])
332 if (d
->as
.i
[i
] < d
->as
.i
[j
])
337 d
->as
.i
[j
] = d
->as
.i
[i
];
344 qsort(d
->as
.i
, d
->nas
, sizeof(d
->as
.i
[0]), &cmp_int
);
345 /* More identical values may become adjacent after sorting. */
346 for (i
= 1, j
= 0; i
< d
->nas
; ++i
)
348 if (d
->as
.i
[i
] != d
->as
.i
[j
])
351 d
->as
.i
[j
] = d
->as
.i
[i
];
360 * See sel_updatefunc() for description of the parameters.
361 * \p data should point to a \c t_methoddata_same.
363 * Calculates which values in \c data->val.i can be found in \c data->as.i
364 * (assumed sorted), and writes the corresponding atoms to output.
365 * If \c data->val is sorted, uses a linear scan of both arrays, otherwise a
366 * binary search of \c data->as is performed for each block of values in
370 evaluate_same_int(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
,
371 gmx_ana_index_t
*g
, gmx_ana_selvalue_t
*out
, void *data
)
373 t_methoddata_same
*d
= (t_methoddata_same
*)data
;
382 /* If we are sorted, we can do a simple linear scan. */
383 while (i
< d
->nas
&& d
->as
.i
[i
] < d
->val
.i
[j
]) ++i
;
387 /* If not, we must do a binary search of all the values. */
394 int itry
= (i1
+ i2
) / 2;
395 if (d
->as
.i
[itry
] <= d
->val
.i
[j
])
404 i
= (d
->as
.i
[i1
] == d
->val
.i
[j
] ? i1
: d
->nas
);
406 /* Check whether the value was found in the as list. */
407 if (i
== d
->nas
|| d
->as
.i
[i
] != d
->val
.i
[j
])
409 /* If not, skip all atoms with the same value. */
410 int tmpval
= d
->val
.i
[j
];
412 while (j
< g
->isize
&& d
->val
.i
[j
] == tmpval
) ++j
;
416 /* Copy all the atoms with this value to the output. */
417 while (j
< g
->isize
&& d
->val
.i
[j
] == d
->as
.i
[i
])
419 out
->u
.g
->index
[out
->u
.g
->isize
++] = g
->index
[j
];
423 if (j
< g
->isize
&& d
->val
.i
[j
] < d
->val
.i
[j
- 1])
432 * Helper function for comparison of two strings.
435 cmp_str(const void *a
, const void *b
)
437 return strcmp(*(char **)a
, *(char **)b
);
441 * \param[in] top Not used.
442 * \param[in] fr Current frame.
443 * \param[in] pbc PBC structure.
444 * \param data Should point to a \ref t_methoddata_same.
445 * \returns 0 on success, a non-zero error code on error.
447 * Sorts the \c data->as.s array and removes identical values for faster and
451 init_frame_same_str(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
, void *data
)
453 t_methoddata_same
*d
= (t_methoddata_same
*)data
;
456 /* Collapse adjacent values.
457 * For strings, it's unlikely that the values would be sorted originally,
458 * so set bSorted always to FALSE. */
460 d
->as_s_sorted
[0] = d
->as
.s
[0];
461 for (i
= 1, j
= 0; i
< d
->nas
; ++i
)
463 if (strcmp(d
->as
.s
[i
], d
->as_s_sorted
[j
]) != 0)
466 d
->as_s_sorted
[j
] = d
->as
.s
[i
];
471 qsort(d
->as_s_sorted
, d
->nas
, sizeof(d
->as_s_sorted
[0]), &cmp_str
);
472 /* More identical values may become adjacent after sorting. */
473 for (i
= 1, j
= 0; i
< d
->nas
; ++i
)
475 if (strcmp(d
->as_s_sorted
[i
], d
->as_s_sorted
[j
]) != 0)
478 d
->as_s_sorted
[j
] = d
->as_s_sorted
[i
];
486 * See sel_updatefunc() for description of the parameters.
487 * \p data should point to a \c t_methoddata_same.
489 * Calculates which strings in \c data->val.s can be found in \c data->as.s
490 * (assumed sorted), and writes the corresponding atoms to output.
491 * A binary search of \c data->as is performed for each block of values in
495 evaluate_same_str(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
,
496 gmx_ana_index_t
*g
, gmx_ana_selvalue_t
*out
, void *data
)
498 t_methoddata_same
*d
= (t_methoddata_same
*)data
;
505 /* Do a binary search of the strings. */
507 ptr
= bsearch(&d
->val
.s
[j
], d
->as_s_sorted
, d
->nas
,
508 sizeof(d
->as_s_sorted
[0]), &cmp_str
);
509 /* Check whether the value was found in the as list. */
512 /* If not, skip all atoms with the same value. */
513 const char *tmpval
= d
->val
.s
[j
];
515 while (j
< g
->isize
&& strcmp(d
->val
.s
[j
], tmpval
) == 0) ++j
;
519 const char *tmpval
= d
->val
.s
[j
];
520 /* Copy all the atoms with this value to the output. */
521 while (j
< g
->isize
&& strcmp(d
->val
.s
[j
], tmpval
) == 0)
523 out
->u
.g
->index
[out
->u
.g
->isize
++] = g
->index
[j
];