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