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) 2012,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.
42 #include "gromacs/legacyheaders/typedefs.h"
43 #include "gromacs/utility/smalloc.h"
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/legacyheaders/macros.h"
46 #include "gromacs/utility/futil.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/fileio/trxio.h"
50 static int calc_ntype(int nft, int *ft, t_idef *idef)
54 for (i = 0; (i < idef->ntypes); i++)
56 for (f = 0; f < nft; f++)
58 if (idef->functype[i] == ft[f])
68 static void fill_ft_ind(int nft, int *ft, t_idef *idef,
69 int ft_ind[], char *grpnames[])
72 int i, f, ftype, ind = 0;
74 /* Loop over all the function types in the topology */
75 for (i = 0; (i < idef->ntypes); i++)
78 /* Check all the selected function types */
79 for (f = 0; f < nft; f++)
82 if (idef->functype[i] == ftype)
88 sprintf(buf, "Theta=%.1f_%.2f", idef->iparams[i].harmonic.rA,
89 idef->iparams[i].harmonic.krA);
92 sprintf(buf, "Cos_th=%.1f_%.2f", idef->iparams[i].harmonic.rA,
93 idef->iparams[i].harmonic.krA);
96 sprintf(buf, "UB_th=%.1f_%.2f2f", idef->iparams[i].u_b.thetaA,
97 idef->iparams[i].u_b.kthetaA);
99 case F_QUARTIC_ANGLES:
100 sprintf(buf, "Q_th=%.1f_%.2f_%.2f", idef->iparams[i].qangle.theta,
101 idef->iparams[i].qangle.c[0], idef->iparams[i].qangle.c[1]);
104 sprintf(buf, "Table=%d_%.2f", idef->iparams[i].tab.table,
105 idef->iparams[i].tab.kA);
108 sprintf(buf, "Phi=%.1f_%d_%.2f", idef->iparams[i].pdihs.phiA,
109 idef->iparams[i].pdihs.mult, idef->iparams[i].pdihs.cpA);
112 sprintf(buf, "Xi=%.1f_%.2f", idef->iparams[i].harmonic.rA,
113 idef->iparams[i].harmonic.krA);
116 sprintf(buf, "RB-A1=%.2f", idef->iparams[i].rbdihs.rbcA[1]);
119 sprintf(buf, "Theta=%.1f_%.2f", idef->iparams[i].harmonic.rA,
120 idef->iparams[i].harmonic.krA);
123 sprintf(buf, "Theta=%.1f_%.2f", idef->iparams[i].harmonic.rA,
124 idef->iparams[i].harmonic.krA);
127 sprintf(buf, "CBT-A1=%.2f", idef->iparams[i].cbtdihs.cbtcA[1]);
131 gmx_fatal(FARGS, "Unsupported function type '%s' selected",
132 interaction_function[ftype].longname);
134 grpnames[ind] = gmx_strdup(buf);
141 static void fill_ang(int nft, int *ft, int fac,
142 int nr[], int *index[], int ft_ind[], t_topology *top,
143 gmx_bool bNoH, real hq)
145 int f, ftype, i, j, indg, nr_fac;
153 atom = top->atoms.atom;
155 for (f = 0; f < nft; f++)
158 ia = idef->il[ftype].iatoms;
159 for (i = 0; (i < idef->il[ftype].nr); )
161 indg = ft_ind[ia[0]];
164 gmx_incons("Routine fill_ang");
169 for (j = 0; j < fac; j++)
171 if (atom[ia[1+j]].m < 1.5)
179 for (j = 0; j < fac; j++)
181 if (atom[ia[1+j]].m < 1.5 && fabs(atom[ia[1+j]].q) < hq)
189 if (nr[indg] % 1000 == 0)
191 srenew(index[indg], fac*(nr[indg]+1000));
193 nr_fac = fac*nr[indg];
194 for (j = 0; (j < fac); j++)
196 index[indg][nr_fac+j] = ia[j+1];
200 ia += interaction_function[ftype].nratoms+1;
201 i += interaction_function[ftype].nratoms+1;
206 static int *select_ftype(const char *opt, int *nft, int *mult)
208 int *ft = NULL, ftype;
213 for (ftype = 0; ftype < F_NRE; ftype++)
215 if ((interaction_function[ftype].flags & IF_ATYPE) ||
216 ftype == F_TABANGLES)
248 int gmx_mk_angndx(int argc, char *argv[])
250 static const char *desc[] = {
251 "[THISMODULE] makes an index file for calculation of",
252 "angle distributions etc. It uses a run input file ([TT].tpx[tt]) for the",
253 "definitions of the angles, dihedrals etc."
255 static const char *opt[] = { NULL, "angle", "dihedral", "improper", "ryckaert-bellemans", NULL };
256 static gmx_bool bH = TRUE;
259 { "-type", FALSE, etENUM, {opt},
261 { "-hyd", FALSE, etBOOL, {&bH},
262 "Include angles with atoms with mass < 1.5" },
263 { "-hq", FALSE, etREAL, {&hq},
264 "Ignore angles with atoms with mass < 1.5 and magnitude of their charge less than this value" }
271 int nft = 0, *ft, mult = 0;
277 { efTPX, NULL, NULL, ffREAD },
278 { efNDX, NULL, "angle", ffWRITE }
280 #define NFILE asize(fnm)
282 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
283 asize(desc), desc, 0, NULL, &oenv))
289 ft = select_ftype(opt[0], &nft, &mult);
291 top = read_top(ftp2fn(efTPX, NFILE, fnm), NULL);
293 ntype = calc_ntype(nft, ft, &(top->idef));
294 snew(grpnames, ntype);
295 snew(ft_ind, top->idef.ntypes);
296 fill_ft_ind(nft, ft, &top->idef, ft_ind, grpnames);
300 fill_ang(nft, ft, mult, nr, index, ft_ind, top, !bH, hq);
302 out = ftp2FILE(efNDX, NFILE, fnm, "w");
303 for (i = 0; (i < ntype); i++)
307 fprintf(out, "[ %s ]\n", grpnames[i]);
308 for (j = 0; (j < nr[i]*mult); j++)
310 fprintf(out, " %5d", index[i][j]+1);