9e6362b579b0244514fac1b488c773133af62da0
[alexxy/gromacs.git] / src / gromacs / selection / selhelp.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2021, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /*! \internal \file
38  * \brief
39  * Implements functions in selhelp.h.
40  *
41  * \author Teemu Murtola <teemu.murtola@gmail.com>
42  * \ingroup module_selection
43  */
44 #include "gmxpre.h"
45
46 #include "selhelp.h"
47
48 #include <memory>
49 #include <set>
50 #include <string>
51 #include <utility>
52 #include <vector>
53
54 #include "gromacs/onlinehelp/helptopic.h"
55 #include "gromacs/onlinehelp/helpwritercontext.h"
56 #include "gromacs/utility/classhelpers.h"
57 #include "gromacs/utility/exceptions.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/stringutil.h"
60 #include "gromacs/utility/textwriter.h"
61
62 #include "selmethod.h"
63 #include "symrec.h"
64
65 namespace
66 {
67
68 struct CommonHelpText
69 {
70     static const char        name[];
71     static const char        title[];
72     static const char* const text[];
73 };
74
75 const char        CommonHelpText::name[]  = "selections";
76 const char        CommonHelpText::title[] = "Selection syntax and usage";
77 const char* const CommonHelpText::text[]  = {
78     "Selections are used to select atoms/molecules/residues for analysis.",
79     "In contrast to traditional index files, selections can be dynamic, i.e.,",
80     "select different atoms for different trajectory frames. The GROMACS",
81     "manual contains a short introductory section to selections in the",
82     "Analysis chapter, including suggestions on how to get familiar with",
83     "selections if you are new to the concept. The subtopics listed below",
84     "provide more details on the technical and syntactic aspects of",
85     "selections.[PAR]",
86
87     "Each analysis tool requires a different number of selections and the",
88     "selections are interpreted differently. The general idea is still the",
89     "same: each selection evaluates to a set of positions, where a position",
90     "can be an atom position or center-of-mass or center-of-geometry of",
91     "a set of atoms. The tool then uses these positions for its analysis to",
92     "allow very flexible processing. Some analysis tools may have limitations",
93     "on the types of selections allowed."
94 };
95
96 struct ArithmeticHelpText
97 {
98     static const char        name[];
99     static const char        title[];
100     static const char* const text[];
101 };
102
103 const char        ArithmeticHelpText::name[]  = "arithmetic";
104 const char        ArithmeticHelpText::title[] = "Arithmetic expressions in selections";
105 const char* const ArithmeticHelpText::text[]  = {
106     "Basic arithmetic evaluation is supported for numeric expressions.",
107     "Supported operations are addition, subtraction, negation, multiplication,",
108     "division, and exponentiation (using ^).",
109     "Result of a division by zero or other illegal operations is undefined.",
110 };
111
112 struct CmdLineHelpText
113 {
114     static const char        name[];
115     static const char        title[];
116     static const char* const text[];
117 };
118
119 const char        CmdLineHelpText::name[]  = "cmdline";
120 const char        CmdLineHelpText::title[] = "Specifying selections from command line";
121 const char* const CmdLineHelpText::text[]  = {
122     "If no selections are provided on the command line, you are prompted to",
123     "type the selections interactively (a pipe can also be used to provide",
124     "the selections in this case for most tools). While this works well for",
125     "testing, it is easier to provide the selections from the command line",
126     "if they are complex or for scripting.[PAR]",
127
128     "Each tool has different command-line arguments for specifying selections",
129     "(see the help for the individual tools).",
130     "You can either pass a single string containing all selections (separated",
131     "by semicolons), or multiple strings, each containing one selection.",
132     "Note that you need to quote the selections to protect them from the",
133     "shell.[PAR]",
134
135     "If you set a selection command-line argument, but do not provide any",
136     "selections, you are prompted to type the selections for that argument",
137     "interactively. This is useful if that selection argument is optional,",
138     "in which case it is not normally prompted for.[PAR]",
139
140     "To provide selections from a file, use [TT]-sf file.dat[tt] in the place",
141     "of the selection for a selection argument (e.g.,",
142     "[TT]-select -sf file.dat[tt]). In general, the [TT]-sf[tt] argument reads",
143     "selections from the provided file and assigns them to selection arguments",
144     "that have been specified up to that point, but for which no selections",
145     "have been provided.",
146     "As a special case, [TT]-sf[tt] provided on its own, without preceding",
147     "selection arguments, assigns the selections to all (yet unset) required",
148     "selections (i.e., those that would be promted interactively if no",
149     "selections are provided on the command line).[PAR]",
150
151     "To use groups from a traditional index file, use argument [TT]-n[tt]",
152     "to provide a file. See the \"syntax\" subtopic for how to use them.",
153     "If this option is not provided, default groups are generated.",
154     "The default groups are generated with the same logic as for",
155     "non-selection tools.",
156     "",
157     "Depending on the tool, two additional command-line arguments may be",
158     "available to control the behavior:",
159     "",
160     "* [TT]-seltype[tt] can be used to specify the default type of",
161     "  positions to calculate for each selection.",
162     "* [TT]-selrpos[tt] can be used to specify the default type of",
163     "  positions used in selecting atoms by coordinates.",
164     "",
165     "See the \"positions\" subtopic for more information on these options.",
166     "",
167     "Tools that take selections apply them to a structure/topology and/or",
168     "a trajectory file. If the tool takes both (typically as [TT]-s[tt]",
169     "for structure/topology and [TT]-f[tt] for trajectory), then the",
170     "trajectory file is only used for coordinate information, and all other",
171     "information, such as atom names and residue information, is read from",
172     "the structure/topology file. If the tool only takes a structure file,",
173     "or if only that input parameter is provided, then also the coordinates",
174     "are taken from that file.",
175     "For example, to select atoms from a [TT].pdb[tt]/[TT].gro[tt] file in",
176     "a tool that provides both options, pass it as [TT]-s[tt] (only).",
177     "There is no warning if the trajectory file specifies, e.g., different",
178     "atom names than the structure file. Only the number of atoms is checked.",
179     "Many selection-enabled tools also provide an [TT]-fgroup[tt] option",
180     "to specify the atom indices that are present in the trajectory for cases",
181     "where the trajectory only has a subset of atoms from the",
182     "topology/structure file."
183 };
184
185 struct EvaluationHelpText
186 {
187     static const char        name[];
188     static const char        title[];
189     static const char* const text[];
190 };
191
192 const char        EvaluationHelpText::name[]  = "evaluation";
193 const char        EvaluationHelpText::title[] = "Selection evaluation and optimization";
194 const char* const EvaluationHelpText::text[]  = {
195     "Boolean evaluation proceeds from left to right and is short-circuiting",
196     "i.e., as soon as it is known whether an atom will be selected, the",
197     "remaining expressions are not evaluated at all.",
198     "This can be used to optimize the selections: you should write the",
199     "most restrictive and/or the most inexpensive expressions first in",
200     "boolean expressions.",
201     "The relative ordering between dynamic and static expressions does not",
202     "matter: all static expressions are evaluated only once, before the first",
203     "frame, and the result becomes the leftmost expression.[PAR]",
204
205     "Another point for optimization is in common subexpressions: they are not",
206     "automatically recognized, but can be manually optimized by the use of",
207     "variables. This can have a big impact on the performance of complex",
208     "selections, in particular if you define several index groups like this::",
209     "",
210     "  rdist = distance from com of resnr 1 to 5;",
211     "  resname RES and rdist < 2;",
212     "  resname RES and rdist < 4;",
213     "  resname RES and rdist < 6;",
214     "",
215     "Without the variable assignment, the distances would be evaluated three",
216     "times, although they are exactly the same within each selection.",
217     "Anything assigned into a variable becomes a common subexpression that",
218     "is evaluated only once during a frame.",
219     "Currently, in some cases the use of variables can actually lead to a small",
220     "performance loss because of the checks necessary to determine for which",
221     "atoms the expression has already been evaluated, but this should not be",
222     "a major problem.",
223 };
224
225 struct ExamplesHelpText
226 {
227     static const char        name[];
228     static const char        title[];
229     static const char* const text[];
230 };
231
232 const char        ExamplesHelpText::name[]  = "examples";
233 const char        ExamplesHelpText::title[] = "Selection examples";
234 const char* const ExamplesHelpText::text[]  = {
235     "Below, examples of different types of selections are given.",
236     "",
237     "* Selection of all water oxygens::",
238     "",
239     "    resname SOL and name OW",
240     "",
241     "* Centers of mass of residues 1 to 5 and 10::",
242     "",
243     "    res_com of resnr 1 to 5 10",
244     "",
245     "* All atoms farther than 1 nm of a fixed position::",
246     "",
247     "    not within 1 of [1.2, 3.1, 2.4]",
248     "",
249     "* All atoms of a residue LIG within 0.5 nm of a protein (with a custom name)::",
250     "",
251     "    \"Close to protein\" resname LIG and within 0.5 of group \"Protein\"",
252     "",
253     "* All protein residues that have at least one atom within 0.5 nm of a residue LIG::",
254     "",
255     "    group \"Protein\" and same residue as within 0.5 of resname LIG",
256     "",
257     "* All RES residues whose COM is between 2 and 4 nm from the COM of all of them::",
258     "",
259     "    rdist = res_com distance from com of resname RES",
260     "    resname RES and rdist >= 2 and rdist <= 4",
261     "",
262     // TODO: Make it possible to use links below.
263     "* Selection like with duplicate atoms like C1 C2 C2 C3 C3 C4 ... C8 C9::",
264     "",
265     "    name \"C[1-8]\" merge name \"C[2-9]\"",
266     "",
267     "  This can be used with [TT]gmx distance[tt] to compute C1-C2, C2-C3 etc.",
268     "  distances.",
269     "",
270     "* Selection with atoms in order C2 C1::",
271     "",
272     "    name C1 C2 permute 2 1",
273     "",
274     "  This can be used with [TT]gmx gangle[tt] to get C2->C1 vectors instead of",
275     "  C1->C2.",
276     "",
277     "* Selection with COMs of two index groups::",
278     "",
279     "    com of group 1 plus com of group 2",
280     "",
281     "  This can be used with [TT]gmx distance[tt] to compute the distance between",
282     "  these two COMs.",
283     "",
284     "* Fixed vector along x (can be used as a reference with [TT]gmx gangle[tt])::",
285     "",
286     "    [0, 0, 0] plus [1, 0, 0]",
287     "",
288     "* The following examples explain the difference between the various",
289     "  position types.  This selection selects a position for each residue",
290     "  where any of the three atoms C[123] has [TT]x < 2[tt].  The positions",
291     "  are computed as the COM of all three atoms.",
292     "  This is the default behavior if you just write [TT]res_com of[tt]. ::",
293     "",
294     "    part_res_com of name C1 C2 C3 and x < 2",
295     "",
296     "  This selection does the same, but the positions are computed as COM",
297     "  positions of whole residues::",
298     "",
299     "    whole_res_com of name C1 C2 C3 and x < 2",
300     "",
301     "  Finally, this selection selects the same residues, but the positions",
302     "  are computed as COM of exactly those atoms atoms that match the",
303     "  [TT]x < 2[tt] criterion::",
304     "",
305     "    dyn_res_com of name C1 C2 C3 and x < 2",
306     "",
307     "* Without the [TT]of[tt] keyword, the default behavior is different from",
308     "  above, but otherwise the rules are the same::",
309     "",
310     "    name C1 C2 C3 and res_com x < 2",
311     "",
312     "  works as if [TT]whole_res_com[tt] was specified, and selects the three",
313     "  atoms from residues whose COM satisfiex [TT]x < 2[tt].",
314     "  Using ::",
315     "",
316     "    name C1 C2 C3 and part_res_com x < 2",
317     "",
318     "  instead selects residues based on the COM computed from the C[123] atoms.",
319 };
320
321 struct KeywordsHelpText
322 {
323     static const char        name[];
324     static const char        title[];
325     static const char* const text[];
326 };
327
328 const char        KeywordsHelpText::name[]  = "keywords";
329 const char        KeywordsHelpText::title[] = "Selection keywords";
330 const char* const KeywordsHelpText::text[]  = {
331     "The following selection keywords are currently available.",
332     "For keywords marked with a plus, additional help is available through",
333     "a subtopic KEYWORD, where KEYWORD is the name of the keyword.",
334 };
335
336 struct LimitationsHelpText
337 {
338     static const char        name[];
339     static const char        title[];
340     static const char* const text[];
341 };
342
343 const char        LimitationsHelpText::name[]  = "limitations";
344 const char        LimitationsHelpText::title[] = "Selection limitations";
345 const char* const LimitationsHelpText::text[]  = {
346     "* Some analysis programs may require a special structure for the input",
347     "  selections (e.g., some options of [TT]gmx gangle[tt] require the index",
348     "  group to be made of groups of three or four atoms).",
349     "  For such programs, it is up to the user to provide a proper selection",
350     "  expression that always returns such positions.",
351     "",
352     "* All selection keywords select atoms in increasing order, i.e., you can",
353     "  consider them as set operations that in the end return the atoms in",
354     "  sorted numerical order.  For example, the following selections select",
355     "  the same atoms in the same order::",
356     "",
357     "    resname RA RB RC",
358     "    resname RB RC RA",
359     "",
360     "  ::",
361     "",
362     "    atomnr 10 11 12 13",
363     "    atomnr 12 13 10 11",
364     "    atomnr 10 to 13",
365     "    atomnr 13 to 10",
366     "",
367     "  If you need atoms/positions in a different order, you can:",
368     "",
369     "  * use external index groups (for some static selections),",
370     "  * use the [TT]permute[tt] keyword to change the final order, or",
371     "  * use the [TT]merge[tt] or [TT]plus[tt] keywords to compose the",
372     "    final selection from multiple distinct selections.",
373     "",
374     "* Due to technical reasons, having a negative value as the first value in",
375     "  expressions like ::",
376     "",
377     "    charge -1 to -0.7",
378     "",
379     "  result in a syntax error. A workaround is to write ::",
380     "",
381     "    charge {-1 to -0.7}",
382     "",
383     "  instead.",
384     "",
385     "* When [TT]name[tt] selection keyword is used together with PDB input",
386     "  files, the behavior may be unintuitive. When GROMACS reads in a PDB",
387     "  file, 4 character atom names that start with a digit are transformed",
388     "  such that, e.g., 1HG2 becomes HG21, and the latter is what is matched",
389     "  by the [TT]name[tt] keyword. Use [TT]pdbname[tt] to match the atom name",
390     "  as it appears in the input PDB file.",
391 };
392
393 struct PositionsHelpText
394 {
395     static const char        name[];
396     static const char        title[];
397     static const char* const text[];
398 };
399
400 const char        PositionsHelpText::name[]  = "positions";
401 const char        PositionsHelpText::title[] = "Specifying positions in selections";
402 const char* const PositionsHelpText::text[]  = {
403     "Possible ways of specifying positions in selections are:",
404     "",
405     "1. A constant position can be defined as [TT][XX, YY, ZZ][tt], where",
406     "   [TT]XX[tt], [TT]YY[tt] and [TT]ZZ[tt] are real numbers.[PAR]",
407     "",
408     "2. [TT]com of ATOM_EXPR [pbc][tt] or [TT]cog of ATOM_EXPR [pbc][tt]",
409     "   calculate the center of mass/geometry of [TT]ATOM_EXPR[tt]. If",
410     "   [TT]pbc[tt] is specified, the center is calculated iteratively to try",
411     "   to deal with cases where [TT]ATOM_EXPR[tt] wraps around periodic",
412     "   boundary conditions.",
413     "",
414     "3. [TT]POSTYPE of ATOM_EXPR[tt] calculates the specified positions for",
415     "   the atoms in [TT]ATOM_EXPR[tt].",
416     "   [TT]POSTYPE[tt] can be [TT]atom[tt], [TT]res_com[tt], [TT]res_cog[tt],",
417     "   [TT]mol_com[tt] or [TT]mol_cog[tt], with an optional prefix [TT]whole_[tt]",
418     "   [TT]part_[tt] or [TT]dyn_[tt].",
419     "   [TT]whole_[tt] calculates the centers for the whole residue/molecule,",
420     "   even if only part of it is selected.",
421     "   [TT]part_[tt] prefix calculates the centers for the selected atoms, but",
422     "   uses always the same atoms for the same residue/molecule. The used atoms",
423     "   are determined from the largest group allowed by the selection.",
424     "   [TT]dyn_[tt] calculates the centers strictly only for the selected atoms.",
425     "   If no prefix is specified, whole selections default to [TT]part_[tt] and",
426     "   other places default to [TT]whole_[tt].",
427     "   The latter is often desirable to select the same molecules in different",
428     "   tools, while the first is a compromise between speed ([TT]dyn_[tt]",
429     "   positions can be slower to evaluate than [TT]part_[tt]) and intuitive",
430     "   behavior.",
431     "",
432     "4. [TT]ATOM_EXPR[tt], when given for whole selections, is handled as 3.",
433     "   above, using the position type from the command-line argument",
434     "   [TT]-seltype[tt].",
435     "",
436     "Selection keywords that select atoms based on their positions, such as",
437     "[TT]dist from[tt], use by default the positions defined by the",
438     "[TT]-selrpos[tt] command-line option.",
439     "This can be overridden by prepending a [TT]POSTYPE[tt] specifier to the",
440     "keyword. For example, [TT]res_com dist from POS[tt] evaluates the",
441     "residue center of mass distances. In the example, all atoms of a residue",
442     "are either selected or not, based on the single distance calculated.",
443 };
444
445 struct SyntaxHelpText
446 {
447     static const char        name[];
448     static const char        title[];
449     static const char* const text[];
450 };
451
452 const char        SyntaxHelpText::name[]  = "syntax";
453 const char        SyntaxHelpText::title[] = "Selection syntax";
454 const char* const SyntaxHelpText::text[]  = {
455     "A set of selections consists of one or more selections, separated by",
456     "semicolons. Each selection defines a set of positions for the analysis.",
457     "Each selection can also be preceded by a string that gives a name for",
458     "the selection for use in, e.g., graph legends.",
459     "If no name is provided, the string used for the selection is used",
460     "automatically as the name.[PAR]",
461
462     "For interactive input, the syntax is slightly altered: line breaks can",
463     "also be used to separate selections. \\ followed by a line break can",
464     "be used to continue a line if necessary.",
465     "Notice that the above only applies to real interactive input,",
466     "not if you provide the selections, e.g., from a pipe.[PAR]",
467
468     "It is possible to use variables to store selection expressions.",
469     "A variable is defined with the following syntax::",
470     "",
471     "  VARNAME = EXPR ;",
472     "",
473     "where [TT]EXPR[tt] is any valid selection expression.",
474     "After this, [TT]VARNAME[tt] can be used anywhere where [TT]EXPR[tt]",
475     "would be valid.[PAR]",
476
477     "Selections are composed of three main types of expressions, those that",
478     "define atoms ([TT]ATOM_EXPR[tt]), those that define positions",
479     "([TT]POS_EXPR[tt]), and those that evaluate to numeric values",
480     "([TT]NUM_EXPR[tt]). Each selection should be a [TT]POS_EXPR[tt]",
481     "or a [TT]ATOM_EXPR[tt] (the latter is automatically converted to",
482     "positions). The basic rules are as follows:",
483     "",
484     "* An expression like [TT]NUM_EXPR1 < NUM_EXPR2[tt] evaluates to an",
485     "  [TT]ATOM_EXPR[tt] that selects all the atoms for which the comparison",
486     "  is true.",
487     "* Atom expressions can be combined with boolean operations such as",
488     "  [TT]not ATOM_EXPR[tt], [TT]ATOM_EXPR and ATOM_EXPR[tt], or",
489     "  [TT]ATOM_EXPR or ATOM_EXPR[tt]. Parentheses can be used to alter the",
490     "  evaluation order.",
491     "* [TT]ATOM_EXPR[tt] expressions can be converted into [TT]POS_EXPR[tt]",
492     "  expressions in various ways, see the \"positions\" subtopic for more",
493     "  details.",
494     "* [TT]POS_EXPR[tt] can be converted into [TT]NUM_EXPR[tt] using syntax",
495     "  like \"[TT]x of POS_EXPR[tt]\". Currently, this is only supported for single",
496     "  positions like in expression \"[TT]x of cog of ATOM_EXPR[tt]\".",
497     "",
498
499     "Some keywords select atoms based on string values such as the atom name.",
500     "For these keywords, it is possible to use wildcards ([TT]name \"C*\"[tt])",
501     "or regular expressions (e.g., [TT]resname \"R[AB]\"[tt]).",
502     "The match type is automatically guessed from the string: if it contains",
503     "other characters than letters, numbers, '*', or '?', it is interpreted",
504     "as a regular expression.",
505     "To force the matching to use literal string matching, use",
506     "[TT]name = \"C*\"[tt] to match a literal C*.",
507     "To force other type of matching, use '?' or '~' in place of '=' to force",
508     "wildcard or regular expression matching, respectively.[PAR]",
509
510     "Strings that contain non-alphanumeric characters should be enclosed in",
511     "double quotes as in the examples. For other strings, the quotes are",
512     "optional, but if the value conflicts with a reserved keyword, a syntax",
513     "error will occur. If your strings contain uppercase letters, this should",
514     "not happen.[PAR]",
515
516     "Index groups provided with the [TT]-n[tt] command-line option or",
517     "generated by default can be accessed with [TT]group NR[tt] or",
518     "[TT]group NAME[tt], where [TT]NR[tt] is a zero-based index of the group",
519     "and [TT]NAME[tt] is part of the name of the desired group.",
520     "The keyword [TT]group[tt] is optional if the whole selection is",
521     "provided from an index group.",
522     "To see a list of available groups in the interactive mode, press enter",
523     "in the beginning of a line.",
524 };
525
526 } // namespace
527
528 namespace gmx
529 {
530
531 namespace
532 {
533
534 /*! \internal \brief
535  * Help topic implementation for an individual selection method.
536  *
537  * \ingroup module_selection
538  */
539 class KeywordDetailsHelpTopic : public AbstractSimpleHelpTopic
540 {
541 public:
542     //! Initialize help topic for the given selection method.
543     KeywordDetailsHelpTopic(const std::string& name, const gmx_ana_selmethod_t& method) :
544         name_(name),
545         method_(method)
546     {
547     }
548
549     const char* name() const override { return name_.c_str(); }
550     const char* title() const override { return method_.help.helpTitle; }
551
552 protected:
553     std::string helpText() const override
554     {
555         return joinStrings(method_.help.help, method_.help.help + method_.help.nlhelp, "\n");
556     }
557
558 private:
559     std::string                name_;
560     const gmx_ana_selmethod_t& method_;
561
562     GMX_DISALLOW_COPY_AND_ASSIGN(KeywordDetailsHelpTopic);
563 };
564
565 /*! \internal \brief
566  * Custom help topic for printing a list of selection keywords.
567  *
568  * \ingroup module_selection
569  */
570 class KeywordsHelpTopic : public CompositeHelpTopic<KeywordsHelpText>
571 {
572 public:
573     KeywordsHelpTopic();
574
575     void writeHelp(const HelpWriterContext& context) const override;
576
577 private:
578     /*! \brief
579      * Container for known selection methods.
580      *
581      * The first item in the pair is the name of the selection method, and
582      * the second points to the static data structure that describes the
583      * method.
584      * The name in the first item may differ from the name of the static
585      * data structure if an alias is defined for that method.
586      */
587     typedef std::vector<std::pair<std::string, const gmx_ana_selmethod_t*>> MethodList;
588
589     /*! \brief
590      * Prints markup for starting a list of keywords.
591      */
592     static void writeKeywordListStart(const HelpWriterContext& context, const char* heading);
593     /*! \brief
594      * Prints markup for ending a list of keywords.
595      */
596     static void writeKeywordListEnd(const HelpWriterContext& context, const char* extraInfo);
597
598     /*! \brief
599      * Prints a brief list of keywords (selection methods) available.
600      *
601      * \param[in] context  Context for printing the help.
602      * \param[in] type     Only methods that return this type are printed.
603      * \param[in] bModifiers  If false, \ref SMETH_MODIFIER methods are
604      *      excluded, otherwise only them are printed.
605      */
606     void printKeywordList(const HelpWriterContext& context, e_selvalue_t type, bool bModifiers) const;
607
608     /*! \brief
609      * Prints the detailed help for keywords for rst export.
610      */
611     void writeKeywordSubTopics(const HelpWriterContext& context) const;
612
613     MethodList methods_;
614 };
615
616 KeywordsHelpTopic::KeywordsHelpTopic()
617 {
618     // TODO: This is not a very elegant way of getting the list of selection
619     // methods, but this needs to be rewritten in any case if/when #652 is
620     // implemented.
621     const std::unique_ptr<SelectionParserSymbolTable> symtab(new SelectionParserSymbolTable);
622     gmx_ana_selmethod_register_defaults(symtab.get());
623
624     SelectionParserSymbolIterator symbol = symtab->beginIterator(SelectionParserSymbol::MethodSymbol);
625     while (symbol != symtab->endIterator())
626     {
627         const std::string&         symname = symbol->name();
628         const gmx_ana_selmethod_t* method  = symbol->methodValue();
629         methods_.push_back(std::make_pair(std::string(symname), method));
630         if (method->help.nlhelp > 0 && method->help.help != nullptr)
631         {
632             addSubTopic(HelpTopicPointer(new KeywordDetailsHelpTopic(symname, *method)));
633         }
634         ++symbol;
635     }
636 }
637
638 void KeywordsHelpTopic::writeHelp(const HelpWriterContext& context) const
639 {
640     context.writeTextBlock(helpText());
641
642     // Print the list of keywords
643     writeKeywordListStart(context, "Keywords that select atoms by an integer property:");
644     printKeywordList(context, INT_VALUE, false);
645     writeKeywordListEnd(context, "(use in expressions or like \"atomnr 1 to 5 7 9\")");
646
647     writeKeywordListStart(context, "Keywords that select atoms by a numeric property:");
648     printKeywordList(context, REAL_VALUE, false);
649     writeKeywordListEnd(context, "(use in expressions or like \"occupancy 0.5 to 1\")");
650
651     writeKeywordListStart(context, "Keywords that select atoms by a string property:");
652     printKeywordList(context, STR_VALUE, false);
653     writeKeywordListEnd(context, "(use like \"name PATTERN [PATTERN] ...\")");
654
655     writeKeywordListStart(context, "Additional keywords that directly select atoms:");
656     printKeywordList(context, GROUP_VALUE, false);
657     writeKeywordListEnd(context, nullptr);
658
659     writeKeywordListStart(context, "Keywords that directly evaluate to positions:");
660     printKeywordList(context, POS_VALUE, false);
661     writeKeywordListEnd(context, "(see also \"positions\" subtopic)");
662
663     writeKeywordListStart(context, "Additional keywords:");
664     printKeywordList(context, POS_VALUE, true);
665     printKeywordList(context, NO_VALUE, true);
666     writeKeywordListEnd(context, nullptr);
667
668     writeKeywordSubTopics(context);
669 }
670
671 void KeywordsHelpTopic::writeKeywordListStart(const HelpWriterContext& context, const char* heading)
672 {
673     context.paragraphBreak();
674     std::string fullHeading("* ");
675     fullHeading.append(heading);
676     context.writeTextBlock(fullHeading);
677     if (context.outputFormat() == eHelpOutputFormat_Rst)
678     {
679         context.paragraphBreak();
680         context.writeTextBlock("  ::");
681         context.paragraphBreak();
682     }
683 }
684
685 void KeywordsHelpTopic::writeKeywordListEnd(const HelpWriterContext& context, const char* extraInfo)
686 {
687     if (context.outputFormat() == eHelpOutputFormat_Rst)
688     {
689         context.paragraphBreak();
690     }
691     if (!isNullOrEmpty(extraInfo))
692     {
693         std::string fullInfo("  ");
694         fullInfo.append(extraInfo);
695         context.writeTextBlock(fullInfo);
696     }
697 }
698
699 void KeywordsHelpTopic::printKeywordList(const HelpWriterContext& context, e_selvalue_t type, bool bModifiers) const
700 {
701     TextWriter&                file = context.outputFile();
702     MethodList::const_iterator iter;
703     for (iter = methods_.begin(); iter != methods_.end(); ++iter)
704     {
705         const gmx_ana_selmethod_t& method      = *iter->second;
706         const bool                 bIsModifier = (method.flags & SMETH_MODIFIER) != 0;
707         if (method.type == type && bModifiers == bIsModifier)
708         {
709             const bool bHasHelp       = (method.help.nlhelp > 0 && method.help.help != nullptr);
710             const bool bPrintHelpMark = bHasHelp && context.outputFormat() == eHelpOutputFormat_Console;
711             file.writeString(formatString("   %c ", bPrintHelpMark ? '+' : ' '));
712             if (method.help.syntax != nullptr)
713             {
714                 file.writeLine(method.help.syntax);
715             }
716             else
717             {
718                 std::string symname = iter->first;
719                 if (symname != method.name)
720                 {
721                     symname.append(formatString(" (synonym for %s)", method.name));
722                 }
723                 file.writeLine(symname);
724             }
725         }
726     }
727 }
728
729 void KeywordsHelpTopic::writeKeywordSubTopics(const HelpWriterContext& context) const
730 {
731     if (context.outputFormat() != eHelpOutputFormat_Rst)
732     {
733         return;
734     }
735     std::set<std::string>      usedSymbols;
736     MethodList::const_iterator iter;
737     for (iter = methods_.begin(); iter != methods_.end(); ++iter)
738     {
739         const gmx_ana_selmethod_t& method = *iter->second;
740         const bool bHasHelp               = (method.help.nlhelp > 0 && method.help.help != nullptr);
741         if (!bHasHelp || usedSymbols.count(iter->first) > 0)
742         {
743             continue;
744         }
745
746         std::string title;
747         if (method.help.helpTitle != nullptr)
748         {
749             title = method.help.helpTitle;
750             title.append(" - ");
751         }
752         title.append(iter->first);
753         MethodList::const_iterator mergeIter = iter;
754         for (++mergeIter; mergeIter != methods_.end(); ++mergeIter)
755         {
756             if (mergeIter->second->help.help == method.help.help)
757             {
758                 title.append(", ");
759                 title.append(mergeIter->first);
760                 usedSymbols.insert(mergeIter->first);
761             }
762         }
763
764         const IHelpTopic* subTopic = findSubTopic(iter->first.c_str());
765         GMX_RELEASE_ASSERT(subTopic != nullptr, "Keyword subtopic no longer exists");
766         HelpWriterContext subContext(context);
767         subContext.enterSubSection(title);
768         subTopic->writeHelp(subContext);
769     }
770 }
771
772 } // namespace
773
774 //! \cond libapi */
775 HelpTopicPointer createSelectionHelpTopic()
776 {
777     CompositeHelpTopicPointer root(new CompositeHelpTopic<CommonHelpText>);
778     root->registerSubTopic<SimpleHelpTopic<CmdLineHelpText>>();
779     root->registerSubTopic<SimpleHelpTopic<SyntaxHelpText>>();
780     root->registerSubTopic<SimpleHelpTopic<PositionsHelpText>>();
781     root->registerSubTopic<SimpleHelpTopic<ArithmeticHelpText>>();
782     root->registerSubTopic<KeywordsHelpTopic>();
783     root->registerSubTopic<SimpleHelpTopic<EvaluationHelpText>>();
784     root->registerSubTopic<SimpleHelpTopic<LimitationsHelpText>>();
785     root->registerSubTopic<SimpleHelpTopic<ExamplesHelpText>>();
786     return HelpTopicPointer(root.release());
787 }
788 //! \endcond
789
790 } // namespace gmx