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