Merge branch 'release-4-6' into release-5-0
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_principal.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include <string.h>
43
44 #include "gromacs/commandline/pargs.h"
45 #include "sysstuff.h"
46 #include "typedefs.h"
47 #include "gromacs/utility/smalloc.h"
48 #include "macros.h"
49 #include "vec.h"
50 #include "pbc.h"
51 #include "gromacs/fileio/futil.h"
52 #include "index.h"
53 #include "mshift.h"
54 #include "xvgr.h"
55 #include "princ.h"
56 #include "rmpbc.h"
57 #include "txtdump.h"
58 #include "gromacs/fileio/tpxio.h"
59 #include "gromacs/fileio/trxio.h"
60 #include "gstat.h"
61 #include "gmx_ana.h"
62
63
64 void
65 calc_principal_axes(t_topology *   top,
66                     rvec *         x,
67                     atom_id *      index,
68                     int            n,
69                     matrix         axes,
70                     rvec           inertia)
71 {
72     rvec   xcm;
73
74     sub_xcm(x, n, index, top->atoms.atom, xcm, FALSE);
75     principal_comp(n, index, top->atoms.atom, x, axes, inertia);
76 }
77
78 int gmx_principal(int argc, char *argv[])
79 {
80     const char     *desc[] = {
81         "[THISMODULE] calculates the three principal axes of inertia for a group",
82         "of atoms. NOTE: Old versions of Gromacs wrote the output data in a",
83         "strange transposed way. As of Gromacs-5.0, the output file paxis1.dat",
84         "contains the x/y/z components of the first (major) principal axis for",
85         "each frame, and similarly for the middle and minor axes in paxis2.dat",
86         "and paxis3.dat."
87     };
88     static gmx_bool foo = FALSE;
89
90     t_pargs         pa[] = {
91         { "-foo",      FALSE, etBOOL, {&foo}, "Dummy option to avoid empty array" }
92     };
93     t_trxstatus    *status;
94     t_topology      top;
95     int             ePBC;
96     real            t;
97     rvec      *     x;
98
99     int             natoms;
100     char           *grpname, title[256];
101     int             i, j, m, gnx, nam, mol;
102     atom_id        *index;
103     rvec            a1, a2, a3, moi;
104     FILE      *     axis1;
105     FILE      *     axis2;
106     FILE      *     axis3;
107     FILE      *     fmoi;
108     matrix          axes, box;
109     output_env_t    oenv;
110     gmx_rmpbc_t     gpbc = NULL;
111     char **         legend;
112
113     t_filenm        fnm[] = {
114         { efTRX, "-f",   NULL,       ffREAD },
115         { efTPS, NULL,   NULL,       ffREAD },
116         { efNDX, NULL,   NULL,       ffOPTRD },
117         { efXVG, "-a1",  "paxis1",   ffWRITE },
118         { efXVG, "-a2",  "paxis2",   ffWRITE },
119         { efXVG, "-a3",  "paxis3",   ffWRITE },
120         { efXVG, "-om",  "moi",      ffWRITE }
121     };
122 #define NFILE asize(fnm)
123
124     if (!parse_common_args(&argc, argv,
125                            PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW | PCA_BE_NICE,
126                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
127     {
128         return 0;
129     }
130
131     snew(legend, DIM);
132     for (i = 0; i < DIM; i++)
133     {
134         snew(legend[i], STRLEN);
135         sprintf(legend[i], "%c component", 'X'+i);
136     }
137
138     axis1 = xvgropen(opt2fn("-a1", NFILE, fnm), "Principal axis 1 (major axis)",
139                      output_env_get_xvgr_tlabel(oenv), "Component (nm)", oenv);
140     xvgr_legend(axis1, DIM, (const char **)legend, oenv);
141
142     axis2 = xvgropen(opt2fn("-a2", NFILE, fnm), "Principal axis 2 (middle axis)",
143                      output_env_get_xvgr_tlabel(oenv), "Component (nm)", oenv);
144     xvgr_legend(axis2, DIM, (const char **)legend, oenv);
145
146     axis3 = xvgropen(opt2fn("-a3", NFILE, fnm), "Principal axis 3 (minor axis)",
147                      output_env_get_xvgr_tlabel(oenv), "Component (nm)", oenv);
148     xvgr_legend(axis3, DIM, (const char **)legend, oenv);
149
150     sprintf(legend[XX], "Axis 1 (major)");
151     sprintf(legend[YY], "Axis 2 (middle)");
152     sprintf(legend[ZZ], "Axis 3 (minor)");
153
154     fmoi  = xvgropen(opt2fn("-om", NFILE, fnm), "Moments of inertia around inertial axes",
155                      output_env_get_xvgr_tlabel(oenv), "I (au nm\\S2\\N)", oenv);
156     xvgr_legend(fmoi, DIM, (const char **)legend, oenv);
157
158     for (i = 0; i < DIM; i++)
159     {
160         sfree(legend[i]);
161     }
162     sfree(legend);
163
164     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, NULL, NULL, box, TRUE);
165
166     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
167
168     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
169
170     gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
171
172     do
173     {
174         gmx_rmpbc(gpbc, natoms, box, x);
175
176         calc_principal_axes(&top, x, index, gnx, axes, moi);
177
178         fprintf(axis1, "%15.10f     %15.10f  %15.10f  %15.10f\n", t, axes[XX][XX], axes[XX][YY], axes[XX][ZZ]);
179         fprintf(axis2, "%15.10f     %15.10f  %15.10f  %15.10f\n", t, axes[YY][XX], axes[YY][YY], axes[YY][ZZ]);
180         fprintf(axis3, "%15.10f     %15.10f  %15.10f  %15.10f\n", t, axes[ZZ][XX], axes[ZZ][YY], axes[ZZ][ZZ]);
181         fprintf(fmoi,  "%15.10f     %15.10f  %15.10f  %15.10f\n", t, moi[XX], moi[YY], moi[ZZ]);
182     }
183     while (read_next_x(oenv, status, &t, x, box));
184
185     gmx_rmpbc_done(gpbc);
186
187     close_trj(status);
188
189     xvgrclose(axis1);
190     xvgrclose(axis2);
191     xvgrclose(axis3);
192     xvgrclose(fmoi);
193
194     return 0;
195 }