2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013 by the GROMACS development team.
5 * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
6 * Copyright (c) 2019,2020, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * Implements the \p merge and \p plus selection modifiers.
41 * \author Teemu Murtola <teemu.murtola@gmail.com>
42 * \ingroup module_selection
46 #include "gromacs/math/vec.h"
47 #include "gromacs/utility/arraysize.h"
48 #include "gromacs/utility/basedefinitions.h"
49 #include "gromacs/utility/exceptions.h"
50 #include "gromacs/utility/smalloc.h"
53 #include "selmethod.h"
54 #include "selmethod_impl.h"
57 * Data structure for the merging selection modifiers.
59 typedef struct methoddata_merge
61 /** Input positions. */
63 /** Other input positions. */
65 /** Stride for merging (\c stride values from \c p1 for each in \c p2). */
69 /** Allocates data for the merging selection modifiers. */
70 static void* init_data_merge(int npar
, gmx_ana_selparam_t
* param
);
72 * Initializes data for the merging selection modifiers.
74 * \param[in] top Not used.
75 * \param[in] npar Not used (should be 2 or 3).
76 * \param[in] param Method parameters (should point to \ref smparams_merge).
77 * \param[in] data Should point to a \p t_methoddata_merge.
78 * \returns 0 if everything is successful, -1 on error.
80 static void init_merge(const gmx_mtop_t
* top
, int npar
, gmx_ana_selparam_t
* param
, void* data
);
81 /** Initializes output for the \p merge selection modifier. */
82 static void init_output_merge(const gmx_mtop_t
* top
, gmx_ana_selvalue_t
* out
, void* data
);
83 /** Initializes output for the \p plus selection modifier. */
84 static void init_output_plus(const gmx_mtop_t
* top
, gmx_ana_selvalue_t
* out
, void* data
);
85 /** Frees the memory allocated for the merging selection modifiers. */
86 static void free_data_merge(void* data
);
88 * Evaluates the \p merge selection modifier.
90 * \param[in] context Not used.
91 * \param[in] p Positions to merge (should point to \p data->p1).
92 * \param[out] out Output data structure (\p out->u.p is used).
93 * \param[in] data Should point to a \p t_methoddata_merge.
95 static void evaluate_merge(const gmx::SelMethodEvalContext
& context
,
97 gmx_ana_selvalue_t
* out
,
100 * Evaluates the \p plus selection modifier.
102 * \param[in] context Not used.
103 * \param[in] p Positions to merge (should point to \p data->p1).
104 * \param[out] out Output data structure (\p out->u.p is used).
105 * \param[in] data Should point to a \p t_methoddata_merge.
107 static void evaluate_plus(const gmx::SelMethodEvalContext
& context
,
109 gmx_ana_selvalue_t
* out
,
112 /** Parameters for the merging selection modifiers. */
113 static gmx_ana_selparam_t smparams_merge
[] = {
114 { nullptr, { POS_VALUE
, -1, { nullptr } }, nullptr, SPAR_DYNAMIC
| SPAR_VARNUM
},
115 { nullptr, { POS_VALUE
, -1, { nullptr } }, nullptr, SPAR_DYNAMIC
| SPAR_VARNUM
},
116 { "stride", { INT_VALUE
, 1, { nullptr } }, nullptr, SPAR_OPTIONAL
},
119 //! Help title for the merging selection modifiers.
120 static const char helptitle_merge
[] = "Merging selections";
121 //! Help text for the merging selection modifiers.
122 static const char* const help_merge
[] = {
125 " POSEXPR merge POSEXPR [stride INT]",
126 " POSEXPR merge POSEXPR [merge POSEXPR ...]",
127 " POSEXPR plus POSEXPR [plus POSEXPR ...]",
129 "Basic selection keywords can only create selections where each atom",
130 "occurs at most once. The [TT]merge[tt] and [TT]plus[tt] selection",
131 "keywords can be used to work around this limitation. Both create",
132 "a selection that contains the positions from all the given position",
133 "expressions, even if they contain duplicates.",
134 "The difference between the two is that [TT]merge[tt] expects two or more",
135 "selections with the same number of positions, and the output contains",
136 "the input positions selected from each expression in turn, i.e.,",
137 "the output is like A1 B1 A2 B2 and so on. It is also possible to merge",
138 "selections of unequal size as long as the size of the first is a",
139 "multiple of the second one. The [TT]stride[tt] parameter can be used",
140 "to explicitly provide this multiplicity.",
141 "[TT]plus[tt] simply concatenates the positions after each other, and",
142 "can work also with selections of different sizes.",
143 "These keywords are valid only at the selection level, not in any",
147 /** Selection method data for the \p plus modifier. */
148 gmx_ana_selmethod_t sm_merge
= {
152 asize(smparams_merge
),
162 { "merge POSEXPR", helptitle_merge
, asize(help_merge
), help_merge
},
165 /** Selection method data for the \p plus modifier. */
166 gmx_ana_selmethod_t sm_plus
= {
170 asize(smparams_merge
) - 1,
180 { "plus POSEXPR", helptitle_merge
, asize(help_merge
), help_merge
},
184 * \param[in] npar Should be 2 for \c plus and 3 for \c merge.
185 * \param[in,out] param Method parameters (should point to a copy of
186 * \ref smparams_merge).
187 * \returns Pointer to the allocated data (\p t_methoddata_merge).
189 * Allocates memory for a \p t_methoddata_merge structure.
191 static void* init_data_merge(int npar
, gmx_ana_selparam_t
* param
)
193 t_methoddata_merge
* data
= new t_methoddata_merge();
195 param
[0].val
.u
.p
= &data
->p1
;
196 param
[1].val
.u
.p
= &data
->p2
;
199 param
[2].val
.u
.i
= &data
->stride
;
204 static void init_merge(const gmx_mtop_t
* /* top */, int /* npar */, gmx_ana_selparam_t
* /* param */, void* data
)
206 t_methoddata_merge
* d
= static_cast<t_methoddata_merge
*>(data
);
210 GMX_THROW(gmx::InvalidInputError("Stride for merging should be positive"));
212 /* If no stride given, deduce it from the input sizes */
215 d
->stride
= d
->p1
.count() / d
->p2
.count();
217 if (d
->p1
.count() != d
->stride
* d
->p2
.count())
219 GMX_THROW(gmx::InconsistentInputError(
220 "The number of positions to be merged are not compatible"));
225 * Does common initialization to all merging modifiers.
227 * \param[in] top Topology data structure.
228 * \param[in,out] out Pointer to output data structure.
229 * \param[in,out] data Should point to \c t_methoddata_merge.
231 static void init_output_common(const gmx_mtop_t
* top
, gmx_ana_selvalue_t
* out
, void* data
)
233 t_methoddata_merge
* d
= static_cast<t_methoddata_merge
*>(data
);
235 GMX_UNUSED_VALUE(top
);
236 if (d
->p1
.m
.type
!= d
->p2
.m
.type
)
238 /* TODO: Maybe we could pick something else here? */
239 out
->u
.p
->m
.type
= INDEX_UNKNOWN
;
243 out
->u
.p
->m
.type
= d
->p1
.m
.type
;
245 gmx_ana_pos_reserve_for_append(out
->u
.p
, d
->p1
.count() + d
->p2
.count(),
246 d
->p1
.m
.b
.nra
+ d
->p2
.m
.b
.nra
, d
->p1
.v
!= nullptr, d
->p1
.f
!= nullptr);
247 gmx_ana_pos_empty_init(out
->u
.p
);
251 * \param[in] top Topology data structure.
252 * \param[in,out] out Pointer to output data structure.
253 * \param[in,out] data Should point to \c t_methoddata_merge.
255 static void init_output_merge(const gmx_mtop_t
* top
, gmx_ana_selvalue_t
* out
, void* data
)
257 t_methoddata_merge
* d
= static_cast<t_methoddata_merge
*>(data
);
260 init_output_common(top
, out
, data
);
261 for (i
= 0; i
< d
->p2
.count(); ++i
)
263 for (j
= 0; j
< d
->stride
; ++j
)
265 gmx_ana_pos_append_init(out
->u
.p
, &d
->p1
, d
->stride
* i
+ j
);
267 gmx_ana_pos_append_init(out
->u
.p
, &d
->p2
, i
);
272 * \param[in] top Topology data structure.
273 * \param[in,out] out Pointer to output data structure.
274 * \param[in,out] data Should point to \c t_methoddata_merge.
276 static void init_output_plus(const gmx_mtop_t
* top
, gmx_ana_selvalue_t
* out
, void* data
)
278 t_methoddata_merge
* d
= static_cast<t_methoddata_merge
*>(data
);
281 init_output_common(top
, out
, data
);
282 for (i
= 0; i
< d
->p1
.count(); ++i
)
284 gmx_ana_pos_append_init(out
->u
.p
, &d
->p1
, i
);
286 for (i
= 0; i
< d
->p2
.count(); ++i
)
288 gmx_ana_pos_append_init(out
->u
.p
, &d
->p2
, i
);
293 * \param data Data to free (should point to a \p t_methoddata_merge).
295 * Frees the memory allocated for \c t_methoddata_merge.
297 static void free_data_merge(void* data
)
299 t_methoddata_merge
* d
= static_cast<t_methoddata_merge
*>(data
);
303 static void evaluate_merge(const gmx::SelMethodEvalContext
& /*context*/,
304 gmx_ana_pos_t
* /* p */,
305 gmx_ana_selvalue_t
* out
,
308 t_methoddata_merge
* d
= static_cast<t_methoddata_merge
*>(data
);
312 if (d
->p1
.count() != d
->stride
* d
->p2
.count())
314 GMX_THROW(gmx::InconsistentInputError(
315 "The number of positions to be merged are not compatible"));
317 gmx_ana_pos_empty(out
->u
.p
);
318 for (i
= 0; i
< d
->p2
.count(); ++i
)
320 for (j
= 0; j
< d
->stride
; ++j
)
322 refid
= d
->p1
.m
.refid
[d
->stride
* i
+ j
];
325 refid
= (d
->stride
+ 1) * (refid
/ d
->stride
) + (refid
% d
->stride
);
327 gmx_ana_pos_append(out
->u
.p
, &d
->p1
, d
->stride
* i
+ j
, refid
);
329 refid
= (d
->stride
+ 1) * d
->p2
.m
.refid
[i
] + d
->stride
;
330 gmx_ana_pos_append(out
->u
.p
, &d
->p2
, i
, refid
);
332 gmx_ana_pos_append_finish(out
->u
.p
);
335 static void evaluate_plus(const gmx::SelMethodEvalContext
& /*context*/,
336 gmx_ana_pos_t
* /* p */,
337 gmx_ana_selvalue_t
* out
,
340 t_methoddata_merge
* d
= static_cast<t_methoddata_merge
*>(data
);
344 gmx_ana_pos_empty(out
->u
.p
);
345 for (i
= 0; i
< d
->p1
.count(); ++i
)
347 refid
= d
->p1
.m
.refid
[i
];
348 gmx_ana_pos_append(out
->u
.p
, &d
->p1
, i
, refid
);
350 for (i
= 0; i
< d
->p2
.count(); ++i
)
352 refid
= d
->p2
.m
.refid
[i
];
355 refid
+= d
->p1
.m
.b
.nr
;
357 gmx_ana_pos_append(out
->u
.p
, &d
->p2
, i
, refid
);
359 gmx_ana_pos_append_finish(out
->u
.p
);