Merge branch 'release-4-6'
[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, 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 "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.",
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
108
109     t_filenm fnm[] = {
110         { efTRX, "-f",   NULL,       ffREAD },
111         { efTPS, NULL,   NULL,       ffREAD },
112         { efNDX, NULL,   NULL,       ffOPTRD },
113         { efDAT, "-a1",  "axis1",    ffWRITE },
114         { efDAT, "-a2",  "axis2",    ffWRITE },
115         { efDAT, "-a3",  "axis3",    ffWRITE },
116         { efDAT, "-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 | PCA_BE_NICE,
122                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
123     {
124         return 0;
125     }
126
127     axis1 = ffopen(opt2fn("-a1", NFILE, fnm), "w");
128     axis2 = ffopen(opt2fn("-a2", NFILE, fnm), "w");
129     axis3 = ffopen(opt2fn("-a3", NFILE, fnm), "w");
130     fmoi  = ffopen(opt2fn("-om", NFILE, fnm), "w");
131
132     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, NULL, NULL, box, TRUE);
133
134     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
135
136     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
137
138     gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
139
140     do
141     {
142         gmx_rmpbc(gpbc, natoms, box, x);
143
144         calc_principal_axes(&top, x, index, gnx, axes, moi);
145
146         fprintf(axis1, "%15.10f     %15.10f  %15.10f  %15.10f\n", t, axes[XX][XX], axes[YY][XX], axes[ZZ][XX]);
147         fprintf(axis2, "%15.10f     %15.10f  %15.10f  %15.10f\n", t, axes[XX][YY], axes[YY][YY], axes[ZZ][YY]);
148         fprintf(axis3, "%15.10f     %15.10f  %15.10f  %15.10f\n", t, axes[XX][ZZ], axes[YY][ZZ], axes[ZZ][ZZ]);
149         fprintf(fmoi,  "%15.10f     %15.10f  %15.10f  %15.10f\n", t, moi[XX], moi[YY], moi[ZZ]);
150     }
151     while (read_next_x(oenv, status, &t, x, box));
152
153     gmx_rmpbc_done(gpbc);
154
155
156     close_trj(status);
157     ffclose(axis1);
158     ffclose(axis2);
159     ffclose(axis3);
160     ffclose(fmoi);
161
162     return 0;
163 }