Removed include of types/ifunc.h from typedefs.h
[gromacs.git] / src / gromacs / gmxana / gmx_mk_angndx.cpp
blob4554e216d19b941e38a2830e2025a77b91985e0b
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) 2012,2013,2014,2015, 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>
41 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/fileio/trxio.h"
43 #include "gromacs/gmxana/gmx_ana.h"
44 #include "gromacs/legacyheaders/typedefs.h"
45 #include "gromacs/legacyheaders/types/ifunc.h"
46 #include "gromacs/utility/arraysize.h"
47 #include "gromacs/utility/cstringutil.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/futil.h"
50 #include "gromacs/utility/gmxassert.h"
51 #include "gromacs/utility/smalloc.h"
53 static int calc_ntype(int nft, int *ft, t_idef *idef)
55 int i, f, nf = 0;
57 for (i = 0; (i < idef->ntypes); i++)
59 for (f = 0; f < nft; f++)
61 if (idef->functype[i] == ft[f])
63 nf++;
68 return nf;
71 static void fill_ft_ind(int nft, int *ft, t_idef *idef,
72 int ft_ind[], char *grpnames[])
74 char buf[125];
75 int i, f, ftype, ind = 0;
77 /* Loop over all the function types in the topology */
78 for (i = 0; (i < idef->ntypes); i++)
80 ft_ind[i] = -1;
81 /* Check all the selected function types */
82 for (f = 0; f < nft; f++)
84 ftype = ft[f];
85 if (idef->functype[i] == ftype)
87 ft_ind[i] = ind;
88 switch (ftype)
90 case F_ANGLES:
91 sprintf(buf, "Theta=%.1f_%.2f", idef->iparams[i].harmonic.rA,
92 idef->iparams[i].harmonic.krA);
93 break;
94 case F_G96ANGLES:
95 sprintf(buf, "Cos_th=%.1f_%.2f", idef->iparams[i].harmonic.rA,
96 idef->iparams[i].harmonic.krA);
97 break;
98 case F_UREY_BRADLEY:
99 sprintf(buf, "UB_th=%.1f_%.2f2f", idef->iparams[i].u_b.thetaA,
100 idef->iparams[i].u_b.kthetaA);
101 break;
102 case F_QUARTIC_ANGLES:
103 sprintf(buf, "Q_th=%.1f_%.2f_%.2f", idef->iparams[i].qangle.theta,
104 idef->iparams[i].qangle.c[0], idef->iparams[i].qangle.c[1]);
105 break;
106 case F_TABANGLES:
107 sprintf(buf, "Table=%d_%.2f", idef->iparams[i].tab.table,
108 idef->iparams[i].tab.kA);
109 break;
110 case F_PDIHS:
111 sprintf(buf, "Phi=%.1f_%d_%.2f", idef->iparams[i].pdihs.phiA,
112 idef->iparams[i].pdihs.mult, idef->iparams[i].pdihs.cpA);
113 break;
114 case F_IDIHS:
115 sprintf(buf, "Xi=%.1f_%.2f", idef->iparams[i].harmonic.rA,
116 idef->iparams[i].harmonic.krA);
117 break;
118 case F_RBDIHS:
119 sprintf(buf, "RB-A1=%.2f", idef->iparams[i].rbdihs.rbcA[1]);
120 break;
121 case F_RESTRANGLES:
122 sprintf(buf, "Theta=%.1f_%.2f", idef->iparams[i].harmonic.rA,
123 idef->iparams[i].harmonic.krA);
124 break;
125 case F_RESTRDIHS:
126 sprintf(buf, "Theta=%.1f_%.2f", idef->iparams[i].harmonic.rA,
127 idef->iparams[i].harmonic.krA);
128 break;
129 case F_CBTDIHS:
130 sprintf(buf, "CBT-A1=%.2f", idef->iparams[i].cbtdihs.cbtcA[1]);
131 break;
133 default:
134 gmx_fatal(FARGS, "Unsupported function type '%s' selected",
135 interaction_function[ftype].longname);
137 grpnames[ind] = gmx_strdup(buf);
138 ind++;
144 static void fill_ang(int nft, int *ft, int fac,
145 int nr[], int *index[], int ft_ind[], t_topology *top,
146 gmx_bool bNoH, real hq)
148 int f, ftype, i, j, indg, nr_fac;
149 gmx_bool bUse;
150 t_idef *idef;
151 t_atom *atom;
152 t_iatom *ia;
155 idef = &top->idef;
156 atom = top->atoms.atom;
158 for (f = 0; f < nft; f++)
160 ftype = ft[f];
161 ia = idef->il[ftype].iatoms;
162 for (i = 0; (i < idef->il[ftype].nr); )
164 indg = ft_ind[ia[0]];
165 if (indg == -1)
167 gmx_incons("Routine fill_ang");
169 bUse = TRUE;
170 if (bNoH)
172 for (j = 0; j < fac; j++)
174 if (atom[ia[1+j]].m < 1.5)
176 bUse = FALSE;
180 if (hq)
182 for (j = 0; j < fac; j++)
184 if (atom[ia[1+j]].m < 1.5 && std::abs(atom[ia[1+j]].q) < hq)
186 bUse = FALSE;
190 if (bUse)
192 if (nr[indg] % 1000 == 0)
194 srenew(index[indg], fac*(nr[indg]+1000));
196 nr_fac = fac*nr[indg];
197 for (j = 0; (j < fac); j++)
199 index[indg][nr_fac+j] = ia[j+1];
201 nr[indg]++;
203 ia += interaction_function[ftype].nratoms+1;
204 i += interaction_function[ftype].nratoms+1;
209 static int *select_ftype(const char *opt, int *nft, int *mult)
211 int *ft = NULL, ftype;
213 if (opt[0] == 'a')
215 *mult = 3;
216 for (ftype = 0; ftype < F_NRE; ftype++)
218 if ((interaction_function[ftype].flags & IF_ATYPE) ||
219 ftype == F_TABANGLES)
221 (*nft)++;
222 srenew(ft, *nft);
223 ft[*nft-1] = ftype;
227 else
229 *mult = 4;
230 *nft = 1;
231 snew(ft, *nft);
232 switch (opt[0])
234 case 'd':
235 ft[0] = F_PDIHS;
236 break;
237 case 'i':
238 ft[0] = F_IDIHS;
239 break;
240 case 'r':
241 ft[0] = F_RBDIHS;
242 break;
243 default:
244 break;
248 return ft;
251 int gmx_mk_angndx(int argc, char *argv[])
253 static const char *desc[] = {
254 "[THISMODULE] makes an index file for calculation of",
255 "angle distributions etc. It uses a run input file ([REF].tpx[ref]) for the",
256 "definitions of the angles, dihedrals etc."
258 static const char *opt[] = { NULL, "angle", "dihedral", "improper", "ryckaert-bellemans", NULL };
259 static gmx_bool bH = TRUE;
260 static real hq = -1;
261 t_pargs pa[] = {
262 { "-type", FALSE, etENUM, {opt},
263 "Type of angle" },
264 { "-hyd", FALSE, etBOOL, {&bH},
265 "Include angles with atoms with mass < 1.5" },
266 { "-hq", FALSE, etREAL, {&hq},
267 "Ignore angles with atoms with mass < 1.5 and magnitude of their charge less than this value" }
270 output_env_t oenv;
271 FILE *out;
272 t_topology *top;
273 int i, j, ntype;
274 int nft = 0, *ft, mult = 0;
275 int **index;
276 int *ft_ind;
277 int *nr;
278 char **grpnames;
279 t_filenm fnm[] = {
280 { efTPR, NULL, NULL, ffREAD },
281 { efNDX, NULL, "angle", ffWRITE }
283 #define NFILE asize(fnm)
285 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
286 asize(desc), desc, 0, NULL, &oenv))
288 return 0;
291 GMX_RELEASE_ASSERT(opt[0] != 0, "Options inconsistency; opt[0] is NULL");
293 ft = select_ftype(opt[0], &nft, &mult);
295 top = read_top(ftp2fn(efTPR, NFILE, fnm), NULL);
297 ntype = calc_ntype(nft, ft, &(top->idef));
298 snew(grpnames, ntype);
299 snew(ft_ind, top->idef.ntypes);
300 fill_ft_ind(nft, ft, &top->idef, ft_ind, grpnames);
302 snew(nr, ntype);
303 snew(index, ntype);
304 fill_ang(nft, ft, mult, nr, index, ft_ind, top, !bH, hq);
306 out = ftp2FILE(efNDX, NFILE, fnm, "w");
307 for (i = 0; (i < ntype); i++)
309 if (nr[i] > 0)
311 fprintf(out, "[ %s ]\n", grpnames[i]);
312 for (j = 0; (j < nr[i]*mult); j++)
314 fprintf(out, " %5d", index[i][j]+1);
315 if ((j % 12) == 11)
317 fprintf(out, "\n");
320 fprintf(out, "\n");
323 gmx_ffclose(out);
325 return 0;