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