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