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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
50 #include "gmx_fatal.h"
52 static int calc_ntype(int nft, int *ft, t_idef *idef)
56 for (i = 0; (i < idef->ntypes); i++)
58 for (f = 0; f < nft; f++)
60 if (idef->functype[i] == ft[f])
70 static void fill_ft_ind(int nft, int *ft, t_idef *idef,
71 int ft_ind[], char *grpnames[])
74 int i, f, ftype, ind = 0;
76 /* Loop over all the function types in the topology */
77 for (i = 0; (i < idef->ntypes); i++)
80 /* Check all the selected function types */
81 for (f = 0; f < nft; f++)
84 if (idef->functype[i] == ftype)
90 sprintf(buf, "Theta=%.1f_%.2f", idef->iparams[i].harmonic.rA,
91 idef->iparams[i].harmonic.krA);
94 sprintf(buf, "Cos_th=%.1f_%.2f", idef->iparams[i].harmonic.rA,
95 idef->iparams[i].harmonic.krA);
98 sprintf(buf, "UB_th=%.1f_%.2f2f", idef->iparams[i].u_b.thetaA,
99 idef->iparams[i].u_b.kthetaA);
101 case F_QUARTIC_ANGLES:
102 sprintf(buf, "Q_th=%.1f_%.2f_%.2f", idef->iparams[i].qangle.theta,
103 idef->iparams[i].qangle.c[0], idef->iparams[i].qangle.c[1]);
106 sprintf(buf, "Table=%d_%.2f", idef->iparams[i].tab.table,
107 idef->iparams[i].tab.kA);
110 sprintf(buf, "Phi=%.1f_%d_%.2f", idef->iparams[i].pdihs.phiA,
111 idef->iparams[i].pdihs.mult, idef->iparams[i].pdihs.cpA);
114 sprintf(buf, "Xi=%.1f_%.2f", idef->iparams[i].harmonic.rA,
115 idef->iparams[i].harmonic.krA);
118 sprintf(buf, "RB-A1=%.2f", idef->iparams[i].rbdihs.rbcA[1]);
121 gmx_fatal(FARGS, "Unsupported function type '%s' selected",
122 interaction_function[ftype].longname);
124 grpnames[ind] = strdup(buf);
131 static void fill_ang(int nft, int *ft, int fac,
132 int nr[], int *index[], int ft_ind[], t_topology *top,
133 gmx_bool bNoH, real hq)
135 int f, ftype, i, j, indg, nr_fac;
143 atom = top->atoms.atom;
145 for (f = 0; f < nft; f++)
148 ia = idef->il[ftype].iatoms;
149 for (i = 0; (i < idef->il[ftype].nr); )
151 indg = ft_ind[ia[0]];
154 gmx_incons("Routine fill_ang");
159 for (j = 0; j < fac; j++)
161 if (atom[ia[1+j]].m < 1.5)
169 for (j = 0; j < fac; j++)
171 if (atom[ia[1+j]].m < 1.5 && fabs(atom[ia[1+j]].q) < hq)
179 if (nr[indg] % 1000 == 0)
181 srenew(index[indg], fac*(nr[indg]+1000));
183 nr_fac = fac*nr[indg];
184 for (j = 0; (j < fac); j++)
186 index[indg][nr_fac+j] = ia[j+1];
190 ia += interaction_function[ftype].nratoms+1;
191 i += interaction_function[ftype].nratoms+1;
196 static int *select_ftype(const char *opt, int *nft, int *mult)
198 int *ft = NULL, ftype;
203 for (ftype = 0; ftype < F_NRE; ftype++)
205 if ((interaction_function[ftype].flags & IF_ATYPE) ||
206 ftype == F_TABANGLES)
238 int gmx_mk_angndx(int argc, char *argv[])
240 static const char *desc[] = {
241 "[TT]mk_angndx[tt] makes an index file for calculation of",
242 "angle distributions etc. It uses a run input file ([TT].tpx[tt]) for the",
243 "definitions of the angles, dihedrals etc."
245 static const char *opt[] = { NULL, "angle", "dihedral", "improper", "ryckaert-bellemans", NULL };
246 static gmx_bool bH = TRUE;
249 { "-type", FALSE, etENUM, {opt},
251 { "-hyd", FALSE, etBOOL, {&bH},
252 "Include angles with atoms with mass < 1.5" },
253 { "-hq", FALSE, etREAL, {&hq},
254 "Ignore angles with atoms with mass < 1.5 and magnitude of their charge less than this value" }
261 int nft = 0, *ft, mult = 0;
267 { efTPX, NULL, NULL, ffREAD },
268 { efNDX, NULL, "angle", ffWRITE }
270 #define NFILE asize(fnm)
272 parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
273 asize(desc), desc, 0, NULL, &oenv);
276 ft = select_ftype(opt[0], &nft, &mult);
278 top = read_top(ftp2fn(efTPX, NFILE, fnm), NULL);
280 ntype = calc_ntype(nft, ft, &(top->idef));
281 snew(grpnames, ntype);
282 snew(ft_ind, top->idef.ntypes);
283 fill_ft_ind(nft, ft, &top->idef, ft_ind, grpnames);
287 fill_ang(nft, ft, mult, nr, index, ft_ind, top, !bH, hq);
289 out = ftp2FILE(efNDX, NFILE, fnm, "w");
290 for (i = 0; (i < ntype); i++)
294 fprintf(out, "[ %s ]\n", grpnames[i]);
295 for (j = 0; (j < nr[i]*mult); j++)
297 fprintf(out, " %5d", index[i][j]+1);