Merge branch release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / selection / sm_merge.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012,2013, by the GROMACS development team, led by
5  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
6  * others, as listed in the AUTHORS file in the top-level source
7  * directory and at http://www.gromacs.org.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief
37  * Implements the \p merge and \p plus selection modifiers.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_selection
41  */
42 #include "gromacs/legacyheaders/macros.h"
43 #include "gromacs/legacyheaders/smalloc.h"
44 #include "gromacs/legacyheaders/vec.h"
45
46 #include "gromacs/selection/position.h"
47 #include "gromacs/selection/selmethod.h"
48 #include "gromacs/utility/exceptions.h"
49
50 /*! \internal \brief
51  * Data structure for the merging selection modifiers.
52  */
53 typedef struct
54 {
55     /** Input positions. */
56     gmx_ana_pos_t    p1;
57     /** Other input positions. */
58     gmx_ana_pos_t    p2;
59     /** Stride for merging (\c stride values from \c p1 for each in \c p2). */
60     int              stride;
61 } t_methoddata_merge;
62
63 /** Allocates data for the merging selection modifiers. */
64 static void *
65 init_data_merge(int npar, gmx_ana_selparam_t *param);
66 /** Initializes data for the merging selection modifiers. */
67 static void
68 init_merge(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
69 /** Initializes output for the \p merge selection modifier. */
70 static void
71 init_output_merge(t_topology *top, gmx_ana_selvalue_t *out, void *data);
72 /** Initializes output for the \p plus selection modifier. */
73 static void
74 init_output_plus(t_topology *top, gmx_ana_selvalue_t *out, void *data);
75 /** Frees the memory allocated for the merging selection modifiers. */
76 static void
77 free_data_merge(void *data);
78 /** Evaluates the \p merge selection modifier. */
79 static void
80 evaluate_merge(t_topology *top, t_trxframe *fr, t_pbc *pbc,
81                gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data);
82 /** Evaluates the \p plus selection modifier. */
83 static void
84 evaluate_plus(t_topology *top, t_trxframe *fr, t_pbc *pbc,
85               gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data);
86
87 /** Parameters for the merging selection modifiers. */
88 static gmx_ana_selparam_t smparams_merge[] = {
89     {NULL,       {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
90     {NULL,       {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
91     {"stride",   {INT_VALUE,  1, {NULL}}, NULL, SPAR_OPTIONAL},
92 };
93
94 /** Help text for the merging selection modifiers. */
95 static const char *help_merge[] = {
96     "MERGING SELECTIONS[PAR]",
97
98     "[TT]POSEXPR merge POSEXPR [stride INT][tt][BR]",
99     "[TT]POSEXPR merge POSEXPR [merge POSEXPR ...][tt][BR]",
100     "[TT]POSEXPR plus POSEXPR [plus POSEXPR ...][tt][PAR]",
101
102     "Basic selection keywords can only create selections where each atom",
103     "occurs at most once. The [TT]merge[tt] and [TT]plus[tt] selection",
104     "keywords can be used to work around this limitation. Both create",
105     "a selection that contains the positions from all the given position",
106     "expressions, even if they contain duplicates.",
107     "The difference between the two is that [TT]merge[tt] expects two or more",
108     "selections with the same number of positions, and the output contains",
109     "the input positions selected from each expression in turn, i.e.,",
110     "the output is like A1 B1 A2 B2 and so on. It is also possible to merge",
111     "selections of unequal size as long as the size of the first is a",
112     "multiple of the second one. The [TT]stride[tt] parameter can be used",
113     "to explicitly provide this multiplicity.",
114     "[TT]plus[tt] simply concatenates the positions after each other, and",
115     "can work also with selections of different sizes.",
116     "These keywords are valid only at the selection level, not in any",
117     "subexpressions.",
118 };
119
120 /** \internal Selection method data for the \p plus modifier. */
121 gmx_ana_selmethod_t sm_merge = {
122     "merge", POS_VALUE, SMETH_MODIFIER,
123     asize(smparams_merge), smparams_merge,
124     &init_data_merge,
125     NULL,
126     &init_merge,
127     &init_output_merge,
128     &free_data_merge,
129     NULL,
130     NULL,
131     &evaluate_merge,
132     {"merge POSEXPR", asize(help_merge), help_merge},
133 };
134
135 /** \internal Selection method data for the \p plus modifier. */
136 gmx_ana_selmethod_t sm_plus = {
137     "plus", POS_VALUE, SMETH_MODIFIER,
138     asize(smparams_merge)-1, smparams_merge,
139     &init_data_merge,
140     NULL,
141     &init_merge,
142     &init_output_plus,
143     &free_data_merge,
144     NULL,
145     NULL,
146     &evaluate_plus,
147     {"plus POSEXPR", asize(help_merge), help_merge},
148 };
149
150 /*!
151  * \param[in]     npar  Should be 2 for \c plus and 3 for \c merge.
152  * \param[in,out] param Method parameters (should point to a copy of
153  *   \ref smparams_merge).
154  * \returns Pointer to the allocated data (\p t_methoddata_merge).
155  *
156  * Allocates memory for a \p t_methoddata_merge structure.
157  */
158 static void *
159 init_data_merge(int npar, gmx_ana_selparam_t *param)
160 {
161     t_methoddata_merge *data = new t_methoddata_merge();
162     data->stride     = 0;
163     param[0].val.u.p = &data->p1;
164     param[1].val.u.p = &data->p2;
165     if (npar > 2)
166     {
167         param[2].val.u.i = &data->stride;
168     }
169     return data;
170 }
171
172 /*!
173  * \param[in] top   Not used.
174  * \param[in] npar  Not used (should be 2 or 3).
175  * \param[in] param Method parameters (should point to \ref smparams_merge).
176  * \param[in] data  Should point to a \p t_methoddata_merge.
177  * \returns   0 if everything is successful, -1 on error.
178  */
179 static void
180 init_merge(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
181 {
182     t_methoddata_merge *d = (t_methoddata_merge *)data;
183
184     if (d->stride < 0)
185     {
186         GMX_THROW(gmx::InvalidInputError("Stride for merging should be positive"));
187     }
188     /* If no stride given, deduce it from the input sizes */
189     if (d->stride == 0)
190     {
191         d->stride = d->p1.count() / d->p2.count();
192     }
193     if (d->p1.count() != d->stride*d->p2.count())
194     {
195         GMX_THROW(gmx::InconsistentInputError("The number of positions to be merged are not compatible"));
196     }
197 }
198
199 /*! \brief
200  * Does common initialization to all merging modifiers.
201  *
202  * \param[in]     top   Topology data structure.
203  * \param[in,out] out   Pointer to output data structure.
204  * \param[in,out] data  Should point to \c t_methoddata_merge.
205  */
206 static void
207 init_output_common(t_topology *top, gmx_ana_selvalue_t *out, void *data)
208 {
209     t_methoddata_merge *d = (t_methoddata_merge *)data;
210
211     if (d->p1.m.type != d->p2.m.type)
212     {
213         /* TODO: Maybe we could pick something else here? */
214         out->u.p->m.type = INDEX_UNKNOWN;
215     }
216     else
217     {
218         out->u.p->m.type = d->p1.m.type;
219     }
220     gmx_ana_pos_reserve_for_append(out->u.p, d->p1.count() + d->p2.count(),
221                                    d->p1.m.b.nra + d->p2.m.b.nra,
222                                    d->p1.v != NULL, d->p1.f != NULL);
223     gmx_ana_pos_empty_init(out->u.p);
224 }
225
226 /*!
227  * \param[in]     top   Topology data structure.
228  * \param[in,out] out   Pointer to output data structure.
229  * \param[in,out] data  Should point to \c t_methoddata_merge.
230  */
231 static void
232 init_output_merge(t_topology *top, gmx_ana_selvalue_t *out, void *data)
233 {
234     t_methoddata_merge *d = (t_methoddata_merge *)data;
235     int                 i, j;
236
237     init_output_common(top, out, data);
238     for (i = 0; i < d->p2.count(); ++i)
239     {
240         for (j = 0; j < d->stride; ++j)
241         {
242             gmx_ana_pos_append_init(out->u.p, &d->p1, d->stride * i + j);
243         }
244         gmx_ana_pos_append_init(out->u.p, &d->p2, i);
245     }
246 }
247
248 /*!
249  * \param[in]     top   Topology data structure.
250  * \param[in,out] out   Pointer to output data structure.
251  * \param[in,out] data  Should point to \c t_methoddata_merge.
252  */
253 static void
254 init_output_plus(t_topology *top, gmx_ana_selvalue_t *out, void *data)
255 {
256     t_methoddata_merge *d = (t_methoddata_merge *)data;
257     int                 i;
258
259     init_output_common(top, out, data);
260     for (i = 0; i < d->p1.count(); ++i)
261     {
262         gmx_ana_pos_append_init(out->u.p, &d->p1, i);
263     }
264     for (i = 0; i < d->p2.count(); ++i)
265     {
266         gmx_ana_pos_append_init(out->u.p, &d->p2, i);
267     }
268 }
269
270 /*!
271  * \param data Data to free (should point to a \p t_methoddata_merge).
272  *
273  * Frees the memory allocated for \c t_methoddata_merge.
274  */
275 static void
276 free_data_merge(void *data)
277 {
278     t_methoddata_merge *d = (t_methoddata_merge *)data;
279     delete d;
280 }
281
282 /*!
283  * \param[in]  top   Not used.
284  * \param[in]  fr    Not used.
285  * \param[in]  pbc   Not used.
286  * \param[in]  p     Positions to merge (should point to \p data->p1).
287  * \param[out] out   Output data structure (\p out->u.p is used).
288  * \param[in]  data  Should point to a \p t_methoddata_merge.
289  */
290 static void
291 evaluate_merge(t_topology *top, t_trxframe *fr, t_pbc *pbc,
292                gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data)
293 {
294     t_methoddata_merge *d = (t_methoddata_merge *)data;
295     int                 i, j;
296     int                 refid;
297
298     if (d->p1.count() != d->stride*d->p2.count())
299     {
300         GMX_THROW(gmx::InconsistentInputError("The number of positions to be merged are not compatible"));
301     }
302     gmx_ana_pos_empty(out->u.p);
303     for (i = 0; i < d->p2.count(); ++i)
304     {
305         for (j = 0; j < d->stride; ++j)
306         {
307             refid = d->p1.m.refid[d->stride*i+j];
308             if (refid != -1)
309             {
310                 refid = (d->stride+1) * (refid / d->stride) + (refid % d->stride);
311             }
312             gmx_ana_pos_append(out->u.p, &d->p1, d->stride*i+j, refid);
313         }
314         refid = (d->stride+1)*d->p2.m.refid[i] + d->stride;
315         gmx_ana_pos_append(out->u.p, &d->p2, i, refid);
316     }
317     gmx_ana_pos_append_finish(out->u.p);
318 }
319
320 /*!
321  * \param[in]  top   Not used.
322  * \param[in]  fr    Not used.
323  * \param[in]  pbc   Not used.
324  * \param[in]  p     Positions to merge (should point to \p data->p1).
325  * \param[out] out   Output data structure (\p out->u.p is used).
326  * \param[in]  data  Should point to a \p t_methoddata_merge.
327  */
328 static void
329 evaluate_plus(t_topology *top, t_trxframe *fr, t_pbc *pbc,
330               gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data)
331 {
332     t_methoddata_merge *d = (t_methoddata_merge *)data;
333     int                 i;
334     int                 refid;
335
336     gmx_ana_pos_empty(out->u.p);
337     for (i = 0; i < d->p1.count(); ++i)
338     {
339         refid = d->p1.m.refid[i];
340         gmx_ana_pos_append(out->u.p, &d->p1, i, refid);
341     }
342     for (i = 0; i < d->p2.count(); ++i)
343     {
344         refid = d->p2.m.refid[i];
345         if (refid != -1)
346         {
347             refid += d->p1.m.b.nr;
348         }
349         gmx_ana_pos_append(out->u.p, &d->p2, i, refid);
350     }
351     gmx_ana_pos_append_finish(out->u.p);
352 }