f85fe269318ba8729cf720ad880d7eae035ff187
[alexxy/gromacs.git] / src / programs / mdrun / mdrun.cpp
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, 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 /*! \defgroup module_mdrun Implementation of mdrun
38  * \ingroup group_mdrun
39  *
40  * \brief This module contains code that implements mdrun.
41  */
42 /*! \internal \file
43  *
44  * \brief This file implements mdrun
45  *
46  * \author Berk Hess <hess@kth.se>
47  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
48  * \author Erik Lindahl <erik@kth.se>
49  * \author Mark Abraham <mark.j.abraham@gmail.com>
50  *
51  * \ingroup module_mdrun
52  */
53 #include "gmxpre.h"
54
55 #include "config.h"
56
57 #include <memory>
58
59 #include "gromacs/commandline/pargs.h"
60 #include "gromacs/domdec/options.h"
61 #include "gromacs/fileio/gmxfio.h"
62 #include "gromacs/mdrun/legacymdrunoptions.h"
63 #include "gromacs/mdrun/runner.h"
64 #include "gromacs/mdrun/simulationcontext.h"
65 #include "gromacs/mdrunutility/handlerestart.h"
66 #include "gromacs/mdrunutility/logging.h"
67 #include "gromacs/mdrunutility/multisim.h"
68 #include "gromacs/utility/arrayref.h"
69 #include "gromacs/utility/smalloc.h"
70
71 #include "mdrun_main.h"
72
73 namespace gmx
74 {
75
76 //! Implements C-style main function for mdrun
77 int gmx_mdrun(int argc, char* argv[])
78 {
79     auto mdModules = std::make_unique<MDModules>();
80
81     std::vector<const char*> desc = {
82         "[THISMODULE] is the main computational chemistry engine",
83         "within GROMACS. Obviously, it performs Molecular Dynamics simulations,",
84         "but it can also perform Stochastic Dynamics, Energy Minimization,",
85         "test particle insertion or (re)calculation of energies.",
86         "Normal mode analysis is another option. In this case [TT]mdrun[tt]",
87         "builds a Hessian matrix from single conformation.",
88         "For usual Normal Modes-like calculations, make sure that",
89         "the structure provided is properly energy-minimized.",
90         "The generated matrix can be diagonalized by [gmx-nmeig].[PAR]",
91         "The [TT]mdrun[tt] program reads the run input file ([TT]-s[tt])",
92         "and distributes the topology over ranks if needed.",
93         "[TT]mdrun[tt] produces at least four output files.",
94         "A single log file ([TT]-g[tt]) is written.",
95         "The trajectory file ([TT]-o[tt]), contains coordinates, velocities and",
96         "optionally forces.",
97         "The structure file ([TT]-c[tt]) contains the coordinates and",
98         "velocities of the last step.",
99         "The energy file ([TT]-e[tt]) contains energies, the temperature,",
100         "pressure, etc, a lot of these things are also printed in the log file.",
101         "Optionally coordinates can be written to a compressed trajectory file",
102         "([TT]-x[tt]).[PAR]",
103         "The option [TT]-dhdl[tt] is only used when free energy calculation is",
104         "turned on.[PAR]",
105         "Running mdrun efficiently in parallel is a complex topic,",
106         "many aspects of which are covered in the online User Guide. You",
107         "should look there for practical advice on using many of the options",
108         "available in mdrun.[PAR]",
109         "ED (essential dynamics) sampling and/or additional flooding potentials",
110         "are switched on by using the [TT]-ei[tt] flag followed by an [REF].edi[ref]",
111         "file. The [REF].edi[ref] file can be produced with the [TT]make_edi[tt] tool",
112         "or by using options in the essdyn menu of the WHAT IF program.",
113         "[TT]mdrun[tt] produces a [REF].xvg[ref] output file that",
114         "contains projections of positions, velocities and forces onto selected",
115         "eigenvectors.[PAR]",
116         "When user-defined potential functions have been selected in the",
117         "[REF].mdp[ref] file the [TT]-table[tt] option is used to pass [TT]mdrun[tt]",
118         "a formatted table with potential functions. The file is read from",
119         "either the current directory or from the [TT]GMXLIB[tt] directory.",
120         "A number of pre-formatted tables are presented in the [TT]GMXLIB[tt] dir,",
121         "for 6-8, 6-9, 6-10, 6-11, 6-12 Lennard-Jones potentials with",
122         "normal Coulomb.",
123         "When pair interactions are present, a separate table for pair interaction",
124         "functions is read using the [TT]-tablep[tt] option.[PAR]",
125         "When tabulated bonded functions are present in the topology,",
126         "interaction functions are read using the [TT]-tableb[tt] option.",
127         "For each different tabulated interaction type used, a table file name must",
128         "be given. For the topology to work, a file name given here must match a",
129         "character sequence before the file extension. That sequence is: an underscore,",
130         "then a 'b' for bonds, an 'a' for angles or a 'd' for dihedrals,",
131         "and finally the matching table number index used in the topology. Note that,",
132         "these options are deprecated, and in future will be available via grompp.[PAR]",
133         "The options [TT]-px[tt] and [TT]-pf[tt] are used for writing pull COM",
134         "coordinates and forces when pulling is selected",
135         "in the [REF].mdp[ref] file.",
136         "[PAR]",
137         "The option [TT]-membed[tt] does what used to be g_membed, i.e. embed",
138         "a protein into a membrane. This module requires a number of settings",
139         "that are provided in a data file that is the argument of this option.",
140         "For more details in membrane embedding, see the documentation in the",
141         "user guide. The options [TT]-mn[tt] and [TT]-mp[tt] are used to provide",
142         "the index and topology files used for the embedding.",
143         "[PAR]",
144         "The option [TT]-pforce[tt] is useful when you suspect a simulation",
145         "crashes due to too large forces. With this option coordinates and",
146         "forces of atoms with a force larger than a certain value will",
147         "be printed to stderr. It will also terminate the run when non-finite",
148         "forces are present.",
149         "[PAR]",
150         "Checkpoints containing the complete state of the system are written",
151         "at regular intervals (option [TT]-cpt[tt]) to the file [TT]-cpo[tt],",
152         "unless option [TT]-cpt[tt] is set to -1.",
153         "The previous checkpoint is backed up to [TT]state_prev.cpt[tt] to",
154         "make sure that a recent state of the system is always available,",
155         "even when the simulation is terminated while writing a checkpoint.",
156         "With [TT]-cpnum[tt] all checkpoint files are kept and appended",
157         "with the step number.",
158         "A simulation can be continued by reading the full state from file",
159         "with option [TT]-cpi[tt]. This option is intelligent in the way that",
160         "if no checkpoint file is found, GROMACS just assumes a normal run and",
161         "starts from the first step of the [REF].tpr[ref] file. By default the output",
162         "will be appending to the existing output files. The checkpoint file",
163         "contains checksums of all output files, such that you will never",
164         "loose data when some output files are modified, corrupt or removed.",
165         "There are three scenarios with [TT]-cpi[tt]:[PAR]",
166         "[TT]*[tt] no files with matching names are present: new output files are written[PAR]",
167         "[TT]*[tt] all files are present with names and checksums matching those stored",
168         "in the checkpoint file: files are appended[PAR]",
169         "[TT]*[tt] otherwise no files are modified and a fatal error is generated[PAR]",
170         "With [TT]-noappend[tt] new output files are opened and the simulation",
171         "part number is added to all output file names.",
172         "Note that in all cases the checkpoint file itself is not renamed",
173         "and will be overwritten, unless its name does not match",
174         "the [TT]-cpo[tt] option.",
175         "[PAR]",
176         "With checkpointing the output is appended to previously written",
177         "output files, unless [TT]-noappend[tt] is used or none of the previous",
178         "output files are present (except for the checkpoint file).",
179         "The integrity of the files to be appended is verified using checksums",
180         "which are stored in the checkpoint file. This ensures that output can",
181         "not be mixed up or corrupted due to file appending. When only some",
182         "of the previous output files are present, a fatal error is generated",
183         "and no old output files are modified and no new output files are opened.",
184         "The result with appending will be the same as from a single run.",
185         "The contents will be binary identical, unless you use a different number",
186         "of ranks or dynamic load balancing or the FFT library uses optimizations",
187         "through timing.",
188         "[PAR]",
189         "With option [TT]-maxh[tt] a simulation is terminated and a checkpoint",
190         "file is written at the first neighbor search step where the run time",
191         "exceeds [TT]-maxh[tt]\\*0.99 hours. This option is particularly useful in",
192         "combination with setting [TT]nsteps[tt] to -1 either in the mdp or using the",
193         "similarly named command line option (although the latter is deprecated).",
194         "This results in an infinite run,",
195         "terminated only when the time limit set by [TT]-maxh[tt] is reached (if any)",
196         "or upon receiving a signal.",
197         "[PAR]",
198         "Interactive molecular dynamics (IMD) can be activated by using at least one",
199         "of the three IMD switches: The [TT]-imdterm[tt] switch allows one to terminate",
200         "the simulation from the molecular viewer (e.g. VMD). With [TT]-imdwait[tt],",
201         "[TT]mdrun[tt] pauses whenever no IMD client is connected. Pulling from the",
202         "IMD remote can be turned on by [TT]-imdpull[tt].",
203         "The port [TT]mdrun[tt] listens to can be altered by [TT]-imdport[tt].The",
204         "file pointed to by [TT]-if[tt] contains atom indices and forces if IMD",
205         "pulling is used."
206     };
207
208     LegacyMdrunOptions options;
209
210     if (options.updateFromCommandLine(argc, argv, desc) == 0)
211     {
212         return 0;
213     }
214
215     ArrayRef<const std::string> multiSimDirectoryNames =
216             opt2fnsIfOptionSet("-multidir", ssize(options.filenames), options.filenames.data());
217
218     // Set up the communicator, where possible (see docs for
219     // SimulationContext).
220     MPI_Comm communicator = GMX_LIB_MPI ? MPI_COMM_WORLD : MPI_COMM_NULL;
221     // The SimulationContext is necessary with gmxapi so that
222     // resources owned by the client code can have suitable
223     // lifetime. The gmx wrapper binary uses the same infrastructure,
224     // but the lifetime is now trivially that of the invocation of the
225     // wrapper binary.
226     SimulationContext simulationContext(communicator, multiSimDirectoryNames);
227
228     StartingBehavior startingBehavior        = StartingBehavior::NewSimulation;
229     LogFilePtr       logFileGuard            = nullptr;
230     gmx_multisim_t*  ms                      = simulationContext.multiSimulation_.get();
231     std::tie(startingBehavior, logFileGuard) = handleRestart(
232             findIsSimulationMasterRank(ms, communicator), communicator, ms,
233             options.mdrunOptions.appendingBehavior, ssize(options.filenames), options.filenames.data());
234
235     /* The named components for the builder exposed here are descriptive of the
236      * state of mdrun at implementation and are not intended to be prescriptive
237      * of future design. (Note the ICommandLineOptions... framework used elsewhere.)
238      * The modules should ultimately take part in composing the Director code
239      * for an extensible Builder.
240      *
241      * In the near term, we assume that resources like domain decomposition and
242      * neighbor lists must be reinitialized between simulation segments.
243      * We would prefer to rebuild resources only as necessary, but we defer such
244      * details to future optimizations.
245      */
246     auto builder = MdrunnerBuilder(std::move(mdModules),
247                                    compat::not_null<SimulationContext*>(&simulationContext));
248     builder.addSimulationMethod(options.mdrunOptions, options.pforce, startingBehavior);
249     builder.addDomainDecomposition(options.domdecOptions);
250     // \todo pass by value
251     builder.addNonBonded(options.nbpu_opt_choices[0]);
252     // \todo pass by value
253     builder.addElectrostatics(options.pme_opt_choices[0], options.pme_fft_opt_choices[0]);
254     builder.addBondedTaskAssignment(options.bonded_opt_choices[0]);
255     builder.addUpdateTaskAssignment(options.update_opt_choices[0]);
256     builder.addNeighborList(options.nstlist_cmdline);
257     builder.addReplicaExchange(options.replExParams);
258     // Need to establish run-time values from various inputs to provide a resource handle to Mdrunner
259     builder.addHardwareOptions(options.hw_opt);
260     // \todo File names are parameters that should be managed modularly through further factoring.
261     builder.addFilenames(options.filenames);
262     builder.addInput(makeSimulationInput(options));
263     // Note: The gmx_output_env_t life time is not managed after the call to parse_common_args.
264     // \todo Implement lifetime management for gmx_output_env_t.
265     // \todo Output environment should be configured outside of Mdrunner and provided as a resource.
266     builder.addOutputEnvironment(options.oenv);
267     builder.addLogFile(logFileGuard.get());
268
269     auto runner = builder.build();
270
271     return runner.mdrunner();
272 }
273
274 } // namespace gmx