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