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