95abef6a57f49e4bcdc9e303cf19ce79ded5701d
[alexxy/gromacs.git] / src / gromacs / mdrun / legacymdrunoptions.h
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) 2011-2020,2021, 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 /*! \libinternal \file
38  *
39  * \brief This file declares helper functionality for legacy option handling for mdrun
40  *
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>
45  *
46  * \ingroup module_mdrun
47  * \inlibraryapi
48  */
49 #ifndef GMX_MDRUN_LEGACYMDRUNOPTIONS_H
50 #define GMX_MDRUN_LEGACYMDRUNOPTIONS_H
51
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"
57
58 #include "replicaexchange.h"
59
60 namespace gmx
61 {
62
63 /*! \libinternal
64  * \brief This class provides the same command-line option
65  * functionality to both CLI and API sessions.
66  *
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.
70  *
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.
75  *
76  * \todo Modules in mdrun should acquire proper option handling so
77  *       that all of these declarations and defaults are local to the
78  *       modules.
79  *
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.
83  */
84 class LegacyMdrunOptions
85 {
86 public:
87     //! Ongoing collection of mdrun options
88     MdrunOptions mdrunOptions;
89     //! Options for the domain decomposition.
90     DomdecOptions domdecOptions;
91     //! Parallelism-related user options.
92     gmx_hw_opt_t hw_opt;
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;
97
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 } } };
131
132     //! Print a warning if any force is larger than this (in kJ/mol nm).
133     real pforce = -1;
134
135     //! The value of the -append option
136     bool appendOption = true;
137
138     /*! \brief Output context for writing text files
139      *
140      * \todo Clarify initialization, ownership, and lifetime. */
141     gmx_output_env_t* oenv = nullptr;
142
143     /*! \brief Command line options, defaults, docs and storage for them to fill. */
144     /*! \{ */
145     rvec        realddxyz                                                           = { 0, 0, 0 };
146     const char* ddrank_opt_choices[static_cast<int>(DdRankOrder::Count) + 1]        = { nullptr,
147                                                                                  "interleave",
148                                                                                  "pp_pme",
149                                                                                  "cartesian",
150                                                                                  nullptr };
151     const char* dddlb_opt_choices[static_cast<int>(DlbOption::Count) + 1]           = { nullptr,
152                                                                               "auto",
153                                                                               "no",
154                                                                               "yes",
155                                                                               nullptr };
156     const char* thread_aff_opt_choices[static_cast<int>(ThreadAffinity::Count) + 1] = { nullptr,
157                                                                                         "auto",
158                                                                                         "on",
159                                                                                         "off",
160                                                                                         nullptr };
161     const char* nbpu_opt_choices[5]    = { nullptr, "auto", "cpu", "gpu", nullptr };
162     const char* pme_opt_choices[5]     = { nullptr, "auto", "cpu", "gpu", nullptr };
163     const char* pme_fft_opt_choices[5] = { nullptr, "auto", "cpu", "gpu", nullptr };
164     const char* bonded_opt_choices[5]  = { nullptr, "auto", "cpu", "gpu", nullptr };
165     const char* update_opt_choices[5]  = { nullptr, "auto", "cpu", "gpu", nullptr };
166     const char* devicesSelectedByUser  = "";
167     const char* userGpuTaskAssignment  = "";
168
169
170     ImdOptions& imdOptions = mdrunOptions.imdOptions;
171
172     t_pargs pa[48] = {
173
174         { "-dd", FALSE, etRVEC, { &realddxyz }, "Domain decomposition grid, 0 is optimize" },
175         { "-ddorder", FALSE, etENUM, { ddrank_opt_choices }, "DD rank order" },
176         { "-npme",
177           FALSE,
178           etINT,
179           { &domdecOptions.numPmeRanks },
180           "Number of separate ranks to be used for PME, -1 is guess" },
181         { "-nt",
182           FALSE,
183           etINT,
184           { &hw_opt.nthreads_tot },
185           "Total number of threads to start (0 is guess)" },
186         { "-ntmpi",
187           FALSE,
188           etINT,
189           { &hw_opt.nthreads_tmpi },
190           "Number of thread-MPI ranks to start (0 is guess)" },
191         { "-ntomp",
192           FALSE,
193           etINT,
194           { &hw_opt.nthreads_omp },
195           "Number of OpenMP threads per MPI rank to start (0 is guess)" },
196         { "-ntomp_pme",
197           FALSE,
198           etINT,
199           { &hw_opt.nthreads_omp_pme },
200           "Number of OpenMP threads per MPI rank to start (0 is -ntomp)" },
201         { "-pin",
202           FALSE,
203           etENUM,
204           { thread_aff_opt_choices },
205           "Whether mdrun should try to set thread affinities" },
206         { "-pinoffset",
207           FALSE,
208           etINT,
209           { &hw_opt.core_pinning_offset },
210           "The lowest logical core number to which mdrun should pin the first thread" },
211         { "-pinstride",
212           FALSE,
213           etINT,
214           { &hw_opt.core_pinning_stride },
215           "Pinning distance in logical cores for threads, use 0 to minimize the number of threads "
216           "per physical core" },
217         { "-gpu_id",
218           FALSE,
219           etSTR,
220           { &devicesSelectedByUser },
221           "List of unique GPU device IDs available to use" },
222         { "-gputasks",
223           FALSE,
224           etSTR,
225           { &userGpuTaskAssignment },
226           "List of GPU device IDs, mapping each PP task on each node to a device" },
227         { "-ddcheck",
228           FALSE,
229           etBOOL,
230           { &domdecOptions.checkBondedInteractions },
231           "Check for all bonded interactions with DD" },
232         { "-ddbondcomm",
233           FALSE,
234           etBOOL,
235           { &domdecOptions.useBondedCommunication },
236           "HIDDENUse special bonded atom communication when [TT]-rdd[tt] > cut-off" },
237         { "-rdd",
238           FALSE,
239           etREAL,
240           { &domdecOptions.minimumCommunicationRange },
241           "The maximum distance for bonded interactions with DD (nm), 0 is determine from initial "
242           "coordinates" },
243         { "-rcon",
244           FALSE,
245           etREAL,
246           { &domdecOptions.constraintCommunicationRange },
247           "Maximum distance for P-LINCS (nm), 0 is estimate" },
248         { "-dlb", FALSE, etENUM, { dddlb_opt_choices }, "Dynamic load balancing (with DD)" },
249         { "-dds",
250           FALSE,
251           etREAL,
252           { &domdecOptions.dlbScaling },
253           "Fraction in (0,1) by whose reciprocal the initial DD cell size will be increased in "
254           "order to "
255           "provide a margin in which dynamic load balancing can act while preserving the minimum "
256           "cell size." },
257         { "-ddcsx",
258           FALSE,
259           etSTR,
260           { &domdecOptions.cellSizeX },
261           "HIDDENA string containing a vector of the relative sizes in the x "
262           "direction of the corresponding DD cells. Only effective with static "
263           "load balancing." },
264         { "-ddcsy",
265           FALSE,
266           etSTR,
267           { &domdecOptions.cellSizeY },
268           "HIDDENA string containing a vector of the relative sizes in the y "
269           "direction of the corresponding DD cells. Only effective with static "
270           "load balancing." },
271         { "-ddcsz",
272           FALSE,
273           etSTR,
274           { &domdecOptions.cellSizeZ },
275           "HIDDENA string containing a vector of the relative sizes in the z "
276           "direction of the corresponding DD cells. Only effective with static "
277           "load balancing." },
278         { "-nb", FALSE, etENUM, { nbpu_opt_choices }, "Calculate non-bonded interactions on" },
279         { "-nstlist",
280           FALSE,
281           etINT,
282           { &nstlist_cmdline },
283           "Set nstlist when using a Verlet buffer tolerance (0 is guess)" },
284         { "-tunepme",
285           FALSE,
286           etBOOL,
287           { &mdrunOptions.tunePme },
288           "Optimize PME load between PP/PME ranks or GPU/CPU" },
289         { "-pme", FALSE, etENUM, { pme_opt_choices }, "Perform PME calculations on" },
290         { "-pmefft", FALSE, etENUM, { pme_fft_opt_choices }, "Perform PME FFT calculations on" },
291         { "-bonded", FALSE, etENUM, { bonded_opt_choices }, "Perform bonded calculations on" },
292         { "-update", FALSE, etENUM, { update_opt_choices }, "Perform update and constraints on" },
293         { "-v", FALSE, etBOOL, { &mdrunOptions.verbose }, "Be loud and noisy" },
294         { "-pforce", FALSE, etREAL, { &pforce }, "Print all forces larger than this (kJ/mol nm)" },
295         { "-reprod",
296           FALSE,
297           etBOOL,
298           { &mdrunOptions.reproducible },
299           "Try to avoid optimizations that affect binary reproducibility" },
300         { "-cpt",
301           FALSE,
302           etREAL,
303           { &mdrunOptions.checkpointOptions.period },
304           "Checkpoint interval (minutes)" },
305         { "-cpnum",
306           FALSE,
307           etBOOL,
308           { &mdrunOptions.checkpointOptions.keepAndNumberCheckpointFiles },
309           "Keep and number checkpoint files" },
310         { "-append",
311           FALSE,
312           etBOOL,
313           { &appendOption },
314           "Append to previous output files when continuing from checkpoint instead of adding the "
315           "simulation part number to all file names" },
316         { "-nsteps",
317           FALSE,
318           etINT64,
319           { &mdrunOptions.numStepsCommandline },
320           "Run this number of steps (-1 means infinite, -2 means use mdp option, smaller is "
321           "invalid)" },
322         { "-maxh",
323           FALSE,
324           etREAL,
325           { &mdrunOptions.maximumHoursToRun },
326           "Terminate after 0.99 times this time (hours)" },
327         { "-replex",
328           FALSE,
329           etINT,
330           { &replExParams.exchangeInterval },
331           "Attempt replica exchange periodically with this period (steps)" },
332         { "-nex",
333           FALSE,
334           etINT,
335           { &replExParams.numExchanges },
336           "Number of random exchanges to carry out each exchange interval (N^3 is one suggestion). "
337           " -nex zero or not specified gives neighbor replica exchange." },
338         { "-reseed",
339           FALSE,
340           etINT,
341           { &replExParams.randomSeed },
342           "Seed for replica exchange, -1 is generate a seed" },
343         { "-imdport", FALSE, etINT, { &imdOptions.port }, "HIDDENIMD listening port" },
344         { "-imdwait",
345           FALSE,
346           etBOOL,
347           { &imdOptions.wait },
348           "HIDDENPause the simulation while no IMD client is connected" },
349         { "-imdterm",
350           FALSE,
351           etBOOL,
352           { &imdOptions.terminatable },
353           "HIDDENAllow termination of the simulation from IMD client" },
354         { "-imdpull",
355           FALSE,
356           etBOOL,
357           { &imdOptions.pull },
358           "HIDDENAllow pulling in the simulation from IMD client" },
359         { "-rerunvsite",
360           FALSE,
361           etBOOL,
362           { &mdrunOptions.rerunConstructVsites },
363           "HIDDENRecalculate virtual site coordinates with [TT]-rerun[tt]" },
364         { "-confout",
365           FALSE,
366           etBOOL,
367           { &mdrunOptions.writeConfout },
368           "HIDDENWrite the last configuration with [TT]-c[tt] and force checkpointing at the last "
369           "step" },
370         { "-stepout",
371           FALSE,
372           etINT,
373           { &mdrunOptions.verboseStepPrintInterval },
374           "HIDDENFrequency of writing the remaining wall clock time for the run" },
375         { "-resetstep",
376           FALSE,
377           etINT,
378           { &mdrunOptions.timingOptions.resetStep },
379           "HIDDENReset cycle counters after these many time steps" },
380         { "-resethway",
381           FALSE,
382           etBOOL,
383           { &mdrunOptions.timingOptions.resetHalfway },
384           "HIDDENReset the cycle counters after half the number of steps or halfway "
385           "[TT]-maxh[tt]" }
386     };
387     /*! \} */
388
389     //! Parses the command-line input and prepares to start mdrun.
390     int updateFromCommandLine(int argc, char** argv, ArrayRef<const char*> desc);
391
392     ~LegacyMdrunOptions();
393 };
394
395 } // end namespace gmx
396
397 #endif