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