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