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