e86c5e0534ef487a8ed29bd2a16434e22465cedf
[alexxy/gromacs.git] / src / programs / legacymodules.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \brief
36  * Registers command-line modules for pre-5.0 binaries.
37  *
38  * \author Teemu Murtola <teemu.murtola@gmail.com>
39  */
40 #include "gmxpre.h"
41
42 #include "legacymodules.h"
43
44 #include <cstdio>
45
46 #include "gromacs/commandline/cmdlinemodule.h"
47 #include "gromacs/commandline/cmdlinemodulemanager.h"
48 #include "gromacs/commandline/cmdlineoptionsmodule.h"
49 #include "gromacs/gmxana/gmx_ana.h"
50 #include "gromacs/gmxpreprocess/editconf.h"
51 #include "gromacs/gmxpreprocess/genconf.h"
52 #include "gromacs/gmxpreprocess/genion.h"
53 #include "gromacs/gmxpreprocess/genrestr.h"
54 #include "gromacs/gmxpreprocess/grompp.h"
55 #include "gromacs/gmxpreprocess/insert_molecules.h"
56 #include "gromacs/gmxpreprocess/pdb2gmx.h"
57 #include "gromacs/gmxpreprocess/solvate.h"
58 #include "gromacs/gmxpreprocess/x2top.h"
59 #include "gromacs/tools/check.h"
60 #include "gromacs/tools/convert_tpr.h"
61 #include "gromacs/tools/dump.h"
62 #include "gromacs/tools/eneconv.h"
63 #include "gromacs/tools/make_ndx.h"
64 #include "gromacs/tools/mk_angndx.h"
65 #include "gromacs/tools/pme_error.h"
66 #include "gromacs/tools/report_methods.h"
67 #include "gromacs/tools/trjcat.h"
68 #include "gromacs/tools/trjconv.h"
69 #include "gromacs/tools/tune_pme.h"
70
71 #include "mdrun/mdrun_main.h"
72 #include "mdrun/nonbonded_bench.h"
73 #include "view/view.h"
74
75 namespace
76 {
77
78 /*! \brief
79  * Command line module that provides information about obsolescence.
80  *
81  * Prints a message directing the user to a wiki page describing replacement
82  * options.
83  */
84 class ObsoleteToolModule : public gmx::ICommandLineModule
85 {
86 public:
87     //! Creates an obsolete tool module for a tool with the given name.
88     explicit ObsoleteToolModule(const char* name) : name_(name) {}
89
90     const char* name() const override { return name_; }
91     const char* shortDescription() const override { return nullptr; }
92
93     void init(gmx::CommandLineModuleSettings* /*settings*/) override {}
94     int  run(int /*argc*/, char* /*argv*/[]) override
95     {
96         printMessage();
97         return 0;
98     }
99     void writeHelp(const gmx::CommandLineHelpContext& /*context*/) const override
100     {
101         printMessage();
102     }
103
104 private:
105     void printMessage() const
106     {
107         std::fprintf(stderr,
108                      "This tool is no longer present in GROMACS. Please see\n"
109                      "  "
110                      "http://jenkins.gromacs.org/job/Documentation_Nightly_master/javadoc/"
111                      "user-guide/cmdline.html#command-changes\n"
112                      "for ideas how to perform the same tasks with the "
113                      "new tools.\n");
114     }
115
116     const char* name_;
117 };
118
119 //! Initializer for a module that defaults to nice level zero.
120 void initSettingsNoNice(gmx::CommandLineModuleSettings* settings)
121 {
122     settings->setDefaultNiceLevel(0);
123 }
124
125 /*! \brief
126  * Convenience function for creating and registering a module.
127  *
128  * \param[in] manager          Module manager to which to register the module.
129  * \param[in] mainFunction     Main function to wrap.
130  * \param[in] name             Name for the new module.
131  * \param[in] shortDescription One-line description for the new module.
132  */
133 void registerModule(gmx::CommandLineModuleManager*               manager,
134                     gmx::CommandLineModuleManager::CMainFunction mainFunction,
135                     const char*                                  name,
136                     const char*                                  shortDescription)
137 {
138     manager->addModuleCMain(name, shortDescription, mainFunction);
139 }
140
141 /*! \brief
142  * Convenience function for creating and registering a module that defaults to
143  * -nice 0.
144  *
145  * \param[in] manager          Module manager to which to register the module.
146  * \param[in] mainFunction     Main function to wrap.
147  * \param[in] name             Name for the new module.
148  * \param[in] shortDescription One-line description for the new module.
149  */
150 void registerModuleNoNice(gmx::CommandLineModuleManager*               manager,
151                           gmx::CommandLineModuleManager::CMainFunction mainFunction,
152                           const char*                                  name,
153                           const char*                                  shortDescription)
154 {
155     manager->addModuleCMainWithSettings(name, shortDescription, mainFunction, &initSettingsNoNice);
156 }
157
158 /*! \brief
159  * Convenience function for registering a module for an obsolete tool.
160  *
161  * \param[in] manager          Module manager to which to register the module.
162  * \param[in] name             Name for the obsolete tool.
163  */
164 void registerObsoleteTool(gmx::CommandLineModuleManager* manager, const char* name)
165 {
166     gmx::CommandLineModulePointer module(new ObsoleteToolModule(name));
167     manager->addModule(std::move(module));
168 }
169
170 } // namespace
171
172 void registerLegacyModules(gmx::CommandLineModuleManager* manager)
173 {
174     registerModule(manager, &gmx_check, "check", "Check and compare files");
175     gmx::ICommandLineOptionsModule::registerModuleFactory(
176             manager, gmx::DumpInfo::name, gmx::DumpInfo::shortDescription, &gmx::DumpInfo::create);
177     registerModule(manager, &gmx_grompp, "grompp", "Make a run input file");
178     registerModule(manager, &gmx_convert_tpr, "convert-tpr", "Make a modifed run-input file");
179     registerObsoleteTool(manager, "tpbconv");
180     registerModule(manager, &gmx_x2top, "x2top", "Generate a primitive topology from coordinates");
181
182     registerModuleNoNice(
183             manager, &gmx::gmx_mdrun, "mdrun",
184             "Perform a simulation, do a normal mode analysis or an energy minimization");
185
186     gmx::ICommandLineOptionsModule::registerModuleFactory(
187             manager, gmx::NonbondedBenchmarkInfo::name,
188             gmx::NonbondedBenchmarkInfo::shortDescription, &gmx::NonbondedBenchmarkInfo::create);
189
190     gmx::ICommandLineOptionsModule::registerModuleFactory(manager, gmx::InsertMoleculesInfo::name(),
191                                                           gmx::InsertMoleculesInfo::shortDescription(),
192                                                           &gmx::InsertMoleculesInfo::create);
193
194     gmx::ICommandLineOptionsModule::registerModuleFactory(manager, gmx::ReportMethodsInfo::name,
195                                                           gmx::ReportMethodsInfo::shortDescription,
196                                                           &gmx::ReportMethodsInfo::create);
197
198     gmx::ICommandLineOptionsModule::registerModuleFactory(manager, gmx::pdb2gmxInfo::name,
199                                                           gmx::pdb2gmxInfo::shortDescription,
200                                                           &gmx::pdb2gmxInfo::create);
201
202     // Modules from gmx_ana.h.
203     registerModule(manager, &gmx_do_dssp, "do_dssp",
204                    "Assign secondary structure and calculate solvent accessible surface area");
205     registerModule(manager, &gmx_editconf, "editconf", "Convert and manipulates structure files");
206     registerModule(manager, &gmx_eneconv, "eneconv", "Convert energy files");
207     registerModule(manager, &gmx_solvate, "solvate", "Solvate a system");
208     registerObsoleteTool(manager, "genbox");
209     registerModule(manager, &gmx_genconf, "genconf",
210                    "Multiply a conformation in 'random' orientations");
211     registerModule(manager, &gmx_genion, "genion",
212                    "Generate monoatomic ions on energetically favorable positions");
213     registerModule(manager, &gmx_genrestr, "genrestr",
214                    "Generate position restraints or distance restraints for index groups");
215     registerModule(manager, &gmx_make_edi, "make_edi",
216                    "Generate input files for essential dynamics sampling");
217     registerModule(manager, &gmx_make_ndx, "make_ndx", "Make index files");
218     registerModule(manager, &gmx_mk_angndx, "mk_angndx", "Generate index files for 'gmx angle'");
219     registerModule(manager, &gmx_trjcat, "trjcat", "Concatenate trajectory files");
220     registerModule(manager, &gmx_trjconv, "trjconv", "Convert and manipulates trajectory files");
221     registerModule(manager, &gmx_trjorder, "trjorder",
222                    "Order molecules according to their distance to a group");
223     registerModule(manager, &gmx_xpm2ps, "xpm2ps",
224                    "Convert XPM (XPixelMap) matrices to postscript or XPM");
225
226     registerModule(manager, &gmx_anaeig, "anaeig", "Analyze eigenvectors/normal modes");
227     registerModule(manager, &gmx_analyze, "analyze", "Analyze data sets");
228     registerModule(manager, &gmx_g_angle, "angle",
229                    "Calculate distributions and correlations for angles and dihedrals");
230     registerModule(manager, &gmx_awh, "awh",
231                    "Extract data from an accelerated weight histogram (AWH) run");
232     registerModule(manager, &gmx_bar, "bar",
233                    "Calculate free energy difference estimates through Bennett's acceptance ratio");
234     registerObsoleteTool(manager, "bond");
235     registerObsoleteTool(manager, "dist");
236     registerObsoleteTool(manager, "sas");
237     registerObsoleteTool(manager, "sgangle");
238
239     registerModule(manager, &gmx_bundle, "bundle", "Analyze bundles of axes, e.g., helices");
240     registerModule(manager, &gmx_chi, "chi",
241                    "Calculate everything you want to know about chi and other dihedrals");
242     registerModule(manager, &gmx_cluster, "cluster", "Cluster structures");
243     registerModule(manager, &gmx_clustsize, "clustsize",
244                    "Calculate size distributions of atomic clusters");
245     registerModule(manager, &gmx_confrms, "confrms", "Fit two structures and calculates the RMSD");
246     registerModule(manager, &gmx_covar, "covar", "Calculate and diagonalize the covariance matrix");
247     registerModule(manager, &gmx_current, "current",
248                    "Calculate dielectric constants and current autocorrelation function");
249     registerModule(manager, &gmx_density, "density", "Calculate the density of the system");
250     registerModule(manager, &gmx_densmap, "densmap",
251                    "Calculate 2D planar or axial-radial density maps");
252     registerModule(manager, &gmx_densorder, "densorder", "Calculate surface fluctuations");
253     registerModule(manager, &gmx_dielectric, "dielectric",
254                    "Calculate frequency dependent dielectric constants");
255     registerModule(manager, &gmx_dipoles, "dipoles", "Compute the total dipole plus fluctuations");
256     registerModule(manager, &gmx_disre, "disre", "Analyze distance restraints");
257     registerModule(manager, &gmx_dos, "dos",
258                    "Analyze density of states and properties based on that");
259     registerModule(manager, &gmx_dyecoupl, "dyecoupl", "Extract dye dynamics from trajectories");
260     registerModule(manager, &gmx_enemat, "enemat", "Extract an energy matrix from an energy file");
261     registerModule(manager, &gmx_energy, "energy",
262                    "Writes energies to xvg files and display averages");
263     registerModule(manager, &gmx_filter, "filter",
264                    "Frequency filter trajectories, useful for making smooth movies");
265     registerModule(manager, &gmx_gyrate, "gyrate", "Calculate the radius of gyration");
266     registerModule(manager, &gmx_h2order, "h2order", "Compute the orientation of water molecules");
267     registerModule(manager, &gmx_hbond, "hbond", "Compute and analyze hydrogen bonds");
268     registerModule(manager, &gmx_helix, "helix", "Calculate basic properties of alpha helices");
269     registerModule(manager, &gmx_helixorient, "helixorient",
270                    "Calculate local pitch/bending/rotation/orientation inside helices");
271     registerModule(manager, &gmx_hydorder, "hydorder",
272                    "Compute tetrahedrality parameters around a given atom");
273     registerModule(manager, &gmx_lie, "lie", "Estimate free energy from linear combinations");
274     registerModule(manager, &gmx_mdmat, "mdmat", "Calculate residue contact maps");
275     registerModule(manager, &gmx_mindist, "mindist",
276                    "Calculate the minimum distance between two groups");
277     registerModule(manager, &gmx_msd, "msd", "Calculates mean square displacements");
278     registerModule(manager, &gmx_nmeig, "nmeig", "Diagonalize the Hessian for normal mode analysis");
279     registerModule(manager, &gmx_nmens, "nmens",
280                    "Generate an ensemble of structures from the normal modes");
281     registerModule(manager, &gmx_nmr, "nmr",
282                    "Analyze nuclear magnetic resonance properties from an energy file");
283     registerModule(manager, &gmx_nmtraj, "nmtraj",
284                    "Generate a virtual oscillating trajectory from an eigenvector");
285     registerModule(manager, &gmx_order, "order",
286                    "Compute the order parameter per atom for carbon tails");
287     registerModule(manager, &gmx_pme_error, "pme_error",
288                    "Estimate the error of using PME with a given input file");
289     registerModule(manager, &gmx_polystat, "polystat", "Calculate static properties of polymers");
290     registerModule(manager, &gmx_potential, "potential",
291                    "Calculate the electrostatic potential across the box");
292     registerModule(manager, &gmx_principal, "principal",
293                    "Calculate principal axes of inertia for a group of atoms");
294     registerModule(manager, &gmx_rama, "rama", "Compute Ramachandran plots");
295     registerModule(manager, &gmx_rms, "rms",
296                    "Calculate RMSDs with a reference structure and RMSD matrices");
297     registerModule(manager, &gmx_rmsdist, "rmsdist",
298                    "Calculate atom pair distances averaged with power -2, -3 or -6");
299     registerModule(manager, &gmx_rmsf, "rmsf", "Calculate atomic fluctuations");
300     registerModule(manager, &gmx_rotacf, "rotacf",
301                    "Calculate the rotational correlation function for molecules");
302     registerModule(manager, &gmx_rotmat, "rotmat",
303                    "Plot the rotation matrix for fitting to a reference structure");
304     registerModule(manager, &gmx_saltbr, "saltbr", "Compute salt bridges");
305     registerModule(manager, &gmx_sans, "sans", "Compute small angle neutron scattering spectra");
306     registerModule(manager, &gmx_saxs, "saxs", "Compute small angle X-ray scattering spectra");
307     registerModule(manager, &gmx_sham, "sham",
308                    "Compute free energies or other histograms from histograms");
309     registerModule(manager, &gmx_sigeps, "sigeps",
310                    "Convert c6/12 or c6/cn combinations to and from sigma/epsilon");
311     registerModule(manager, &gmx_sorient, "sorient", "Analyze solvent orientation around solutes");
312     registerModule(manager, &gmx_spatial, "spatial", "Calculate the spatial distribution function");
313     registerModule(manager, &gmx_spol, "spol",
314                    "Analyze solvent dipole orientation and polarization around solutes");
315     registerModule(manager, &gmx_tcaf, "tcaf", "Calculate viscosities of liquids");
316     registerModule(manager, &gmx_traj, "traj",
317                    "Plot x, v, f, box, temperature and rotational energy from trajectories");
318     registerModule(manager, &gmx_tune_pme, "tune_pme",
319                    "Time mdrun as a function of PME ranks to optimize settings");
320     registerModule(manager, &gmx_vanhove, "vanhove",
321                    "Compute Van Hove displacement and correlation functions");
322     registerModule(manager, &gmx_velacc, "velacc", "Calculate velocity autocorrelation functions");
323     registerModule(manager, &gmx_wham, "wham",
324                    "Perform weighted histogram analysis after umbrella sampling");
325     registerModule(manager, &gmx_wheel, "wheel", "Plot helical wheels");
326     registerModuleNoNice(manager, &gmx_view, "view", "View a trajectory on an X-Windows terminal");
327
328     {
329         gmx::CommandLineModuleGroup group =
330                 manager->addModuleGroup("Generating topologies and coordinates");
331         group.addModuleWithDescription("editconf", "Edit the box and write subgroups");
332         group.addModule("x2top");
333         group.addModule("solvate");
334         group.addModule("insert-molecules");
335         group.addModule("genconf");
336         group.addModule("genion");
337         group.addModule("genrestr");
338         group.addModule("pdb2gmx");
339     }
340     {
341         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Running a simulation");
342         group.addModule("grompp");
343         group.addModule("mdrun");
344         group.addModule("convert-tpr");
345     }
346     {
347         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Viewing trajectories");
348         group.addModule("nmtraj");
349         group.addModule("view");
350     }
351     {
352         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Processing energies");
353         group.addModule("enemat");
354         group.addModule("energy");
355         group.addModuleWithDescription("mdrun",
356                                        "(Re)calculate energies for trajectory frames with -rerun");
357     }
358     {
359         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Converting files");
360         group.addModule("editconf");
361         group.addModule("eneconv");
362         group.addModule("sigeps");
363         group.addModule("trjcat");
364         group.addModule("trjconv");
365         group.addModule("xpm2ps");
366     }
367     {
368         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Tools");
369         group.addModule("analyze");
370         group.addModule("awh");
371         group.addModule("filter");
372         group.addModule("lie");
373         group.addModule("pme_error");
374         group.addModule("sham");
375         group.addModule("spatial");
376         group.addModule("traj");
377         group.addModule("tune_pme");
378         group.addModule("wham");
379         group.addModule("check");
380         group.addModule("dump");
381         group.addModule("make_ndx");
382         group.addModule("mk_angndx");
383         group.addModule("trjorder");
384         group.addModule("xpm2ps");
385         group.addModule("report-methods");
386     }
387     {
388         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Distances between structures");
389         group.addModule("cluster");
390         group.addModule("confrms");
391         group.addModule("rms");
392         group.addModule("rmsf");
393     }
394     {
395         gmx::CommandLineModuleGroup group =
396                 manager->addModuleGroup("Distances in structures over time");
397         group.addModule("mindist");
398         group.addModule("mdmat");
399         group.addModule("polystat");
400         group.addModule("rmsdist");
401     }
402     {
403         gmx::CommandLineModuleGroup group =
404                 manager->addModuleGroup("Mass distribution properties over time");
405         group.addModule("gyrate");
406         group.addModule("msd");
407         group.addModule("polystat");
408         group.addModule("rdf");
409         group.addModule("rotacf");
410         group.addModule("rotmat");
411         group.addModule("sans");
412         group.addModule("saxs");
413         group.addModule("traj");
414         group.addModule("vanhove");
415     }
416     {
417         gmx::CommandLineModuleGroup group =
418                 manager->addModuleGroup("Analyzing bonded interactions");
419         group.addModule("angle");
420         group.addModule("mk_angndx");
421     }
422     {
423         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Structural properties");
424         group.addModule("bundle");
425         group.addModule("clustsize");
426         group.addModule("disre");
427         group.addModule("hbond");
428         group.addModule("order");
429         group.addModule("principal");
430         group.addModule("rdf");
431         group.addModule("saltbr");
432         group.addModule("sorient");
433         group.addModule("spol");
434     }
435     {
436         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Kinetic properties");
437         group.addModule("bar");
438         group.addModule("current");
439         group.addModule("dos");
440         group.addModule("dyecoupl");
441         group.addModule("principal");
442         group.addModule("tcaf");
443         group.addModule("traj");
444         group.addModule("vanhove");
445         group.addModule("velacc");
446     }
447     {
448         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Electrostatic properties");
449         group.addModule("current");
450         group.addModule("dielectric");
451         group.addModule("dipoles");
452         group.addModule("potential");
453         group.addModule("spol");
454         group.addModule("genion");
455     }
456     {
457         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Protein-specific analysis");
458         group.addModule("do_dssp");
459         group.addModule("chi");
460         group.addModule("helix");
461         group.addModule("helixorient");
462         group.addModule("rama");
463         group.addModule("wheel");
464     }
465     {
466         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Interfaces");
467         group.addModule("bundle");
468         group.addModule("density");
469         group.addModule("densmap");
470         group.addModule("densorder");
471         group.addModule("h2order");
472         group.addModule("hydorder");
473         group.addModule("order");
474         group.addModule("potential");
475     }
476     {
477         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Covariance analysis");
478         group.addModuleWithDescription("anaeig", "Analyze the eigenvectors");
479         group.addModule("covar");
480         group.addModule("make_edi");
481     }
482     {
483         gmx::CommandLineModuleGroup group = manager->addModuleGroup("Normal modes");
484         group.addModuleWithDescription("anaeig", "Analyze the normal modes");
485         group.addModule("nmeig");
486         group.addModule("nmtraj");
487         group.addModule("nmens");
488         group.addModule("grompp");
489         group.addModuleWithDescription("mdrun",
490                                        "Find a potential energy minimum and calculate the Hessian");
491     }
492 }