Add gmx convert-trj
[gromacs.git] / src / gromacs / gmxana / gmx_principal.cpp
blob65bfa3a2026bcf496fd5ec970a443e8db6ea2f5f
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,2015,2016,2017,2018, 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 <cmath>
40 #include <cstring>
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/math/vec.h"
50 #include "gromacs/pbcutil/rmpbc.h"
51 #include "gromacs/topology/index.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/arraysize.h"
54 #include "gromacs/utility/cstringutil.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/utility/smalloc.h"
59 static void
60 calc_principal_axes(const t_topology *top,
61 rvec *x,
62 int *index,
63 int n,
64 matrix axes,
65 rvec inertia)
67 rvec xcm;
69 sub_xcm(x, n, index, top->atoms.atom, xcm, FALSE);
70 principal_comp(n, index, top->atoms.atom, x, axes, inertia);
73 int gmx_principal(int argc, char *argv[])
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."
83 static gmx_bool foo = FALSE;
85 t_pargs pa[] = {
86 { "-foo", FALSE, etBOOL, {&foo}, "Dummy option to avoid empty array" }
88 t_trxstatus *status;
89 t_topology top;
90 int ePBC;
91 real t;
92 rvec * x;
94 int natoms;
95 char *grpname;
96 int i, gnx;
97 int *index;
98 rvec moi;
99 FILE * axis1;
100 FILE * axis2;
101 FILE * axis3;
102 FILE * fmoi;
103 matrix axes, box;
104 gmx_output_env_t *oenv;
105 gmx_rmpbc_t gpbc = nullptr;
106 char ** legend;
108 t_filenm fnm[] = {
109 { efTRX, "-f", nullptr, ffREAD },
110 { efTPS, nullptr, nullptr, ffREAD },
111 { efNDX, nullptr, nullptr, ffOPTRD },
112 { efXVG, "-a1", "paxis1", ffWRITE },
113 { efXVG, "-a2", "paxis2", ffWRITE },
114 { efXVG, "-a3", "paxis3", ffWRITE },
115 { efXVG, "-om", "moi", ffWRITE }
117 #define NFILE asize(fnm)
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, nullptr, &oenv))
123 return 0;
126 snew(legend, DIM);
127 for (i = 0; i < DIM; i++)
129 snew(legend[i], STRLEN);
130 sprintf(legend[i], "%c component", 'X'+i);
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, legend, oenv);
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, legend, oenv);
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, legend, oenv);
145 sprintf(legend[XX], "Axis 1 (major)");
146 sprintf(legend[YY], "Axis 2 (middle)");
147 sprintf(legend[ZZ], "Axis 3 (minor)");
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, legend, oenv);
153 for (i = 0; i < DIM; i++)
155 sfree(legend[i]);
157 sfree(legend);
159 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, nullptr, nullptr, box, TRUE);
161 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
163 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
165 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
169 gmx_rmpbc(gpbc, natoms, box, x);
171 calc_principal_axes(&top, x, index, gnx, axes, moi);
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]);
178 while (read_next_x(oenv, status, &t, x, box));
180 gmx_rmpbc_done(gpbc);
182 close_trx(status);
184 xvgrclose(axis1);
185 xvgrclose(axis2);
186 xvgrclose(axis3);
187 xvgrclose(fmoi);
189 return 0;