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