Split g96 I/O routines from confio.cpp
[gromacs.git] / src / gromacs / selection / sm_position.cpp
blobad174aa1591118fbec7f170d7de7348663886637
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015, 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.
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.
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.
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.
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.
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.
35 /*! \internal \file
36 * \brief
37 * Implements position evaluation selection methods.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
42 #include "gmxpre.h"
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/selection/indexutil.h"
46 #include "gromacs/selection/position.h"
47 #include "gromacs/utility/cstringutil.h"
48 #include "gromacs/utility/smalloc.h"
50 #include "keywords.h"
51 #include "poscalc.h"
52 #include "selelem.h"
53 #include "selmethod.h"
55 /*! \internal \brief
56 * Data structure for position keyword evaluation.
58 typedef struct
60 /** Position calculation collection to use. */
61 gmx::PositionCalculationCollection *pcc;
62 /** Index group for which the center should be evaluated. */
63 gmx_ana_index_t g;
64 /** Position evaluation data structure. */
65 gmx_ana_poscalc_t *pc;
66 /** true if periodic boundary conditions should be used. */
67 bool bPBC;
68 /** Type of positions to calculate. */
69 char *type;
70 /** Flags for the position calculation. */
71 int flags;
72 } t_methoddata_pos;
74 /** Allocates data for position evaluation selection methods. */
75 static void *
76 init_data_pos(int npar, gmx_ana_selparam_t *param);
77 /** Sets the position calculation collection for position evaluation selection methods. */
78 static void
79 set_poscoll_pos(gmx::PositionCalculationCollection *pcc, void *data);
80 /*! \brief
81 * Initializes position evaluation keywords.
83 * \param[in] top Not used.
84 * \param[in] npar Not used.
85 * \param[in] param Not used.
86 * \param[in,out] data Should point to \c t_methoddata_pos.
87 * \returns 0 on success, a non-zero error code on error.
89 * The \c t_methoddata_pos::type field should have been initialized
90 * externally using _gmx_selelem_set_kwpos_type().
92 static void
93 init_kwpos(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
94 /*! \brief
95 * Initializes the \p cog selection method.
97 * \param[in] top Topology data structure.
98 * \param[in] npar Not used.
99 * \param[in] param Not used.
100 * \param[in,out] data Should point to \c t_methoddata_pos.
101 * \returns 0 on success, a non-zero error code on error.
103 static void
104 init_cog(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
105 /*! \brief
106 * Initializes the \p cog selection method.
108 * \param[in] top Topology data structure.
109 * \param[in] npar Not used.
110 * \param[in] param Not used.
111 * \param[in,out] data Should point to \c t_methoddata_pos.
112 * \returns 0 on success, a non-zero error code on error.
114 static void
115 init_com(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
116 /*! \brief
117 * Initializes output for position evaluation selection methods.
119 * \param[in] top Topology data structure.
120 * \param[in,out] out Pointer to output data structure.
121 * \param[in,out] data Should point to \c t_methoddata_pos.
122 * \returns 0 for success.
124 static void
125 init_output_pos(t_topology *top, gmx_ana_selvalue_t *out, void *data);
126 /** Frees the data allocated for position evaluation selection methods. */
127 static void
128 free_data_pos(void *data);
129 /** Evaluates position evaluation selection methods. */
130 static void
131 evaluate_pos(t_topology * /* top */, t_trxframe *fr, t_pbc *pbc, gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void *data);
133 /** Parameters for position keyword evaluation. */
134 static gmx_ana_selparam_t smparams_keyword_pos[] = {
135 {NULL, {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
138 /** Parameters for the \p cog and \p com selection methods. */
139 static gmx_ana_selparam_t smparams_com[] = {
140 {"of", {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
141 {"pbc", {NO_VALUE, 0, {NULL}}, NULL, 0},
144 /** Selection method data for position keyword evaluation. */
145 gmx_ana_selmethod_t sm_keyword_pos = {
146 "kw_pos", POS_VALUE, SMETH_DYNAMIC | SMETH_VARNUMVAL | SMETH_ALLOW_UNSORTED,
147 asize(smparams_keyword_pos), smparams_keyword_pos,
148 &init_data_pos,
149 &set_poscoll_pos,
150 &init_kwpos,
151 &init_output_pos,
152 &free_data_pos,
153 NULL,
154 &evaluate_pos,
155 NULL,
156 {NULL, NULL, 0, NULL},
159 /** Selection method data for the \p cog method. */
160 gmx_ana_selmethod_t sm_cog = {
161 "cog", POS_VALUE, SMETH_DYNAMIC | SMETH_SINGLEVAL,
162 asize(smparams_com), smparams_com,
163 &init_data_pos,
164 &set_poscoll_pos,
165 &init_cog,
166 &init_output_pos,
167 &free_data_pos,
168 NULL,
169 &evaluate_pos,
170 NULL,
171 {"cog of ATOM_EXPR [pbc]", NULL, 0, NULL},
174 /** Selection method data for the \p com method. */
175 gmx_ana_selmethod_t sm_com = {
176 "com", POS_VALUE, SMETH_REQTOP | SMETH_DYNAMIC | SMETH_SINGLEVAL,
177 asize(smparams_com), smparams_com,
178 &init_data_pos,
179 &set_poscoll_pos,
180 &init_com,
181 &init_output_pos,
182 &free_data_pos,
183 NULL,
184 &evaluate_pos,
185 NULL,
186 {"com of ATOM_EXPR [pbc]", NULL, 0, NULL},
190 * \param[in] npar Should be 1 or 2.
191 * \param[in,out] param Method parameters (should point to
192 * \ref smparams_keyword_pos or \ref smparams_com).
193 * \returns Pointer to the allocated data (\c t_methoddata_pos).
195 * Allocates memory for a \c t_methoddata_pos structure and initializes
196 * the first parameter to define the value for \c t_methoddata_pos::g.
197 * If a second parameter is present, it is used for setting the
198 * \c t_methoddata_pos::bPBC flag.
200 static void *
201 init_data_pos(int npar, gmx_ana_selparam_t *param)
203 t_methoddata_pos *data;
205 snew(data, 1);
206 param[0].val.u.g = &data->g;
207 if (npar > 1)
209 param[1].val.u.b = &data->bPBC;
211 data->pc = NULL;
212 data->bPBC = false;
213 data->type = NULL;
214 data->flags = -1;
215 return data;
219 * \param[in] pcc Position calculation collection to use.
220 * \param[in,out] data Should point to \c t_methoddata_pos.
222 static void
223 set_poscoll_pos(gmx::PositionCalculationCollection *pcc, void *data)
225 ((t_methoddata_pos *)data)->pcc = pcc;
228 bool
229 _gmx_selelem_is_default_kwpos(const gmx::SelectionTreeElement &sel)
231 if (sel.type != SEL_EXPRESSION || !sel.u.expr.method
232 || sel.u.expr.method->name != sm_keyword_pos.name)
234 return false;
237 t_methoddata_pos *d = static_cast<t_methoddata_pos *>(sel.u.expr.mdata);
238 return d->type == NULL;
242 * \param[in,out] sel Selection element to initialize.
243 * \param[in] type One of the enum values acceptable for
244 * PositionCalculationCollection::typeFromEnum().
246 * Initializes the reference position type for position evaluation.
247 * If called multiple times, the first setting takes effect, and later calls
248 * are neglected.
250 void
251 _gmx_selelem_set_kwpos_type(gmx::SelectionTreeElement *sel, const char *type)
253 t_methoddata_pos *d = (t_methoddata_pos *)sel->u.expr.mdata;
255 if (sel->type != SEL_EXPRESSION || !sel->u.expr.method
256 || sel->u.expr.method->name != sm_keyword_pos.name)
258 return;
260 if (!d->type && type)
262 d->type = gmx_strdup(type);
263 /* FIXME: It would be better not to have the string here hardcoded. */
264 if (type[0] != 'a')
266 sel->u.expr.method->flags |= SMETH_REQTOP;
272 * \param[in,out] sel Selection element to initialize.
273 * \param[in] flags Default completion flags
274 * (see PositionCalculationCollection::typeFromEnum()).
276 * Initializes the flags for position evaluation.
277 * If called multiple times, the first setting takes effect, and later calls
278 * are neglected.
280 void
281 _gmx_selelem_set_kwpos_flags(gmx::SelectionTreeElement *sel, int flags)
283 t_methoddata_pos *d = (t_methoddata_pos *)sel->u.expr.mdata;
285 if (sel->type != SEL_EXPRESSION || !sel->u.expr.method
286 || sel->u.expr.method->name != sm_keyword_pos.name)
288 return;
290 if (d->flags == -1)
292 d->flags = flags;
296 static void
297 init_kwpos(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
299 t_methoddata_pos *d = (t_methoddata_pos *)data;
301 if (!(param[0].flags & SPAR_DYNAMIC))
303 d->flags &= ~(POS_DYNAMIC | POS_MASKONLY);
305 else if (!(d->flags & POS_MASKONLY))
307 d->flags |= POS_DYNAMIC;
309 d->pc = d->pcc->createCalculationFromEnum(d->type, d->flags);
310 gmx_ana_poscalc_set_maxindex(d->pc, &d->g);
313 static void
314 init_cog(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
316 t_methoddata_pos *d = (t_methoddata_pos *)data;
318 d->flags = (param[0].flags & SPAR_DYNAMIC) ? POS_DYNAMIC : 0;
319 d->pc = d->pcc->createCalculation(d->bPBC ? POS_ALL_PBC : POS_ALL, d->flags);
320 gmx_ana_poscalc_set_maxindex(d->pc, &d->g);
323 static void
324 init_com(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
326 t_methoddata_pos *d = (t_methoddata_pos *)data;
328 d->flags = (param[0].flags & SPAR_DYNAMIC) ? POS_DYNAMIC : 0;
329 d->flags |= POS_MASS;
330 d->pc = d->pcc->createCalculation(d->bPBC ? POS_ALL_PBC : POS_ALL, d->flags);
331 gmx_ana_poscalc_set_maxindex(d->pc, &d->g);
334 static void
335 init_output_pos(t_topology * /* top */, gmx_ana_selvalue_t *out, void *data)
337 t_methoddata_pos *d = (t_methoddata_pos *)data;
339 gmx_ana_poscalc_init_pos(d->pc, out->u.p);
343 * \param data Data to free (should point to a \c t_methoddata_pos).
345 * Frees the memory allocated for \c t_methoddata_pos::g and
346 * \c t_methoddata_pos::pc.
348 static void
349 free_data_pos(void *data)
351 t_methoddata_pos *d = (t_methoddata_pos *)data;
353 sfree(d->type);
354 gmx_ana_poscalc_free(d->pc);
355 sfree(d);
359 * See sel_updatefunc() for description of the parameters.
360 * \p data should point to a \c t_methoddata_pos.
362 * Calculates the positions using \c t_methoddata_pos::pc for the index group
363 * in \c t_methoddata_pos::g and stores the results in \p out->u.p.
365 static void
366 evaluate_pos(t_topology * /* top */, t_trxframe *fr, t_pbc *pbc,
367 gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void *data)
369 t_methoddata_pos *d = (t_methoddata_pos *)data;
371 gmx_ana_poscalc_update(d->pc, out->u.p, &d->g, fr, pbc);