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