Move read_tps_conf() to confio.h
[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,2015, 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 #include "gmxpre.h"
38
39 #include <math.h>
40 #include <string.h>
41
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/fileio/confio.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/gmx_ana.h"
47 #include "gromacs/gmxana/gstat.h"
48 #include "gromacs/gmxana/princ.h"
49 #include "gromacs/legacyheaders/macros.h"
50 #include "gromacs/legacyheaders/txtdump.h"
51 #include "gromacs/legacyheaders/typedefs.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/pbcutil/rmpbc.h"
54 #include "gromacs/topology/index.h"
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/utility/futil.h"
57 #include "gromacs/utility/smalloc.h"
58
59
60 void
61 calc_principal_axes(t_topology *   top,
62                     rvec *         x,
63                     atom_id *      index,
64                     int            n,
65                     matrix         axes,
66                     rvec           inertia)
67 {
68     rvec   xcm;
69
70     sub_xcm(x, n, index, top->atoms.atom, xcm, FALSE);
71     principal_comp(n, index, top->atoms.atom, x, axes, inertia);
72 }
73
74 int gmx_principal(int argc, char *argv[])
75 {
76     const char     *desc[] = {
77         "[THISMODULE] calculates the three principal axes of inertia for a group",
78         "of atoms. NOTE: Old versions of GROMACS wrote the output data in a",
79         "strange transposed way. As of GROMACS 5.0, the output file paxis1.dat",
80         "contains the x/y/z components of the first (major) principal axis for",
81         "each frame, and similarly for the middle and minor axes in paxis2.dat",
82         "and paxis3.dat."
83     };
84     static gmx_bool foo = FALSE;
85
86     t_pargs         pa[] = {
87         { "-foo",      FALSE, etBOOL, {&foo}, "Dummy option to avoid empty array" }
88     };
89     t_trxstatus    *status;
90     t_topology      top;
91     int             ePBC;
92     real            t;
93     rvec      *     x;
94
95     int             natoms;
96     char           *grpname, title[256];
97     int             i, j, m, gnx, nam, mol;
98     atom_id        *index;
99     rvec            a1, a2, a3, moi;
100     FILE      *     axis1;
101     FILE      *     axis2;
102     FILE      *     axis3;
103     FILE      *     fmoi;
104     matrix          axes, box;
105     output_env_t    oenv;
106     gmx_rmpbc_t     gpbc = NULL;
107     char **         legend;
108
109     t_filenm        fnm[] = {
110         { efTRX, "-f",   NULL,       ffREAD },
111         { efTPS, NULL,   NULL,       ffREAD },
112         { efNDX, NULL,   NULL,       ffOPTRD },
113         { efXVG, "-a1",  "paxis1",   ffWRITE },
114         { efXVG, "-a2",  "paxis2",   ffWRITE },
115         { efXVG, "-a3",  "paxis3",   ffWRITE },
116         { efXVG, "-om",  "moi",      ffWRITE }
117     };
118 #define NFILE asize(fnm)
119
120     if (!parse_common_args(&argc, argv,
121                            PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW,
122                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
123     {
124         return 0;
125     }
126
127     snew(legend, DIM);
128     for (i = 0; i < DIM; i++)
129     {
130         snew(legend[i], STRLEN);
131         sprintf(legend[i], "%c component", 'X'+i);
132     }
133
134     axis1 = xvgropen(opt2fn("-a1", NFILE, fnm), "Principal axis 1 (major axis)",
135                      output_env_get_xvgr_tlabel(oenv), "Component (nm)", oenv);
136     xvgr_legend(axis1, DIM, (const char **)legend, oenv);
137
138     axis2 = xvgropen(opt2fn("-a2", NFILE, fnm), "Principal axis 2 (middle axis)",
139                      output_env_get_xvgr_tlabel(oenv), "Component (nm)", oenv);
140     xvgr_legend(axis2, DIM, (const char **)legend, oenv);
141
142     axis3 = xvgropen(opt2fn("-a3", NFILE, fnm), "Principal axis 3 (minor axis)",
143                      output_env_get_xvgr_tlabel(oenv), "Component (nm)", oenv);
144     xvgr_legend(axis3, DIM, (const char **)legend, oenv);
145
146     sprintf(legend[XX], "Axis 1 (major)");
147     sprintf(legend[YY], "Axis 2 (middle)");
148     sprintf(legend[ZZ], "Axis 3 (minor)");
149
150     fmoi  = xvgropen(opt2fn("-om", NFILE, fnm), "Moments of inertia around inertial axes",
151                      output_env_get_xvgr_tlabel(oenv), "I (au nm\\S2\\N)", oenv);
152     xvgr_legend(fmoi, DIM, (const char **)legend, oenv);
153
154     for (i = 0; i < DIM; i++)
155     {
156         sfree(legend[i]);
157     }
158     sfree(legend);
159
160     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, NULL, NULL, box, TRUE);
161
162     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
163
164     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
165
166     gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
167
168     do
169     {
170         gmx_rmpbc(gpbc, natoms, box, x);
171
172         calc_principal_axes(&top, x, index, gnx, axes, moi);
173
174         fprintf(axis1, "%15.10f     %15.10f  %15.10f  %15.10f\n", t, axes[XX][XX], axes[XX][YY], axes[XX][ZZ]);
175         fprintf(axis2, "%15.10f     %15.10f  %15.10f  %15.10f\n", t, axes[YY][XX], axes[YY][YY], axes[YY][ZZ]);
176         fprintf(axis3, "%15.10f     %15.10f  %15.10f  %15.10f\n", t, axes[ZZ][XX], axes[ZZ][YY], axes[ZZ][ZZ]);
177         fprintf(fmoi,  "%15.10f     %15.10f  %15.10f  %15.10f\n", t, moi[XX], moi[YY], moi[ZZ]);
178     }
179     while (read_next_x(oenv, status, &t, x, box));
180
181     gmx_rmpbc_done(gpbc);
182
183     close_trj(status);
184
185     xvgrclose(axis1);
186     xvgrclose(axis2);
187     xvgrclose(axis3);
188     xvgrclose(fmoi);
189
190     return 0;
191 }