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