Merge branch release-4-6 into release-5-0
[alexxy/gromacs.git] / src / gromacs / selection / sm_permute.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012,2013,2014, 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 permute selection modifier.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_selection
41  */
42 #include "gromacs/legacyheaders/macros.h"
43 #include "gromacs/legacyheaders/vec.h"
44
45 #include "gromacs/selection/position.h"
46 #include "gromacs/selection/selmethod.h"
47 #include "gromacs/utility/exceptions.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/utility/stringutil.h"
50
51 /*! \internal \brief
52  * Data structure for the \p permute selection modifier.
53  */
54 typedef struct
55 {
56     /** Positions to permute. */
57     gmx_ana_pos_t    p;
58     /** Number of elements in the permutation. */
59     int              n;
60     /** Array describing the permutation. */
61     int             *perm;
62     /** Array that has the permutation reversed. */
63     int             *rperm;
64 } t_methoddata_permute;
65
66 /*! \brief
67  * Allocates data for the \p permute selection modifier.
68  *
69  * \param[in]     npar  Not used (should be 2).
70  * \param[in,out] param Method parameters (should point to a copy of
71  *   \ref smparams_permute).
72  * \returns Pointer to the allocated data (\p t_methoddata_permute).
73  *
74  * Allocates memory for a \p t_methoddata_permute structure.
75  */
76 static void *
77 init_data_permute(int npar, gmx_ana_selparam_t *param);
78 /*! \brief
79  * Initializes data for the \p permute selection modifier.
80  *
81  * \param[in] top   Not used.
82  * \param[in] npar  Not used (should be 2).
83  * \param[in] param Method parameters (should point to \ref smparams_permute).
84  * \param[in] data  Should point to a \p t_methoddata_permute.
85  * \returns   0 if the input permutation is valid, -1 on error.
86  */
87 static void
88 init_permute(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
89 /*! \brief
90  * Initializes output for the \p permute selection modifier.
91  *
92  * \param[in]     top   Topology data structure.
93  * \param[in,out] out   Pointer to output data structure.
94  * \param[in,out] data  Should point to \c t_methoddata_permute.
95  */
96 static void
97 init_output_permute(t_topology *top, gmx_ana_selvalue_t *out, void *data);
98 /** Frees the memory allocated for the \p permute selection modifier. */
99 static void
100 free_data_permute(void *data);
101 static void
102 /*! \brief
103  * Evaluates the \p permute selection modifier.
104  *
105  * \param[in]  top   Not used.
106  * \param[in]  fr    Not used.
107  * \param[in]  pbc   Not used.
108  * \param[in]  p     Positions to permute (should point to \p data->p).
109  * \param[out] out   Output data structure (\p out->u.p is used).
110  * \param[in]  data  Should point to a \p t_methoddata_permute.
111  * \returns    0 if \p p could be permuted, -1 on error.
112  *
113  * Returns -1 if the size of \p p is not divisible by the number of
114  * elements in the permutation.
115  */
116 evaluate_permute(t_topology *top, t_trxframe *fr, t_pbc *pbc,
117                  gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data);
118
119 /** Parameters for the \p permute selection modifier. */
120 static gmx_ana_selparam_t smparams_permute[] = {
121     {NULL,       {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
122     {NULL,       {INT_VALUE, -1, {NULL}}, NULL, SPAR_VARNUM},
123 };
124
125 /** Help text for the \p permute selection modifier. */
126 static const char *help_permute[] = {
127     "PERMUTING SELECTIONS[PAR]",
128
129     "[TT]permute P1 ... PN[tt][PAR]",
130
131     "By default, all selections are evaluated such that the atom indices are",
132     "returned in ascending order. This can be changed by appending",
133     "[TT]permute P1 P2 ... PN[tt] to an expression.",
134     "The [TT]Pi[tt] should form a permutation of the numbers 1 to N.",
135     "This keyword permutes each N-position block in the selection such that",
136     "the i'th position in the block becomes Pi'th.",
137     "Note that it is the positions that are permuted, not individual atoms.",
138     "A fatal error occurs if the size of the selection is not a multiple of n.",
139     "It is only possible to permute the whole selection expression, not any",
140     "subexpressions, i.e., the [TT]permute[tt] keyword should appear last in",
141     "a selection.",
142 };
143
144 /** Selection method data for the \p permute modifier. */
145 gmx_ana_selmethod_t sm_permute = {
146     "permute", POS_VALUE, SMETH_MODIFIER,
147     asize(smparams_permute), smparams_permute,
148     &init_data_permute,
149     NULL,
150     &init_permute,
151     &init_output_permute,
152     &free_data_permute,
153     NULL,
154     NULL,
155     &evaluate_permute,
156     {"permute P1 ... PN", asize(help_permute), help_permute},
157 };
158
159 static void *
160 init_data_permute(int /* npar */, gmx_ana_selparam_t *param)
161 {
162     t_methoddata_permute *data = new t_methoddata_permute();
163     data->n          = 0;
164     data->perm       = NULL;
165     data->rperm      = NULL;
166     param[0].val.u.p = &data->p;
167     return data;
168 }
169
170 static void
171 init_permute(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
172 {
173     t_methoddata_permute *d = (t_methoddata_permute *)data;
174     int                   i;
175
176     d->n    = param[1].val.nr;
177     d->perm = param[1].val.u.i;
178     if (d->p.count() % d->n != 0)
179     {
180         GMX_THROW(gmx::InconsistentInputError(
181                           gmx::formatString("The number of positions to be permuted is not divisible by %d", d->n)));
182     }
183     snew(d->rperm, d->n);
184     for (i = 0; i < d->n; ++i)
185     {
186         d->rperm[i] = -1;
187     }
188     for (i = 0; i < d->n; ++i)
189     {
190         d->perm[i]--;
191         if (d->perm[i] < 0 || d->perm[i] >= d->n)
192         {
193             GMX_THROW(gmx::InvalidInputError("Invalid permutation"));
194         }
195         if (d->rperm[d->perm[i]] >= 0)
196         {
197             GMX_THROW(gmx::InvalidInputError("Invalid permutation"));
198         }
199         d->rperm[d->perm[i]] = i;
200     }
201 }
202
203 static void
204 init_output_permute(t_topology * /* top */, gmx_ana_selvalue_t *out, void *data)
205 {
206     t_methoddata_permute *d = (t_methoddata_permute *)data;
207     int                   i, j, b;
208
209     out->u.p->m.type = d->p.m.type;
210     gmx_ana_pos_reserve_for_append(out->u.p, d->p.count(), d->p.m.b.nra,
211                                    d->p.v != NULL, d->p.f != NULL);
212     gmx_ana_pos_empty_init(out->u.p);
213     for (i = 0; i < d->p.count(); i += d->n)
214     {
215         for (j = 0; j < d->n; ++j)
216         {
217             b = i + d->rperm[j];
218             gmx_ana_pos_append_init(out->u.p, &d->p, b);
219         }
220     }
221 }
222
223 /*!
224  * \param data Data to free (should point to a \p t_methoddata_permute).
225  *
226  * Frees the memory allocated for \c t_methoddata_permute.
227  */
228 static void
229 free_data_permute(void *data)
230 {
231     t_methoddata_permute *d = (t_methoddata_permute *)data;
232
233     sfree(d->rperm);
234     delete d;
235 }
236
237 static void
238 evaluate_permute(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
239                  gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data)
240 {
241     t_methoddata_permute *d = (t_methoddata_permute *)data;
242     int                   i, j, b;
243     int                   refid;
244
245     if (d->p.count() % d->n != 0)
246     {
247         GMX_THROW(gmx::InconsistentInputError(
248                           gmx::formatString("The number of positions to be permuted is not divisible by %d", d->n)));
249     }
250     gmx_ana_pos_empty(out->u.p);
251     for (i = 0; i < d->p.count(); i += d->n)
252     {
253         for (j = 0; j < d->n; ++j)
254         {
255             b     = i + d->rperm[j];
256             refid = d->p.m.refid[b];
257             if (refid != -1)
258             {
259                 /* De-permute the reference ID */
260                 refid = refid - (refid % d->n) + d->perm[refid % d->n];
261             }
262             gmx_ana_pos_append(out->u.p, p, b, refid);
263         }
264     }
265     gmx_ana_pos_append_finish(out->u.p);
266 }