553c10fc543eaf2f59444aa01f8315003255cacb
[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,2014,2015,2016,2017,2018,2019, 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.
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 "gmxpre.h"
43
44 #include "gromacs/math/vec.h"
45 #include "gromacs/utility/arraysize.h"
46 #include "gromacs/utility/basedefinitions.h"
47 #include "gromacs/utility/exceptions.h"
48 #include "gromacs/utility/smalloc.h"
49
50 #include "position.h"
51 #include "selmethod.h"
52 #include "selmethod_impl.h"
53
54 /*! \internal \brief
55  * Data structure for the merging selection modifiers.
56  */
57 typedef struct
58 {
59     /** Input positions. */
60     gmx_ana_pos_t p1;
61     /** Other input positions. */
62     gmx_ana_pos_t p2;
63     /** Stride for merging (\c stride values from \c p1 for each in \c p2). */
64     int stride;
65 } t_methoddata_merge;
66
67 /** Allocates data for the merging selection modifiers. */
68 static void* init_data_merge(int npar, gmx_ana_selparam_t* param);
69 /*! \brief
70  * Initializes data for the merging selection modifiers.
71  *
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.
77  */
78 static void init_merge(const gmx_mtop_t* top, int npar, gmx_ana_selparam_t* param, void* data);
79 /** Initializes output for the \p merge selection modifier. */
80 static void init_output_merge(const gmx_mtop_t* top, gmx_ana_selvalue_t* out, void* data);
81 /** Initializes output for the \p plus selection modifier. */
82 static void init_output_plus(const gmx_mtop_t* top, gmx_ana_selvalue_t* out, void* data);
83 /** Frees the memory allocated for the merging selection modifiers. */
84 static void free_data_merge(void* data);
85 /*! \brief
86  * Evaluates the \p merge selection modifier.
87  *
88  * \param[in]  context Not used.
89  * \param[in]  p     Positions to merge (should point to \p data->p1).
90  * \param[out] out   Output data structure (\p out->u.p is used).
91  * \param[in]  data  Should point to a \p t_methoddata_merge.
92  */
93 static void evaluate_merge(const gmx::SelMethodEvalContext& context,
94                            gmx_ana_pos_t*                   p,
95                            gmx_ana_selvalue_t*              out,
96                            void*                            data);
97 /*! \brief
98  * Evaluates the \p plus selection modifier.
99  *
100  * \param[in]  context Not used.
101  * \param[in]  p     Positions to merge (should point to \p data->p1).
102  * \param[out] out   Output data structure (\p out->u.p is used).
103  * \param[in]  data  Should point to a \p t_methoddata_merge.
104  */
105 static void evaluate_plus(const gmx::SelMethodEvalContext& context,
106                           gmx_ana_pos_t*                   p,
107                           gmx_ana_selvalue_t*              out,
108                           void*                            data);
109
110 /** Parameters for the merging selection modifiers. */
111 static gmx_ana_selparam_t smparams_merge[] = {
112     { nullptr, { POS_VALUE, -1, { nullptr } }, nullptr, SPAR_DYNAMIC | SPAR_VARNUM },
113     { nullptr, { POS_VALUE, -1, { nullptr } }, nullptr, SPAR_DYNAMIC | SPAR_VARNUM },
114     { "stride", { INT_VALUE, 1, { nullptr } }, nullptr, SPAR_OPTIONAL },
115 };
116
117 //! Help title for the merging selection modifiers.
118 static const char helptitle_merge[] = "Merging selections";
119 //! Help text for the merging selection modifiers.
120 static const char* const help_merge[] = {
121     "::",
122     "",
123     "  POSEXPR merge POSEXPR [stride INT]",
124     "  POSEXPR merge POSEXPR [merge POSEXPR ...]",
125     "  POSEXPR plus POSEXPR [plus POSEXPR ...]",
126     "",
127     "Basic selection keywords can only create selections where each atom",
128     "occurs at most once. The [TT]merge[tt] and [TT]plus[tt] selection",
129     "keywords can be used to work around this limitation. Both create",
130     "a selection that contains the positions from all the given position",
131     "expressions, even if they contain duplicates.",
132     "The difference between the two is that [TT]merge[tt] expects two or more",
133     "selections with the same number of positions, and the output contains",
134     "the input positions selected from each expression in turn, i.e.,",
135     "the output is like A1 B1 A2 B2 and so on. It is also possible to merge",
136     "selections of unequal size as long as the size of the first is a",
137     "multiple of the second one. The [TT]stride[tt] parameter can be used",
138     "to explicitly provide this multiplicity.",
139     "[TT]plus[tt] simply concatenates the positions after each other, and",
140     "can work also with selections of different sizes.",
141     "These keywords are valid only at the selection level, not in any",
142     "subexpressions.",
143 };
144
145 /** Selection method data for the \p plus modifier. */
146 gmx_ana_selmethod_t sm_merge = {
147     "merge",
148     POS_VALUE,
149     SMETH_MODIFIER,
150     asize(smparams_merge),
151     smparams_merge,
152     &init_data_merge,
153     nullptr,
154     &init_merge,
155     &init_output_merge,
156     &free_data_merge,
157     nullptr,
158     nullptr,
159     &evaluate_merge,
160     { "merge POSEXPR", helptitle_merge, asize(help_merge), help_merge },
161 };
162
163 /** Selection method data for the \p plus modifier. */
164 gmx_ana_selmethod_t sm_plus = {
165     "plus",
166     POS_VALUE,
167     SMETH_MODIFIER,
168     asize(smparams_merge) - 1,
169     smparams_merge,
170     &init_data_merge,
171     nullptr,
172     &init_merge,
173     &init_output_plus,
174     &free_data_merge,
175     nullptr,
176     nullptr,
177     &evaluate_plus,
178     { "plus POSEXPR", helptitle_merge, asize(help_merge), help_merge },
179 };
180
181 /*!
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).
186  *
187  * Allocates memory for a \p t_methoddata_merge structure.
188  */
189 static void* init_data_merge(int npar, gmx_ana_selparam_t* param)
190 {
191     t_methoddata_merge* data = new t_methoddata_merge();
192     data->stride             = 0;
193     param[0].val.u.p         = &data->p1;
194     param[1].val.u.p         = &data->p2;
195     if (npar > 2)
196     {
197         param[2].val.u.i = &data->stride;
198     }
199     return data;
200 }
201
202 static void init_merge(const gmx_mtop_t* /* top */, int /* npar */, gmx_ana_selparam_t* /* param */, void* data)
203 {
204     t_methoddata_merge* d = static_cast<t_methoddata_merge*>(data);
205
206     if (d->stride < 0)
207     {
208         GMX_THROW(gmx::InvalidInputError("Stride for merging should be positive"));
209     }
210     /* If no stride given, deduce it from the input sizes */
211     if (d->stride == 0)
212     {
213         d->stride = d->p1.count() / d->p2.count();
214     }
215     if (d->p1.count() != d->stride * d->p2.count())
216     {
217         GMX_THROW(gmx::InconsistentInputError(
218                 "The number of positions to be merged are not compatible"));
219     }
220 }
221
222 /*! \brief
223  * Does common initialization to all merging modifiers.
224  *
225  * \param[in]     top   Topology data structure.
226  * \param[in,out] out   Pointer to output data structure.
227  * \param[in,out] data  Should point to \c t_methoddata_merge.
228  */
229 static void init_output_common(const gmx_mtop_t* top, gmx_ana_selvalue_t* out, void* data)
230 {
231     t_methoddata_merge* d = static_cast<t_methoddata_merge*>(data);
232
233     GMX_UNUSED_VALUE(top);
234     if (d->p1.m.type != d->p2.m.type)
235     {
236         /* TODO: Maybe we could pick something else here? */
237         out->u.p->m.type = INDEX_UNKNOWN;
238     }
239     else
240     {
241         out->u.p->m.type = d->p1.m.type;
242     }
243     gmx_ana_pos_reserve_for_append(out->u.p, d->p1.count() + d->p2.count(),
244                                    d->p1.m.b.nra + d->p2.m.b.nra, d->p1.v != nullptr, d->p1.f != nullptr);
245     gmx_ana_pos_empty_init(out->u.p);
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 init_output_merge(const gmx_mtop_t* top, gmx_ana_selvalue_t* out, void* data)
254 {
255     t_methoddata_merge* d = static_cast<t_methoddata_merge*>(data);
256     int                 i, j;
257
258     init_output_common(top, out, data);
259     for (i = 0; i < d->p2.count(); ++i)
260     {
261         for (j = 0; j < d->stride; ++j)
262         {
263             gmx_ana_pos_append_init(out->u.p, &d->p1, d->stride * i + j);
264         }
265         gmx_ana_pos_append_init(out->u.p, &d->p2, i);
266     }
267 }
268
269 /*!
270  * \param[in]     top   Topology data structure.
271  * \param[in,out] out   Pointer to output data structure.
272  * \param[in,out] data  Should point to \c t_methoddata_merge.
273  */
274 static void init_output_plus(const gmx_mtop_t* top, gmx_ana_selvalue_t* out, void* data)
275 {
276     t_methoddata_merge* d = static_cast<t_methoddata_merge*>(data);
277     int                 i;
278
279     init_output_common(top, out, data);
280     for (i = 0; i < d->p1.count(); ++i)
281     {
282         gmx_ana_pos_append_init(out->u.p, &d->p1, i);
283     }
284     for (i = 0; i < d->p2.count(); ++i)
285     {
286         gmx_ana_pos_append_init(out->u.p, &d->p2, i);
287     }
288 }
289
290 /*!
291  * \param data Data to free (should point to a \p t_methoddata_merge).
292  *
293  * Frees the memory allocated for \c t_methoddata_merge.
294  */
295 static void free_data_merge(void* data)
296 {
297     t_methoddata_merge* d = static_cast<t_methoddata_merge*>(data);
298     delete d;
299 }
300
301 static void evaluate_merge(const gmx::SelMethodEvalContext& /*context*/,
302                            gmx_ana_pos_t* /* p */,
303                            gmx_ana_selvalue_t* out,
304                            void*               data)
305 {
306     t_methoddata_merge* d = static_cast<t_methoddata_merge*>(data);
307     int                 i, j;
308     int                 refid;
309
310     if (d->p1.count() != d->stride * d->p2.count())
311     {
312         GMX_THROW(gmx::InconsistentInputError(
313                 "The number of positions to be merged are not compatible"));
314     }
315     gmx_ana_pos_empty(out->u.p);
316     for (i = 0; i < d->p2.count(); ++i)
317     {
318         for (j = 0; j < d->stride; ++j)
319         {
320             refid = d->p1.m.refid[d->stride * i + j];
321             if (refid != -1)
322             {
323                 refid = (d->stride + 1) * (refid / d->stride) + (refid % d->stride);
324             }
325             gmx_ana_pos_append(out->u.p, &d->p1, d->stride * i + j, refid);
326         }
327         refid = (d->stride + 1) * d->p2.m.refid[i] + d->stride;
328         gmx_ana_pos_append(out->u.p, &d->p2, i, refid);
329     }
330     gmx_ana_pos_append_finish(out->u.p);
331 }
332
333 static void evaluate_plus(const gmx::SelMethodEvalContext& /*context*/,
334                           gmx_ana_pos_t* /* p */,
335                           gmx_ana_selvalue_t* out,
336                           void*               data)
337 {
338     t_methoddata_merge* d = static_cast<t_methoddata_merge*>(data);
339     int                 i;
340     int                 refid;
341
342     gmx_ana_pos_empty(out->u.p);
343     for (i = 0; i < d->p1.count(); ++i)
344     {
345         refid = d->p1.m.refid[i];
346         gmx_ana_pos_append(out->u.p, &d->p1, i, refid);
347     }
348     for (i = 0; i < d->p2.count(); ++i)
349     {
350         refid = d->p2.m.refid[i];
351         if (refid != -1)
352         {
353             refid += d->p1.m.b.nr;
354         }
355         gmx_ana_pos_append(out->u.p, &d->p2, i, refid);
356     }
357     gmx_ana_pos_append_finish(out->u.p);
358 }