Preparation for de-duplicating mdrun setup
[gromacs.git] / src / gromacs / mdrun / legacyoptions.h
blob742e04a873545f4c24c637c95df09af3ee0da03c
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) 2011,2012,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 /*! \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
47 * \inlibraryapi
49 #include "gromacs/commandline/filenm.h"
50 #include "gromacs/commandline/pargs.h"
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/hardware/hw_info.h"
53 #include "gromacs/mdlib/mdrun.h"
54 #include "gromacs/mdrun/logging.h"
56 #include "replicaexchange.h"
58 struct gmx_multisim_t;
60 namespace gmx
63 /*! \libinternal
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 * \todo Modules in mdrun should acquire proper option handling so
72 * that all of these declarations and defaults are local to the
73 * modules.
75 * \todo Contextual aspects, such as working directory, MPI
76 * environment, and environment variable handling are more properly
77 * the role of SimulationContext, and should be moved there */
78 class LegacyMdrunOptions
80 public:
81 //! Ongoing collection of mdrun options
82 MdrunOptions mdrunOptions;
83 //! Options for the domain decomposition.
84 DomdecOptions domdecOptions;
85 //! Parallelism-related user options.
86 gmx_hw_opt_t hw_opt;
87 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
88 int nstlist_cmdline = 0;
89 //! Parameters for replica-exchange simulations.
90 ReplicaExchangeParameters replExParams;
92 //! Filename options to fill from command-line argument values.
93 std::vector<t_filenm> filenames =
94 {{{ efTPR, nullptr, nullptr, ffREAD },
95 { efTRN, "-o", nullptr, ffWRITE },
96 { efCOMPRESSED, "-x", nullptr, ffOPTWR },
97 { efCPT, "-cpi", nullptr, ffOPTRD | ffALLOW_MISSING },
98 { efCPT, "-cpo", nullptr, ffOPTWR },
99 { efSTO, "-c", "confout", ffWRITE },
100 { efEDR, "-e", "ener", ffWRITE },
101 { efLOG, "-g", "md", ffWRITE },
102 { efXVG, "-dhdl", "dhdl", ffOPTWR },
103 { efXVG, "-field", "field", ffOPTWR },
104 { efXVG, "-table", "table", ffOPTRD },
105 { efXVG, "-tablep", "tablep", ffOPTRD },
106 { efXVG, "-tableb", "table", ffOPTRDMULT },
107 { efTRX, "-rerun", "rerun", ffOPTRD },
108 { efXVG, "-tpi", "tpi", ffOPTWR },
109 { efXVG, "-tpid", "tpidist", ffOPTWR },
110 { efEDI, "-ei", "sam", ffOPTRD },
111 { efXVG, "-eo", "edsam", ffOPTWR },
112 { efXVG, "-devout", "deviatie", ffOPTWR },
113 { efXVG, "-runav", "runaver", ffOPTWR },
114 { efXVG, "-px", "pullx", ffOPTWR },
115 { efXVG, "-pf", "pullf", ffOPTWR },
116 { efXVG, "-ro", "rotation", ffOPTWR },
117 { efLOG, "-ra", "rotangles", ffOPTWR },
118 { efLOG, "-rs", "rotslabs", ffOPTWR },
119 { efLOG, "-rt", "rottorque", ffOPTWR },
120 { efMTX, "-mtx", "nm", ffOPTWR },
121 { efRND, "-multidir", nullptr, ffOPTRDMULT},
122 { efXVG, "-awh", "awhinit", ffOPTRD },
123 { efDAT, "-membed", "membed", ffOPTRD },
124 { efTOP, "-mp", "membed", ffOPTRD },
125 { efNDX, "-mn", "membed", ffOPTRD },
126 { efXVG, "-if", "imdforces", ffOPTWR },
127 { efXVG, "-swap", "swapions", ffOPTWR }}};
129 //! Print a warning if any force is larger than this (in kJ/mol nm).
130 real pforce = -1;
132 /*! \brief Output context for writing text files
134 * \todo Clarify initialization, ownership, and lifetime. */
135 gmx_output_env_t *oenv = nullptr;
137 //! Handle to file used for logging.
138 LogFilePtr logFileGuard = nullptr;
140 /*! \brief Command line options, defaults, docs and storage for them to fill. */
141 /*! \{ */
142 rvec realddxyz = {0, 0, 0};
143 const char *ddrank_opt_choices[static_cast<int>(DdRankOrder::nr)+1] =
144 { nullptr, "interleave", "pp_pme", "cartesian", nullptr };
145 const char *dddlb_opt_choices[static_cast<int>(DlbOption::nr)+1] =
146 { nullptr, "auto", "no", "yes", nullptr };
147 const char *thread_aff_opt_choices[threadaffNR+1] =
148 { nullptr, "auto", "on", "off", nullptr };
149 const char *nbpu_opt_choices[5] =
150 { nullptr, "auto", "cpu", "gpu", nullptr };
151 const char *pme_opt_choices[5] =
152 { nullptr, "auto", "cpu", "gpu", nullptr };
153 const char *pme_fft_opt_choices[5] =
154 { nullptr, "auto", "cpu", "gpu", nullptr };
155 gmx_bool bTryToAppendFiles = TRUE;
156 const char *gpuIdsAvailable = "";
157 const char *userGpuTaskAssignment = "";
159 ImdOptions &imdOptions = mdrunOptions.imdOptions;
161 t_pargs pa[47] = {
163 { "-dd", FALSE, etRVEC, {&realddxyz},
164 "Domain decomposition grid, 0 is optimize" },
165 { "-ddorder", FALSE, etENUM, {ddrank_opt_choices},
166 "DD rank order" },
167 { "-npme", FALSE, etINT, {&domdecOptions.numPmeRanks},
168 "Number of separate ranks to be used for PME, -1 is guess" },
169 { "-nt", FALSE, etINT, {&hw_opt.nthreads_tot},
170 "Total number of threads to start (0 is guess)" },
171 { "-ntmpi", FALSE, etINT, {&hw_opt.nthreads_tmpi},
172 "Number of thread-MPI ranks to start (0 is guess)" },
173 { "-ntomp", FALSE, etINT, {&hw_opt.nthreads_omp},
174 "Number of OpenMP threads per MPI rank to start (0 is guess)" },
175 { "-ntomp_pme", FALSE, etINT, {&hw_opt.nthreads_omp_pme},
176 "Number of OpenMP threads per MPI rank to start (0 is -ntomp)" },
177 { "-pin", FALSE, etENUM, {thread_aff_opt_choices},
178 "Whether mdrun should try to set thread affinities" },
179 { "-pinoffset", FALSE, etINT, {&hw_opt.core_pinning_offset},
180 "The lowest logical core number to which mdrun should pin the first thread" },
181 { "-pinstride", FALSE, etINT, {&hw_opt.core_pinning_stride},
182 "Pinning distance in logical cores for threads, use 0 to minimize the number of threads per physical core" },
183 { "-gpu_id", FALSE, etSTR, {&gpuIdsAvailable},
184 "List of unique GPU device IDs available to use" },
185 { "-gputasks", FALSE, etSTR, {&userGpuTaskAssignment},
186 "List of GPU device IDs, mapping each PP task on each node to a device" },
187 { "-ddcheck", FALSE, etBOOL, {&domdecOptions.checkBondedInteractions},
188 "Check for all bonded interactions with DD" },
189 { "-ddbondcomm", FALSE, etBOOL, {&domdecOptions.useBondedCommunication},
190 "HIDDENUse special bonded atom communication when [TT]-rdd[tt] > cut-off" },
191 { "-rdd", FALSE, etREAL, {&domdecOptions.minimumCommunicationRange},
192 "The maximum distance for bonded interactions with DD (nm), 0 is determine from initial coordinates" },
193 { "-rcon", FALSE, etREAL, {&domdecOptions.constraintCommunicationRange},
194 "Maximum distance for P-LINCS (nm), 0 is estimate" },
195 { "-dlb", FALSE, etENUM, {dddlb_opt_choices},
196 "Dynamic load balancing (with DD)" },
197 { "-dds", FALSE, etREAL, {&domdecOptions.dlbScaling},
198 "Fraction in (0,1) by whose reciprocal the initial DD cell size will be increased in order to "
199 "provide a margin in which dynamic load balancing can act while preserving the minimum cell size." },
200 { "-ddcsx", FALSE, etSTR, {&domdecOptions.cellSizeX},
201 "HIDDENA string containing a vector of the relative sizes in the x "
202 "direction of the corresponding DD cells. Only effective with static "
203 "load balancing." },
204 { "-ddcsy", FALSE, etSTR, {&domdecOptions.cellSizeY},
205 "HIDDENA string containing a vector of the relative sizes in the y "
206 "direction of the corresponding DD cells. Only effective with static "
207 "load balancing." },
208 { "-ddcsz", FALSE, etSTR, {&domdecOptions.cellSizeZ},
209 "HIDDENA string containing a vector of the relative sizes in the z "
210 "direction of the corresponding DD cells. Only effective with static "
211 "load balancing." },
212 { "-gcom", FALSE, etINT, {&mdrunOptions.globalCommunicationInterval},
213 "Global communication frequency" },
214 { "-nb", FALSE, etENUM, {nbpu_opt_choices},
215 "Calculate non-bonded interactions on" },
216 { "-nstlist", FALSE, etINT, {&nstlist_cmdline},
217 "Set nstlist when using a Verlet buffer tolerance (0 is guess)" },
218 { "-tunepme", FALSE, etBOOL, {&mdrunOptions.tunePme},
219 "Optimize PME load between PP/PME ranks or GPU/CPU (only with the Verlet cut-off scheme)" },
220 { "-pme", FALSE, etENUM, {pme_opt_choices},
221 "Perform PME calculations on" },
222 { "-pmefft", FALSE, etENUM, {pme_fft_opt_choices},
223 "Perform PME FFT calculations on" },
224 { "-v", FALSE, etBOOL, {&mdrunOptions.verbose},
225 "Be loud and noisy" },
226 { "-pforce", FALSE, etREAL, {&pforce},
227 "Print all forces larger than this (kJ/mol nm)" },
228 { "-reprod", FALSE, etBOOL, {&mdrunOptions.reproducible},
229 "Try to avoid optimizations that affect binary reproducibility" },
230 { "-cpt", FALSE, etREAL, {&mdrunOptions.checkpointOptions.period},
231 "Checkpoint interval (minutes)" },
232 { "-cpnum", FALSE, etBOOL, {&mdrunOptions.checkpointOptions.keepAndNumberCheckpointFiles},
233 "Keep and number checkpoint files" },
234 { "-append", FALSE, etBOOL, {&bTryToAppendFiles},
235 "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names" },
236 { "-nsteps", FALSE, etINT64, {&mdrunOptions.numStepsCommandline},
237 "Run this number of steps, overrides .mdp file option (-1 means infinite, -2 means use mdp option, smaller is invalid)" },
238 { "-maxh", FALSE, etREAL, {&mdrunOptions.maximumHoursToRun},
239 "Terminate after 0.99 times this time (hours)" },
240 { "-replex", FALSE, etINT, {&replExParams.exchangeInterval},
241 "Attempt replica exchange periodically with this period (steps)" },
242 { "-nex", FALSE, etINT, {&replExParams.numExchanges},
243 "Number of random exchanges to carry out each exchange interval (N^3 is one suggestion). -nex zero or not specified gives neighbor replica exchange." },
244 { "-reseed", FALSE, etINT, {&replExParams.randomSeed},
245 "Seed for replica exchange, -1 is generate a seed" },
246 { "-imdport", FALSE, etINT, {&imdOptions.port},
247 "HIDDENIMD listening port" },
248 { "-imdwait", FALSE, etBOOL, {&imdOptions.wait},
249 "HIDDENPause the simulation while no IMD client is connected" },
250 { "-imdterm", FALSE, etBOOL, {&imdOptions.terminatable},
251 "HIDDENAllow termination of the simulation from IMD client" },
252 { "-imdpull", FALSE, etBOOL, {&imdOptions.pull},
253 "HIDDENAllow pulling in the simulation from IMD client" },
254 { "-rerunvsite", FALSE, etBOOL, {&mdrunOptions.rerunConstructVsites},
255 "HIDDENRecalculate virtual site coordinates with [TT]-rerun[tt]" },
256 { "-confout", FALSE, etBOOL, {&mdrunOptions.writeConfout},
257 "HIDDENWrite the last configuration with [TT]-c[tt] and force checkpointing at the last step" },
258 { "-stepout", FALSE, etINT, {&mdrunOptions.verboseStepPrintInterval},
259 "HIDDENFrequency of writing the remaining wall clock time for the run" },
260 { "-resetstep", FALSE, etINT, {&mdrunOptions.timingOptions.resetStep},
261 "HIDDENReset cycle counters after these many time steps" },
262 { "-resethway", FALSE, etBOOL, {&mdrunOptions.timingOptions.resetHalfway},
263 "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt]" }
265 /*! \} */
267 //! Handle to communication object.
268 t_commrec *cr = nullptr;
269 //! Multi-simulation object.
270 gmx_multisim_t *ms = nullptr;
272 //! Parses the command-line input and prepares to start mdrun.
273 int updateFromCommandLine(int argc, char **argv, ArrayRef<const char *> desc);
275 ~LegacyMdrunOptions();
278 } // end namespace gmx