Big rework of selection compilation/evaluation.
[gromacs/qmmm-gamess-us.git] / src / gmxlib / selection / sm_permute.c
bloba318912eb2cf77685ae4de4d177544127eb1a328
1 /*
3 * This source code is part of
5 * G R O M A C S
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
29 * For more info, check our website at http://www.gromacs.org
31 /*! \internal \file
32 * \brief Implementation of the \p permute selection modifier.
34 #ifdef HAVE_CONFIG_H
35 #include <config.h>
36 #endif
38 #include <macros.h>
39 #include <smalloc.h>
40 #include <vec.h>
42 #include <position.h>
43 #include <selmethod.h>
45 /*! \internal \brief
46 * Data structure for the \p permute selection modifier.
48 typedef struct
50 /** Positions to permute. */
51 gmx_ana_pos_t p;
52 /** Group to receive the output permutation. */
53 gmx_ana_index_t g;
54 /** Number of elements in the permutation. */
55 int n;
56 /** Array describing the permutation. */
57 int *perm;
58 /** Array that has the permutation reversed. */
59 int *rperm;
60 } t_methoddata_permute;
62 /** Allocates data for the \p permute selection modifier. */
63 static void *
64 init_data_permute(int npar, gmx_ana_selparam_t *param);
65 /** Initializes data for the \p permute selection modifier. */
66 static int
67 init_permute(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
68 /** Initializes output for the \p permute selection modifier. */
69 static int
70 init_output_permute(t_topology *top, gmx_ana_selvalue_t *out, void *data);
71 /** Frees the memory allocated for the \p permute selection modifier. */
72 static void
73 free_data_permute(void *data);
74 /** Evaluates the \p permute selection modifier. */
75 static int
76 evaluate_permute(t_topology *top, t_trxframe *fr, t_pbc *pbc,
77 gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data);
79 /** Parameters for the \p permute selection modifier. */
80 static gmx_ana_selparam_t smparams_permute[] = {
81 {NULL, {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
82 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_VARNUM},
85 /** Help text for the \p permute selection modifier. */
86 static const char *help_permute[] = {
87 "PERMUTING SELECTIONS[PAR]",
89 "[TT]permute P1 ... PN[tt][PAR]",
91 "By default, all selections are evaluated such that the atom indices are",
92 "returned in ascending order. This can be changed by appending",
93 "[TT]permute P1 P2 ... PN[tt] to an expression.",
94 "The [TT]Pi[tt] should form a permutation of the numbers 1 to N.",
95 "This keyword permutes each N-position block in the selection such that",
96 "the i'th position in the block becomes Pi'th.",
97 "Note that it is the positions that are permuted, not individual atoms.",
98 "A fatal error occurs if the size of the selection is not a multiple of n.",
99 "It is only possible to permute the whole selection expression, not any",
100 "subexpressions, i.e., the [TT]permute[tt] keyword should appear last in",
101 "a selection.",
104 /** \internal Selection method data for the \p permute modifier. */
105 gmx_ana_selmethod_t sm_permute = {
106 "permute", POS_VALUE, SMETH_MODIFIER,
107 asize(smparams_permute), smparams_permute,
108 &init_data_permute,
109 NULL,
110 &init_permute,
111 &init_output_permute,
112 &free_data_permute,
113 NULL,
114 NULL,
115 &evaluate_permute,
116 {"permute P1 ... PN", asize(help_permute), help_permute},
120 * \param[in] npar Not used (should be 2).
121 * \param[in,out] param Method parameters (should point to a copy of
122 * \ref smparams_permute).
123 * \returns Pointer to the allocated data (\p t_methoddata_permute).
125 * Allocates memory for a \p t_methoddata_permute structure.
127 static void *
128 init_data_permute(int npar, gmx_ana_selparam_t *param)
130 t_methoddata_permute *data;
132 snew(data, 1);
133 param[0].val.u.p = &data->p;
134 return data;
138 * \param[in] top Not used.
139 * \param[in] npar Not used (should be 2).
140 * \param[in] param Method parameters (should point to \ref smparams_permute).
141 * \param[in] data Should point to a \p t_methoddata_permute.
142 * \returns 0 if the input permutation is valid, -1 on error.
144 static int
145 init_permute(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
147 t_methoddata_permute *d = (t_methoddata_permute *)data;
148 int i;
150 gmx_ana_index_reserve(&d->g, d->p.g->isize);
151 d->n = param[1].val.nr;
152 d->perm = param[1].val.u.i;
153 if (d->p.nr % d->n != 0)
155 fprintf(stderr, "error: the number of positions to be permuted is not divisible by %d\n",
156 d->n);
157 return -1;
159 snew(d->rperm, d->n);
160 for (i = 0; i < d->n; ++i)
162 d->rperm[i] = -1;
164 for (i = 0; i < d->n; ++i)
166 d->perm[i]--;
167 if (d->perm[i] < 0 || d->perm[i] >= d->n)
169 fprintf(stderr, "invalid permutation");
170 return -1;
172 if (d->rperm[d->perm[i]] >= 0)
174 fprintf(stderr, "invalid permutation");
175 return -1;
177 d->rperm[d->perm[i]] = i;
179 return 0;
183 * \param[in] top Topology data structure.
184 * \param[in,out] out Pointer to output data structure.
185 * \param[in,out] data Should point to \c t_methoddata_permute.
186 * \returns 0 for success.
188 static int
189 init_output_permute(t_topology *top, gmx_ana_selvalue_t *out, void *data)
191 t_methoddata_permute *d = (t_methoddata_permute *)data;
192 int i, j, b, k;
194 gmx_ana_pos_copy(out->u.p, &d->p, TRUE);
195 gmx_ana_pos_set_evalgrp(out->u.p, &d->g);
196 d->g.isize = 0;
197 gmx_ana_pos_empty_init(out->u.p);
198 for (i = 0; i < d->p.nr; i += d->n)
200 for (j = 0; j < d->n; ++j)
202 b = i + d->rperm[j];
203 gmx_ana_pos_append_init(out->u.p, &d->g, &d->p, b);
206 return 0;
210 * \param data Data to free (should point to a \p t_methoddata_permute).
212 * Frees the memory allocated for \c t_methoddata_permute.
214 static void
215 free_data_permute(void *data)
217 t_methoddata_permute *d = (t_methoddata_permute *)data;
219 gmx_ana_index_deinit(&d->g);
220 sfree(d->rperm);
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.
232 * Returns -1 if the size of \p p is not divisible by the number of
233 * elements in the permutation.
235 static int
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)
239 t_methoddata_permute *d = (t_methoddata_permute *)data;
240 int i, j, b, k;
241 int refid;
243 if (d->p.nr % d->n != 0)
245 fprintf(stderr, "error: the number of positions to be permuted is not divisible by %d\n",
246 d->n);
247 return -1;
249 d->g.isize = 0;
250 gmx_ana_pos_empty(out->u.p);
251 for (i = 0; i < d->p.nr; i += d->n)
253 for (j = 0; j < d->n; ++j)
255 b = i + d->rperm[j];
256 refid = d->p.m.refid[b];
257 if (refid != -1)
259 /* De-permute the reference ID */
260 refid = refid - (refid % d->n) + d->perm[refid % d->n];
262 gmx_ana_pos_append(out->u.p, &d->g, p, b, refid);
265 gmx_ana_pos_append_finish(out->u.p);
266 return 0;