Merge remote-tracking branch 'origin/release-4-6' into HEAD
[alexxy/gromacs.git] / src / gromacs / selection / sm_distance.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 distance-based selection methods.
34  *
35  * This file implements the \p distance, \p mindistance and \p within
36  * selection methods.
37  *
38  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
39  * \ingroup module_selection
40  */
41 #ifdef HAVE_CONFIG_H
42 #include <config.h>
43 #endif
44
45 #include "macros.h"
46 #include "pbc.h"
47 #include "smalloc.h"
48 #include "vec.h"
49
50 #include "gromacs/selection/nbsearch.h"
51 #include "gromacs/selection/position.h"
52 #include "gromacs/selection/selmethod.h"
53 #include "gromacs/utility/exceptions.h"
54
55 /*! \internal \brief
56  * Data structure for distance-based selection method.
57  *
58  * The same data structure is used by all the distance-based methods.
59  */
60 typedef struct
61 {
62     /** Cutoff distance. */
63     real                cutoff;
64     /** Positions of the reference points. */
65     gmx_ana_pos_t       p;
66     /** Neighborhood search data. */
67     gmx_ana_nbsearch_t *nb;
68 } t_methoddata_distance;
69
70 /** Allocates data for distance-based selection methods. */
71 static void *
72 init_data_common(int npar, gmx_ana_selparam_t *param);
73 /** Initializes a distance-based selection method. */
74 static void
75 init_common(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
76 /** Frees the data allocated for a distance-based selection method. */
77 static void
78 free_data_common(void *data);
79 /** Initializes the evaluation of a distance-based within selection method for a frame. */
80 static void
81 init_frame_common(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
82 /** Evaluates the \p distance selection method. */
83 static void
84 evaluate_distance(t_topology *top, t_trxframe *fr, t_pbc *pbc,
85                   gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
86 /** Evaluates the \p within selection method. */
87 static void
88 evaluate_within(t_topology *top, t_trxframe *fr, t_pbc *pbc,
89                 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
90
91 /** Parameters for the \p distance selection method. */
92 static gmx_ana_selparam_t smparams_distance[] = {
93     {"cutoff", {REAL_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
94     {"from",   {POS_VALUE,  1, {NULL}}, NULL, SPAR_DYNAMIC},
95 };
96
97 /** Parameters for the \p mindistance selection method. */
98 static gmx_ana_selparam_t smparams_mindistance[] = {
99     {"cutoff", {REAL_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
100     {"from",   {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
101 };
102
103 /** Parameters for the \p within selection method. */
104 static gmx_ana_selparam_t smparams_within[] = {
105     {NULL, {REAL_VALUE,  1, {NULL}}, NULL, 0},
106     {"of", {POS_VALUE,  -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
107 };
108
109 /** Help text for the distance selection methods. */
110 static const char *help_distance[] = {
111     "DISTANCE-BASED SELECTION KEYWORDS[PAR]",
112
113     "[TT]distance from POS [cutoff REAL][tt][BR]",
114     "[TT]mindistance from POS_EXPR [cutoff REAL][tt][BR]",
115     "[TT]within REAL of POS_EXPR[tt][PAR]",
116
117     "[TT]distance[tt] and [TT]mindistance[tt] calculate the distance from the",
118     "given position(s), the only difference being in that [TT]distance[tt]",
119     "only accepts a single position, while any number of positions can be",
120     "given for [TT]mindistance[tt], which then calculates the distance to the",
121     "closest position.",
122     "[TT]within[tt] directly selects atoms that are within [TT]REAL[tt] of",
123     "[TT]POS_EXPR[tt].[PAR]",
124
125     "For the first two keywords, it is possible to specify a cutoff to speed",
126     "up the evaluation: all distances above the specified cutoff are",
127     "returned as equal to the cutoff.",
128 };
129
130 /** \internal Selection method data for the \p distance method. */
131 gmx_ana_selmethod_t sm_distance = {
132     "distance", REAL_VALUE, SMETH_DYNAMIC,
133     asize(smparams_distance), smparams_distance,
134     &init_data_common,
135     NULL,
136     &init_common,
137     NULL,
138     &free_data_common,
139     &init_frame_common,
140     NULL,
141     &evaluate_distance,
142     {"distance from POS [cutoff REAL]", asize(help_distance), help_distance},
143 };
144
145 /** \internal Selection method data for the \p distance method. */
146 gmx_ana_selmethod_t sm_mindistance = {
147     "mindistance", REAL_VALUE, SMETH_DYNAMIC,
148     asize(smparams_mindistance), smparams_mindistance,
149     &init_data_common,
150     NULL,
151     &init_common,
152     NULL,
153     &free_data_common,
154     &init_frame_common,
155     NULL,
156     &evaluate_distance,
157     {"mindistance from POS_EXPR [cutoff REAL]", asize(help_distance), help_distance},
158 };
159
160 /** \internal Selection method data for the \p within method. */
161 gmx_ana_selmethod_t sm_within = {
162     "within", GROUP_VALUE, SMETH_DYNAMIC,
163     asize(smparams_within), smparams_within,
164     &init_data_common,
165     NULL,
166     &init_common,
167     NULL,
168     &free_data_common,
169     &init_frame_common,
170     NULL,
171     &evaluate_within,
172     {"within REAL of POS_EXPR", asize(help_distance), help_distance},
173 };
174
175 /*!
176  * \param[in]     npar  Not used (should be 2).
177  * \param[in,out] param Method parameters (should point to one of the distance
178  *   parameter arrays).
179  * \returns       Pointer to the allocated data (\c t_methoddata_distance).
180  *
181  * Allocates memory for a \c t_methoddata_distance structure and
182  * initializes the parameter as follows:
183  *  - the first parameter defines the value for
184  *    \c t_methoddata_distance::cutoff.
185  *  - the second parameter defines the reference positions and the value is
186  *    stored in \c t_methoddata_distance::p.
187  */
188 static void *
189 init_data_common(int npar, gmx_ana_selparam_t *param)
190 {
191     t_methoddata_distance *data;
192
193     snew(data, 1);
194     data->cutoff     = -1;
195     param[0].val.u.r = &data->cutoff;
196     param[1].val.u.p = &data->p;
197     return data;
198 }
199
200 /*!
201  * \param   top   Not used.
202  * \param   npar  Not used (should be 2).
203  * \param   param Method parameters (should point to one of the distance
204  *   parameter arrays).
205  * \param   data  Pointer to \c t_methoddata_distance to initialize.
206  * \returns 0 on success, a non-zero error code on failure.
207  *
208  * Initializes the neighborhood search data structure
209  * (\c t_methoddata_distance::nb).
210  * Also checks that the cutoff is valid.
211  */
212 static void
213 init_common(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
214 {
215     t_methoddata_distance *d = (t_methoddata_distance *)data;
216
217     if ((param[0].flags & SPAR_SET) && d->cutoff <= 0)
218     {
219         GMX_THROW(gmx::InvalidInputError("Distance cutoff should be > 0"));
220     }
221     d->nb = gmx_ana_nbsearch_create(d->cutoff, d->p.nr);
222 }
223
224 /*!
225  * \param data Data to free (should point to a \c t_methoddata_distance).
226  *
227  * Frees the memory allocated for \c t_methoddata_distance::xref and
228  * \c t_methoddata_distance::nb.
229  */
230 static void
231 free_data_common(void *data)
232 {
233     t_methoddata_distance *d = (t_methoddata_distance *)data;
234
235     if (d->nb)
236     {
237         gmx_ana_nbsearch_free(d->nb);
238     }
239 }
240
241 /*!
242  * \param[in]  top  Not used.
243  * \param[in]  fr   Current frame.
244  * \param[in]  pbc  PBC structure.
245  * \param      data Should point to a \c t_methoddata_distance.
246  * \returns    0 on success, a non-zero error code on error.
247  *
248  * Initializes the neighborhood search for the current frame.
249  */
250 static void
251 init_frame_common(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
252 {
253     t_methoddata_distance *d = (t_methoddata_distance *)data;
254
255     gmx_ana_nbsearch_pos_init(d->nb, pbc, &d->p);
256 }
257
258 /*!
259  * See sel_updatefunc_pos() for description of the parameters.
260  * \p data should point to a \c t_methoddata_distance.
261  *
262  * Calculates the distance of each position from \c t_methoddata_distance::p
263  * and puts them in \p out->u.r.
264  */
265 static void
266 evaluate_distance(t_topology *top, t_trxframe *fr, t_pbc *pbc,
267                   gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
268 {
269     t_methoddata_distance *d = (t_methoddata_distance *)data;
270     int  b, i;
271     real n;
272
273     out->nr = pos->g->isize;
274     for (b = 0; b < pos->nr; ++b)
275     {
276         n = gmx_ana_nbsearch_pos_mindist(d->nb, pos, b);
277         for (i = pos->m.mapb.index[b]; i < pos->m.mapb.index[b+1]; ++i)
278         {
279             out->u.r[i] = n;
280         }
281     }
282 }
283
284 /*!
285  * See sel_updatefunc() for description of the parameters.
286  * \p data should point to a \c t_methoddata_distance.
287  *
288  * Finds the atoms that are closer than the defined cutoff to
289  * \c t_methoddata_distance::xref and puts them in \p out.g.
290  */
291 static void
292 evaluate_within(t_topology *top, t_trxframe *fr, t_pbc *pbc,
293                 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
294 {
295     t_methoddata_distance *d = (t_methoddata_distance *)data;
296     int                    b;
297
298     out->u.g->isize = 0;
299     for (b = 0; b < pos->nr; ++b)
300     {
301         if (gmx_ana_nbsearch_pos_is_within(d->nb, pos, b))
302         {
303             gmx_ana_pos_append(NULL, out->u.g, pos, b, 0);
304         }
305     }
306 }