Improve interactive selection input
[alexxy/gromacs.git] / src / gromacs / selection / selectioncollection.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014, 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 gmx::SelectionCollection.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_selection
41  */
42 #include "selectioncollection.h"
43
44 #include <cstdio>
45
46 #include <boost/shared_ptr.hpp>
47
48 #include "gromacs/legacyheaders/oenv.h"
49 #include "gromacs/legacyheaders/smalloc.h"
50 #include "gromacs/legacyheaders/xvgr.h"
51
52 #include "gromacs/options/basicoptions.h"
53 #include "gromacs/options/options.h"
54 #include "gromacs/selection/selection.h"
55 #include "gromacs/utility/exceptions.h"
56 #include "gromacs/utility/file.h"
57 #include "gromacs/utility/gmxassert.h"
58 #include "gromacs/utility/messagestringcollector.h"
59 #include "gromacs/utility/stringutil.h"
60
61 #include "compiler.h"
62 #include "mempool.h"
63 #include "parser.h"
64 #include "poscalc.h"
65 #include "scanner.h"
66 #include "selection.h"
67 #include "selectioncollection-impl.h"
68 #include "selelem.h"
69 #include "selhelp.h"
70 #include "selmethod.h"
71 #include "symrec.h"
72
73 namespace gmx
74 {
75
76 /********************************************************************
77  * SelectionCollection::Impl
78  */
79
80 SelectionCollection::Impl::Impl()
81     : debugLevel_(0), bExternalGroupsSet_(false), grps_(NULL)
82 {
83     sc_.nvars     = 0;
84     sc_.varstrs   = NULL;
85     sc_.top       = NULL;
86     gmx_ana_index_clear(&sc_.gall);
87     sc_.mempool   = NULL;
88     sc_.symtab.reset(new SelectionParserSymbolTable);
89     gmx_ana_selmethod_register_defaults(sc_.symtab.get());
90 }
91
92
93 SelectionCollection::Impl::~Impl()
94 {
95     clearSymbolTable();
96     // The tree must be freed before the SelectionData objects, since the
97     // tree may hold references to the position data in SelectionData.
98     sc_.root.reset();
99     sc_.sel.clear();
100     for (int i = 0; i < sc_.nvars; ++i)
101     {
102         sfree(sc_.varstrs[i]);
103     }
104     sfree(sc_.varstrs);
105     gmx_ana_index_deinit(&sc_.gall);
106     if (sc_.mempool)
107     {
108         _gmx_sel_mempool_destroy(sc_.mempool);
109     }
110 }
111
112
113 void
114 SelectionCollection::Impl::clearSymbolTable()
115 {
116     sc_.symtab.reset();
117 }
118
119
120 namespace
121 {
122
123 /*! \brief
124  * Reads a single selection line from stdin.
125  *
126  * \param[in]  infile        File to read from (typically File::standardInput()).
127  * \param[in]  bInteractive  Whether to print interactive prompts.
128  * \param[out] line          The read line in stored here.
129  * \returns true if something was read, false if at end of input.
130  *
131  * Handles line continuation, reading also the continuing line(s) in one call.
132  */
133 bool promptLine(File *infile, bool bInteractive, std::string *line)
134 {
135     if (bInteractive)
136     {
137         fprintf(stderr, "> ");
138     }
139     if (!infile->readLineWithTrailingSpace(line))
140     {
141         return false;
142     }
143     while (endsWith(*line, "\\\n"))
144     {
145         line->resize(line->length() - 2);
146         if (bInteractive)
147         {
148             fprintf(stderr, "... ");
149         }
150         std::string buffer;
151         // Return value ignored, buffer remains empty and works correctly
152         // if there is nothing to read.
153         infile->readLineWithTrailingSpace(&buffer);
154         line->append(buffer);
155     }
156     if (endsWith(*line, "\n"))
157     {
158         line->resize(line->length() - 1);
159     }
160     else if (bInteractive)
161     {
162         fprintf(stderr, "\n");
163     }
164     return true;
165 }
166
167 /*! \brief
168  * Helper function for tokenizing the input and pushing them to the parser.
169  *
170  * \param     scanner       Tokenizer data structure.
171  * \param     parserState   Parser data structure.
172  * \param[in] bInteractive  Whether to operate in interactive mode.
173  *
174  * Repeatedly reads tokens using \p scanner and pushes them to the parser with
175  * \p parserState until there is no more input, or until enough input is given
176  * (only in interactive mode).
177  */
178 int runParserLoop(yyscan_t scanner, _gmx_sel_yypstate *parserState,
179                   bool bInteractive)
180 {
181     int status    = YYPUSH_MORE;
182     int prevToken = 0;
183     do
184     {
185         YYSTYPE value;
186         int     token = _gmx_sel_yylex(&value, scanner);
187         if (bInteractive)
188         {
189             if (token == 0)
190             {
191                 break;
192             }
193             // Empty commands cause the interactive parser to print out
194             // status information. This avoids producing those unnecessarily,
195             // e.g., from "resname RA;;".
196             if (prevToken == CMD_SEP && token == CMD_SEP)
197             {
198                 continue;
199             }
200             prevToken = token;
201         }
202         status = _gmx_sel_yypush_parse(parserState, token, &value, scanner);
203     }
204     while (status == YYPUSH_MORE);
205     _gmx_sel_lexer_rethrow_exception_if_occurred(scanner);
206     return status;
207 }
208
209 /*! \brief
210  * Print current status in response to empty line in interactive input.
211  *
212  * \param[in] sc             Selection collection data structure.
213  * \param[in] grps           Available index groups.
214  * \param[in] firstSelection Index of first selection from this interactive
215  *     session.
216  * \param[in] maxCount       Maximum number of selections.
217  * \param[in] context        Context to print for what the selections are for.
218  * \param[in] bFirst         Whether this is the header that is printed before
219  *     any user input.
220  *
221  * Prints the available index groups and currently provided selections.
222  */
223 void printCurrentStatus(gmx_ana_selcollection_t *sc, gmx_ana_indexgrps_t *grps,
224                         size_t firstSelection, int maxCount,
225                         const std::string &context, bool bFirst)
226 {
227     if (grps != NULL)
228     {
229         std::fprintf(stderr, "Available static index groups:\n");
230         gmx_ana_indexgrps_print(stderr, grps, 0);
231     }
232     std::fprintf(stderr, "Specify ");
233     if (maxCount < 0)
234     {
235         std::fprintf(stderr, "any number of selections");
236     }
237     else if (maxCount == 1)
238     {
239         std::fprintf(stderr, "a selection");
240     }
241     else
242     {
243         std::fprintf(stderr, "%d selections", maxCount);
244     }
245     std::fprintf(stderr, "%s%s:\n",
246                  context.empty() ? "" : " ", context.c_str());
247     std::fprintf(stderr,
248                  "(one per line, <enter> for status/groups, 'help' for help%s)\n",
249                  maxCount < 0 ? ", Ctrl-D to end" : "");
250     if (!bFirst && (sc->nvars > 0 || sc->sel.size() > firstSelection))
251     {
252         std::fprintf(stderr, "Currently provided selections:\n");
253         for (int i = 0; i < sc->nvars; ++i)
254         {
255             std::fprintf(stderr, "     %s\n", sc->varstrs[i]);
256         }
257         for (size_t i = firstSelection; i < sc->sel.size(); ++i)
258         {
259             std::fprintf(stderr, " %2d. %s\n",
260                          static_cast<int>(i - firstSelection + 1),
261                          sc->sel[i]->selectionText());
262         }
263         if (maxCount > 0)
264         {
265             const int remaining
266                 = maxCount - static_cast<int>(sc->sel.size() - firstSelection);
267             std::fprintf(stderr, "(%d more selection%s required)\n",
268                          remaining, remaining > 1 ? "s" : "");
269         }
270     }
271 }
272
273 /*! \brief
274  * Helper function that runs the parser once the tokenizer has been
275  * initialized.
276  *
277  * \param[in,out] scanner Scanner data structure.
278  * \param[in]     bStdIn  Whether to use a line-based reading
279  *      algorithm designed for interactive input.
280  * \param[in]     maxnr   Maximum number of selections to parse
281  *      (if -1, parse as many as provided by the user).
282  * \param[in]     context Context to print for what the selections are for.
283  * \returns       Vector of parsed selections.
284  * \throws        std::bad_alloc if out of memory.
285  * \throws        InvalidInputError if there is a parsing error.
286  *
287  * Used internally to implement parseFromStdin(), parseFromFile() and
288  * parseFromString().
289  */
290 SelectionList runParser(yyscan_t scanner, bool bStdIn, int maxnr,
291                         const std::string &context)
292 {
293     boost::shared_ptr<void>  scannerGuard(scanner, &_gmx_sel_free_lexer);
294     gmx_ana_selcollection_t *sc   = _gmx_sel_lexer_selcollection(scanner);
295     gmx_ana_indexgrps_t     *grps = _gmx_sel_lexer_indexgrps(scanner);
296
297     MessageStringCollector   errors;
298     _gmx_sel_set_lexer_error_reporter(scanner, &errors);
299
300     size_t oldCount = sc->sel.size();
301     bool   bOk      = false;
302     {
303         boost::shared_ptr<_gmx_sel_yypstate> parserState(
304                 _gmx_sel_yypstate_new(), &_gmx_sel_yypstate_delete);
305         if (bStdIn)
306         {
307             File       &stdinFile(File::standardInput());
308             const bool  bInteractive = _gmx_sel_is_lexer_interactive(scanner);
309             if (bInteractive)
310             {
311                 printCurrentStatus(sc, grps, oldCount, maxnr, context, true);
312             }
313             std::string line;
314             int         status;
315             while (promptLine(&stdinFile, bInteractive, &line))
316             {
317                 if (bInteractive)
318                 {
319                     if (line.empty())
320                     {
321                         printCurrentStatus(sc, grps, oldCount, maxnr, context, false);
322                         continue;
323                     }
324                 }
325                 line.append("\n");
326                 _gmx_sel_set_lex_input_str(scanner, line.c_str());
327                 status = runParserLoop(scanner, parserState.get(), true);
328                 if (status != YYPUSH_MORE)
329                 {
330                     // TODO: Check if there is more input, and issue an
331                     // error/warning if some input was ignored.
332                     goto early_termination;
333                 }
334                 if (!errors.isEmpty() && bInteractive)
335                 {
336                     fprintf(stderr, "%s", errors.toString().c_str());
337                     errors.clear();
338                 }
339             }
340             status = _gmx_sel_yypush_parse(parserState.get(), 0, NULL,
341                                            scanner);
342             _gmx_sel_lexer_rethrow_exception_if_occurred(scanner);
343 early_termination:
344             bOk = (status == 0);
345         }
346         else
347         {
348             int status = runParserLoop(scanner, parserState.get(), false);
349             bOk = (status == 0);
350         }
351     }
352     scannerGuard.reset();
353     int nr = sc->sel.size() - oldCount;
354     if (maxnr > 0 && nr != maxnr)
355     {
356         bOk = false;
357         errors.append("Too few selections provided");
358     }
359
360     // TODO: Remove added selections from the collection if parsing failed?
361     if (!bOk || !errors.isEmpty())
362     {
363         GMX_ASSERT(!bOk && !errors.isEmpty(), "Inconsistent error reporting");
364         GMX_THROW(InvalidInputError(errors.toString()));
365     }
366
367     SelectionList                     result;
368     SelectionDataList::const_iterator i;
369     result.reserve(nr);
370     for (i = sc->sel.begin() + oldCount; i != sc->sel.end(); ++i)
371     {
372         result.push_back(Selection(i->get()));
373     }
374     return result;
375 }
376
377 }   // namespace
378
379
380 void SelectionCollection::Impl::resolveExternalGroups(
381         const SelectionTreeElementPointer &root,
382         ExceptionInitializer              *errors)
383 {
384
385     if (root->type == SEL_GROUPREF)
386     {
387         try
388         {
389             root->resolveIndexGroupReference(grps_);
390         }
391         catch (const UserInputError &)
392         {
393             errors->addCurrentExceptionAsNested();
394         }
395     }
396
397     SelectionTreeElementPointer child = root->child;
398     while (child)
399     {
400         resolveExternalGroups(child, errors);
401         child = child->next;
402     }
403 }
404
405
406 /********************************************************************
407  * SelectionCollection
408  */
409
410 SelectionCollection::SelectionCollection()
411     : impl_(new Impl)
412 {
413 }
414
415
416 SelectionCollection::~SelectionCollection()
417 {
418 }
419
420
421 void
422 SelectionCollection::initOptions(Options *options)
423 {
424     const char * const debug_levels[]
425         = { "no", "basic", "compile", "eval", "full" };
426
427     bool bAllowNonAtomOutput = false;
428     SelectionDataList::const_iterator iter;
429     for (iter = impl_->sc_.sel.begin(); iter != impl_->sc_.sel.end(); ++iter)
430     {
431         const internal::SelectionData &sel = **iter;
432         if (!sel.hasFlag(efSelection_OnlyAtoms))
433         {
434             bAllowNonAtomOutput = true;
435         }
436     }
437
438     const char *const *postypes = PositionCalculationCollection::typeEnumValues;
439     options->addOption(StringOption("selrpos")
440                            .enumValueFromNullTerminatedArray(postypes)
441                            .store(&impl_->rpost_).defaultValue(postypes[0])
442                            .description("Selection reference positions"));
443     if (bAllowNonAtomOutput)
444     {
445         options->addOption(StringOption("seltype")
446                                .enumValueFromNullTerminatedArray(postypes)
447                                .store(&impl_->spost_).defaultValue(postypes[0])
448                                .description("Default selection output positions"));
449     }
450     else
451     {
452         impl_->spost_ = postypes[0];
453     }
454     GMX_RELEASE_ASSERT(impl_->debugLevel_ >= 0 && impl_->debugLevel_ <= 4,
455                        "Debug level out of range");
456     options->addOption(StringOption("seldebug").hidden(impl_->debugLevel_ == 0)
457                            .enumValue(debug_levels)
458                            .defaultValue(debug_levels[impl_->debugLevel_])
459                            .storeEnumIndex(&impl_->debugLevel_)
460                            .description("Print out selection trees for debugging"));
461 }
462
463
464 void
465 SelectionCollection::setReferencePosType(const char *type)
466 {
467     GMX_RELEASE_ASSERT(type != NULL, "Cannot assign NULL position type");
468     // Check that the type is valid, throw if it is not.
469     e_poscalc_t  dummytype;
470     int          dummyflags;
471     PositionCalculationCollection::typeFromEnum(type, &dummytype, &dummyflags);
472     impl_->rpost_ = type;
473 }
474
475
476 void
477 SelectionCollection::setOutputPosType(const char *type)
478 {
479     GMX_RELEASE_ASSERT(type != NULL, "Cannot assign NULL position type");
480     // Check that the type is valid, throw if it is not.
481     e_poscalc_t  dummytype;
482     int          dummyflags;
483     PositionCalculationCollection::typeFromEnum(type, &dummytype, &dummyflags);
484     impl_->spost_ = type;
485 }
486
487
488 void
489 SelectionCollection::setDebugLevel(int debugLevel)
490 {
491     impl_->debugLevel_ = debugLevel;
492 }
493
494
495 void
496 SelectionCollection::setTopology(t_topology *top, int natoms)
497 {
498     GMX_RELEASE_ASSERT(natoms > 0 || top != NULL,
499                        "The number of atoms must be given if there is no topology");
500     // Get the number of atoms from the topology if it is not given.
501     if (natoms <= 0)
502     {
503         natoms = top->atoms.nr;
504     }
505     gmx_ana_selcollection_t *sc = &impl_->sc_;
506     // Do this first, as it allocates memory, while the others don't throw.
507     gmx_ana_index_init_simple(&sc->gall, natoms);
508     sc->pcc.setTopology(top);
509     sc->top = top;
510 }
511
512
513 void
514 SelectionCollection::setIndexGroups(gmx_ana_indexgrps_t *grps)
515 {
516     GMX_RELEASE_ASSERT(grps == NULL || !impl_->bExternalGroupsSet_,
517                        "Can only set external groups once or clear them afterwards");
518     impl_->grps_               = grps;
519     impl_->bExternalGroupsSet_ = true;
520
521     ExceptionInitializer        errors("Unknown index group references encountered");
522     SelectionTreeElementPointer root = impl_->sc_.root;
523     while (root)
524     {
525         impl_->resolveExternalGroups(root, &errors);
526         root = root->next;
527     }
528     if (errors.hasNestedExceptions())
529     {
530         GMX_THROW(InconsistentInputError(errors));
531     }
532     for (size_t i = 0; i < impl_->sc_.sel.size(); ++i)
533     {
534         impl_->sc_.sel[i]->refreshName();
535     }
536 }
537
538
539 bool
540 SelectionCollection::requiresTopology() const
541 {
542     e_poscalc_t  type;
543     int          flags;
544
545     if (!impl_->rpost_.empty())
546     {
547         flags = 0;
548         // Should not throw, because has been checked earlier.
549         PositionCalculationCollection::typeFromEnum(impl_->rpost_.c_str(),
550                                                     &type, &flags);
551         if (type != POS_ATOM)
552         {
553             return true;
554         }
555     }
556     if (!impl_->spost_.empty())
557     {
558         flags = 0;
559         // Should not throw, because has been checked earlier.
560         PositionCalculationCollection::typeFromEnum(impl_->spost_.c_str(),
561                                                     &type, &flags);
562         if (type != POS_ATOM)
563         {
564             return true;
565         }
566     }
567
568     SelectionTreeElementPointer sel = impl_->sc_.root;
569     while (sel)
570     {
571         if (_gmx_selelem_requires_top(*sel))
572         {
573             return true;
574         }
575         sel = sel->next;
576     }
577     return false;
578 }
579
580
581 SelectionList
582 SelectionCollection::parseFromStdin(int nr, bool bInteractive,
583                                     const std::string &context)
584 {
585     yyscan_t scanner;
586
587     _gmx_sel_init_lexer(&scanner, &impl_->sc_, bInteractive, nr,
588                         impl_->bExternalGroupsSet_,
589                         impl_->grps_);
590     return runParser(scanner, true, nr, context);
591 }
592
593
594 SelectionList
595 SelectionCollection::parseFromFile(const std::string &filename)
596 {
597
598     try
599     {
600         yyscan_t scanner;
601         File     file(filename, "r");
602         // TODO: Exception-safe way of using the lexer.
603         _gmx_sel_init_lexer(&scanner, &impl_->sc_, false, -1,
604                             impl_->bExternalGroupsSet_,
605                             impl_->grps_);
606         _gmx_sel_set_lex_input_file(scanner, file.handle());
607         return runParser(scanner, false, -1, std::string());
608     }
609     catch (GromacsException &ex)
610     {
611         ex.prependContext(formatString(
612                                   "Error in parsing selections from file '%s'",
613                                   filename.c_str()));
614         throw;
615     }
616 }
617
618
619 SelectionList
620 SelectionCollection::parseFromString(const std::string &str)
621 {
622     yyscan_t scanner;
623
624     _gmx_sel_init_lexer(&scanner, &impl_->sc_, false, -1,
625                         impl_->bExternalGroupsSet_,
626                         impl_->grps_);
627     _gmx_sel_set_lex_input_str(scanner, str.c_str());
628     return runParser(scanner, false, -1, std::string());
629 }
630
631
632 void
633 SelectionCollection::compile()
634 {
635     if (impl_->sc_.top == NULL && requiresTopology())
636     {
637         GMX_THROW(InconsistentInputError("Selection requires topology information, but none provided"));
638     }
639     if (!impl_->bExternalGroupsSet_)
640     {
641         setIndexGroups(NULL);
642     }
643     if (impl_->debugLevel_ >= 1)
644     {
645         printTree(stderr, false);
646     }
647
648     SelectionCompiler compiler;
649     compiler.compile(this);
650
651     if (impl_->debugLevel_ >= 1)
652     {
653         std::fprintf(stderr, "\n");
654         printTree(stderr, false);
655         std::fprintf(stderr, "\n");
656         impl_->sc_.pcc.printTree(stderr);
657         std::fprintf(stderr, "\n");
658     }
659     impl_->sc_.pcc.initEvaluation();
660     if (impl_->debugLevel_ >= 1)
661     {
662         impl_->sc_.pcc.printTree(stderr);
663         std::fprintf(stderr, "\n");
664     }
665
666     // TODO: It would be nicer to associate the name of the selection option
667     // (if available) to the error message.
668     SelectionDataList::const_iterator iter;
669     for (iter = impl_->sc_.sel.begin(); iter != impl_->sc_.sel.end(); ++iter)
670     {
671         const internal::SelectionData &sel = **iter;
672         if (sel.hasFlag(efSelection_OnlyAtoms))
673         {
674             if (sel.type() != INDEX_ATOM)
675             {
676                 std::string message = formatString(
677                             "Selection '%s' does not evaluate to individual atoms. "
678                             "This is not allowed in this context.",
679                             sel.selectionText());
680                 GMX_THROW(InvalidInputError(message));
681             }
682         }
683         if (sel.hasFlag(efSelection_DisallowEmpty))
684         {
685             if (sel.posCount() == 0)
686             {
687                 std::string message = formatString(
688                             "Selection '%s' never matches any atoms.",
689                             sel.selectionText());
690                 GMX_THROW(InvalidInputError(message));
691             }
692         }
693     }
694 }
695
696
697 void
698 SelectionCollection::evaluate(t_trxframe *fr, t_pbc *pbc)
699 {
700     impl_->sc_.pcc.initFrame();
701
702     SelectionEvaluator evaluator;
703     evaluator.evaluate(this, fr, pbc);
704
705     if (impl_->debugLevel_ >= 3)
706     {
707         std::fprintf(stderr, "\n");
708         printTree(stderr, true);
709     }
710 }
711
712
713 void
714 SelectionCollection::evaluateFinal(int nframes)
715 {
716     SelectionEvaluator evaluator;
717     evaluator.evaluateFinal(this, nframes);
718 }
719
720
721 void
722 SelectionCollection::printTree(FILE *fp, bool bValues) const
723 {
724     SelectionTreeElementPointer sel = impl_->sc_.root;
725     while (sel)
726     {
727         _gmx_selelem_print_tree(fp, *sel, bValues, 0);
728         sel = sel->next;
729     }
730 }
731
732
733 void
734 SelectionCollection::printXvgrInfo(FILE *out, output_env_t oenv) const
735 {
736     if (output_env_get_xvg_format(oenv) != exvgNONE)
737     {
738         const gmx_ana_selcollection_t &sc = impl_->sc_;
739         std::fprintf(out, "# Selections:\n");
740         for (int i = 0; i < sc.nvars; ++i)
741         {
742             std::fprintf(out, "#   %s\n", sc.varstrs[i]);
743         }
744         for (size_t i = 0; i < sc.sel.size(); ++i)
745         {
746             std::fprintf(out, "#   %s\n", sc.sel[i]->selectionText());
747         }
748         std::fprintf(out, "#\n");
749     }
750 }
751
752 // static
753 HelpTopicPointer
754 SelectionCollection::createDefaultHelpTopic()
755 {
756     return createSelectionHelpTopic();
757 }
758
759 } // namespace gmx