31f3c0a8b2dd9c6f84182d5255be263680eecf43
[alexxy/gromacs.git] / src / gromacs / selection / sm_position.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012,2013,2014, 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 position evaluation selection methods.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_selection
41  */
42 #include "gromacs/legacyheaders/macros.h"
43
44 #include "gromacs/selection/indexutil.h"
45 #include "gromacs/selection/poscalc.h"
46 #include "gromacs/selection/position.h"
47 #include "gromacs/selection/selmethod.h"
48 #include "gromacs/utility/cstringutil.h"
49 #include "gromacs/utility/smalloc.h"
50
51 #include "keywords.h"
52 #include "selelem.h"
53
54 /*! \internal \brief
55  * Data structure for position keyword evaluation.
56  */
57 typedef struct
58 {
59     /** Position calculation collection to use. */
60     gmx::PositionCalculationCollection *pcc;
61     /** Index group for which the center should be evaluated. */
62     gmx_ana_index_t                     g;
63     /** Position evaluation data structure. */
64     gmx_ana_poscalc_t                  *pc;
65     /** true if periodic boundary conditions should be used. */
66     bool                                bPBC;
67     /** Type of positions to calculate. */
68     char                               *type;
69     /** Flags for the position calculation. */
70     int                                 flags;
71 } t_methoddata_pos;
72
73 /** Allocates data for position evaluation selection methods. */
74 static void *
75 init_data_pos(int npar, gmx_ana_selparam_t *param);
76 /** Sets the position calculation collection for position evaluation selection methods. */
77 static void
78 set_poscoll_pos(gmx::PositionCalculationCollection *pcc, void *data);
79 /*! \brief
80  * Initializes position evaluation keywords.
81  *
82  * \param[in] top   Not used.
83  * \param[in] npar  Not used.
84  * \param[in] param Not used.
85  * \param[in,out] data  Should point to \c t_methoddata_pos.
86  * \returns       0 on success, a non-zero error code on error.
87  *
88  * The \c t_methoddata_pos::type field should have been initialized
89  * externally using _gmx_selelem_set_kwpos_type().
90  */
91 static void
92 init_kwpos(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
93 /*! \brief
94  * Initializes the \p cog selection method.
95  *
96  * \param[in]     top   Topology data structure.
97  * \param[in]     npar  Not used.
98  * \param[in]     param Not used.
99  * \param[in,out] data  Should point to \c t_methoddata_pos.
100  * \returns       0 on success, a non-zero error code on error.
101  */
102 static void
103 init_cog(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
104 /*! \brief
105  * Initializes the \p cog selection method.
106  *
107  * \param[in]     top   Topology data structure.
108  * \param[in]     npar  Not used.
109  * \param[in]     param Not used.
110  * \param[in,out] data  Should point to \c t_methoddata_pos.
111  * \returns       0 on success, a non-zero error code on error.
112  */
113 static void
114 init_com(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
115 /*! \brief
116  * Initializes output for position evaluation selection methods.
117  *
118  * \param[in]     top   Topology data structure.
119  * \param[in,out] out   Pointer to output data structure.
120  * \param[in,out] data  Should point to \c t_methoddata_pos.
121  * \returns       0 for success.
122  */
123 static void
124 init_output_pos(t_topology *top, gmx_ana_selvalue_t *out, void *data);
125 /** Frees the data allocated for position evaluation selection methods. */
126 static void
127 free_data_pos(void *data);
128 /** Evaluates position evaluation selection methods. */
129 static void
130 evaluate_pos(t_topology * /* top */, t_trxframe *fr, t_pbc *pbc, gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void *data);
131
132 /** Parameters for position keyword evaluation. */
133 static gmx_ana_selparam_t smparams_keyword_pos[] = {
134     {NULL,   {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
135 };
136
137 /** Parameters for the \p cog and \p com selection methods. */
138 static gmx_ana_selparam_t smparams_com[] = {
139     {"of",   {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
140     {"pbc",  {NO_VALUE,    0, {NULL}}, NULL, 0},
141 };
142
143 /** Selection method data for position keyword evaluation. */
144 gmx_ana_selmethod_t sm_keyword_pos = {
145     "kw_pos", POS_VALUE, SMETH_DYNAMIC | SMETH_VARNUMVAL | SMETH_ALLOW_UNSORTED,
146     asize(smparams_keyword_pos), smparams_keyword_pos,
147     &init_data_pos,
148     &set_poscoll_pos,
149     &init_kwpos,
150     &init_output_pos,
151     &free_data_pos,
152     NULL,
153     &evaluate_pos,
154     NULL,
155     {NULL, 0, NULL},
156 };
157
158 /** Selection method data for the \p cog method. */
159 gmx_ana_selmethod_t sm_cog = {
160     "cog", POS_VALUE, SMETH_DYNAMIC | SMETH_SINGLEVAL,
161     asize(smparams_com), smparams_com,
162     &init_data_pos,
163     &set_poscoll_pos,
164     &init_cog,
165     &init_output_pos,
166     &free_data_pos,
167     NULL,
168     &evaluate_pos,
169     NULL,
170     {"cog of ATOM_EXPR [pbc]", 0, NULL},
171 };
172
173 /** Selection method data for the \p com method. */
174 gmx_ana_selmethod_t sm_com = {
175     "com", POS_VALUE, SMETH_REQTOP | SMETH_DYNAMIC | SMETH_SINGLEVAL,
176     asize(smparams_com), smparams_com,
177     &init_data_pos,
178     &set_poscoll_pos,
179     &init_com,
180     &init_output_pos,
181     &free_data_pos,
182     NULL,
183     &evaluate_pos,
184     NULL,
185     {"com of ATOM_EXPR [pbc]", 0, NULL},
186 };
187
188 /*!
189  * \param[in]     npar  Should be 1 or 2.
190  * \param[in,out] param Method parameters (should point to
191  *   \ref smparams_keyword_pos or \ref smparams_com).
192  * \returns       Pointer to the allocated data (\c t_methoddata_pos).
193  *
194  * Allocates memory for a \c t_methoddata_pos structure and initializes
195  * the first parameter to define the value for \c t_methoddata_pos::g.
196  * If a second parameter is present, it is used for setting the
197  * \c t_methoddata_pos::bPBC flag.
198  */
199 static void *
200 init_data_pos(int npar, gmx_ana_selparam_t *param)
201 {
202     t_methoddata_pos *data;
203
204     snew(data, 1);
205     param[0].val.u.g = &data->g;
206     if (npar > 1)
207     {
208         param[1].val.u.b = &data->bPBC;
209     }
210     data->pc       = NULL;
211     data->bPBC     = false;
212     data->type     = NULL;
213     data->flags    = -1;
214     return data;
215 }
216
217 /*!
218  * \param[in]     pcc   Position calculation collection to use.
219  * \param[in,out] data  Should point to \c t_methoddata_pos.
220  */
221 static void
222 set_poscoll_pos(gmx::PositionCalculationCollection *pcc, void *data)
223 {
224     ((t_methoddata_pos *)data)->pcc = pcc;
225 }
226
227 /*!
228  * \param[in,out] sel   Selection element to initialize.
229  * \param[in]     type  One of the enum values acceptable for
230  *     PositionCalculationCollection::typeFromEnum().
231  *
232  * Initializes the reference position type for position evaluation.
233  * If called multiple times, the first setting takes effect, and later calls
234  * are neglected.
235  */
236 void
237 _gmx_selelem_set_kwpos_type(gmx::SelectionTreeElement *sel, const char *type)
238 {
239     t_methoddata_pos *d = (t_methoddata_pos *)sel->u.expr.mdata;
240
241     if (sel->type != SEL_EXPRESSION || !sel->u.expr.method
242         || sel->u.expr.method->name != sm_keyword_pos.name)
243     {
244         return;
245     }
246     if (!d->type && type)
247     {
248         d->type  = gmx_strdup(type);
249         /* FIXME: It would be better not to have the string here hardcoded. */
250         if (type[0] != 'a')
251         {
252             sel->u.expr.method->flags |= SMETH_REQTOP;
253         }
254     }
255 }
256
257 /*!
258  * \param[in,out] sel   Selection element to initialize.
259  * \param[in]     flags Default completion flags
260  *     (see PositionCalculationCollection::typeFromEnum()).
261  *
262  * Initializes the flags for position evaluation.
263  * If called multiple times, the first setting takes effect, and later calls
264  * are neglected.
265  */
266 void
267 _gmx_selelem_set_kwpos_flags(gmx::SelectionTreeElement *sel, int flags)
268 {
269     t_methoddata_pos *d = (t_methoddata_pos *)sel->u.expr.mdata;
270
271     if (sel->type != SEL_EXPRESSION || !sel->u.expr.method
272         || sel->u.expr.method->name != sm_keyword_pos.name)
273     {
274         return;
275     }
276     if (d->flags == -1)
277     {
278         d->flags = flags;
279     }
280 }
281
282 static void
283 init_kwpos(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
284 {
285     t_methoddata_pos *d = (t_methoddata_pos *)data;
286
287     if (!(param[0].flags & SPAR_DYNAMIC))
288     {
289         d->flags &= ~(POS_DYNAMIC | POS_MASKONLY);
290     }
291     else if (!(d->flags & POS_MASKONLY))
292     {
293         d->flags |= POS_DYNAMIC;
294     }
295     d->pc = d->pcc->createCalculationFromEnum(d->type, d->flags);
296     gmx_ana_poscalc_set_maxindex(d->pc, &d->g);
297 }
298
299 static void
300 init_cog(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
301 {
302     t_methoddata_pos *d = (t_methoddata_pos *)data;
303
304     d->flags = (param[0].flags & SPAR_DYNAMIC) ? POS_DYNAMIC : 0;
305     d->pc    = d->pcc->createCalculation(d->bPBC ? POS_ALL_PBC : POS_ALL, d->flags);
306     gmx_ana_poscalc_set_maxindex(d->pc, &d->g);
307 }
308
309 static void
310 init_com(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
311 {
312     t_methoddata_pos *d = (t_methoddata_pos *)data;
313
314     d->flags  = (param[0].flags & SPAR_DYNAMIC) ? POS_DYNAMIC : 0;
315     d->flags |= POS_MASS;
316     d->pc     = d->pcc->createCalculation(d->bPBC ? POS_ALL_PBC : POS_ALL, d->flags);
317     gmx_ana_poscalc_set_maxindex(d->pc, &d->g);
318 }
319
320 static void
321 init_output_pos(t_topology * /* top */, gmx_ana_selvalue_t *out, void *data)
322 {
323     t_methoddata_pos *d = (t_methoddata_pos *)data;
324
325     gmx_ana_poscalc_init_pos(d->pc, out->u.p);
326 }
327
328 /*!
329  * \param data Data to free (should point to a \c t_methoddata_pos).
330  *
331  * Frees the memory allocated for \c t_methoddata_pos::g and
332  * \c t_methoddata_pos::pc.
333  */
334 static void
335 free_data_pos(void *data)
336 {
337     t_methoddata_pos *d = (t_methoddata_pos *)data;
338
339     sfree(d->type);
340     gmx_ana_poscalc_free(d->pc);
341     sfree(d);
342 }
343
344 /*!
345  * See sel_updatefunc() for description of the parameters.
346  * \p data should point to a \c t_methoddata_pos.
347  *
348  * Calculates the positions using \c t_methoddata_pos::pc for the index group
349  * in \c t_methoddata_pos::g and stores the results in \p out->u.p.
350  */
351 static void
352 evaluate_pos(t_topology * /* top */, t_trxframe *fr, t_pbc *pbc,
353              gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void *data)
354 {
355     t_methoddata_pos *d = (t_methoddata_pos *)data;
356
357     gmx_ana_poscalc_update(d->pc, out->u.p, &d->g, fr, pbc);
358 }