Merge branch 'release-4-6' into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_mk_angndx.c
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  * 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <math.h>
43 #include "typedefs.h"
44 #include "smalloc.h"
45 #include "copyrite.h"
46 #include "statutil.h"
47 #include "macros.h"
48 #include "string2.h"
49 #include "futil.h"
50 #include "gmx_fatal.h"
51
52 static int calc_ntype(int nft, int *ft, t_idef *idef)
53 {
54     int  i, f, nf = 0;
55
56     for (i = 0; (i < idef->ntypes); i++)
57     {
58         for (f = 0; f < nft; f++)
59         {
60             if (idef->functype[i] == ft[f])
61             {
62                 nf++;
63             }
64         }
65     }
66
67     return nf;
68 }
69
70 static void fill_ft_ind(int nft, int *ft, t_idef *idef,
71                         int ft_ind[], char *grpnames[])
72 {
73     char buf[125];
74     int  i, f, ftype, ind = 0;
75
76     /* Loop over all the function types in the topology */
77     for (i = 0; (i < idef->ntypes); i++)
78     {
79         ft_ind[i] = -1;
80         /* Check all the selected function types */
81         for (f = 0; f < nft; f++)
82         {
83             ftype = ft[f];
84             if (idef->functype[i] == ftype)
85             {
86                 ft_ind[i] = ind;
87                 switch (ftype)
88                 {
89                     case F_ANGLES:
90                         sprintf(buf, "Theta=%.1f_%.2f", idef->iparams[i].harmonic.rA,
91                                 idef->iparams[i].harmonic.krA);
92                         break;
93                     case F_G96ANGLES:
94                         sprintf(buf, "Cos_th=%.1f_%.2f", idef->iparams[i].harmonic.rA,
95                                 idef->iparams[i].harmonic.krA);
96                         break;
97                     case F_UREY_BRADLEY:
98                         sprintf(buf, "UB_th=%.1f_%.2f2f", idef->iparams[i].u_b.thetaA,
99                                 idef->iparams[i].u_b.kthetaA);
100                         break;
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]);
104                         break;
105                     case F_TABANGLES:
106                         sprintf(buf, "Table=%d_%.2f", idef->iparams[i].tab.table,
107                                 idef->iparams[i].tab.kA);
108                         break;
109                     case F_PDIHS:
110                         sprintf(buf, "Phi=%.1f_%d_%.2f", idef->iparams[i].pdihs.phiA,
111                                 idef->iparams[i].pdihs.mult, idef->iparams[i].pdihs.cpA);
112                         break;
113                     case F_IDIHS:
114                         sprintf(buf, "Xi=%.1f_%.2f", idef->iparams[i].harmonic.rA,
115                                 idef->iparams[i].harmonic.krA);
116                         break;
117                     case F_RBDIHS:
118                         sprintf(buf, "RB-A1=%.2f", idef->iparams[i].rbdihs.rbcA[1]);
119                         break;
120                     default:
121                         gmx_fatal(FARGS, "Unsupported function type '%s' selected",
122                                   interaction_function[ftype].longname);
123                 }
124                 grpnames[ind] = strdup(buf);
125                 ind++;
126             }
127         }
128     }
129 }
130
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)
134 {
135     int         f, ftype, i, j, indg, nr_fac;
136     gmx_bool    bUse;
137     t_idef     *idef;
138     t_atom     *atom;
139     t_iatom    *ia;
140
141
142     idef = &top->idef;
143     atom = top->atoms.atom;
144
145     for (f = 0; f < nft; f++)
146     {
147         ftype = ft[f];
148         ia    = idef->il[ftype].iatoms;
149         for (i = 0; (i < idef->il[ftype].nr); )
150         {
151             indg = ft_ind[ia[0]];
152             if (indg == -1)
153             {
154                 gmx_incons("Routine fill_ang");
155             }
156             bUse = TRUE;
157             if (bNoH)
158             {
159                 for (j = 0; j < fac; j++)
160                 {
161                     if (atom[ia[1+j]].m < 1.5)
162                     {
163                         bUse = FALSE;
164                     }
165                 }
166             }
167             if (hq)
168             {
169                 for (j = 0; j < fac; j++)
170                 {
171                     if (atom[ia[1+j]].m < 1.5 && fabs(atom[ia[1+j]].q) < hq)
172                     {
173                         bUse = FALSE;
174                     }
175                 }
176             }
177             if (bUse)
178             {
179                 if (nr[indg] % 1000 == 0)
180                 {
181                     srenew(index[indg], fac*(nr[indg]+1000));
182                 }
183                 nr_fac = fac*nr[indg];
184                 for (j = 0; (j < fac); j++)
185                 {
186                     index[indg][nr_fac+j] = ia[j+1];
187                 }
188                 nr[indg]++;
189             }
190             ia += interaction_function[ftype].nratoms+1;
191             i  += interaction_function[ftype].nratoms+1;
192         }
193     }
194 }
195
196 static int *select_ftype(const char *opt, int *nft, int *mult)
197 {
198     int *ft = NULL, ftype;
199
200     if (opt[0] == 'a')
201     {
202         *mult = 3;
203         for (ftype = 0; ftype < F_NRE; ftype++)
204         {
205             if ((interaction_function[ftype].flags & IF_ATYPE) ||
206                 ftype == F_TABANGLES)
207             {
208                 (*nft)++;
209                 srenew(ft, *nft);
210                 ft[*nft-1] = ftype;
211             }
212         }
213     }
214     else
215     {
216         *mult = 4;
217         *nft  = 1;
218         snew(ft, *nft);
219         switch (opt[0])
220         {
221             case 'd':
222                 ft[0] = F_PDIHS;
223                 break;
224             case 'i':
225                 ft[0] = F_IDIHS;
226                 break;
227             case 'r':
228                 ft[0] = F_RBDIHS;
229                 break;
230             default:
231                 break;
232         }
233     }
234
235     return ft;
236 }
237
238 int gmx_mk_angndx(int argc, char *argv[])
239 {
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."
244     };
245     static const char *opt[] = { NULL, "angle", "dihedral", "improper", "ryckaert-bellemans", NULL };
246     static gmx_bool    bH    = TRUE;
247     static real        hq    = -1;
248     t_pargs            pa[]  = {
249         { "-type", FALSE, etENUM, {opt},
250           "Type of angle" },
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" }
255     };
256
257     output_env_t       oenv;
258     FILE              *out;
259     t_topology        *top;
260     int                i, j, ntype;
261     int                nft = 0, *ft, mult = 0;
262     int              **index;
263     int               *ft_ind;
264     int               *nr;
265     char             **grpnames;
266     t_filenm           fnm[] = {
267         { efTPX, NULL, NULL, ffREAD  },
268         { efNDX, NULL, "angle", ffWRITE }
269     };
270 #define NFILE asize(fnm)
271
272     parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
273                       asize(desc), desc, 0, NULL, &oenv);
274
275
276     ft = select_ftype(opt[0], &nft, &mult);
277
278     top = read_top(ftp2fn(efTPX, NFILE, fnm), NULL);
279
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);
284
285     snew(nr, ntype);
286     snew(index, ntype);
287     fill_ang(nft, ft, mult, nr, index, ft_ind, top, !bH, hq);
288
289     out = ftp2FILE(efNDX, NFILE, fnm, "w");
290     for (i = 0; (i < ntype); i++)
291     {
292         if (nr[i] > 0)
293         {
294             fprintf(out, "[ %s ]\n", grpnames[i]);
295             for (j = 0; (j < nr[i]*mult); j++)
296             {
297                 fprintf(out, " %5d", index[i][j]+1);
298                 if ((j % 12) == 11)
299                 {
300                     fprintf(out, "\n");
301                 }
302             }
303             fprintf(out, "\n");
304         }
305     }
306     ffclose(out);
307
308     thanx(stderr);
309
310     return 0;
311 }