Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_principal.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <math.h>
40 #include <string.h>
41
42 #include "gromacs/commandline/pargs.h"
43 #include "sysstuff.h"
44 #include "typedefs.h"
45 #include "smalloc.h"
46 #include "macros.h"
47 #include "vec.h"
48 #include "pbc.h"
49 #include "gromacs/fileio/futil.h"
50 #include "index.h"
51 #include "mshift.h"
52 #include "xvgr.h"
53 #include "princ.h"
54 #include "rmpbc.h"
55 #include "txtdump.h"
56 #include "gromacs/fileio/tpxio.h"
57 #include "gromacs/fileio/trxio.h"
58 #include "gstat.h"
59 #include "gmx_ana.h"
60
61
62 void
63 calc_principal_axes(t_topology *   top,
64                     rvec *         x,
65                     atom_id *      index,
66                     int            n,
67                     matrix         axes,
68                     rvec           inertia)
69 {
70     rvec   xcm;
71
72     sub_xcm(x, n, index, top->atoms.atom, xcm, FALSE);
73     principal_comp(n, index, top->atoms.atom, x, axes, inertia);
74 }
75
76 int gmx_principal(int argc, char *argv[])
77 {
78     const char     *desc[] = {
79         "[THISMODULE] calculates the three principal axes of inertia for a group",
80         "of atoms.",
81     };
82     static gmx_bool foo = FALSE;
83
84     t_pargs         pa[] = {
85         { "-foo",      FALSE, etBOOL, {&foo}, "Dummy option to avoid empty array" }
86     };
87     t_trxstatus    *status;
88     t_topology      top;
89     int             ePBC;
90     real            t;
91     rvec      *     x;
92
93     int             natoms;
94     char           *grpname, title[256];
95     int             i, j, m, gnx, nam, mol;
96     atom_id        *index;
97     rvec            a1, a2, a3, moi;
98     FILE      *     axis1;
99     FILE      *     axis2;
100     FILE      *     axis3;
101     FILE      *     fmoi;
102     matrix          axes, box;
103     output_env_t    oenv;
104     gmx_rmpbc_t     gpbc = NULL;
105
106
107     t_filenm fnm[] = {
108         { efTRX, "-f",   NULL,       ffREAD },
109         { efTPS, NULL,   NULL,       ffREAD },
110         { efNDX, NULL,   NULL,       ffOPTRD },
111         { efDAT, "-a1",  "axis1",    ffWRITE },
112         { efDAT, "-a2",  "axis2",    ffWRITE },
113         { efDAT, "-a3",  "axis3",    ffWRITE },
114         { efDAT, "-om",  "moi",      ffWRITE }
115     };
116 #define NFILE asize(fnm)
117
118     if (!parse_common_args(&argc, argv,
119                            PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW | PCA_BE_NICE,
120                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
121     {
122         return 0;
123     }
124
125     axis1 = ffopen(opt2fn("-a1", NFILE, fnm), "w");
126     axis2 = ffopen(opt2fn("-a2", NFILE, fnm), "w");
127     axis3 = ffopen(opt2fn("-a3", NFILE, fnm), "w");
128     fmoi  = ffopen(opt2fn("-om", NFILE, fnm), "w");
129
130     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, NULL, NULL, box, TRUE);
131
132     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
133
134     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
135
136     gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
137
138     do
139     {
140         gmx_rmpbc(gpbc, natoms, box, x);
141
142         calc_principal_axes(&top, x, index, gnx, axes, moi);
143
144         fprintf(axis1, "%15.10f     %15.10f  %15.10f  %15.10f\n", t, axes[XX][XX], axes[YY][XX], axes[ZZ][XX]);
145         fprintf(axis2, "%15.10f     %15.10f  %15.10f  %15.10f\n", t, axes[XX][YY], axes[YY][YY], axes[ZZ][YY]);
146         fprintf(axis3, "%15.10f     %15.10f  %15.10f  %15.10f\n", t, axes[XX][ZZ], axes[YY][ZZ], axes[ZZ][ZZ]);
147         fprintf(fmoi,  "%15.10f     %15.10f  %15.10f  %15.10f\n", t, moi[XX], moi[YY], moi[ZZ]);
148     }
149     while (read_next_x(oenv, status, &t, x, box));
150
151     gmx_rmpbc_done(gpbc);
152
153
154     close_trj(status);
155     ffclose(axis1);
156     ffclose(axis2);
157     ffclose(axis3);
158     ffclose(fmoi);
159
160     return 0;
161 }