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