Add conserved quantity for Berendsen P-couple
[gromacs.git] / src / gromacs / commandline / pargs.cpp
blob0ceb6a777fc23af35234e002baca4843b0fc5793
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, 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/cstringutil.h"
62 #include "gromacs/utility/exceptions.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/gmxassert.h"
65 #include "gromacs/utility/path.h"
66 #include "gromacs/utility/programcontext.h"
67 #include "gromacs/utility/smalloc.h"
68 #include "gromacs/utility/stringutil.h"
70 /* The source code in this file should be thread-safe.
71 Please keep it that way. */
73 int nenum(const char *const enumc[])
75 int i;
77 i = 1;
78 /* we *can* compare pointers directly here! */
79 while (enumc[i] && enumc[0] != enumc[i])
81 i++;
84 return i;
87 int opt2parg_int(const char *option, int nparg, t_pargs pa[])
89 int i;
91 for (i = 0; (i < nparg); i++)
93 if (strcmp(pa[i].option, option) == 0)
95 return *pa[i].u.i;
99 gmx_fatal(FARGS, "No integer option %s in pargs", option);
101 return 0;
104 gmx_bool opt2parg_bool(const char *option, int nparg, t_pargs pa[])
106 int i;
108 for (i = 0; (i < nparg); i++)
110 if (strcmp(pa[i].option, option) == 0)
112 return *pa[i].u.b;
116 gmx_fatal(FARGS, "No boolean option %s in pargs", option);
118 return FALSE;
121 real opt2parg_real(const char *option, int nparg, t_pargs pa[])
123 int i;
125 for (i = 0; (i < nparg); i++)
127 if (strcmp(pa[i].option, option) == 0)
129 return *pa[i].u.r;
133 gmx_fatal(FARGS, "No real option %s in pargs", option);
135 return 0.0;
138 const char *opt2parg_str(const char *option, int nparg, t_pargs pa[])
140 int i;
142 for (i = 0; (i < nparg); i++)
144 if (strcmp(pa[i].option, option) == 0)
146 return *(pa[i].u.c);
150 gmx_fatal(FARGS, "No string option %s in pargs", option);
152 return nullptr;
155 gmx_bool opt2parg_bSet(const char *option, int nparg, t_pargs pa[])
157 int i;
159 for (i = 0; (i < nparg); i++)
161 if (strcmp(pa[i].option, option) == 0)
163 return pa[i].bSet;
167 gmx_fatal(FARGS, "No such option %s in pargs", option);
169 return FALSE; /* Too make some compilers happy */
172 const char *opt2parg_enum(const char *option, int nparg, t_pargs pa[])
174 int i;
176 for (i = 0; (i < nparg); i++)
178 if (strcmp(pa[i].option, option) == 0)
180 return pa[i].u.c[0];
184 gmx_fatal(FARGS, "No such option %s in pargs", option);
186 return nullptr;
189 /********************************************************************
190 * parse_common_args()
193 namespace gmx
196 namespace
199 /*! \brief
200 * Returns the index of the default xvg format.
202 * \ingroup module_commandline
204 int getDefaultXvgFormat(gmx::ConstArrayRef<const char *> xvgFormats)
206 const char *const select = getenv("GMX_VIEW_XVG");
207 if (select != nullptr)
209 ConstArrayRef<const char *>::const_iterator i =
210 std::find(xvgFormats.begin(), xvgFormats.end(), std::string(select));
211 if (i != xvgFormats.end())
213 return std::distance(xvgFormats.begin(), i);
215 else
217 return exvgNONE - 1;
220 /* The default is the first option */
221 return 0;
224 /*! \brief
225 * Conversion helper between t_pargs/t_filenm and Options.
227 * This class holds the necessary mapping between the old C structures and
228 * the new C++ options to allow copying values back after parsing for cases
229 * where the C++ options do not directly provide the type of value required for
230 * the C structures.
232 * \ingroup module_commandline
234 class OptionsAdapter
236 public:
237 /*! \brief
238 * Initializes the adapter to convert from a specified command line.
240 * The command line is required, because t_pargs wants to return
241 * strings by reference to the original command line.
242 * OptionsAdapter creates a copy of the `argv` array (but not the
243 * strings) to make this possible, even if the parser removes
244 * options it has recognized.
246 OptionsAdapter(int argc, const char *const argv[])
247 : argv_(argv, argv + argc)
251 /*! \brief
252 * Converts a t_filenm option into an Options option.
254 * \param options Options object to add the new option to.
255 * \param fnm t_filenm option to convert.
257 void filenmToOptions(Options *options, t_filenm *fnm);
258 /*! \brief
259 * Converts a t_pargs option into an Options option.
261 * \param options Options object to add the new option to.
262 * \param pa t_pargs option to convert.
264 void pargsToOptions(Options *options, t_pargs *pa);
266 /*! \brief
267 * Copies values back from options to t_pargs/t_filenm.
269 void copyValues(bool bReadNode);
271 private:
272 struct FileNameData
274 //! Creates a conversion helper for a given `t_filenm` struct.
275 explicit FileNameData(t_filenm *fnm) : fnm(fnm), optionInfo(nullptr)
279 //! t_filenm structure to receive the final values.
280 t_filenm *fnm;
281 //! Option info object for the created FileNameOption.
282 FileNameOptionInfo *optionInfo;
283 //! Value storage for the created FileNameOption.
284 std::vector<std::string> values;
286 struct ProgramArgData
288 //! Creates a conversion helper for a given `t_pargs` struct.
289 explicit ProgramArgData(t_pargs *pa)
290 : pa(pa), optionInfo(nullptr), enumIndex(0), boolValue(false)
294 //! t_pargs structure to receive the final values.
295 t_pargs *pa;
296 //! Option info object for the created option.
297 OptionInfo *optionInfo;
298 //! Value storage for a non-enum StringOption (unused for other types).
299 std::string stringValue;
300 //! Value storage for an enum option (unused for other types).
301 int enumIndex;
302 //! Value storage for a BooleanOption (unused for other types).
303 bool boolValue;
306 std::vector<const char *> argv_;
307 // These are lists instead of vectors to avoid relocating existing
308 // objects in case the container is reallocated (the Options object
309 // contains pointes to members of the objects, which would get
310 // invalidated).
311 std::list<FileNameData> fileNameOptions_;
312 std::list<ProgramArgData> programArgs_;
314 GMX_DISALLOW_COPY_AND_ASSIGN(OptionsAdapter);
317 void OptionsAdapter::filenmToOptions(Options *options, t_filenm *fnm)
319 if (fnm->opt == nullptr)
321 // Existing code may use opt2fn() instead of ftp2fn() for
322 // options that use the default option name, so we need to
323 // keep the old behavior instead of fixing opt2fn().
324 // TODO: Check that this is not the case, remove this, and make
325 // opt2*() work even if fnm->opt is NULL for some options.
326 fnm->opt = ftp2defopt(fnm->ftp);
328 const bool bRead = ((fnm->flag & ffREAD) != 0);
329 const bool bWrite = ((fnm->flag & ffWRITE) != 0);
330 const bool bOptional = ((fnm->flag & ffOPT) != 0);
331 const bool bLibrary = ((fnm->flag & ffLIB) != 0);
332 const bool bMultiple = ((fnm->flag & ffMULT) != 0);
333 const bool bMissing = ((fnm->flag & ffALLOW_MISSING) != 0);
334 const char *const name = &fnm->opt[1];
335 const char * defName = fnm->fn;
336 int defType = -1;
337 if (defName == nullptr)
339 defName = ftp2defnm(fnm->ftp);
341 else if (Path::hasExtension(defName))
343 defType = fn2ftp(defName);
344 GMX_RELEASE_ASSERT(defType != efNR,
345 "File name option specifies an invalid extension");
347 fileNameOptions_.emplace_back(fnm);
348 FileNameData &data = fileNameOptions_.back();
349 data.optionInfo = options->addOption(
350 FileNameOption(name).storeVector(&data.values)
351 .defaultBasename(defName).defaultType(defType)
352 .legacyType(fnm->ftp).legacyOptionalBehavior()
353 .readWriteFlags(bRead, bWrite).required(!bOptional)
354 .libraryFile(bLibrary).multiValue(bMultiple)
355 .allowMissing(bMissing)
356 .description(ftp2desc(fnm->ftp)));
359 void OptionsAdapter::pargsToOptions(Options *options, t_pargs *pa)
361 const bool bHidden = startsWith(pa->desc, "HIDDEN");
362 const char *const name = &pa->option[1];
363 const char *const desc = (bHidden ? &pa->desc[6] : pa->desc);
364 programArgs_.emplace_back(pa);
365 ProgramArgData &data = programArgs_.back();
366 switch (pa->type)
368 case etINT:
369 data.optionInfo = options->addOption(
370 IntegerOption(name).store(pa->u.i)
371 .description(desc).hidden(bHidden));
372 return;
373 case etINT64:
374 data.optionInfo = options->addOption(
375 Int64Option(name).store(pa->u.is)
376 .description(desc).hidden(bHidden));
377 return;
378 case etREAL:
379 data.optionInfo = options->addOption(
380 RealOption(name).store(pa->u.r)
381 .description(desc).hidden(bHidden));
382 return;
383 case etTIME:
384 data.optionInfo = options->addOption(
385 RealOption(name).store(pa->u.r).timeValue()
386 .description(desc).hidden(bHidden));
387 return;
388 case etSTR:
390 const char *const defValue = (*pa->u.c != nullptr ? *pa->u.c : "");
391 data.optionInfo = options->addOption(
392 StringOption(name).store(&data.stringValue)
393 .defaultValue(defValue)
394 .description(desc).hidden(bHidden));
395 return;
397 case etBOOL:
398 data.optionInfo = options->addOption(
399 BooleanOption(name).store(&data.boolValue)
400 .defaultValue(*pa->u.b)
401 .description(desc).hidden(bHidden));
402 return;
403 case etRVEC:
404 data.optionInfo = options->addOption(
405 RealOption(name).store(*pa->u.rv).vector()
406 .description(desc).hidden(bHidden));
407 return;
408 case etENUM:
410 const int defaultIndex = (pa->u.c[0] != nullptr ? nenum(pa->u.c) - 1 : 0);
411 data.optionInfo = options->addOption(
412 EnumIntOption(name).store(&data.enumIndex)
413 .defaultValue(defaultIndex)
414 .enumValueFromNullTerminatedArray(pa->u.c + 1)
415 .description(desc).hidden(bHidden));
416 return;
419 GMX_THROW(NotImplementedError("Argument type not implemented"));
422 void OptionsAdapter::copyValues(bool bReadNode)
424 std::list<FileNameData>::const_iterator file;
425 for (file = fileNameOptions_.begin(); file != fileNameOptions_.end(); ++file)
427 if (!bReadNode && (file->fnm->flag & ffREAD))
429 continue;
431 if (file->optionInfo->isSet())
433 file->fnm->flag |= ffSET;
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)
439 file->fnm->fns[i] = gmx_strdup(file->values[i].c_str());
442 std::list<ProgramArgData>::const_iterator arg;
443 for (arg = programArgs_.begin(); arg != programArgs_.end(); ++arg)
445 arg->pa->bSet = arg->optionInfo->isSet();
446 switch (arg->pa->type)
448 case etSTR:
450 if (arg->pa->bSet)
452 std::vector<const char *>::const_iterator pos =
453 std::find(argv_.begin(), argv_.end(), arg->stringValue);
454 GMX_RELEASE_ASSERT(pos != argv_.end(),
455 "String argument got a value not in argv");
456 *arg->pa->u.c = *pos;
458 break;
460 case etBOOL:
461 *arg->pa->u.b = arg->boolValue;
462 break;
463 case etENUM:
464 *arg->pa->u.c = arg->pa->u.c[arg->enumIndex + 1];
465 break;
466 default:
467 // For other types, there is nothing type-specific to do.
468 break;
473 } // namespace
475 } // namespace gmx
477 gmx_bool parse_common_args(int *argc, char *argv[], unsigned long Flags,
478 int nfile, t_filenm fnm[], int npargs, t_pargs *pa,
479 int ndesc, const char **desc,
480 int nbugs, const char **bugs,
481 gmx_output_env_t **oenv)
483 /* This array should match the order of the enum in oenv.h */
484 const char *const xvg_formats[] = { "xmgrace", "xmgr", "none" };
486 // Handle the flags argument, which is a bit field
487 // The FF macro returns whether or not the bit is set
488 #define FF(arg) ((Flags & arg) == arg)
492 double tbegin = 0.0, tend = 0.0, tdelta = 0.0;
493 bool bBeginTimeSet = false, bEndTimeSet = false, bDtSet = false;
494 bool bView = false;
495 int xvgFormat = 0;
496 gmx::OptionsAdapter adapter(*argc, argv);
497 gmx::Options options;
498 gmx::OptionsBehaviorCollection behaviors(&options);
499 gmx::FileNameOptionManager fileOptManager;
501 fileOptManager.disableInputOptionChecking(
502 FF(PCA_NOT_READ_NODE) || FF(PCA_DISABLE_INPUT_FILE_CHECKING));
503 options.addManager(&fileOptManager);
505 if (FF(PCA_CAN_SET_DEFFNM))
507 fileOptManager.addDefaultFileNameOption(&options, "deffnm");
509 if (FF(PCA_CAN_BEGIN))
511 options.addOption(
512 gmx::DoubleOption("b")
513 .store(&tbegin).storeIsSet(&bBeginTimeSet).timeValue()
514 .description("First frame (%t) to read from trajectory"));
516 if (FF(PCA_CAN_END))
518 options.addOption(
519 gmx::DoubleOption("e")
520 .store(&tend).storeIsSet(&bEndTimeSet).timeValue()
521 .description("Last frame (%t) to read from trajectory"));
523 if (FF(PCA_CAN_DT))
525 options.addOption(
526 gmx::DoubleOption("dt")
527 .store(&tdelta).storeIsSet(&bDtSet).timeValue()
528 .description("Only use frame when t MOD dt = first time (%t)"));
530 gmx::TimeUnit timeUnit = gmx::TimeUnit_Default;
531 if (FF(PCA_TIME_UNIT))
533 std::shared_ptr<gmx::TimeUnitBehavior> timeUnitBehavior(
534 new gmx::TimeUnitBehavior());
535 timeUnitBehavior->setTimeUnitStore(&timeUnit);
536 timeUnitBehavior->setTimeUnitFromEnvironment();
537 timeUnitBehavior->addTimeUnitOption(&options, "tu");
538 behaviors.addBehavior(timeUnitBehavior);
540 if (FF(PCA_CAN_VIEW))
542 options.addOption(
543 gmx::BooleanOption("w").store(&bView)
544 .description("View output [REF].xvg[ref], [REF].xpm[ref], "
545 "[REF].eps[ref] and [REF].pdb[ref] files"));
548 bool bXvgr = false;
549 for (int i = 0; i < nfile; i++)
551 bXvgr = bXvgr || (fnm[i].ftp == efXVG);
553 xvgFormat = gmx::getDefaultXvgFormat(xvg_formats);
554 if (bXvgr)
556 options.addOption(
557 gmx::EnumIntOption("xvg").enumValue(xvg_formats)
558 .store(&xvgFormat)
559 .description("xvg plot formatting"));
562 /* Now append the program specific arguments */
563 for (int i = 0; i < nfile; i++)
565 adapter.filenmToOptions(&options, &fnm[i]);
567 for (int i = 0; i < npargs; i++)
569 adapter.pargsToOptions(&options, &pa[i]);
572 const gmx::CommandLineHelpContext *context =
573 gmx::GlobalCommandLineHelpContext::get();
574 if (context != nullptr)
576 GMX_RELEASE_ASSERT(gmx_node_rank() == 0,
577 "Help output should be handled higher up and "
578 "only get called only on the master rank");
579 gmx::CommandLineHelpWriter(options)
580 .setHelpText(gmx::constArrayRefFromArray<const char *>(desc, ndesc))
581 .setKnownIssues(gmx::constArrayRefFromArray(bugs, nbugs))
582 .writeHelp(*context);
583 return FALSE;
586 /* Now parse all the command-line options */
587 gmx::CommandLineParser(&options).skipUnknown(FF(PCA_NOEXIT_ON_ARGS))
588 .parse(argc, argv);
589 behaviors.optionsFinishing();
590 options.finish();
592 /* set program name, command line, and default values for output options */
593 output_env_init(oenv, gmx::getProgramContext(),
594 (time_unit_t)(timeUnit + 1), bView,
595 (xvg_format_t)(xvgFormat + 1), 0);
597 /* Extract Time info from arguments */
598 if (bBeginTimeSet)
600 setTimeValue(TBEGIN, tbegin);
602 if (bEndTimeSet)
604 setTimeValue(TEND, tend);
606 if (bDtSet)
608 setTimeValue(TDELTA, tdelta);
611 adapter.copyValues(!FF(PCA_NOT_READ_NODE));
613 return TRUE;
615 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
616 #undef FF