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