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