52eaafb075d8bfaf93d352c91d3fdf2ce597b4e1
[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, 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
49 #include "gromacs/gmxana/gmx_ana.h"
50 #include "gromacs/gmxpreprocess/genconf.h"
51 #include "gromacs/gmxpreprocess/grompp.h"
52 #include "gromacs/gmxpreprocess/insert-molecules.h"
53 #include "gromacs/gmxpreprocess/pdb2gmx.h"
54 #include "gromacs/gmxpreprocess/protonate.h"
55 #include "gromacs/gmxpreprocess/solvate.h"
56 #include "gromacs/gmxpreprocess/x2top.h"
57 #include "gromacs/tools/check.h"
58 #include "gromacs/tools/convert_tpr.h"
59 #include "gromacs/tools/dump.h"
60
61 #include "mdrun/mdrun_main.h"
62 #include "view/view.h"
63
64 namespace
65 {
66
67 /*! \brief
68  * Command line module that provides information about obsolescence.
69  *
70  * Prints a message directing the user to a wiki page describing replacement
71  * options.
72  */
73 class ObsoleteToolModule : public gmx::CommandLineModuleInterface
74 {
75     public:
76         //! Creates an obsolete tool module for a tool with the given name.
77         explicit ObsoleteToolModule(const char *name)
78             : name_(name)
79         {
80         }
81
82         virtual const char *name() const
83         {
84             return name_;
85         }
86         virtual const char *shortDescription() const
87         {
88             return NULL;
89         }
90
91         virtual void init(gmx::CommandLineModuleSettings * /*settings*/)
92         {
93         }
94         virtual int run(int /*argc*/, char * /*argv*/[])
95         {
96             printMessage();
97             return 0;
98         }
99         virtual void writeHelp(const gmx::CommandLineHelpContext & /*context*/) const
100         {
101             printMessage();
102         }
103
104     private:
105         void printMessage() const
106         {
107             std::fprintf(stderr,
108                          "This tool has been removed from Gromacs 5.0. Please see\n"
109                          "  http://www.gromacs.org/Documentation/How-tos/Tool_Changes_for_5.0\n"
110                          "for ideas how to perform the same tasks with the "
111                          "new tools.\n");
112         }
113
114         const char             *name_;
115 };
116
117 // TODO: Consider removing duplication with CMainCommandLineModule from
118 // cmdlinemodulemanager.cpp.
119 class NoNiceModule : public gmx::CommandLineModuleInterface
120 {
121     public:
122         //! \copydoc gmx::CommandLineModuleManager::CMainFunction
123         typedef gmx::CommandLineModuleManager::CMainFunction CMainFunction;
124
125         /*! \brief
126          * Creates a wrapper module for the given main function.
127          *
128          * \param[in] name             Name for the module.
129          * \param[in] shortDescription One-line description for the module.
130          * \param[in] mainFunction     Main function to wrap.
131          *
132          * Does not throw.
133          */
134         NoNiceModule(const char *name, const char *shortDescription,
135                      CMainFunction mainFunction)
136             : name_(name), shortDescription_(shortDescription),
137               mainFunction_(mainFunction)
138         {
139         }
140
141         virtual const char *name() const
142         {
143             return name_;
144         }
145         virtual const char *shortDescription() const
146         {
147             return shortDescription_;
148         }
149
150         virtual void init(gmx::CommandLineModuleSettings *settings)
151         {
152             settings->setDefaultNiceLevel(0);
153         }
154         virtual int run(int argc, char *argv[])
155         {
156             return mainFunction_(argc, argv);
157         }
158         virtual void writeHelp(const gmx::CommandLineHelpContext &context) const
159         {
160             writeCommandLineHelpCMain(context, name_, mainFunction_);
161         }
162
163     private:
164         const char             *name_;
165         const char             *shortDescription_;
166         CMainFunction           mainFunction_;
167 };
168
169 /*! \brief
170  * Convenience function for creating and registering a module.
171  *
172  * \param[in] manager          Module manager to which to register the module.
173  * \param[in] mainFunction     Main function to wrap.
174  * \param[in] name             Name for the new module.
175  * \param[in] shortDescription One-line description for the new module.
176  */
177 void registerModule(gmx::CommandLineModuleManager                *manager,
178                     gmx::CommandLineModuleManager::CMainFunction  mainFunction,
179                     const char *name, const char *shortDescription)
180 {
181     manager->addModuleCMain(name, shortDescription, mainFunction);
182 }
183
184 /*! \brief
185  * Convenience function for creating and registering a module that defaults to
186  * -nice 0.
187  *
188  * \param[in] manager          Module manager to which to register the module.
189  * \param[in] mainFunction     Main function to wrap.
190  * \param[in] name             Name for the new module.
191  * \param[in] shortDescription One-line description for the new module.
192  */
193 void registerModuleNoNice(gmx::CommandLineModuleManager                *manager,
194                           gmx::CommandLineModuleManager::CMainFunction  mainFunction,
195                           const char *name, const char *shortDescription)
196 {
197     gmx::CommandLineModulePointer module(
198             new NoNiceModule(name, shortDescription, mainFunction));
199     manager->addModule(move(module));
200 }
201
202 /*! \brief
203  * Convenience function for registering a module for an obsolete tool.
204  *
205  * \param[in] manager          Module manager to which to register the module.
206  * \param[in] name             Name for the obsolete tool.
207  */
208 void registerObsoleteTool(gmx::CommandLineModuleManager *manager,
209                           const char                    *name)
210 {
211     gmx::CommandLineModulePointer module(new ObsoleteToolModule(name));
212     manager->addModule(move(module));
213 }
214
215 } // namespace
216
217 void registerLegacyModules(gmx::CommandLineModuleManager *manager)
218 {
219     // Modules from this directory (were in src/kernel/).
220     registerModule(manager, &gmx_check, "check",
221                    "Check and compare files");
222     registerModule(manager, &gmx_dump, "dump",
223                    "Make binary files human readable");
224     registerModule(manager, &gmx_grompp, "grompp",
225                    "Make a run input file");
226     registerModule(manager, &gmx_pdb2gmx, "pdb2gmx",
227                    "Convert coordinate files to topology and FF-compliant coordinate files");
228     registerModule(manager, &gmx_convert_tpr, "convert-tpr",
229                    "Make a modifed run-input file");
230     registerObsoleteTool(manager, "tpbconv");
231
232     registerModule(manager, &gmx_protonate, "protonate",
233                    "Protonate structures");
234     registerModule(manager, &gmx_x2top, "x2top",
235                    "Generate a primitive topology from coordinates");
236
237     registerModuleNoNice(manager, &gmx_mdrun, "mdrun",
238                          "Perform a simulation, do a normal mode analysis or an energy minimization");
239
240     // Modules from gmx_ana.h.
241     registerModule(manager, &gmx_do_dssp, "do_dssp",
242                    "Assign secondary structure and calculate solvent accessible surface area");
243     registerModule(manager, &gmx_editconf, "editconf",
244                    "Convert and manipulates structure files");
245     registerModule(manager, &gmx_eneconv, "eneconv",
246                    "Convert energy files");
247     registerModule(manager, &gmx_solvate, "solvate",
248                    "Solvate a system");
249     registerModule(manager, &gmx_insert_molecules, "insert-molecules",
250                    "Insert molecules into existing vacancies");
251     registerObsoleteTool(manager, "genbox");
252     registerModule(manager, &gmx_genconf, "genconf",
253                    "Multiply a conformation in 'random' orientations");
254     registerModule(manager, &gmx_genion, "genion",
255                    "Generate monoatomic ions on energetically favorable positions");
256     registerModule(manager, &gmx_genpr, "genrestr",
257                    "Generate position restraints or distance restraints for index groups");
258     registerModule(manager, &gmx_make_edi, "make_edi",
259                    "Generate input files for essential dynamics sampling");
260     registerModule(manager, &gmx_make_ndx, "make_ndx",
261                    "Make index files");
262     registerModule(manager, &gmx_mk_angndx, "mk_angndx",
263                    "Generate index files for 'gmx angle'");
264     registerModule(manager, &gmx_trjcat, "trjcat",
265                    "Concatenate trajectory files");
266     registerModule(manager, &gmx_trjconv, "trjconv",
267                    "Convert and manipulates trajectory files");
268     registerModule(manager, &gmx_trjorder, "trjorder",
269                    "Order molecules according to their distance to a group");
270     registerModule(manager, &gmx_xpm2ps, "xpm2ps",
271                    "Convert XPM (XPixelMap) matrices to postscript or XPM");
272
273     registerModule(manager, &gmx_anadock, "anadock",
274                    "Cluster structures from Autodock runs");
275     registerModule(manager, &gmx_anaeig, "anaeig",
276                    "Analyze eigenvectors/normal modes");
277     registerModule(manager, &gmx_analyze, "analyze",
278                    "Analyze data sets");
279     registerModule(manager, &gmx_g_angle, "angle",
280                    "Calculate distributions and correlations for angles and dihedrals");
281     registerModule(manager, &gmx_bar, "bar",
282                    "Calculate free energy difference estimates through Bennett's acceptance ratio");
283     registerObsoleteTool(manager, "bond");
284     registerObsoleteTool(manager, "dist");
285     registerObsoleteTool(manager, "sas");
286     registerObsoleteTool(manager, "sgangle");
287
288     registerModule(manager, &gmx_bundle, "bundle",
289                    "Analyze bundles of axes, e.g., helices");
290     registerModule(manager, &gmx_chi, "chi",
291                    "Calculate everything you want to know about chi and other dihedrals");
292     registerModule(manager, &gmx_cluster, "cluster",
293                    "Cluster structures");
294     registerModule(manager, &gmx_clustsize, "clustsize",
295                    "Calculate size distributions of atomic clusters");
296     registerModule(manager, &gmx_confrms, "confrms",
297                    "Fit two structures and calculates the RMSD");
298     registerModule(manager, &gmx_covar, "covar",
299                    "Calculate and diagonalize the covariance matrix");
300     registerModule(manager, &gmx_current, "current",
301                    "Calculate dielectric constants and current autocorrelation function");
302     registerModule(manager, &gmx_density, "density",
303                    "Calculate the density of the system");
304     registerModule(manager, &gmx_densmap, "densmap",
305                    "Calculate 2D planar or axial-radial density maps");
306     registerModule(manager, &gmx_densorder, "densorder",
307                    "Calculate surface fluctuations");
308     registerModule(manager, &gmx_dielectric, "dielectric",
309                    "Calculate frequency dependent dielectric constants");
310     registerModule(manager, &gmx_dipoles, "dipoles",
311                    "Compute the total dipole plus fluctuations");
312     registerModule(manager, &gmx_disre, "disre",
313                    "Analyze distance restraints");
314     registerModule(manager, &gmx_dos, "dos",
315                    "Analyze density of states and properties based on that");
316     registerModule(manager, &gmx_dyecoupl, "dyecoupl",
317                    "Extract dye dynamics from trajectories");
318     registerModule(manager, &gmx_dyndom, "dyndom",
319                    "Interpolate and extrapolate structure rotations");
320     registerModule(manager, &gmx_enemat, "enemat",
321                    "Extract an energy matrix from an energy file");
322     registerModule(manager, &gmx_energy, "energy",
323                    "Writes energies to xvg files and display averages");
324     registerModule(manager, &gmx_filter, "filter",
325                    "Frequency filter trajectories, useful for making smooth movies");
326     registerModule(manager, &gmx_gyrate, "gyrate",
327                    "Calculate the radius of gyration");
328     registerModule(manager, &gmx_h2order, "h2order",
329                    "Compute the orientation of water molecules");
330     registerModule(manager, &gmx_hbond, "hbond",
331                    "Compute and analyze hydrogen bonds");
332     registerModule(manager, &gmx_helix, "helix",
333                    "Calculate basic properties of alpha helices");
334     registerModule(manager, &gmx_helixorient, "helixorient",
335                    "Calculate local pitch/bending/rotation/orientation inside helices");
336     registerModule(manager, &gmx_hydorder, "hydorder",
337                    "Compute tetrahedrality parameters around a given atom");
338     registerModule(manager, &gmx_lie, "lie",
339                    "Estimate free energy from linear combinations");
340     registerModule(manager, &gmx_mdmat, "mdmat",
341                    "Calculate residue contact maps");
342     registerModule(manager, &gmx_mindist, "mindist",
343                    "Calculate the minimum distance between two groups");
344     registerModule(manager, &gmx_morph, "morph",
345                    "Interpolate linearly between conformations");
346     registerModule(manager, &gmx_msd, "msd",
347                    "Calculates mean square displacements");
348     registerModule(manager, &gmx_nmeig, "nmeig",
349                    "Diagonalize the Hessian for normal mode analysis");
350     registerModule(manager, &gmx_nmens, "nmens",
351                    "Generate an ensemble of structures from the normal modes");
352     registerModule(manager, &gmx_nmtraj, "nmtraj",
353                    "Generate a virtual oscillating trajectory from an eigenvector");
354     registerModule(manager, &gmx_options, "options", NULL);
355     registerModule(manager, &gmx_order, "order",
356                    "Compute the order parameter per atom for carbon tails");
357     registerModule(manager, &gmx_pme_error, "pme_error",
358                    "Estimate the error of using PME with a given input file");
359     registerModule(manager, &gmx_polystat, "polystat",
360                    "Calculate static properties of polymers");
361     registerModule(manager, &gmx_potential, "potential",
362                    "Calculate the electrostatic potential across the box");
363     registerModule(manager, &gmx_principal, "principal",
364                    "Calculate principal axes of inertia for a group of atoms");
365     registerModule(manager, &gmx_rama, "rama",
366                    "Compute Ramachandran plots");
367     registerModule(manager, &gmx_rdf, "rdf",
368                    "Calculate radial distribution functions");
369     registerModule(manager, &gmx_rms, "rms",
370                    "Calculate RMSDs with a reference structure and RMSD matrices");
371     registerModule(manager, &gmx_rmsdist, "rmsdist",
372                    "Calculate atom pair distances averaged with power -2, -3 or -6");
373     registerModule(manager, &gmx_rmsf, "rmsf",
374                    "Calculate atomic fluctuations");
375     registerModule(manager, &gmx_rotacf, "rotacf",
376                    "Calculate the rotational correlation function for molecules");
377     registerModule(manager, &gmx_rotmat, "rotmat",
378                    "Plot the rotation matrix for fitting to a reference structure");
379     registerModule(manager, &gmx_saltbr, "saltbr",
380                    "Compute salt bridges");
381     registerModule(manager, &gmx_sans, "sans",
382                    "Compute small angle neutron scattering spectra");
383     registerModule(manager, &gmx_saxs, "saxs",
384                    "Compute small angle X-ray scattering spectra");
385     registerModule(manager, &gmx_sham, "sham",
386                    "Compute free energies or other histograms from histograms");
387     registerModule(manager, &gmx_sigeps, "sigeps",
388                    "Convert c6/12 or c6/cn combinations to and from sigma/epsilon");
389     registerModule(manager, &gmx_sorient, "sorient",
390                    "Analyze solvent orientation around solutes");
391     registerModule(manager, &gmx_spatial, "spatial",
392                    "Calculate the spatial distribution function");
393     registerModule(manager, &gmx_spol, "spol",
394                    "Analyze solvent dipole orientation and polarization around solutes");
395     registerModule(manager, &gmx_tcaf, "tcaf",
396                    "Calculate viscosities of liquids");
397     registerModule(manager, &gmx_traj, "traj",
398                    "Plot x, v, f, box, temperature and rotational energy from trajectories");
399     registerModule(manager, &gmx_tune_pme, "tune_pme",
400                    "Time mdrun as a function of PME ranks to optimize settings");
401     registerModule(manager, &gmx_vanhove, "vanhove",
402                    "Compute Van Hove displacement and correlation functions");
403     registerModule(manager, &gmx_velacc, "velacc",
404                    "Calculate velocity autocorrelation functions");
405     registerModule(manager, &gmx_wham, "wham",
406                    "Perform weighted histogram analysis after umbrella sampling");
407     registerModule(manager, &gmx_wheel, "wheel",
408                    "Plot helical wheels");
409     registerModuleNoNice(manager, &gmx_view, "view",
410                          "View a trajectory on an X-Windows terminal");
411
412     {
413         gmx::CommandLineModuleGroup group =
414             manager->addModuleGroup("Generating topologies and coordinates");
415         group.addModuleWithDescription("editconf", "Edit the box and write subgroups");
416         group.addModule("protonate");
417         group.addModule("x2top");
418         group.addModule("solvate");
419         group.addModule("insert-molecules");
420         group.addModule("genconf");
421         group.addModule("genion");
422         group.addModule("genrestr");
423         group.addModule("pdb2gmx");
424     }
425     {
426         gmx::CommandLineModuleGroup group =
427             manager->addModuleGroup("Running a simulation");
428         group.addModule("grompp");
429         group.addModule("mdrun");
430         group.addModule("convert-tpr");
431     }
432     {
433         gmx::CommandLineModuleGroup group =
434             manager->addModuleGroup("Viewing trajectories");
435         group.addModule("nmtraj");
436         group.addModule("view");
437     }
438     {
439         gmx::CommandLineModuleGroup group =
440             manager->addModuleGroup("Processing energies");
441         group.addModule("enemat");
442         group.addModule("energy");
443         group.addModuleWithDescription("mdrun", "(Re)calculate energies for trajectory frames with -rerun");
444     }
445     {
446         gmx::CommandLineModuleGroup group =
447             manager->addModuleGroup("Converting files");
448         group.addModule("editconf");
449         group.addModule("eneconv");
450         group.addModule("sigeps");
451         group.addModule("trjcat");
452         group.addModule("trjconv");
453         group.addModule("xpm2ps");
454     }
455     {
456         gmx::CommandLineModuleGroup group =
457             manager->addModuleGroup("Tools");
458         group.addModule("analyze");
459         group.addModule("dyndom");
460         group.addModule("filter");
461         group.addModule("lie");
462         group.addModule("morph");
463         group.addModule("pme_error");
464         group.addModule("sham");
465         group.addModule("spatial");
466         group.addModule("traj");
467         group.addModule("tune_pme");
468         group.addModule("wham");
469         group.addModule("check");
470         group.addModule("dump");
471         group.addModule("make_ndx");
472         group.addModule("mk_angndx");
473         group.addModule("trjorder");
474         group.addModule("xpm2ps");
475     }
476     {
477         gmx::CommandLineModuleGroup group =
478             manager->addModuleGroup("Distances between structures");
479         group.addModule("cluster");
480         group.addModule("confrms");
481         group.addModule("rms");
482         group.addModule("rmsf");
483     }
484     {
485         gmx::CommandLineModuleGroup group =
486             manager->addModuleGroup("Distances in structures over time");
487         group.addModule("mindist");
488         group.addModule("mdmat");
489         group.addModule("polystat");
490         group.addModule("rmsdist");
491     }
492     {
493         gmx::CommandLineModuleGroup group =
494             manager->addModuleGroup("Mass distribution properties over time");
495         group.addModule("gyrate");
496         group.addModule("msd");
497         group.addModule("polystat");
498         group.addModule("rdf");
499         group.addModule("rotacf");
500         group.addModule("rotmat");
501         group.addModule("sans");
502         group.addModule("saxs");
503         group.addModule("traj");
504         group.addModule("vanhove");
505     }
506     {
507         gmx::CommandLineModuleGroup group =
508             manager->addModuleGroup("Analyzing bonded interactions");
509         group.addModule("angle");
510         group.addModule("mk_angndx");
511     }
512     {
513         gmx::CommandLineModuleGroup group =
514             manager->addModuleGroup("Structural properties");
515         group.addModule("anadock");
516         group.addModule("bundle");
517         group.addModule("clustsize");
518         group.addModule("disre");
519         group.addModule("hbond");
520         group.addModule("order");
521         group.addModule("principal");
522         group.addModule("rdf");
523         group.addModule("saltbr");
524         group.addModule("sorient");
525         group.addModule("spol");
526     }
527     {
528         gmx::CommandLineModuleGroup group =
529             manager->addModuleGroup("Kinetic properties");
530         group.addModule("bar");
531         group.addModule("current");
532         group.addModule("dos");
533         group.addModule("dyecoupl");
534         group.addModule("principal");
535         group.addModule("tcaf");
536         group.addModule("traj");
537         group.addModule("vanhove");
538         group.addModule("velacc");
539     }
540     {
541         gmx::CommandLineModuleGroup group =
542             manager->addModuleGroup("Electrostatic properties");
543         group.addModule("current");
544         group.addModule("dielectric");
545         group.addModule("dipoles");
546         group.addModule("potential");
547         group.addModule("spol");
548         group.addModule("genion");
549     }
550     {
551         gmx::CommandLineModuleGroup group =
552             manager->addModuleGroup("Protein-specific analysis");
553         group.addModule("do_dssp");
554         group.addModule("chi");
555         group.addModule("helix");
556         group.addModule("helixorient");
557         group.addModule("rama");
558         group.addModule("wheel");
559     }
560     {
561         gmx::CommandLineModuleGroup group =
562             manager->addModuleGroup("Interfaces");
563         group.addModule("bundle");
564         group.addModule("density");
565         group.addModule("densmap");
566         group.addModule("densorder");
567         group.addModule("h2order");
568         group.addModule("hydorder");
569         group.addModule("order");
570         group.addModule("potential");
571     }
572     {
573         gmx::CommandLineModuleGroup group =
574             manager->addModuleGroup("Covariance analysis");
575         group.addModuleWithDescription("anaeig", "Analyze the eigenvectors");
576         group.addModule("covar");
577         group.addModule("make_edi");
578     }
579     {
580         gmx::CommandLineModuleGroup group =
581             manager->addModuleGroup("Normal modes");
582         group.addModuleWithDescription("anaeig", "Analyze the normal modes");
583         group.addModule("nmeig");
584         group.addModule("nmtraj");
585         group.addModule("nmens");
586         group.addModule("grompp");
587         group.addModuleWithDescription("mdrun", "Find a potential energy minimum and calculate the Hessian");
588     }
589 }