Move -debug to cmdlinemodulemanager.cpp
[alexxy/gromacs.git] / src / gromacs / commandline / pargs.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, 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 /* This file is completely threadsafe - keep it that way! */
38 #include "gromacs/commandline/pargs.h"
39
40 #ifdef HAVE_CONFIG_H
41 #include <config.h>
42 #endif
43
44 #include <cstdlib>
45 #include <cstring>
46
47 #include <algorithm>
48 #include <list>
49
50 #ifdef HAVE_UNISTD_H
51 #include <unistd.h>
52 #endif
53
54 #include "thread_mpi/threads.h"
55
56 #include "gromacs/commandline/cmdlinehelpcontext.h"
57 #include "gromacs/commandline/cmdlinehelpwriter.h"
58 #include "gromacs/commandline/cmdlineparser.h"
59 #include "gromacs/fileio/timecontrol.h"
60 #include "gromacs/options/basicoptions.h"
61 #include "gromacs/options/filenameoption.h"
62 #include "gromacs/options/filenameoptionmanager.h"
63 #include "gromacs/options/options.h"
64 #include "gromacs/options/timeunitmanager.h"
65 #include "gromacs/utility/arrayref.h"
66 #include "gromacs/utility/common.h"
67 #include "gromacs/utility/cstringutil.h"
68 #include "gromacs/utility/exceptions.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/gmxassert.h"
71 #include "gromacs/utility/programcontext.h"
72 #include "gromacs/utility/smalloc.h"
73 #include "gromacs/utility/stringutil.h"
74
75 /* The source code in this file should be thread-safe.
76       Please keep it that way. */
77
78 int nenum(const char *const enumc[])
79 {
80     int i;
81
82     i = 1;
83     /* we *can* compare pointers directly here! */
84     while (enumc[i] && enumc[0] != enumc[i])
85     {
86         i++;
87     }
88
89     return i;
90 }
91
92 int opt2parg_int(const char *option, int nparg, t_pargs pa[])
93 {
94     int i;
95
96     for (i = 0; (i < nparg); i++)
97     {
98         if (strcmp(pa[i].option, option) == 0)
99         {
100             return *pa[i].u.i;
101         }
102     }
103
104     gmx_fatal(FARGS, "No integer option %s in pargs", option);
105
106     return 0;
107 }
108
109 gmx_bool opt2parg_bool(const char *option, int nparg, t_pargs pa[])
110 {
111     int i;
112
113     for (i = 0; (i < nparg); i++)
114     {
115         if (strcmp(pa[i].option, option) == 0)
116         {
117             return *pa[i].u.b;
118         }
119     }
120
121     gmx_fatal(FARGS, "No boolean option %s in pargs", option);
122
123     return FALSE;
124 }
125
126 real opt2parg_real(const char *option, int nparg, t_pargs pa[])
127 {
128     int i;
129
130     for (i = 0; (i < nparg); i++)
131     {
132         if (strcmp(pa[i].option, option) == 0)
133         {
134             return *pa[i].u.r;
135         }
136     }
137
138     gmx_fatal(FARGS, "No real option %s in pargs", option);
139
140     return 0.0;
141 }
142
143 const char *opt2parg_str(const char *option, int nparg, t_pargs pa[])
144 {
145     int i;
146
147     for (i = 0; (i < nparg); i++)
148     {
149         if (strcmp(pa[i].option, option) == 0)
150         {
151             return *(pa[i].u.c);
152         }
153     }
154
155     gmx_fatal(FARGS, "No string option %s in pargs", option);
156
157     return NULL;
158 }
159
160 gmx_bool opt2parg_bSet(const char *option, int nparg, t_pargs pa[])
161 {
162     int i;
163
164     for (i = 0; (i < nparg); i++)
165     {
166         if (strcmp(pa[i].option, option) == 0)
167         {
168             return pa[i].bSet;
169         }
170     }
171
172     gmx_fatal(FARGS, "No such option %s in pargs", option);
173
174     return FALSE; /* Too make some compilers happy */
175 }
176
177 const char *opt2parg_enum(const char *option, int nparg, t_pargs pa[])
178 {
179     int i;
180
181     for (i = 0; (i < nparg); i++)
182     {
183         if (strcmp(pa[i].option, option) == 0)
184         {
185             return pa[i].u.c[0];
186         }
187     }
188
189     gmx_fatal(FARGS, "No such option %s in pargs", option);
190
191     return NULL;
192 }
193
194 /********************************************************************
195  * parse_common_args()
196  */
197
198 namespace gmx
199 {
200
201 namespace
202 {
203
204 /*! \brief
205  * Returns the index of the default xvg format.
206  *
207  * \ingroup module_commandline
208  */
209 int getDefaultXvgFormat(gmx::ConstArrayRef<const char *> xvgFormats)
210 {
211     const char *const select = getenv("GMX_VIEW_XVG");
212     if (select != NULL)
213     {
214         ConstArrayRef<const char *>::const_iterator i =
215             std::find(xvgFormats.begin(), xvgFormats.end(), std::string(select));
216         if (i != xvgFormats.end())
217         {
218             return i - xvgFormats.begin();
219         }
220         else
221         {
222             return exvgNONE - 1;
223         }
224     }
225     /* The default is the first option */
226     return 0;
227 }
228
229 /*! \brief
230  * Conversion helper between t_pargs/t_filenm and Options.
231  *
232  * This class holds the necessary mapping between the old C structures and
233  * the new C++ options to allow copying values back after parsing for cases
234  * where the C++ options do not directly provide the type of value required for
235  * the C structures.
236  *
237  * \ingroup module_commandline
238  */
239 class OptionsAdapter
240 {
241     public:
242         /*! \brief
243          * Initializes the adapter to convert from a specified command line.
244          *
245          * The command line is required, because t_pargs wants to return
246          * strings by reference to the original command line.
247          * OptionsAdapter creates a copy of the `argv` array (but not the
248          * strings) to make this possible, even if the parser removes
249          * options it has recognized.
250          */
251         OptionsAdapter(int argc, const char *const argv[])
252             : argv_(argv, argv + argc)
253         {
254         }
255
256         /*! \brief
257          * Converts a t_filenm option into an Options option.
258          *
259          * \param options Options object to add the new option to.
260          * \param fnm     t_filenm option to convert.
261          */
262         void filenmToOptions(Options *options, t_filenm *fnm);
263         /*! \brief
264          * Converts a t_pargs option into an Options option.
265          *
266          * \param     options Options object to add the new option to.
267          * \param     pa      t_pargs option to convert.
268          */
269         void pargsToOptions(Options *options, t_pargs *pa);
270
271         /*! \brief
272          * Copies values back from options to t_pargs/t_filenm.
273          */
274         void copyValues(bool bReadNode);
275
276     private:
277         struct FileNameData
278         {
279             //! Creates a conversion helper for a given `t_filenm` struct.
280             explicit FileNameData(t_filenm *fnm) : fnm(fnm), optionInfo(NULL)
281             {
282             }
283
284             //! t_filenm structure to receive the final values.
285             t_filenm                 *fnm;
286             //! Option info object for the created FileNameOption.
287             FileNameOptionInfo       *optionInfo;
288             //! Value storage for the created FileNameOption.
289             std::vector<std::string>  values;
290         };
291         struct ProgramArgData
292         {
293             //! Creates a conversion helper for a given `t_pargs` struct.
294             explicit ProgramArgData(t_pargs *pa)
295                 : pa(pa), optionInfo(NULL), enumIndex(0), boolValue(false)
296             {
297             }
298
299             //! t_pargs structure to receive the final values.
300             t_pargs                 *pa;
301             //! Option info object for the created option.
302             OptionInfo              *optionInfo;
303             //! Value storage for a non-enum StringOption (unused for other types).
304             std::string              stringValue;
305             //! Value storage for an enum option (unused for other types).
306             int                      enumIndex;
307             //! Value storage for a BooleanOption (unused for other types).
308             bool                     boolValue;
309         };
310
311         std::vector<const char *>    argv_;
312         // These are lists instead of vectors to avoid relocating existing
313         // objects in case the container is reallocated (the Options object
314         // contains pointes to members of the objects, which would get
315         // invalidated).
316         std::list<FileNameData>      fileNameOptions_;
317         std::list<ProgramArgData>    programArgs_;
318
319         GMX_DISALLOW_COPY_AND_ASSIGN(OptionsAdapter);
320 };
321
322 void OptionsAdapter::filenmToOptions(Options *options, t_filenm *fnm)
323 {
324     if (fnm->opt == NULL)
325     {
326         // Existing code may use opt2fn() instead of ftp2fn() for
327         // options that use the default option name, so we need to
328         // keep the old behavior instead of fixing opt2fn().
329         // TODO: Check that this is not the case, remove this, and make
330         // opt2*() work even if fnm->opt is NULL for some options.
331         fnm->opt = ftp2defopt(fnm->ftp);
332     }
333     const bool        bRead     = ((fnm->flag & ffREAD)  != 0);
334     const bool        bWrite    = ((fnm->flag & ffWRITE) != 0);
335     const bool        bOptional = ((fnm->flag & ffOPT)   != 0);
336     const bool        bLibrary  = ((fnm->flag & ffLIB)   != 0);
337     const bool        bMultiple = ((fnm->flag & ffMULT)  != 0);
338     const char *const name      = &fnm->opt[1];
339     const char *      defName   = fnm->fn;
340     if (defName == NULL)
341     {
342         defName = ftp2defnm(fnm->ftp);
343     }
344     fileNameOptions_.push_back(FileNameData(fnm));
345     FileNameData &data = fileNameOptions_.back();
346     data.optionInfo = options->addOption(
347                 FileNameOption(name).storeVector(&data.values)
348                     .defaultBasename(defName).legacyType(fnm->ftp)
349                     .legacyOptionalBehavior()
350                     .readWriteFlags(bRead, bWrite).required(!bOptional)
351                     .libraryFile(bLibrary).multiValue(bMultiple)
352                     .description(ftp2desc(fnm->ftp)));
353 }
354
355 void OptionsAdapter::pargsToOptions(Options *options, t_pargs *pa)
356 {
357     const bool        bHidden = startsWith(pa->desc, "HIDDEN");
358     const char *const name    = &pa->option[1];
359     const char *const desc    = (bHidden ? &pa->desc[6] : pa->desc);
360     programArgs_.push_back(ProgramArgData(pa));
361     ProgramArgData   &data = programArgs_.back();
362     switch (pa->type)
363     {
364         case etINT:
365             data.optionInfo = options->addOption(
366                         IntegerOption(name).store(pa->u.i)
367                             .description(desc).hidden(bHidden));
368             return;
369         case etINT64:
370             data.optionInfo = options->addOption(
371                         Int64Option(name).store(pa->u.is)
372                             .description(desc).hidden(bHidden));
373             return;
374         case etREAL:
375             data.optionInfo = options->addOption(
376                         RealOption(name).store(pa->u.r)
377                             .description(desc).hidden(bHidden));
378             return;
379         case etTIME:
380             data.optionInfo = options->addOption(
381                         RealOption(name).store(pa->u.r).timeValue()
382                             .description(desc).hidden(bHidden));
383             return;
384         case etSTR:
385         {
386             const char *const defValue = (*pa->u.c != NULL ? *pa->u.c : "");
387             data.optionInfo = options->addOption(
388                         StringOption(name).store(&data.stringValue)
389                             .defaultValue(defValue)
390                             .description(desc).hidden(bHidden));
391             return;
392         }
393         case etBOOL:
394             data.optionInfo = options->addOption(
395                         BooleanOption(name).store(&data.boolValue)
396                             .defaultValue(*pa->u.b)
397                             .description(desc).hidden(bHidden));
398             return;
399         case etRVEC:
400             data.optionInfo = options->addOption(
401                         RealOption(name).store(*pa->u.rv).vector()
402                             .description(desc).hidden(bHidden));
403             return;
404         case etENUM:
405         {
406             const int defaultIndex = (pa->u.c[0] != NULL ? nenum(pa->u.c) - 1 : 0);
407             data.optionInfo = options->addOption(
408                         StringOption(name).storeEnumIndex(&data.enumIndex)
409                             .defaultEnumIndex(defaultIndex)
410                             .enumValueFromNullTerminatedArray(pa->u.c + 1)
411                             .description(desc).hidden(bHidden));
412             return;
413         }
414     }
415     GMX_THROW(NotImplementedError("Argument type not implemented"));
416 }
417
418 void OptionsAdapter::copyValues(bool bReadNode)
419 {
420     std::list<FileNameData>::const_iterator file;
421     for (file = fileNameOptions_.begin(); file != fileNameOptions_.end(); ++file)
422     {
423         // FIXME: FF_NOT_READ_NODE should also skip all fexist() calls in
424         // FileNameOption.  However, it is not currently used, and other
425         // commented out code for using it is also outdated, so left this for
426         // later.
427         if (!bReadNode && (file->fnm->flag & ffREAD))
428         {
429             continue;
430         }
431         if (file->optionInfo->isSet())
432         {
433             file->fnm->flag |= ffSET;
434         }
435         file->fnm->nfiles = file->values.size();
436         snew(file->fnm->fns, file->fnm->nfiles);
437         for (int i = 0; i < file->fnm->nfiles; ++i)
438         {
439             // TODO: Check for out-of-memory.
440             file->fnm->fns[i] = strdup(file->values[i].c_str());
441         }
442     }
443     std::list<ProgramArgData>::const_iterator arg;
444     for (arg = programArgs_.begin(); arg != programArgs_.end(); ++arg)
445     {
446         arg->pa->bSet = arg->optionInfo->isSet();
447         switch (arg->pa->type)
448         {
449             case etSTR:
450             {
451                 if (arg->pa->bSet)
452                 {
453                     std::vector<const char *>::const_iterator pos =
454                         std::find(argv_.begin(), argv_.end(), arg->stringValue);
455                     GMX_RELEASE_ASSERT(pos != argv_.end(),
456                                        "String argument got a value not in argv");
457                     *arg->pa->u.c = *pos;
458                 }
459                 break;
460             }
461             case etBOOL:
462                 *arg->pa->u.b = arg->boolValue;
463                 break;
464             case etENUM:
465                 *arg->pa->u.c = arg->pa->u.c[arg->enumIndex + 1];
466                 break;
467             default:
468                 // For other types, there is nothing type-specific to do.
469                 break;
470         }
471     }
472 }
473
474 } // namespace
475
476 } // namespace gmx
477
478 gmx_bool parse_common_args(int *argc, char *argv[], unsigned long Flags,
479                            int nfile, t_filenm fnm[], int npargs, t_pargs *pa,
480                            int ndesc, const char **desc,
481                            int nbugs, const char **bugs,
482                            output_env_t *oenv)
483 {
484     /* This array should match the order of the enum in oenv.h */
485     const char *const xvg_formats[] = { "xmgrace", "xmgr", "none" };
486
487     // Handle the flags argument, which is a bit field
488     // The FF macro returns whether or not the bit is set
489 #define FF(arg) ((Flags & arg) == arg)
490
491     try
492     {
493         int                        nicelevel = 0;
494         double                     tbegin    = 0.0, tend = 0.0, tdelta = 0.0;
495         bool                       bView     = false;
496         int                        xvgFormat = 0;
497         gmx::TimeUnitManager       timeUnitManager;
498         gmx::OptionsAdapter        adapter(*argc, argv);
499         gmx::Options               options(NULL, NULL);
500         gmx::FileNameOptionManager fileOptManager;
501
502         options.addManager(&fileOptManager);
503         options.setDescription(gmx::ConstArrayRef<const char *>(desc, ndesc));
504
505         options.addOption(
506                 gmx::IntegerOption("nice").store(&nicelevel)
507                     .defaultValue(FF(PCA_BE_NICE) ? 19 : 0)
508                     .description("Set the nicelevel"));
509
510         if (FF(PCA_CAN_SET_DEFFNM))
511         {
512             fileOptManager.addDefaultFileNameOption(&options, "deffnm");
513         }
514         if (FF(PCA_CAN_BEGIN))
515         {
516             options.addOption(
517                     gmx::DoubleOption("b").store(&tbegin).timeValue()
518                         .description("First frame (%t) to read from trajectory"));
519         }
520         if (FF(PCA_CAN_END))
521         {
522             options.addOption(
523                     gmx::DoubleOption("e").store(&tend).timeValue()
524                         .description("Last frame (%t) to read from trajectory"));
525         }
526         if (FF(PCA_CAN_DT))
527         {
528             options.addOption(
529                     gmx::DoubleOption("dt").store(&tdelta).timeValue()
530                         .description("Only use frame when t MOD dt = first time (%t)"));
531         }
532         if (FF(PCA_TIME_UNIT))
533         {
534             timeUnitManager.setTimeUnitFromEnvironment();
535             timeUnitManager.addTimeUnitOption(&options, "tu");
536         }
537         if (FF(PCA_CAN_VIEW))
538         {
539             options.addOption(
540                     gmx::BooleanOption("w").store(&bView)
541                         .description("View output [TT].xvg[tt], [TT].xpm[tt], "
542                                      "[TT].eps[tt] and [TT].pdb[tt] files"));
543         }
544
545         bool bXvgr = false;
546         for (int i = 0; i < nfile; i++)
547         {
548             bXvgr = bXvgr || (fnm[i].ftp == efXVG);
549         }
550         xvgFormat = gmx::getDefaultXvgFormat(xvg_formats);
551         if (bXvgr)
552         {
553             options.addOption(
554                     gmx::StringOption("xvg").enumValue(xvg_formats)
555                         .storeEnumIndex(&xvgFormat)
556                         .description("xvg plot formatting"));
557         }
558
559         /* Now append the program specific arguments */
560         for (int i = 0; i < nfile; i++)
561         {
562             adapter.filenmToOptions(&options, &fnm[i]);
563         }
564         for (int i = 0; i < npargs; i++)
565         {
566             adapter.pargsToOptions(&options, &pa[i]);
567         }
568
569         const gmx::CommandLineHelpContext *context =
570             gmx::GlobalCommandLineHelpContext::get();
571         if (context != NULL)
572         {
573             if (!(FF(PCA_QUIET)))
574             {
575                 gmx::CommandLineHelpWriter(options)
576                     .setShowDescriptions(true)
577                     .setTimeUnitString(timeUnitManager.timeUnitAsString())
578                     .setKnownIssues(gmx::ConstArrayRef<const char *>(bugs, nbugs))
579                     .writeHelp(*context);
580             }
581             return FALSE;
582         }
583
584         /* Now parse all the command-line options */
585         gmx::CommandLineParser(&options).skipUnknown(FF(PCA_NOEXIT_ON_ARGS))
586             .parse(argc, argv);
587         options.finish();
588
589         /* set program name, command line, and default values for output options */
590         output_env_init(oenv, gmx::getProgramContext(),
591                         (time_unit_t)(timeUnitManager.timeUnit() + 1), bView,
592                         (xvg_format_t)(xvgFormat + 1), 0);
593
594         /* Set the nice level */
595 #ifdef HAVE_UNISTD_H
596 #ifndef GMX_NO_NICE
597         /* The some system, e.g. the catamount kernel on cray xt3 do not have nice(2). */
598         if (nicelevel != 0)
599         {
600             static gmx_bool            nice_set   = FALSE; /* only set it once */
601             static tMPI_Thread_mutex_t init_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
602             tMPI_Thread_mutex_lock(&init_mutex);
603             if (!nice_set)
604             {
605                 if (nice(nicelevel) == -1)
606                 {
607                     /* Do nothing, but use the return value to avoid warnings. */
608                 }
609                 nice_set = TRUE;
610             }
611             tMPI_Thread_mutex_unlock(&init_mutex);
612         }
613 #endif
614 #endif
615
616         timeUnitManager.scaleTimeOptions(&options);
617
618         /* Extract Time info from arguments */
619         // TODO: Use OptionInfo objects instead of string constants
620         if (FF(PCA_CAN_BEGIN) && options.isSet("b"))
621         {
622             setTimeValue(TBEGIN, tbegin);
623         }
624         if (FF(PCA_CAN_END) && options.isSet("e"))
625         {
626             setTimeValue(TEND, tend);
627         }
628         if (FF(PCA_CAN_DT) && options.isSet("dt"))
629         {
630             setTimeValue(TDELTA, tdelta);
631         }
632
633         adapter.copyValues(!FF(PCA_NOT_READ_NODE));
634
635         return TRUE;
636     }
637     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
638 #undef FF
639 }