Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / commandline / pargs.cpp
blob0804bbb7badad0389c071001e4bfedf02b4847f8
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include "pargs.h"
41 #include <cstdlib>
42 #include <cstring>
44 #include <algorithm>
45 #include <list>
47 #include "gromacs/commandline/cmdlinehelpcontext.h"
48 #include "gromacs/commandline/cmdlinehelpwriter.h"
49 #include "gromacs/commandline/cmdlineparser.h"
50 #include "gromacs/fileio/oenv.h"
51 #include "gromacs/fileio/timecontrol.h"
52 #include "gromacs/options/basicoptions.h"
53 #include "gromacs/options/behaviorcollection.h"
54 #include "gromacs/options/filenameoption.h"
55 #include "gromacs/options/filenameoptionmanager.h"
56 #include "gromacs/options/options.h"
57 #include "gromacs/options/timeunitmanager.h"
58 #include "gromacs/utility/arrayref.h"
59 #include "gromacs/utility/basenetwork.h"
60 #include "gromacs/utility/classhelpers.h"
61 #include "gromacs/utility/exceptions.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/gmxassert.h"
64 #include "gromacs/utility/path.h"
65 #include "gromacs/utility/programcontext.h"
66 #include "gromacs/utility/stringutil.h"
68 /* The source code in this file should be thread-safe.
69 Please keep it that way. */
71 int nenum(const char *const enumc[])
73 int i;
75 i = 1;
76 /* we *can* compare pointers directly here! */
77 while (enumc[i] && enumc[0] != enumc[i])
79 i++;
82 return i;
85 int opt2parg_int(const char *option, int nparg, t_pargs pa[])
87 int i;
89 for (i = 0; (i < nparg); i++)
91 if (strcmp(pa[i].option, option) == 0)
93 return *pa[i].u.i;
97 gmx_fatal(FARGS, "No integer option %s in pargs", option);
100 gmx_bool opt2parg_bool(const char *option, int nparg, t_pargs pa[])
102 int i;
104 for (i = 0; (i < nparg); i++)
106 if (strcmp(pa[i].option, option) == 0)
108 return *pa[i].u.b;
112 gmx_fatal(FARGS, "No boolean option %s in pargs", option);
114 return FALSE;
117 real opt2parg_real(const char *option, int nparg, t_pargs pa[])
119 int i;
121 for (i = 0; (i < nparg); i++)
123 if (strcmp(pa[i].option, option) == 0)
125 return *pa[i].u.r;
129 gmx_fatal(FARGS, "No real option %s in pargs", option);
132 const char *opt2parg_str(const char *option, int nparg, t_pargs pa[])
134 int i;
136 for (i = 0; (i < nparg); i++)
138 if (strcmp(pa[i].option, option) == 0)
140 return *(pa[i].u.c);
144 gmx_fatal(FARGS, "No string option %s in pargs", option);
147 gmx_bool opt2parg_bSet(const char *option, int nparg, const t_pargs *pa)
149 int i;
151 for (i = 0; (i < nparg); i++)
153 if (strcmp(pa[i].option, option) == 0)
155 return pa[i].bSet;
159 gmx_fatal(FARGS, "No such option %s in pargs", option);
161 return FALSE; /* Too make some compilers happy */
164 const char *opt2parg_enum(const char *option, int nparg, t_pargs pa[])
166 int i;
168 for (i = 0; (i < nparg); i++)
170 if (strcmp(pa[i].option, option) == 0)
172 return pa[i].u.c[0];
176 gmx_fatal(FARGS, "No such option %s in pargs", option);
179 /********************************************************************
180 * parse_common_args()
183 namespace gmx
186 namespace
189 /*! \brief
190 * Returns the index of the default xvg format.
192 * \ingroup module_commandline
194 int getDefaultXvgFormat(gmx::ArrayRef<const char *const> xvgFormats)
196 const char *const select = getenv("GMX_VIEW_XVG");
197 if (select != nullptr)
199 ArrayRef<const char *const>::const_iterator i =
200 std::find(xvgFormats.begin(), xvgFormats.end(), std::string(select));
201 if (i != xvgFormats.end())
203 return std::distance(xvgFormats.begin(), i);
205 else
207 return exvgNONE - 1;
210 /* The default is the first option */
211 return 0;
214 /*! \brief
215 * Conversion helper between t_pargs/t_filenm and Options.
217 * This class holds the necessary mapping between the old C structures and
218 * the new C++ options to allow copying values back after parsing for cases
219 * where the C++ options do not directly provide the type of value required for
220 * the C structures.
222 * \ingroup module_commandline
224 class OptionsAdapter
226 public:
227 /*! \brief
228 * Initializes the adapter to convert from a specified command line.
230 * The command line is required, because t_pargs wants to return
231 * strings by reference to the original command line.
232 * OptionsAdapter creates a copy of the `argv` array (but not the
233 * strings) to make this possible, even if the parser removes
234 * options it has recognized.
236 OptionsAdapter(int argc, const char *const argv[])
237 : argv_(argv, argv + argc)
241 /*! \brief
242 * Converts a t_filenm option into an Options option.
244 * \param options Options object to add the new option to.
245 * \param fnm t_filenm option to convert.
247 void filenmToOptions(Options *options, t_filenm *fnm);
248 /*! \brief
249 * Converts a t_pargs option into an Options option.
251 * \param options Options object to add the new option to.
252 * \param pa t_pargs option to convert.
254 void pargsToOptions(Options *options, t_pargs *pa);
256 /*! \brief
257 * Copies values back from options to t_pargs/t_filenm.
259 void copyValues();
261 private:
262 struct FileNameData
264 //! Creates a conversion helper for a given `t_filenm` struct.
265 explicit FileNameData(t_filenm *fnm) : fnm(fnm), optionInfo(nullptr)
269 //! t_filenm structure to receive the final values.
270 t_filenm *fnm;
271 //! Option info object for the created FileNameOption.
272 FileNameOptionInfo *optionInfo;
273 //! Value storage for the created FileNameOption.
274 std::vector<std::string> values;
276 struct ProgramArgData
278 //! Creates a conversion helper for a given `t_pargs` struct.
279 explicit ProgramArgData(t_pargs *pa)
280 : pa(pa), optionInfo(nullptr), enumIndex(0), boolValue(false)
284 //! t_pargs structure to receive the final values.
285 t_pargs *pa;
286 //! Option info object for the created option.
287 OptionInfo *optionInfo;
288 //! Value storage for a non-enum StringOption (unused for other types).
289 std::string stringValue;
290 //! Value storage for an enum option (unused for other types).
291 int enumIndex;
292 //! Value storage for a BooleanOption (unused for other types).
293 bool boolValue;
296 std::vector<const char *> argv_;
297 // These are lists instead of vectors to avoid relocating existing
298 // objects in case the container is reallocated (the Options object
299 // contains pointes to members of the objects, which would get
300 // invalidated).
301 std::list<FileNameData> fileNameOptions_;
302 std::list<ProgramArgData> programArgs_;
304 GMX_DISALLOW_COPY_AND_ASSIGN(OptionsAdapter);
307 void OptionsAdapter::filenmToOptions(Options *options, t_filenm *fnm)
309 const bool bRead = ((fnm->flag & ffREAD) != 0);
310 const bool bWrite = ((fnm->flag & ffWRITE) != 0);
311 const bool bOptional = ((fnm->flag & ffOPT) != 0);
312 const bool bLibrary = ((fnm->flag & ffLIB) != 0);
313 const bool bMultiple = ((fnm->flag & ffMULT) != 0);
314 const bool bMissing = ((fnm->flag & ffALLOW_MISSING) != 0);
315 const char *const name = (fnm->opt ? &fnm->opt[1] : &ftp2defopt(fnm->ftp)[1]);
316 const char * defName = fnm->fn;
317 int defType = -1;
318 if (defName == nullptr)
320 defName = ftp2defnm(fnm->ftp);
322 else if (Path::hasExtension(defName))
324 defType = fn2ftp(defName);
325 GMX_RELEASE_ASSERT(defType != efNR,
326 "File name option specifies an invalid extension");
328 fileNameOptions_.emplace_back(fnm);
329 FileNameData &data = fileNameOptions_.back();
330 data.optionInfo = options->addOption(
331 FileNameOption(name).storeVector(&data.values)
332 .defaultBasename(defName).defaultType(defType)
333 .legacyType(fnm->ftp).legacyOptionalBehavior()
334 .readWriteFlags(bRead, bWrite).required(!bOptional)
335 .libraryFile(bLibrary).multiValue(bMultiple)
336 .allowMissing(bMissing)
337 .description(ftp2desc(fnm->ftp)));
340 void OptionsAdapter::pargsToOptions(Options *options, t_pargs *pa)
342 const bool bHidden = startsWith(pa->desc, "HIDDEN");
343 const char *const name = &pa->option[1];
344 const char *const desc = (bHidden ? &pa->desc[6] : pa->desc);
345 programArgs_.emplace_back(pa);
346 ProgramArgData &data = programArgs_.back();
347 switch (pa->type)
349 case etINT:
350 data.optionInfo = options->addOption(
351 IntegerOption(name).store(pa->u.i)
352 .description(desc).hidden(bHidden));
353 return;
354 case etINT64:
355 data.optionInfo = options->addOption(
356 Int64Option(name).store(pa->u.is)
357 .description(desc).hidden(bHidden));
358 return;
359 case etREAL:
360 data.optionInfo = options->addOption(
361 RealOption(name).store(pa->u.r)
362 .description(desc).hidden(bHidden));
363 return;
364 case etTIME:
365 data.optionInfo = options->addOption(
366 RealOption(name).store(pa->u.r).timeValue()
367 .description(desc).hidden(bHidden));
368 return;
369 case etSTR:
371 const char *const defValue = (*pa->u.c != nullptr ? *pa->u.c : "");
372 data.optionInfo = options->addOption(
373 StringOption(name).store(&data.stringValue)
374 .defaultValue(defValue)
375 .description(desc).hidden(bHidden));
376 return;
378 case etBOOL:
379 data.optionInfo = options->addOption(
380 BooleanOption(name).store(&data.boolValue)
381 .defaultValue(*pa->u.b)
382 .description(desc).hidden(bHidden));
383 return;
384 case etRVEC:
385 data.optionInfo = options->addOption(
386 RealOption(name).store(*pa->u.rv).vector()
387 .description(desc).hidden(bHidden));
388 return;
389 case etENUM:
391 const int defaultIndex = (pa->u.c[0] != nullptr ? nenum(pa->u.c) - 1 : 0);
392 data.optionInfo = options->addOption(
393 EnumIntOption(name).store(&data.enumIndex)
394 .defaultValue(defaultIndex)
395 .enumValueFromNullTerminatedArray(pa->u.c + 1)
396 .description(desc).hidden(bHidden));
397 return;
400 GMX_THROW(NotImplementedError("Argument type not implemented"));
403 void OptionsAdapter::copyValues()
405 std::list<FileNameData>::const_iterator file;
406 for (file = fileNameOptions_.begin(); file != fileNameOptions_.end(); ++file)
408 if (file->optionInfo->isSet())
410 file->fnm->flag |= ffSET;
412 file->fnm->filenames = file->values;
414 std::list<ProgramArgData>::const_iterator arg;
415 for (arg = programArgs_.begin(); arg != programArgs_.end(); ++arg)
417 arg->pa->bSet = arg->optionInfo->isSet();
418 switch (arg->pa->type)
420 case etSTR:
422 if (arg->pa->bSet)
424 std::vector<const char *>::const_iterator pos =
425 std::find(argv_.begin(), argv_.end(), arg->stringValue);
426 GMX_RELEASE_ASSERT(pos != argv_.end(),
427 "String argument got a value not in argv");
428 *arg->pa->u.c = *pos;
430 break;
432 case etBOOL:
433 *arg->pa->u.b = arg->boolValue;
434 break;
435 case etENUM:
436 *arg->pa->u.c = arg->pa->u.c[arg->enumIndex + 1];
437 break;
438 default:
439 // For other types, there is nothing type-specific to do.
440 break;
445 } // namespace
447 } // namespace gmx
449 gmx_bool parse_common_args(int *argc, char *argv[], unsigned long Flags,
450 int nfile, t_filenm fnm[], int npargs, t_pargs *pa,
451 int ndesc, const char **desc,
452 int nbugs, const char **bugs,
453 gmx_output_env_t **oenv)
455 /* This array should match the order of the enum in oenv.h */
456 const char *const xvg_formats[] = { "xmgrace", "xmgr", "none" };
458 // Lambda function to test the (local) Flags parameter against a bit mask.
459 auto isFlagSet = [Flags](unsigned long bits) {
460 return (Flags & bits) == bits;
465 double tbegin = 0.0, tend = 0.0, tdelta = 0.0;
466 bool bBeginTimeSet = false, bEndTimeSet = false, bDtSet = false;
467 bool bView = false;
468 int xvgFormat = 0;
469 gmx::OptionsAdapter adapter(*argc, argv);
470 gmx::Options options;
471 gmx::OptionsBehaviorCollection behaviors(&options);
472 gmx::FileNameOptionManager fileOptManager;
474 fileOptManager.disableInputOptionChecking(
475 isFlagSet(PCA_DISABLE_INPUT_FILE_CHECKING));
476 options.addManager(&fileOptManager);
478 if (isFlagSet(PCA_CAN_SET_DEFFNM))
480 fileOptManager.addDefaultFileNameOption(&options, "deffnm");
482 if (isFlagSet(PCA_CAN_BEGIN))
484 options.addOption(
485 gmx::DoubleOption("b")
486 .store(&tbegin).storeIsSet(&bBeginTimeSet).timeValue()
487 .description("Time of first frame to read from trajectory (default unit %t)"));
489 if (isFlagSet(PCA_CAN_END))
491 options.addOption(
492 gmx::DoubleOption("e")
493 .store(&tend).storeIsSet(&bEndTimeSet).timeValue()
494 .description("Time of last frame to read from trajectory (default unit %t)"));
496 if (isFlagSet(PCA_CAN_DT))
498 options.addOption(
499 gmx::DoubleOption("dt")
500 .store(&tdelta).storeIsSet(&bDtSet).timeValue()
501 .description("Only use frame when t MOD dt = first time (default unit %t)"));
503 gmx::TimeUnit timeUnit = gmx::TimeUnit_Default;
504 if (isFlagSet(PCA_TIME_UNIT))
506 std::shared_ptr<gmx::TimeUnitBehavior> timeUnitBehavior(
507 new gmx::TimeUnitBehavior());
508 timeUnitBehavior->setTimeUnitStore(&timeUnit);
509 timeUnitBehavior->setTimeUnitFromEnvironment();
510 timeUnitBehavior->addTimeUnitOption(&options, "tu");
511 behaviors.addBehavior(timeUnitBehavior);
513 if (isFlagSet(PCA_CAN_VIEW))
515 options.addOption(
516 gmx::BooleanOption("w").store(&bView)
517 .description("View output [REF].xvg[ref], [REF].xpm[ref], "
518 "[REF].eps[ref] and [REF].pdb[ref] files"));
521 bool bXvgr = false;
522 for (int i = 0; i < nfile; i++)
524 bXvgr = bXvgr || (fnm[i].ftp == efXVG);
526 xvgFormat = gmx::getDefaultXvgFormat(xvg_formats);
527 if (bXvgr)
529 options.addOption(
530 gmx::EnumIntOption("xvg").enumValue(xvg_formats)
531 .store(&xvgFormat)
532 .description("xvg plot formatting"));
535 /* Now append the program specific arguments */
536 for (int i = 0; i < nfile; i++)
538 adapter.filenmToOptions(&options, &fnm[i]);
540 for (int i = 0; i < npargs; i++)
542 adapter.pargsToOptions(&options, &pa[i]);
545 const gmx::CommandLineHelpContext *context =
546 gmx::GlobalCommandLineHelpContext::get();
547 if (context != nullptr)
549 GMX_RELEASE_ASSERT(gmx_node_rank() == 0,
550 "Help output should be handled higher up and "
551 "only get called only on the master rank");
552 gmx::CommandLineHelpWriter(options)
553 .setHelpText(gmx::constArrayRefFromArray<const char *>(desc, ndesc))
554 .setKnownIssues(gmx::constArrayRefFromArray(bugs, nbugs))
555 .writeHelp(*context);
556 return FALSE;
559 /* Now parse all the command-line options */
560 gmx::CommandLineParser(&options)
561 .skipUnknown(isFlagSet(PCA_NOEXIT_ON_ARGS))
562 .allowPositionalArguments(isFlagSet(PCA_NOEXIT_ON_ARGS))
563 .parse(argc, argv);
564 behaviors.optionsFinishing();
565 options.finish();
567 /* set program name, command line, and default values for output options */
568 output_env_init(oenv, gmx::getProgramContext(),
569 static_cast<time_unit_t>(timeUnit + 1), bView, // NOLINT(bugprone-misplaced-widening-cast)
570 static_cast<xvg_format_t>(xvgFormat + 1), 0);
572 /* Extract Time info from arguments */
573 if (bBeginTimeSet)
575 setTimeValue(TBEGIN, tbegin);
577 if (bEndTimeSet)
579 setTimeValue(TEND, tend);
581 if (bDtSet)
583 setTimeValue(TDELTA, tdelta);
586 adapter.copyValues();
588 return TRUE;
590 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;