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) 2011-2020, 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 /*! \libinternal \file
39 * \brief This file declares helper functionality for legacy option handling for mdrun
41 * \author Berk Hess <hess@kth.se>
42 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
43 * \author Erik Lindahl <erik@kth.se>
44 * \author Mark Abraham <mark.j.abraham@gmail.com>
46 * \ingroup module_mdrun
49 #ifndef GMX_MDRUN_LEGACYMDRUNOPTIONS_H
50 #define GMX_MDRUN_LEGACYMDRUNOPTIONS_H
52 #include "gromacs/commandline/filenm.h"
53 #include "gromacs/commandline/pargs.h"
54 #include "gromacs/domdec/options.h"
55 #include "gromacs/hardware/hw_info.h"
56 #include "gromacs/mdtypes/mdrunoptions.h"
58 #include "replicaexchange.h"
64 * \brief This class provides the same command-line option
65 * functionality to both CLI and API sessions.
67 * This class should not exist, but is necessary now to introduce
68 * support for the CLI and API without duplicating code. It should be
69 * eliminated following the TODOs below.
71 * \warning Instances provide lifetime scope for members that do not have
72 * effective lifetime management or which are frequently accessed unsafely.
73 * The caller is responsible for keeping a LegacyMdrunOptions object alive
74 * for as long as any consumers, direct or transitive.
76 * \todo Modules in mdrun should acquire proper option handling so
77 * that all of these declarations and defaults are local to the
80 * \todo Contextual aspects, such as working directory
81 * and environment variable handling are more properly
82 * the role of SimulationContext, and should be moved there.
84 class LegacyMdrunOptions
87 //! Ongoing collection of mdrun options
88 MdrunOptions mdrunOptions
;
89 //! Options for the domain decomposition.
90 DomdecOptions domdecOptions
;
91 //! Parallelism-related user options.
93 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
94 int nstlist_cmdline
= 0;
95 //! Parameters for replica-exchange simulations.
96 ReplicaExchangeParameters replExParams
;
98 //! Filename options to fill from command-line argument values.
99 std::vector
<t_filenm
> filenames
= { { { efTPR
, nullptr, nullptr, ffREAD
},
100 { efTRN
, "-o", nullptr, ffWRITE
},
101 { efCOMPRESSED
, "-x", nullptr, ffOPTWR
},
102 { efCPT
, "-cpi", nullptr, ffOPTRD
| ffALLOW_MISSING
},
103 { efCPT
, "-cpo", nullptr, ffOPTWR
},
104 { efSTO
, "-c", "confout", ffWRITE
},
105 { efEDR
, "-e", "ener", ffWRITE
},
106 { efLOG
, "-g", "md", ffWRITE
},
107 { efXVG
, "-dhdl", "dhdl", ffOPTWR
},
108 { efXVG
, "-field", "field", ffOPTWR
},
109 { efXVG
, "-table", "table", ffOPTRD
},
110 { efXVG
, "-tablep", "tablep", ffOPTRD
},
111 { efXVG
, "-tableb", "table", ffOPTRDMULT
},
112 { efTRX
, "-rerun", "rerun", ffOPTRD
},
113 { efXVG
, "-tpi", "tpi", ffOPTWR
},
114 { efXVG
, "-tpid", "tpidist", ffOPTWR
},
115 { efEDI
, "-ei", "sam", ffOPTRD
},
116 { efXVG
, "-eo", "edsam", ffOPTWR
},
117 { efXVG
, "-px", "pullx", ffOPTWR
},
118 { efXVG
, "-pf", "pullf", ffOPTWR
},
119 { efXVG
, "-ro", "rotation", ffOPTWR
},
120 { efLOG
, "-ra", "rotangles", ffOPTWR
},
121 { efLOG
, "-rs", "rotslabs", ffOPTWR
},
122 { efLOG
, "-rt", "rottorque", ffOPTWR
},
123 { efMTX
, "-mtx", "nm", ffOPTWR
},
124 { efRND
, "-multidir", nullptr, ffOPTRDMULT
},
125 { efXVG
, "-awh", "awhinit", ffOPTRD
},
126 { efDAT
, "-membed", "membed", ffOPTRD
},
127 { efTOP
, "-mp", "membed", ffOPTRD
},
128 { efNDX
, "-mn", "membed", ffOPTRD
},
129 { efXVG
, "-if", "imdforces", ffOPTWR
},
130 { efXVG
, "-swap", "swapions", ffOPTWR
} } };
132 //! Print a warning if any force is larger than this (in kJ/mol nm).
135 //! The value of the -append option
136 bool appendOption
= true;
138 /*! \brief Output context for writing text files
140 * \todo Clarify initialization, ownership, and lifetime. */
141 gmx_output_env_t
* oenv
= nullptr;
143 /*! \brief Command line options, defaults, docs and storage for them to fill. */
145 rvec realddxyz
= { 0, 0, 0 };
146 const char* ddrank_opt_choices
[static_cast<int>(DdRankOrder::Count
) + 1] = {
147 nullptr, "interleave", "pp_pme", "cartesian", nullptr
149 const char* dddlb_opt_choices
[static_cast<int>(DlbOption::Count
) + 1] = { nullptr, "auto", "no",
151 const char* thread_aff_opt_choices
[static_cast<int>(ThreadAffinity::Count
) + 1] = {
152 nullptr, "auto", "on", "off", nullptr
154 const char* nbpu_opt_choices
[5] = { nullptr, "auto", "cpu", "gpu", nullptr };
155 const char* pme_opt_choices
[5] = { nullptr, "auto", "cpu", "gpu", nullptr };
156 const char* pme_fft_opt_choices
[5] = { nullptr, "auto", "cpu", "gpu", nullptr };
157 const char* bonded_opt_choices
[5] = { nullptr, "auto", "cpu", "gpu", nullptr };
158 const char* update_opt_choices
[5] = { nullptr, "auto", "cpu", "gpu", nullptr };
159 const char* gpuIdsAvailable
= "";
160 const char* userGpuTaskAssignment
= "";
163 ImdOptions
& imdOptions
= mdrunOptions
.imdOptions
;
167 { "-dd", FALSE
, etRVEC
, { &realddxyz
}, "Domain decomposition grid, 0 is optimize" },
168 { "-ddorder", FALSE
, etENUM
, { ddrank_opt_choices
}, "DD rank order" },
172 { &domdecOptions
.numPmeRanks
},
173 "Number of separate ranks to be used for PME, -1 is guess" },
177 { &hw_opt
.nthreads_tot
},
178 "Total number of threads to start (0 is guess)" },
182 { &hw_opt
.nthreads_tmpi
},
183 "Number of thread-MPI ranks to start (0 is guess)" },
187 { &hw_opt
.nthreads_omp
},
188 "Number of OpenMP threads per MPI rank to start (0 is guess)" },
192 { &hw_opt
.nthreads_omp_pme
},
193 "Number of OpenMP threads per MPI rank to start (0 is -ntomp)" },
197 { thread_aff_opt_choices
},
198 "Whether mdrun should try to set thread affinities" },
202 { &hw_opt
.core_pinning_offset
},
203 "The lowest logical core number to which mdrun should pin the first thread" },
207 { &hw_opt
.core_pinning_stride
},
208 "Pinning distance in logical cores for threads, use 0 to minimize the number of threads "
209 "per physical core" },
213 { &gpuIdsAvailable
},
214 "List of unique GPU device IDs available to use" },
218 { &userGpuTaskAssignment
},
219 "List of GPU device IDs, mapping each PP task on each node to a device" },
223 { &domdecOptions
.checkBondedInteractions
},
224 "Check for all bonded interactions with DD" },
228 { &domdecOptions
.useBondedCommunication
},
229 "HIDDENUse special bonded atom communication when [TT]-rdd[tt] > cut-off" },
233 { &domdecOptions
.minimumCommunicationRange
},
234 "The maximum distance for bonded interactions with DD (nm), 0 is determine from initial "
239 { &domdecOptions
.constraintCommunicationRange
},
240 "Maximum distance for P-LINCS (nm), 0 is estimate" },
241 { "-dlb", FALSE
, etENUM
, { dddlb_opt_choices
}, "Dynamic load balancing (with DD)" },
245 { &domdecOptions
.dlbScaling
},
246 "Fraction in (0,1) by whose reciprocal the initial DD cell size will be increased in "
248 "provide a margin in which dynamic load balancing can act while preserving the minimum "
253 { &domdecOptions
.cellSizeX
},
254 "HIDDENA string containing a vector of the relative sizes in the x "
255 "direction of the corresponding DD cells. Only effective with static "
260 { &domdecOptions
.cellSizeY
},
261 "HIDDENA string containing a vector of the relative sizes in the y "
262 "direction of the corresponding DD cells. Only effective with static "
267 { &domdecOptions
.cellSizeZ
},
268 "HIDDENA string containing a vector of the relative sizes in the z "
269 "direction of the corresponding DD cells. Only effective with static "
271 { "-nb", FALSE
, etENUM
, { nbpu_opt_choices
}, "Calculate non-bonded interactions on" },
275 { &nstlist_cmdline
},
276 "Set nstlist when using a Verlet buffer tolerance (0 is guess)" },
280 { &mdrunOptions
.tunePme
},
281 "Optimize PME load between PP/PME ranks or GPU/CPU" },
282 { "-pme", FALSE
, etENUM
, { pme_opt_choices
}, "Perform PME calculations on" },
283 { "-pmefft", FALSE
, etENUM
, { pme_fft_opt_choices
}, "Perform PME FFT calculations on" },
284 { "-bonded", FALSE
, etENUM
, { bonded_opt_choices
}, "Perform bonded calculations on" },
285 { "-update", FALSE
, etENUM
, { update_opt_choices
}, "Perform update and constraints on" },
286 { "-v", FALSE
, etBOOL
, { &mdrunOptions
.verbose
}, "Be loud and noisy" },
287 { "-pforce", FALSE
, etREAL
, { &pforce
}, "Print all forces larger than this (kJ/mol nm)" },
291 { &mdrunOptions
.reproducible
},
292 "Try to avoid optimizations that affect binary reproducibility" },
296 { &mdrunOptions
.checkpointOptions
.period
},
297 "Checkpoint interval (minutes)" },
301 { &mdrunOptions
.checkpointOptions
.keepAndNumberCheckpointFiles
},
302 "Keep and number checkpoint files" },
307 "Append to previous output files when continuing from checkpoint instead of adding the "
308 "simulation part number to all file names" },
312 { &mdrunOptions
.numStepsCommandline
},
313 "Run this number of steps (-1 means infinite, -2 means use mdp option, smaller is "
318 { &mdrunOptions
.maximumHoursToRun
},
319 "Terminate after 0.99 times this time (hours)" },
323 { &replExParams
.exchangeInterval
},
324 "Attempt replica exchange periodically with this period (steps)" },
328 { &replExParams
.numExchanges
},
329 "Number of random exchanges to carry out each exchange interval (N^3 is one suggestion). "
330 " -nex zero or not specified gives neighbor replica exchange." },
334 { &replExParams
.randomSeed
},
335 "Seed for replica exchange, -1 is generate a seed" },
336 { "-imdport", FALSE
, etINT
, { &imdOptions
.port
}, "HIDDENIMD listening port" },
340 { &imdOptions
.wait
},
341 "HIDDENPause the simulation while no IMD client is connected" },
345 { &imdOptions
.terminatable
},
346 "HIDDENAllow termination of the simulation from IMD client" },
350 { &imdOptions
.pull
},
351 "HIDDENAllow pulling in the simulation from IMD client" },
355 { &mdrunOptions
.rerunConstructVsites
},
356 "HIDDENRecalculate virtual site coordinates with [TT]-rerun[tt]" },
360 { &mdrunOptions
.writeConfout
},
361 "HIDDENWrite the last configuration with [TT]-c[tt] and force checkpointing at the last "
366 { &mdrunOptions
.verboseStepPrintInterval
},
367 "HIDDENFrequency of writing the remaining wall clock time for the run" },
371 { &mdrunOptions
.timingOptions
.resetStep
},
372 "HIDDENReset cycle counters after these many time steps" },
376 { &mdrunOptions
.timingOptions
.resetHalfway
},
377 "HIDDENReset the cycle counters after half the number of steps or halfway "
382 //! Parses the command-line input and prepares to start mdrun.
383 int updateFromCommandLine(int argc
, char** argv
, ArrayRef
<const char*> desc
);
385 ~LegacyMdrunOptions();
388 } // end namespace gmx