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