Fixing copyright issues and code contributors
[alexxy/gromacs.git] / src / gmxlib / selection / sm_position.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2009, The GROMACS development team,
6  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012,2013, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /*! \internal \file
39  * \brief Implementation of position evaluation selection methods.
40  */
41 #ifdef HAVE_CONFIG_H
42 #include <config.h>
43 #endif
44
45 #include <macros.h>
46 #include <smalloc.h>
47 #include <string2.h>
48
49 #include <indexutil.h>
50 #include <poscalc.h>
51 #include <position.h>
52 #include <selmethod.h>
53
54 #include "keywords.h"
55 #include "selelem.h"
56
57 /*! \internal \brief
58  * Data structure for position keyword evaluation.
59  */
60 typedef struct
61 {
62     /** Position calculation collection to use. */
63     gmx_ana_poscalc_coll_t *pcc;
64     /** Index group for which the center should be evaluated. */
65     gmx_ana_index_t    g;
66     /** Position evaluation data structure. */
67     gmx_ana_poscalc_t *pc;
68     /** TRUE if periodic boundary conditions should be used. */
69     gmx_bool               bPBC;
70     /** Type of positions to calculate. */
71     char              *type;
72     /** Flags for the position calculation. */
73     int                flags;
74 } t_methoddata_pos;
75
76 /** Allocates data for position evaluation selection methods. */
77 static void *
78 init_data_pos(int npar, gmx_ana_selparam_t *param);
79 /** Sets the position calculation collection for position evaluation selection methods. */
80 static void
81 set_poscoll_pos(gmx_ana_poscalc_coll_t *pcc, void *data);
82 /** Initializes position evaluation keywords. */
83 static int
84 init_kwpos(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
85 /** Initializes the \p cog selection method. */
86 static int
87 init_cog(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
88 /** Initializes the \p cog selection method. */
89 static int
90 init_com(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
91 /** Initializes output for position evaluation selection methods. */
92 static int
93 init_output_pos(t_topology *top, gmx_ana_selvalue_t *out, void *data);
94 /** Frees the data allocated for position evaluation selection methods. */
95 static void
96 free_data_pos(void *data);
97 /** Evaluates position evaluation selection methods. */
98 static int
99 evaluate_pos(t_topology *top, t_trxframe *fr, t_pbc *pbc,
100              gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
101
102 /** Parameters for position keyword evaluation. */
103 static gmx_ana_selparam_t smparams_keyword_pos[] = {
104     {NULL,   {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
105 };
106
107 /** Parameters for the \p cog and \p com selection methods. */
108 static gmx_ana_selparam_t smparams_com[] = {
109     {"of",   {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
110     {"pbc",  {NO_VALUE,    0, {NULL}}, NULL, 0},
111 };
112
113 /** \internal Selection method data for position keyword evaluation. */
114 gmx_ana_selmethod_t sm_keyword_pos = {
115     "kw_pos", POS_VALUE, SMETH_DYNAMIC | SMETH_VARNUMVAL,
116     asize(smparams_keyword_pos), smparams_keyword_pos,
117     &init_data_pos,
118     &set_poscoll_pos,
119     &init_kwpos,
120     &init_output_pos,
121     &free_data_pos,
122      NULL,
123     &evaluate_pos,
124      NULL,
125     {NULL, 0, NULL},
126 };
127
128 /** \internal Selection method data for the \p cog method. */
129 gmx_ana_selmethod_t sm_cog = {
130     "cog", POS_VALUE, SMETH_DYNAMIC | SMETH_SINGLEVAL,
131     asize(smparams_com), smparams_com,
132     &init_data_pos,
133     &set_poscoll_pos,
134     &init_cog,
135     &init_output_pos,
136     &free_data_pos,
137      NULL,
138     &evaluate_pos,
139      NULL,
140     {"cog of ATOM_EXPR [pbc]", 0, NULL},
141 };
142
143 /** \internal Selection method data for the \p com method. */
144 gmx_ana_selmethod_t sm_com = {
145     "com", POS_VALUE, SMETH_REQTOP | SMETH_DYNAMIC | SMETH_SINGLEVAL,
146     asize(smparams_com), smparams_com,
147     &init_data_pos,
148     &set_poscoll_pos,
149     &init_com,
150     &init_output_pos,
151     &free_data_pos,
152      NULL,
153     &evaluate_pos,
154      NULL,
155     {"com of ATOM_EXPR [pbc]", 0, NULL},
156 };
157
158 /*!
159  * \param[in]     npar  Should be 1 or 2.
160  * \param[in,out] param Method parameters (should point to
161  *   \ref smparams_keyword_pos or \ref smparams_com).
162  * \returns       Pointer to the allocated data (\c t_methoddata_pos).
163  *
164  * Allocates memory for a \c t_methoddata_pos structure and initializes
165  * the first parameter to define the value for \c t_methoddata_pos::g.
166  * If a second parameter is present, it is used for setting the
167  * \c t_methoddata_pos::bPBC flag.
168  */
169 static void *
170 init_data_pos(int npar, gmx_ana_selparam_t *param)
171 {
172     t_methoddata_pos *data;
173
174     snew(data, 1);
175     param[0].val.u.g = &data->g;
176     if (npar > 1)
177     {
178         param[1].val.u.b = &data->bPBC;
179     }
180     data->pc       = NULL;
181     data->bPBC     = FALSE;
182     data->type     = NULL;
183     data->flags    = -1;
184     return data;
185 }
186
187 /*!
188  * \param[in]     pcc   Position calculation collection to use.
189  * \param[in,out] data  Should point to \c t_methoddata_pos.
190  */
191 static void
192 set_poscoll_pos(gmx_ana_poscalc_coll_t *pcc, void *data)
193 {
194     ((t_methoddata_pos *)data)->pcc = pcc;
195 }
196
197 /*!
198  * \param[in,out] sel   Selection element to initialize.
199  * \param[in]     type  One of the enum values acceptable for
200  *   gmx_ana_poscalc_type_from_enum().
201  *
202  * Initializes the reference position type for position evaluation.
203  * If called multiple times, the first setting takes effect, and later calls
204  * are neglected.
205  */
206 void
207 _gmx_selelem_set_kwpos_type(t_selelem *sel, const char *type)
208 {
209     t_methoddata_pos *d = (t_methoddata_pos *)sel->u.expr.mdata;
210
211     if (sel->type != SEL_EXPRESSION || !sel->u.expr.method
212         || sel->u.expr.method->name != sm_keyword_pos.name)
213     {
214         return;
215     }
216     if (!d->type && type)
217     {
218         d->type  = strdup(type);
219         /* FIXME: It would be better not to have the string here hardcoded. */
220         if (type[0] != 'a')
221         {
222             sel->u.expr.method->flags |= SMETH_REQTOP;
223         }
224     }
225 }
226
227 /*!
228  * \param[in,out] sel   Selection element to initialize.
229  * \param[in]     flags Default completion flags
230  *   (see gmx_ana_poscalc_type_from_enum()).
231  *
232  * Initializes the flags 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_flags(t_selelem *sel, int flags)
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->flags == -1)
247     {
248         d->flags = flags;
249     }
250 }
251
252 /*!
253  * \param[in] top   Not used.
254  * \param[in] npar  Not used.
255  * \param[in] param Not used.
256  * \param[in,out] data  Should point to \c t_methoddata_pos.
257  * \returns       0 on success, a non-zero error code on error.
258  *
259  * The \c t_methoddata_pos::type field should have been initialized
260  * externally using _gmx_selelem_set_kwpos_type().
261  */
262 static int
263 init_kwpos(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
264 {
265     t_methoddata_pos *d = (t_methoddata_pos *)data;
266     int               rc;
267
268     if (!(param[0].flags & SPAR_DYNAMIC))
269     {
270         d->flags &= ~(POS_DYNAMIC | POS_MASKONLY);
271     }
272     else if (!(d->flags & POS_MASKONLY))
273     {
274         d->flags |= POS_DYNAMIC;
275     }
276     rc = gmx_ana_poscalc_create_enum(&d->pc, d->pcc, d->type, d->flags);
277     if (rc != 0)
278     {
279         return rc;
280     }
281     gmx_ana_poscalc_set_maxindex(d->pc, &d->g);
282     return 0;
283 }
284
285 /*!
286  * \param[in]     top   Topology data structure.
287  * \param[in]     npar  Not used.
288  * \param[in]     param Not used.
289  * \param[in,out] data  Should point to \c t_methoddata_pos.
290  * \returns       0 on success, a non-zero error code on error.
291  */
292 static int
293 init_cog(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
294 {
295     t_methoddata_pos *d = (t_methoddata_pos *)data;
296     int               rc;
297
298     d->flags = (param[0].flags & SPAR_DYNAMIC) ? POS_DYNAMIC : 0;
299     rc = gmx_ana_poscalc_create(&d->pc, d->pcc, d->bPBC ? POS_ALL_PBC : POS_ALL,
300                                 d->flags);
301     if (rc != 0)
302     {
303         return rc;
304     }
305     gmx_ana_poscalc_set_maxindex(d->pc, &d->g);
306     return 0;
307 }
308
309 /*!
310  * \param[in]     top   Topology data structure.
311  * \param[in]     npar  Not used.
312  * \param[in]     param Not used.
313  * \param[in,out] data  Should point to \c t_methoddata_pos.
314  * \returns       0 on success, a non-zero error code on error.
315  */
316 static int
317 init_com(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
318 {
319     t_methoddata_pos *d = (t_methoddata_pos *)data;
320     int               rc;
321
322     d->flags  = (param[0].flags & SPAR_DYNAMIC) ? POS_DYNAMIC : 0;
323     d->flags |= POS_MASS;
324     rc = gmx_ana_poscalc_create(&d->pc, d->pcc, d->bPBC ? POS_ALL_PBC : POS_ALL,
325                                 d->flags);
326     if (rc != 0)
327     {
328         return rc;
329     }
330     gmx_ana_poscalc_set_maxindex(d->pc, &d->g);
331     return 0;
332 }
333
334 /*!
335  * \param[in]     top   Topology data structure.
336  * \param[in,out] out   Pointer to output data structure.
337  * \param[in,out] data  Should point to \c t_methoddata_pos.
338  * \returns       0 for success.
339  */
340 static int
341 init_output_pos(t_topology *top, gmx_ana_selvalue_t *out, void *data)
342 {
343     t_methoddata_pos *d = (t_methoddata_pos *)data;
344
345     gmx_ana_poscalc_init_pos(d->pc, out->u.p);
346     gmx_ana_pos_set_evalgrp(out->u.p, &d->g);
347     return 0;
348 }
349
350 /*!
351  * \param data Data to free (should point to a \c t_methoddata_pos).
352  *
353  * Frees the memory allocated for \c t_methoddata_pos::g and
354  * \c t_methoddata_pos::pc.
355  */
356 static void
357 free_data_pos(void *data)
358 {
359     t_methoddata_pos *d = (t_methoddata_pos *)data;
360
361     sfree(d->type);
362     gmx_ana_poscalc_free(d->pc);
363 }
364
365 /*!
366  * See sel_updatefunc() for description of the parameters.
367  * \p data should point to a \c t_methoddata_pos.
368  *
369  * Calculates the positions using \c t_methoddata_pos::pc for the index group
370  * in \c t_methoddata_pos::g and stores the results in \p out->u.p.
371  */
372 static int
373 evaluate_pos(t_topology *top, t_trxframe *fr, t_pbc *pbc,
374              gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
375 {
376     t_methoddata_pos *d = (t_methoddata_pos *)data;
377
378     gmx_ana_poscalc_update(d->pc, out->u.p, &d->g, fr, pbc);
379     return 0;
380 }