Reduce preprocessor dependency in resourcedivision.cpp
[gromacs.git] / src / gromacs / selection / sm_permute.cpp
blob0b942e54451e6da51f57129047bab2b3711185e3
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016,2017,2018, 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.
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.
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.
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.
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.
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.
35 /*! \internal \file
36 * \brief
37 * Implements the \p permute selection modifier.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
42 #include "gmxpre.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/selection/position.h"
46 #include "gromacs/utility/arraysize.h"
47 #include "gromacs/utility/exceptions.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/utility/stringutil.h"
51 #include "selmethod.h"
52 #include "selmethod-impl.h"
54 /*! \internal \brief
55 * Data structure for the \p permute selection modifier.
57 typedef struct
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;
69 /*! \brief
70 * Allocates data for the \p permute selection modifier.
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).
77 * Allocates memory for a \p t_methoddata_permute structure.
79 static void *
80 init_data_permute(int npar, gmx_ana_selparam_t *param);
81 /*! \brief
82 * Initializes data for the \p permute selection modifier.
84 * \param[in] top Not used.
85 * \param[in] npar Not used (should be 2).
86 * \param[in] param Method parameters (should point to \ref smparams_permute).
87 * \param[in] data Should point to a \p t_methoddata_permute.
88 * \returns 0 if the input permutation is valid, -1 on error.
90 static void
91 init_permute(const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data);
92 /*! \brief
93 * Initializes output for the \p permute selection modifier.
95 * \param[in] top Topology data structure.
96 * \param[in,out] out Pointer to output data structure.
97 * \param[in,out] data Should point to \c t_methoddata_permute.
99 static void
100 init_output_permute(const gmx_mtop_t *top, gmx_ana_selvalue_t *out, void *data);
101 /** Frees the memory allocated for the \p permute selection modifier. */
102 static void
103 free_data_permute(void *data);
104 /*! \brief
105 * Evaluates the \p permute selection modifier.
107 * \param[in] context 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.
112 * Throws if the size of \p p is not divisible by the number of
113 * elements in the permutation.
115 static void
116 evaluate_permute(const gmx::SelMethodEvalContext &context,
117 gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data);
119 /** Parameters for the \p permute selection modifier. */
120 static gmx_ana_selparam_t smparams_permute[] = {
121 {nullptr, {POS_VALUE, -1, {nullptr}}, nullptr, SPAR_DYNAMIC | SPAR_VARNUM},
122 {nullptr, {INT_VALUE, -1, {nullptr}}, nullptr, SPAR_VARNUM},
125 /** Help text for the \p permute selection modifier. */
126 static const char *const help_permute[] = {
127 "::",
129 " permute P1 ... PN",
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.",
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 nullptr,
150 &init_permute,
151 &init_output_permute,
152 &free_data_permute,
153 nullptr,
154 nullptr,
155 &evaluate_permute,
156 {"POSEXPR permute P1 ... PN",
157 "Permuting selections", asize(help_permute), help_permute},
160 static void *
161 init_data_permute(int /* npar */, gmx_ana_selparam_t *param)
163 t_methoddata_permute *data = new t_methoddata_permute();
164 data->n = 0;
165 data->perm = nullptr;
166 data->rperm = nullptr;
167 param[0].val.u.p = &data->p;
168 return data;
171 static void
172 init_permute(const gmx_mtop_t * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
174 t_methoddata_permute *d = static_cast<t_methoddata_permute *>(data);
175 int i;
177 d->n = param[1].val.nr;
178 d->perm = param[1].val.u.i;
179 if (d->p.count() % d->n != 0)
181 GMX_THROW(gmx::InconsistentInputError(
182 gmx::formatString("The number of positions to be permuted is not divisible by %d", d->n)));
184 snew(d->rperm, d->n);
185 for (i = 0; i < d->n; ++i)
187 d->rperm[i] = -1;
189 for (i = 0; i < d->n; ++i)
191 d->perm[i]--;
192 if (d->perm[i] < 0 || d->perm[i] >= d->n)
194 GMX_THROW(gmx::InvalidInputError("Invalid permutation"));
196 if (d->rperm[d->perm[i]] >= 0)
198 GMX_THROW(gmx::InvalidInputError("Invalid permutation"));
200 d->rperm[d->perm[i]] = i;
204 static void
205 init_output_permute(const gmx_mtop_t * /* top */, gmx_ana_selvalue_t *out, void *data)
207 t_methoddata_permute *d = static_cast<t_methoddata_permute *>(data);
208 int i, j, b;
210 out->u.p->m.type = d->p.m.type;
211 gmx_ana_pos_reserve_for_append(out->u.p, d->p.count(), d->p.m.b.nra,
212 d->p.v != nullptr, d->p.f != nullptr);
213 gmx_ana_pos_empty_init(out->u.p);
214 for (i = 0; i < d->p.count(); i += d->n)
216 for (j = 0; j < d->n; ++j)
218 b = i + d->rperm[j];
219 gmx_ana_pos_append_init(out->u.p, &d->p, b);
225 * \param data Data to free (should point to a \p t_methoddata_permute).
227 * Frees the memory allocated for \c t_methoddata_permute.
229 static void
230 free_data_permute(void *data)
232 t_methoddata_permute *d = static_cast<t_methoddata_permute *>(data);
234 sfree(d->rperm);
235 delete d;
238 static void
239 evaluate_permute(const gmx::SelMethodEvalContext & /*context*/,
240 gmx_ana_pos_t * /*p*/, gmx_ana_selvalue_t *out, void *data)
242 t_methoddata_permute *d = static_cast<t_methoddata_permute *>(data);
243 int i, j, b;
244 int refid;
246 if (d->p.count() % d->n != 0)
248 GMX_THROW(gmx::InconsistentInputError(
249 gmx::formatString("The number of positions to be permuted is not divisible by %d", d->n)));
251 gmx_ana_pos_empty(out->u.p);
252 for (i = 0; i < d->p.count(); i += d->n)
254 for (j = 0; j < d->n; ++j)
256 b = i + d->rperm[j];
257 refid = d->p.m.refid[b];
258 if (refid != -1)
260 /* De-permute the reference ID */
261 refid = refid - (refid % d->n) + d->perm[refid % d->n];
263 gmx_ana_pos_append(out->u.p, &d->p, b, refid);
266 gmx_ana_pos_append_finish(out->u.p);