Remove mols from gmx_mtop_t
[gromacs.git] / src / gromacs / selection / sm_simple.cpp
blob619a19f2c2af15881bacd696beeeaf08f5abcc38
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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 simple keyword selection methods.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
42 #include "gmxpre.h"
44 #include <cctype>
46 #include "gromacs/selection/position.h"
47 #include "gromacs/topology/mtop_lookup.h"
48 #include "gromacs/topology/topology.h"
49 #include "gromacs/utility/arraysize.h"
50 #include "gromacs/utility/exceptions.h"
51 #include "gromacs/utility/gmxassert.h"
53 #include "selmethod.h"
55 /** Evaluates the \p all selection keyword. */
56 static void
57 evaluate_all(const gmx::SelMethodEvalContext &context,
58 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
59 /** Evaluates the \p none selection keyword. */
60 static void
61 evaluate_none(const gmx::SelMethodEvalContext &context,
62 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
63 /** Evaluates the \p atomnr selection keyword. */
64 static void
65 evaluate_atomnr(const gmx::SelMethodEvalContext &context,
66 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
67 /** Evaluates the \p resnr selection keyword. */
68 static void
69 evaluate_resnr(const gmx::SelMethodEvalContext &context,
70 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
71 /** Evaluates the \p resindex selection keyword. */
72 static void
73 evaluate_resindex(const gmx::SelMethodEvalContext &context,
74 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
75 /*! \brief
76 * Checks whether molecule information is present in the topology.
78 * \param[in] top Topology structure.
79 * \param npar Not used.
80 * \param param Not used.
81 * \param data Not used.
82 * \returns 0 if molecule info is present in the topology, -1 otherwise.
84 * If molecule information is not found, also prints an error message.
86 static void
87 check_molecules(const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data);
88 /** Evaluates the \p molindex selection keyword. */
89 static void
90 evaluate_molindex(const gmx::SelMethodEvalContext &context,
91 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
92 /** Evaluates the \p atomname selection keyword. */
93 static void
94 evaluate_atomname(const gmx::SelMethodEvalContext &context,
95 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
96 /** Evaluates the \p pdbatomname selection keyword. */
97 static void
98 evaluate_pdbatomname(const gmx::SelMethodEvalContext &context,
99 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
100 /*! \brief
101 * Checks whether atom types are present in the topology.
103 * \param[in] top Topology structure.
104 * \param npar Not used.
105 * \param param Not used.
106 * \param data Not used.
108 static void
109 check_atomtype(const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data);
110 /** Evaluates the \p atomtype selection keyword. */
111 static void
112 evaluate_atomtype(const gmx::SelMethodEvalContext &context,
113 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
114 /** Evaluates the \p insertcode selection keyword. */
115 static void
116 evaluate_insertcode(const gmx::SelMethodEvalContext &context,
117 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
118 /** Evaluates the \p chain selection keyword. */
119 static void
120 evaluate_chain(const gmx::SelMethodEvalContext &context,
121 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
122 /** Evaluates the \p mass selection keyword. */
123 static void
124 evaluate_mass(const gmx::SelMethodEvalContext &context,
125 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
126 /*! \brief
127 * Checks whether charges are present in the topology.
129 * \param[in] top Topology structure.
130 * \param npar Not used.
131 * \param param Not used.
132 * \param data Not used.
134 static void
135 check_charge(const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data);
136 /** Evaluates the \p charge selection keyword. */
137 static void
138 evaluate_charge(const gmx::SelMethodEvalContext &context,
139 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
140 /*! \brief
141 * Checks whether PDB info is present in the topology.
143 * \param[in] top Topology structure.
144 * \param npar Not used.
145 * \param param Not used.
146 * \param data Not used.
147 * \returns 0 if PDB info is present in the topology, -1 otherwise.
149 * If PDB info is not found, also prints an error message.
151 static void
152 check_pdbinfo(const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data);
153 /** Evaluates the \p altloc selection keyword. */
154 static void
155 evaluate_altloc(const gmx::SelMethodEvalContext &context,
156 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
157 /** Evaluates the \p occupancy selection keyword. */
158 static void
159 evaluate_occupancy(const gmx::SelMethodEvalContext &context,
160 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
161 /** Evaluates the \p betafactor selection keyword. */
162 static void
163 evaluate_betafactor(const gmx::SelMethodEvalContext &context,
164 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
165 /** Evaluates the \p resname selection keyword. */
166 static void
167 evaluate_resname(const gmx::SelMethodEvalContext &context,
168 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
170 /** Evaluates the \p x selection keyword. */
171 static void
172 evaluate_x(const gmx::SelMethodEvalContext &context,
173 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
174 /** Evaluates the \p y selection keyword. */
175 static void
176 evaluate_y(const gmx::SelMethodEvalContext &context,
177 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
178 /** Evaluates the \p z selection keyword. */
179 static void
180 evaluate_z(const gmx::SelMethodEvalContext &context,
181 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
183 //! Help title for atom name selection keywords.
184 static const char helptitle_atomname[] = "Selecting atoms by name";
185 //! Help text for atom name selection keywords.
186 static const char *const help_atomname[] = {
187 "::",
189 " name",
190 " pdbname",
191 " atomname",
192 " pdbatomname",
194 "These keywords select atoms by name. [TT]name[tt] selects atoms using",
195 "the GROMACS atom naming convention.",
196 "For input formats other than PDB, the atom names are matched exactly",
197 "as they appear in the input file. For PDB files, 4 character atom names",
198 "that start with a digit are matched after moving the digit to the end",
199 "(e.g., to match 3HG2 from a PDB file, use [TT]name HG23[tt]).",
200 "[TT]pdbname[tt] can only be used with a PDB input file, and selects",
201 "atoms based on the exact name given in the input file, without the",
202 "transformation described above.[PAR]",
204 "[TT]atomname[tt] and [TT]pdbatomname[tt] are synonyms for the above two",
205 "keywords."
208 //! Help title for residue index selection keywords.
209 static const char helptitle_resindex[] = "Selecting atoms by residue number";
210 //! Help text for residue index selection keywords.
211 static const char *const help_resindex[] = {
212 "::",
214 " resnr",
215 " resid",
216 " resindex",
217 " residue",
219 "[TT]resnr[tt] selects atoms using the residue numbering in the input",
220 "file. [TT]resid[tt] is synonym for this keyword for VMD compatibility.",
222 "[TT]resindex N[tt] selects the [TT]N[tt]th residue starting from the",
223 "beginning of the input file. This is useful for uniquely identifying",
224 "residues if there are duplicate numbers in the input file (e.g., in",
225 "multiple chains).",
226 "[TT]residue[tt] is a synonym for [TT]resindex[tt]. This allows",
227 "[TT]same residue as[tt] to work as expected."
230 /** Selection method data for \p all selection keyword. */
231 gmx_ana_selmethod_t sm_all = {
232 "all", GROUP_VALUE, 0,
233 0, nullptr,
234 nullptr,
235 nullptr,
236 nullptr,
237 nullptr,
238 nullptr,
239 nullptr,
240 &evaluate_all,
241 nullptr,
244 /** Selection method data for \p none selection keyword. */
245 gmx_ana_selmethod_t sm_none = {
246 "none", GROUP_VALUE, 0,
247 0, nullptr,
248 nullptr,
249 nullptr,
250 nullptr,
251 nullptr,
252 nullptr,
253 nullptr,
254 &evaluate_none,
255 nullptr,
258 /** Selection method data for \p atomnr selection keyword. */
259 gmx_ana_selmethod_t sm_atomnr = {
260 "atomnr", INT_VALUE, 0,
261 0, nullptr,
262 nullptr,
263 nullptr,
264 nullptr,
265 nullptr,
266 nullptr,
267 nullptr,
268 &evaluate_atomnr,
269 nullptr,
272 /** Selection method data for \p resnr selection keyword. */
273 gmx_ana_selmethod_t sm_resnr = {
274 "resnr", INT_VALUE, SMETH_REQTOP,
275 0, nullptr,
276 nullptr,
277 nullptr,
278 nullptr,
279 nullptr,
280 nullptr,
281 nullptr,
282 &evaluate_resnr,
283 nullptr,
284 {nullptr, helptitle_resindex, asize(help_resindex), help_resindex}
287 /** Selection method data for \p resindex selection keyword. */
288 gmx_ana_selmethod_t sm_resindex = {
289 "resindex", INT_VALUE, SMETH_REQTOP,
290 0, nullptr,
291 nullptr,
292 nullptr,
293 nullptr,
294 nullptr,
295 nullptr,
296 nullptr,
297 &evaluate_resindex,
298 nullptr,
299 {nullptr, helptitle_resindex, asize(help_resindex), help_resindex}
302 /** Selection method data for \p molindex selection keyword. */
303 gmx_ana_selmethod_t sm_molindex = {
304 "molindex", INT_VALUE, SMETH_REQTOP,
305 0, nullptr,
306 nullptr,
307 nullptr,
308 &check_molecules,
309 nullptr,
310 nullptr,
311 nullptr,
312 &evaluate_molindex,
313 nullptr,
316 /** Selection method data for \p atomname selection keyword. */
317 gmx_ana_selmethod_t sm_atomname = {
318 "atomname", STR_VALUE, SMETH_REQTOP,
319 0, nullptr,
320 nullptr,
321 nullptr,
322 nullptr,
323 nullptr,
324 nullptr,
325 nullptr,
326 &evaluate_atomname,
327 nullptr,
328 {nullptr, helptitle_atomname, asize(help_atomname), help_atomname}
331 /** Selection method data for \p pdbatomname selection keyword. */
332 gmx_ana_selmethod_t sm_pdbatomname = {
333 "pdbatomname", STR_VALUE, SMETH_REQTOP,
334 0, nullptr,
335 nullptr,
336 nullptr,
337 &check_pdbinfo,
338 nullptr,
339 nullptr,
340 nullptr,
341 &evaluate_pdbatomname,
342 nullptr,
343 {nullptr, helptitle_atomname, asize(help_atomname), help_atomname}
346 /** Selection method data for \p atomtype selection keyword. */
347 gmx_ana_selmethod_t sm_atomtype = {
348 "atomtype", STR_VALUE, SMETH_REQTOP,
349 0, nullptr,
350 nullptr,
351 nullptr,
352 &check_atomtype,
353 nullptr,
354 nullptr,
355 nullptr,
356 &evaluate_atomtype,
357 nullptr,
360 /** Selection method data for \p resname selection keyword. */
361 gmx_ana_selmethod_t sm_resname = {
362 "resname", STR_VALUE, SMETH_REQTOP,
363 0, nullptr,
364 nullptr,
365 nullptr,
366 nullptr,
367 nullptr,
368 nullptr,
369 nullptr,
370 &evaluate_resname,
371 nullptr,
374 /** Selection method data for \p chain selection keyword. */
375 gmx_ana_selmethod_t sm_insertcode = {
376 "insertcode", STR_VALUE, SMETH_REQTOP | SMETH_CHARVAL,
377 0, nullptr,
378 nullptr,
379 nullptr,
380 nullptr,
381 nullptr,
382 nullptr,
383 nullptr,
384 &evaluate_insertcode,
385 nullptr,
388 /** Selection method data for \p chain selection keyword. */
389 gmx_ana_selmethod_t sm_chain = {
390 "chain", STR_VALUE, SMETH_REQTOP | SMETH_CHARVAL,
391 0, nullptr,
392 nullptr,
393 nullptr,
394 nullptr,
395 nullptr,
396 nullptr,
397 nullptr,
398 &evaluate_chain,
399 nullptr,
402 /** Selection method data for \p mass selection keyword. */
403 gmx_ana_selmethod_t sm_mass = {
404 "mass", REAL_VALUE, SMETH_REQMASS,
405 0, nullptr,
406 nullptr,
407 nullptr,
408 nullptr,
409 nullptr,
410 nullptr,
411 nullptr,
412 &evaluate_mass,
413 nullptr,
416 /** Selection method data for \p charge selection keyword. */
417 gmx_ana_selmethod_t sm_charge = {
418 "charge", REAL_VALUE, SMETH_REQTOP,
419 0, nullptr,
420 nullptr,
421 nullptr,
422 &check_charge,
423 nullptr,
424 nullptr,
425 nullptr,
426 &evaluate_charge,
427 nullptr,
430 /** Selection method data for \p chain selection keyword. */
431 gmx_ana_selmethod_t sm_altloc = {
432 "altloc", STR_VALUE, SMETH_REQTOP | SMETH_CHARVAL,
433 0, nullptr,
434 nullptr,
435 nullptr,
436 &check_pdbinfo,
437 nullptr,
438 nullptr,
439 nullptr,
440 &evaluate_altloc,
441 nullptr,
444 /** Selection method data for \p occupancy selection keyword. */
445 gmx_ana_selmethod_t sm_occupancy = {
446 "occupancy", REAL_VALUE, SMETH_REQTOP,
447 0, nullptr,
448 nullptr,
449 nullptr,
450 &check_pdbinfo,
451 nullptr,
452 nullptr,
453 nullptr,
454 &evaluate_occupancy,
455 nullptr,
458 /** Selection method data for \p betafactor selection keyword. */
459 gmx_ana_selmethod_t sm_betafactor = {
460 "betafactor", REAL_VALUE, SMETH_REQTOP,
461 0, nullptr,
462 nullptr,
463 nullptr,
464 &check_pdbinfo,
465 nullptr,
466 nullptr,
467 nullptr,
468 &evaluate_betafactor,
469 nullptr,
472 /** Selection method data for \p x selection keyword. */
473 gmx_ana_selmethod_t sm_x = {
474 "x", REAL_VALUE, SMETH_DYNAMIC,
475 0, nullptr,
476 nullptr,
477 nullptr,
478 nullptr,
479 nullptr,
480 nullptr,
481 nullptr,
482 nullptr,
483 &evaluate_x,
486 /** Selection method data for \p y selection keyword. */
487 gmx_ana_selmethod_t sm_y = {
488 "y", REAL_VALUE, SMETH_DYNAMIC,
489 0, nullptr,
490 nullptr,
491 nullptr,
492 nullptr,
493 nullptr,
494 nullptr,
495 nullptr,
496 nullptr,
497 &evaluate_y,
500 /** Selection method data for \p z selection keyword. */
501 gmx_ana_selmethod_t sm_z = {
502 "z", REAL_VALUE, SMETH_DYNAMIC,
503 0, nullptr,
504 nullptr,
505 nullptr,
506 nullptr,
507 nullptr,
508 nullptr,
509 nullptr,
510 nullptr,
511 &evaluate_z,
515 * See sel_updatefunc() for description of the parameters.
516 * \p data is not used.
518 * Copies \p g to \p out->u.g.
520 static void
521 evaluate_all(const gmx::SelMethodEvalContext & /*context*/,
522 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
524 gmx_ana_index_copy(out->u.g, g, false);
528 * See sel_updatefunc() for description of the parameters.
529 * \p data is not used.
531 * Returns an empty \p out->u.g.
533 static void
534 evaluate_none(const gmx::SelMethodEvalContext & /*context*/,
535 gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void * /* data */)
537 out->u.g->isize = 0;
541 * See sel_updatefunc() for description of the parameters.
542 * \p data is not used.
544 * Returns the indices for each atom in \p out->u.i.
546 static void
547 evaluate_atomnr(const gmx::SelMethodEvalContext & /*context*/,
548 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
550 int i;
552 out->nr = g->isize;
553 for (i = 0; i < g->isize; ++i)
555 out->u.i[i] = g->index[i] + 1;
560 * See sel_updatefunc() for description of the parameters.
561 * \p data is not used.
563 * Returns the residue numbers for each atom in \p out->u.i.
565 static void
566 evaluate_resnr(const gmx::SelMethodEvalContext &context,
567 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
569 out->nr = g->isize;
570 int molb = 0;
571 for (int i = 0; i < g->isize; ++i)
573 mtopGetAtomAndResidueName(context.top, g->index[i], &molb,
574 nullptr, &out->u.i[i], nullptr, nullptr);
579 * See sel_updatefunc() for description of the parameters.
580 * \p data is not used.
582 * Returns the residue indices for each atom in \p out->u.i.
584 static void
585 evaluate_resindex(const gmx::SelMethodEvalContext &context,
586 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
588 out->nr = g->isize;
589 int molb = 0;
590 for (int i = 0; i < g->isize; ++i)
592 int resind;
593 mtopGetAtomAndResidueName(context.top, g->index[i], &molb,
594 nullptr, nullptr, nullptr, &resind);
595 out->u.i[i] = resind + 1;
599 static void
600 check_molecules(const gmx_mtop_t *top, int /* npar */, gmx_ana_selparam_t * /* param */, void * /* data */)
602 bool bOk;
604 bOk = (top != nullptr && top->haveMoleculeIndices);
605 if (!bOk)
607 GMX_THROW(gmx::InconsistentInputError("Molecule information not available in topology"));
612 * See sel_updatefunc() for description of the parameters.
613 * \p data is not used.
615 * Returns the molecule indices for each atom in \p out->u.i.
617 static void
618 evaluate_molindex(const gmx::SelMethodEvalContext &context,
619 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
621 out->nr = g->isize;
622 int molb = 0;
623 for (int i = 0; i < g->isize; ++i)
625 out->u.i[i] = mtopGetMoleculeIndex(context.top, g->index[i], &molb) + 1;
630 * See sel_updatefunc() for description of the parameters.
631 * \p data is not used.
633 * Returns the atom name for each atom in \p out->u.s.
635 static void
636 evaluate_atomname(const gmx::SelMethodEvalContext &context,
637 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
639 out->nr = g->isize;
640 int molb = 0;
641 for (int i = 0; i < g->isize; ++i)
643 const char *atom_name;
644 mtopGetAtomAndResidueName(context.top, g->index[i], &molb,
645 &atom_name, nullptr, nullptr, nullptr);
646 out->u.s[i] = const_cast<char *>(atom_name);
651 * See sel_updatefunc() for description of the parameters.
652 * \p data is not used.
654 * Returns the PDB atom name for each atom in \p out->u.s.
656 static void
657 evaluate_pdbatomname(const gmx::SelMethodEvalContext &context,
658 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
660 out->nr = g->isize;
661 int molb = 0;
662 for (int i = 0; i < g->isize; ++i)
664 const char *s = mtopGetAtomPdbInfo(context.top, g->index[i], &molb).atomnm;
665 while (std::isspace(*s))
667 ++s;
669 out->u.s[i] = const_cast<char *>(s);
673 static void
674 check_atomtype(const gmx_mtop_t *top, int /* npar */, gmx_ana_selparam_t * /* param */, void * /* data */)
676 if (!gmx_mtop_has_atomtypes(top))
678 GMX_THROW(gmx::InconsistentInputError("Atom types not available in topology"));
683 * See sel_updatefunc() for description of the parameters.
684 * \p data is not used.
686 * Returns the atom type for each atom in \p out->u.s.
687 * Segfaults if atom types are not found in the topology.
689 static void
690 evaluate_atomtype(const gmx::SelMethodEvalContext &context,
691 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
693 out->nr = g->isize;
694 int molb = 0;
695 for (int i = 0; i < g->isize; ++i)
697 int atomIndexInMolecule;
698 mtopGetMolblockIndex(context.top, g->index[i], &molb,
699 nullptr, &atomIndexInMolecule);
700 const gmx_moltype_t &moltype = context.top->moltype[context.top->molblock[molb].type];
701 out->u.s[i] = *moltype.atoms.atomtype[atomIndexInMolecule];
706 * See sel_updatefunc() for description of the parameters.
707 * \p data is not used.
709 * Returns the residue name for each atom in \p out->u.s.
711 static void
712 evaluate_resname(const gmx::SelMethodEvalContext &context,
713 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
715 out->nr = g->isize;
716 int molb = 0;
717 for (int i = 0; i < g->isize; ++i)
719 out->u.s[i] = *mtopGetResidueInfo(context.top, g->index[i], &molb).name;
724 * See sel_updatefunc() for description of the parameters.
725 * \p data is not used.
727 * Returns the insertion code for each atom in \p out->u.s.
729 static void
730 evaluate_insertcode(const gmx::SelMethodEvalContext &context,
731 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
733 out->nr = g->isize;
734 int molb = 0;
735 for (int i = 0; i < g->isize; ++i)
737 out->u.s[i][0] = mtopGetResidueInfo(context.top, g->index[i], &molb).ic;
742 * See sel_updatefunc() for description of the parameters.
743 * \p data is not used.
745 * Returns the chain for each atom in \p out->u.s.
747 static void
748 evaluate_chain(const gmx::SelMethodEvalContext &context,
749 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
751 out->nr = g->isize;
752 int molb = 0;
753 for (int i = 0; i < g->isize; ++i)
755 out->u.s[i][0] = mtopGetResidueInfo(context.top, g->index[i], &molb).chainid;
760 * See sel_updatefunc() for description of the parameters.
761 * \p data is not used.
763 * Returns the mass for each atom in \p out->u.r.
765 static void
766 evaluate_mass(const gmx::SelMethodEvalContext &context,
767 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
769 GMX_RELEASE_ASSERT(gmx_mtop_has_masses(context.top),
770 "Masses not available for evaluation");
771 out->nr = g->isize;
772 int molb = 0;
773 for (int i = 0; i < g->isize; ++i)
775 out->u.r[i] = mtopGetAtomMass(context.top, g->index[i], &molb);
780 static void
781 check_charge(const gmx_mtop_t *top, int /* npar */, gmx_ana_selparam_t * /* param */, void * /* data */)
783 if (!gmx_mtop_has_charges(top))
785 GMX_THROW(gmx::InconsistentInputError("Charges not available in topology"));
790 * See sel_updatefunc() for description of the parameters.
791 * \p data is not used.
793 * Returns the charge for each atom in \p out->u.r.
795 static void
796 evaluate_charge(const gmx::SelMethodEvalContext &context,
797 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
799 out->nr = g->isize;
800 int molb = 0;
801 for (int i = 0; i < g->isize; ++i)
803 out->u.r[i] = mtopGetAtomParameters(context.top, g->index[i], &molb).q;
807 static void
808 check_pdbinfo(const gmx_mtop_t *top, int /* npar */, gmx_ana_selparam_t * /* param */, void * /* data */)
810 if (!gmx_mtop_has_pdbinfo(top))
812 GMX_THROW(gmx::InconsistentInputError("PDB info not available in topology"));
817 * See sel_updatefunc() for description of the parameters.
818 * \p data is not used.
820 * Returns the alternate location identifier for each atom in \p out->u.s.
822 static void
823 evaluate_altloc(const gmx::SelMethodEvalContext &context,
824 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
826 out->nr = g->isize;
827 int molb = 0;
828 for (int i = 0; i < g->isize; ++i)
830 out->u.s[i][0] = mtopGetAtomPdbInfo(context.top, g->index[i], &molb).altloc;
835 * See sel_updatefunc() for description of the parameters.
836 * \p data is not used.
838 * Returns the occupancy numbers for each atom in \p out->u.r.
839 * Segfaults if PDB info is not found in the topology.
841 static void
842 evaluate_occupancy(const gmx::SelMethodEvalContext &context,
843 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
845 out->nr = g->isize;
846 int molb = 0;
847 for (int i = 0; i < g->isize; ++i)
849 out->u.r[i] = mtopGetAtomPdbInfo(context.top, g->index[i], &molb).occup;
854 * See sel_updatefunc() for description of the parameters.
855 * \p data is not used.
857 * Returns the B-factors for each atom in \p out->u.r.
858 * Segfaults if PDB info is not found in the topology.
860 static void
861 evaluate_betafactor(const gmx::SelMethodEvalContext &context,
862 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void * /* data */)
864 out->nr = g->isize;
865 int molb = 0;
866 for (int i = 0; i < g->isize; ++i)
868 out->u.r[i] = mtopGetAtomPdbInfo(context.top, g->index[i], &molb).bfac;
872 /*! \brief
873 * Internal utility function for position keyword evaluation.
875 * \param[out] out Output array.
876 * \param[in] pos Position data to use instead of atomic coordinates.
877 * \param[in] d Coordinate index to evaluate (\p XX, \p YY or \p ZZ).
879 * This function is used internally by evaluate_x(), evaluate_y() and
880 * evaluate_z() to do the actual evaluation.
882 static void
883 evaluate_coord(real out[], gmx_ana_pos_t *pos, int d)
885 for (int i = 0; i < pos->count(); ++i)
887 out[i] = pos->x[i][d];
889 // TODO: Make this more efficient by directly extracting the coordinates
890 // from the frame coordinates for atomic positions instead of going through
891 // a position calculation.
895 * See sel_updatefunc_pos() for description of the parameters.
896 * \p data is not used.
898 * Returns the \p x coordinate for each position in \p out->u.r.
900 static void
901 evaluate_x(const gmx::SelMethodEvalContext & /*context*/,
902 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void * /*data*/)
904 out->nr = pos->count();
905 evaluate_coord(out->u.r, pos, XX);
909 * See sel_updatefunc() for description of the parameters.
910 * \p data is not used.
912 * Returns the \p y coordinate for each position in \p out->u.r.
914 static void
915 evaluate_y(const gmx::SelMethodEvalContext & /*context*/,
916 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void * /*data*/)
918 out->nr = pos->count();
919 evaluate_coord(out->u.r, pos, YY);
923 * See sel_updatefunc() for description of the parameters.
924 * \p data is not used.
926 * Returns the \p z coordinate for each position in \p out->u.r.
928 static void
929 evaluate_z(const gmx::SelMethodEvalContext & /*context*/,
930 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void * /*data*/)
932 out->nr = pos->count();
933 evaluate_coord(out->u.r, pos, ZZ);