2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015, 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 the \p merge and \p plus selection modifiers.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/selection/position.h"
47 #include "gromacs/utility/basedefinitions.h"
48 #include "gromacs/utility/exceptions.h"
49 #include "gromacs/utility/smalloc.h"
51 #include "selmethod.h"
54 * Data structure for the merging selection modifiers.
58 /** Input positions. */
60 /** Other input positions. */
62 /** Stride for merging (\c stride values from \c p1 for each in \c p2). */
66 /** Allocates data for the merging selection modifiers. */
68 init_data_merge(int npar, gmx_ana_selparam_t *param);
70 * Initializes data for the merging selection modifiers.
72 * \param[in] top Not used.
73 * \param[in] npar Not used (should be 2 or 3).
74 * \param[in] param Method parameters (should point to \ref smparams_merge).
75 * \param[in] data Should point to a \p t_methoddata_merge.
76 * \returns 0 if everything is successful, -1 on error.
79 init_merge(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
80 /** Initializes output for the \p merge selection modifier. */
82 init_output_merge(t_topology *top, gmx_ana_selvalue_t *out, void *data);
83 /** Initializes output for the \p plus selection modifier. */
85 init_output_plus(t_topology *top, gmx_ana_selvalue_t *out, void *data);
86 /** Frees the memory allocated for the merging selection modifiers. */
88 free_data_merge(void *data);
90 * Evaluates the \p merge selection modifier.
92 * \param[in] top Not used.
93 * \param[in] fr Not used.
94 * \param[in] pbc Not used.
95 * \param[in] p Positions to merge (should point to \p data->p1).
96 * \param[out] out Output data structure (\p out->u.p is used).
97 * \param[in] data Should point to a \p t_methoddata_merge.
100 evaluate_merge(t_topology *top, t_trxframe *fr, t_pbc *pbc,
101 gmx_ana_pos_t * p, gmx_ana_selvalue_t *out, void *data);
103 * Evaluates the \p plus selection modifier.
105 * \param[in] top Not used.
106 * \param[in] fr Not used.
107 * \param[in] pbc Not used.
108 * \param[in] p Positions to merge (should point to \p data->p1).
109 * \param[out] out Output data structure (\p out->u.p is used).
110 * \param[in] data Should point to a \p t_methoddata_merge.
113 evaluate_plus(t_topology *top, t_trxframe *fr, t_pbc *pbc,
114 gmx_ana_pos_t * p, gmx_ana_selvalue_t *out, void *data);
116 /** Parameters for the merging selection modifiers. */
117 static gmx_ana_selparam_t smparams_merge[] = {
118 {NULL, {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
119 {NULL, {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
120 {"stride", {INT_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
123 /** Help text for the merging selection modifiers. */
124 static const char *help_merge[] = {
125 "MERGING SELECTIONS[PAR]",
129 " POSEXPR merge POSEXPR [stride INT]",
130 " POSEXPR merge POSEXPR [merge POSEXPR ...]",
131 " POSEXPR plus POSEXPR [plus POSEXPR ...]",
133 "Basic selection keywords can only create selections where each atom",
134 "occurs at most once. The [TT]merge[tt] and [TT]plus[tt] selection",
135 "keywords can be used to work around this limitation. Both create",
136 "a selection that contains the positions from all the given position",
137 "expressions, even if they contain duplicates.",
138 "The difference between the two is that [TT]merge[tt] expects two or more",
139 "selections with the same number of positions, and the output contains",
140 "the input positions selected from each expression in turn, i.e.,",
141 "the output is like A1 B1 A2 B2 and so on. It is also possible to merge",
142 "selections of unequal size as long as the size of the first is a",
143 "multiple of the second one. The [TT]stride[tt] parameter can be used",
144 "to explicitly provide this multiplicity.",
145 "[TT]plus[tt] simply concatenates the positions after each other, and",
146 "can work also with selections of different sizes.",
147 "These keywords are valid only at the selection level, not in any",
151 /** Selection method data for the \p plus modifier. */
152 gmx_ana_selmethod_t sm_merge = {
153 "merge", POS_VALUE, SMETH_MODIFIER,
154 asize(smparams_merge), smparams_merge,
163 {"merge POSEXPR", asize(help_merge), help_merge},
166 /** Selection method data for the \p plus modifier. */
167 gmx_ana_selmethod_t sm_plus = {
168 "plus", POS_VALUE, SMETH_MODIFIER,
169 asize(smparams_merge)-1, smparams_merge,
178 {"plus POSEXPR", asize(help_merge), help_merge},
182 * \param[in] npar Should be 2 for \c plus and 3 for \c merge.
183 * \param[in,out] param Method parameters (should point to a copy of
184 * \ref smparams_merge).
185 * \returns Pointer to the allocated data (\p t_methoddata_merge).
187 * Allocates memory for a \p t_methoddata_merge structure.
190 init_data_merge(int npar, gmx_ana_selparam_t *param)
192 t_methoddata_merge *data = new t_methoddata_merge();
194 param[0].val.u.p = &data->p1;
195 param[1].val.u.p = &data->p2;
198 param[2].val.u.i = &data->stride;
204 init_merge(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t * /* param */, void *data)
206 t_methoddata_merge *d = (t_methoddata_merge *)data;
210 GMX_THROW(gmx::InvalidInputError("Stride for merging should be positive"));
212 /* If no stride given, deduce it from the input sizes */
215 d->stride = d->p1.count() / d->p2.count();
217 if (d->p1.count() != d->stride*d->p2.count())
219 GMX_THROW(gmx::InconsistentInputError("The number of positions to be merged are not compatible"));
224 * Does common initialization to all merging modifiers.
226 * \param[in] top Topology data structure.
227 * \param[in,out] out Pointer to output data structure.
228 * \param[in,out] data Should point to \c t_methoddata_merge.
231 init_output_common(t_topology *top, gmx_ana_selvalue_t *out, void *data)
233 t_methoddata_merge *d = (t_methoddata_merge *)data;
235 GMX_UNUSED_VALUE(top);
236 if (d->p1.m.type != d->p2.m.type)
238 /* TODO: Maybe we could pick something else here? */
239 out->u.p->m.type = INDEX_UNKNOWN;
243 out->u.p->m.type = d->p1.m.type;
245 gmx_ana_pos_reserve_for_append(out->u.p, d->p1.count() + d->p2.count(),
246 d->p1.m.b.nra + d->p2.m.b.nra,
247 d->p1.v != NULL, d->p1.f != NULL);
248 gmx_ana_pos_empty_init(out->u.p);
252 * \param[in] top Topology data structure.
253 * \param[in,out] out Pointer to output data structure.
254 * \param[in,out] data Should point to \c t_methoddata_merge.
257 init_output_merge(t_topology *top, gmx_ana_selvalue_t *out, void *data)
259 t_methoddata_merge *d = (t_methoddata_merge *)data;
262 init_output_common(top, out, data);
263 for (i = 0; i < d->p2.count(); ++i)
265 for (j = 0; j < d->stride; ++j)
267 gmx_ana_pos_append_init(out->u.p, &d->p1, d->stride * i + j);
269 gmx_ana_pos_append_init(out->u.p, &d->p2, i);
274 * \param[in] top Topology data structure.
275 * \param[in,out] out Pointer to output data structure.
276 * \param[in,out] data Should point to \c t_methoddata_merge.
279 init_output_plus(t_topology *top, gmx_ana_selvalue_t *out, void *data)
281 t_methoddata_merge *d = (t_methoddata_merge *)data;
284 init_output_common(top, out, data);
285 for (i = 0; i < d->p1.count(); ++i)
287 gmx_ana_pos_append_init(out->u.p, &d->p1, i);
289 for (i = 0; i < d->p2.count(); ++i)
291 gmx_ana_pos_append_init(out->u.p, &d->p2, i);
296 * \param data Data to free (should point to a \p t_methoddata_merge).
298 * Frees the memory allocated for \c t_methoddata_merge.
301 free_data_merge(void *data)
303 t_methoddata_merge *d = (t_methoddata_merge *)data;
308 evaluate_merge(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
309 gmx_ana_pos_t * /* p */, gmx_ana_selvalue_t *out, void *data)
311 t_methoddata_merge *d = (t_methoddata_merge *)data;
315 if (d->p1.count() != d->stride*d->p2.count())
317 GMX_THROW(gmx::InconsistentInputError("The number of positions to be merged are not compatible"));
319 gmx_ana_pos_empty(out->u.p);
320 for (i = 0; i < d->p2.count(); ++i)
322 for (j = 0; j < d->stride; ++j)
324 refid = d->p1.m.refid[d->stride*i+j];
327 refid = (d->stride+1) * (refid / d->stride) + (refid % d->stride);
329 gmx_ana_pos_append(out->u.p, &d->p1, d->stride*i+j, refid);
331 refid = (d->stride+1)*d->p2.m.refid[i] + d->stride;
332 gmx_ana_pos_append(out->u.p, &d->p2, i, refid);
334 gmx_ana_pos_append_finish(out->u.p);
338 evaluate_plus(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
339 gmx_ana_pos_t * /* p */, gmx_ana_selvalue_t *out, void *data)
341 t_methoddata_merge *d = (t_methoddata_merge *)data;
345 gmx_ana_pos_empty(out->u.p);
346 for (i = 0; i < d->p1.count(); ++i)
348 refid = d->p1.m.refid[i];
349 gmx_ana_pos_append(out->u.p, &d->p1, i, refid);
351 for (i = 0; i < d->p2.count(); ++i)
353 refid = d->p2.m.refid[i];
356 refid += d->p1.m.b.nr;
358 gmx_ana_pos_append(out->u.p, &d->p2, i, refid);
360 gmx_ana_pos_append_finish(out->u.p);