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