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