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
33 * Implements the \p merge and \p plus selection modifiers.
35 * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36 * \ingroup module_selection
46 #include "gromacs/selection/position.h"
47 #include "gromacs/selection/selmethod.h"
48 #include "gromacs/utility/exceptions.h"
51 * Data structure for the merging selection modifiers.
55 /** Input positions. */
57 /** Other input positions. */
59 /** Group to store the output atom indices. */
61 /** Stride for merging (\c stride values from \c p1 for each in \c p2). */
65 /** Allocates data for the merging selection modifiers. */
67 init_data_merge(int npar, gmx_ana_selparam_t *param);
68 /** Initializes data for the merging selection modifiers. */
70 init_merge(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
71 /** Initializes output for the \p merge selection modifier. */
73 init_output_merge(t_topology *top, gmx_ana_selvalue_t *out, void *data);
74 /** Initializes output for the \p plus selection modifier. */
76 init_output_plus(t_topology *top, gmx_ana_selvalue_t *out, void *data);
77 /** Frees the memory allocated for the merging selection modifiers. */
79 free_data_merge(void *data);
80 /** Evaluates the \p merge selection modifier. */
82 evaluate_merge(t_topology *top, t_trxframe *fr, t_pbc *pbc,
83 gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data);
84 /** Evaluates the \p plus selection modifier. */
86 evaluate_plus(t_topology *top, t_trxframe *fr, t_pbc *pbc,
87 gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data);
89 /** Parameters for the merging selection modifiers. */
90 static gmx_ana_selparam_t smparams_merge[] = {
91 {NULL, {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
92 {NULL, {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
93 {"stride", {INT_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
96 /** Help text for the merging selection modifiers. */
97 static const char *help_merge[] = {
98 "MERGING SELECTIONS[PAR]",
100 "[TT]POSEXPR merge POSEXPR [stride INT][tt][BR]",
101 "[TT]POSEXPR merge POSEXPR [merge POSEXPR ...][tt][BR]",
102 "[TT]POSEXPR plus POSEXPR [plus POSEXPR ...][tt][PAR]",
104 "Basic selection keywords can only create selections where each atom",
105 "occurs at most once. The [TT]merge[tt] and [TT]plus[tt] selection",
106 "keywords can be used to work around this limitation. Both create",
107 "a selection that contains the positions from all the given position",
108 "expressions, even if they contain duplicates.",
109 "The difference between the two is that [TT]merge[tt] expects two or more",
110 "selections with the same number of positions, and the output contains",
111 "the input positions selected from each expression in turn, i.e.,",
112 "the output is like A1 B1 A2 B2 and so on. It is also possible to merge",
113 "selections of unequal size as long as the size of the first is a",
114 "multiple of the second one. The [TT]stride[tt] parameter can be used",
115 "to explicitly provide this multiplicity.",
116 "[TT]plus[tt] simply concatenates the positions after each other, and",
117 "can work also with selections of different sizes.",
118 "These keywords are valid only at the selection level, not in any",
122 /** \internal Selection method data for the \p plus modifier. */
123 gmx_ana_selmethod_t sm_merge = {
124 "merge", POS_VALUE, SMETH_MODIFIER,
125 asize(smparams_merge), smparams_merge,
134 {"merge POSEXPR", asize(help_merge), help_merge},
137 /** \internal Selection method data for the \p plus modifier. */
138 gmx_ana_selmethod_t sm_plus = {
139 "plus", POS_VALUE, SMETH_MODIFIER,
140 asize(smparams_merge)-1, smparams_merge,
149 {"plus POSEXPR", asize(help_merge), help_merge},
153 * \param[in] npar Should be 2 for \c plus and 3 for \c merge.
154 * \param[in,out] param Method parameters (should point to a copy of
155 * \ref smparams_merge).
156 * \returns Pointer to the allocated data (\p t_methoddata_merge).
158 * Allocates memory for a \p t_methoddata_merge structure.
161 init_data_merge(int npar, gmx_ana_selparam_t *param)
163 t_methoddata_merge *data;
167 param[0].val.u.p = &data->p1;
168 param[1].val.u.p = &data->p2;
171 param[2].val.u.i = &data->stride;
177 * \param[in] top Not used.
178 * \param[in] npar Not used (should be 2 or 3).
179 * \param[in] param Method parameters (should point to \ref smparams_merge).
180 * \param[in] data Should point to a \p t_methoddata_merge.
181 * \returns 0 if everything is successful, -1 on error.
184 init_merge(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
186 t_methoddata_merge *d = (t_methoddata_merge *)data;
190 GMX_THROW(gmx::InvalidInputError("Stride for merging should be positive"));
192 /* If no stride given, deduce it from the input sizes */
195 d->stride = d->p1.nr / d->p2.nr;
197 if (d->p1.nr != d->stride*d->p2.nr)
199 GMX_THROW(gmx::InconsistentInputError("The number of positions to be merged are not compatible"));
201 /* We access the m.b.nra field instead of g->isize in the position
202 * data structures to handle cases where g is NULL
203 * (this occurs with constant positions. */
204 gmx_ana_index_reserve(&d->g, d->p1.m.b.nra + d->p2.m.b.nra);
205 d->g.isize = d->p1.m.b.nra + d->p2.m.b.nra;
209 * Does common initialization to all merging modifiers.
211 * \param[in] top Topology data structure.
212 * \param[in,out] out Pointer to output data structure.
213 * \param[in,out] data Should point to \c t_methoddata_merge.
216 init_output_common(t_topology *top, gmx_ana_selvalue_t *out, void *data)
218 t_methoddata_merge *d = (t_methoddata_merge *)data;
220 if (d->p1.m.type != d->p2.m.type)
222 /* TODO: Maybe we could pick something else here? */
223 out->u.p->m.type = INDEX_UNKNOWN;
227 out->u.p->m.type = d->p1.m.type;
229 gmx_ana_pos_reserve(out->u.p, d->p1.nr + d->p2.nr, d->g.isize);
232 gmx_ana_pos_reserve_velocities(out->u.p);
236 gmx_ana_pos_reserve_forces(out->u.p);
238 gmx_ana_pos_set_evalgrp(out->u.p, &d->g);
239 gmx_ana_pos_empty_init(out->u.p);
244 * \param[in] top Topology data structure.
245 * \param[in,out] out Pointer to output data structure.
246 * \param[in,out] data Should point to \c t_methoddata_merge.
249 init_output_merge(t_topology *top, gmx_ana_selvalue_t *out, void *data)
251 t_methoddata_merge *d = (t_methoddata_merge *)data;
254 init_output_common(top, out, data);
255 for (i = 0; i < d->p2.nr; ++i)
257 for (j = 0; j < d->stride; ++j)
259 gmx_ana_pos_append_init(out->u.p, &d->g, &d->p1, d->stride*i+j);
261 gmx_ana_pos_append_init(out->u.p, &d->g, &d->p2, i);
266 * \param[in] top Topology data structure.
267 * \param[in,out] out Pointer to output data structure.
268 * \param[in,out] data Should point to \c t_methoddata_merge.
271 init_output_plus(t_topology *top, gmx_ana_selvalue_t *out, void *data)
273 t_methoddata_merge *d = (t_methoddata_merge *)data;
276 init_output_common(top, out, data);
277 for (i = 0; i < d->p1.nr; ++i)
279 gmx_ana_pos_append_init(out->u.p, &d->g, &d->p1, i);
281 for (i = 0; i < d->p2.nr; ++i)
283 gmx_ana_pos_append_init(out->u.p, &d->g, &d->p2, i);
288 * \param data Data to free (should point to a \p t_methoddata_merge).
290 * Frees the memory allocated for \c t_methoddata_merge.
293 free_data_merge(void *data)
295 t_methoddata_merge *d = (t_methoddata_merge *)data;
297 gmx_ana_index_deinit(&d->g);
301 * \param[in] top Not used.
302 * \param[in] fr Not used.
303 * \param[in] pbc Not used.
304 * \param[in] p Positions to merge (should point to \p data->p1).
305 * \param[out] out Output data structure (\p out->u.p is used).
306 * \param[in] data Should point to a \p t_methoddata_merge.
309 evaluate_merge(t_topology *top, t_trxframe *fr, t_pbc *pbc,
310 gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data)
312 t_methoddata_merge *d = (t_methoddata_merge *)data;
316 if (d->p1.nr != d->stride*d->p2.nr)
318 GMX_THROW(gmx::InconsistentInputError("The number of positions to be merged are not compatible"));
321 gmx_ana_pos_empty(out->u.p);
322 for (i = 0; i < d->p2.nr; ++i)
324 for (j = 0; j < d->stride; ++j)
326 refid = d->p1.m.refid[d->stride*i+j];
329 refid = (d->stride+1) * (refid / d->stride) + (refid % d->stride);
331 gmx_ana_pos_append(out->u.p, &d->g, &d->p1, d->stride*i+j, refid);
333 refid = (d->stride+1)*d->p2.m.refid[i]+d->stride;
334 gmx_ana_pos_append(out->u.p, &d->g, &d->p2, i, refid);
336 gmx_ana_pos_append_finish(out->u.p);
340 * \param[in] top Not used.
341 * \param[in] fr Not used.
342 * \param[in] pbc Not used.
343 * \param[in] p Positions to merge (should point to \p data->p1).
344 * \param[out] out Output data structure (\p out->u.p is used).
345 * \param[in] data Should point to a \p t_methoddata_merge.
348 evaluate_plus(t_topology *top, t_trxframe *fr, t_pbc *pbc,
349 gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data)
351 t_methoddata_merge *d = (t_methoddata_merge *)data;
356 gmx_ana_pos_empty(out->u.p);
357 for (i = 0; i < d->p1.nr; ++i)
359 refid = d->p1.m.refid[i];
360 gmx_ana_pos_append(out->u.p, &d->g, &d->p1, i, refid);
362 for (i = 0; i < d->p2.nr; ++i)
364 refid = d->p2.m.refid[i];
367 refid += d->p1.m.b.nr;
369 gmx_ana_pos_append(out->u.p, &d->g, &d->p2, i, refid);
371 gmx_ana_pos_append_finish(out->u.p);