Apply clang-format to source tree
[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, int ft_ind[], char* grpnames[])
73 {
74     char buf[125];
75     int  i, f, ftype, ind = 0;
76
77     /* Loop over all the function types in the topology */
78     for (i = 0; (i < idef->ntypes); i++)
79     {
80         ft_ind[i] = -1;
81         /* Check all the selected function types */
82         for (f = 0; f < nft; f++)
83         {
84             ftype = ft[f];
85             if (idef->functype[i] == ftype)
86             {
87                 ft_ind[i] = ind;
88                 switch (ftype)
89                 {
90                     case F_ANGLES:
91                         sprintf(buf, "Theta=%.1f_%.2f", idef->iparams[i].harmonic.rA,
92                                 idef->iparams[i].harmonic.krA);
93                         break;
94                     case F_G96ANGLES:
95                         sprintf(buf, "Cos_th=%.1f_%.2f", idef->iparams[i].harmonic.rA,
96                                 idef->iparams[i].harmonic.krA);
97                         break;
98                     case F_UREY_BRADLEY:
99                         sprintf(buf, "UB_th=%.1f_%.2f2f", idef->iparams[i].u_b.thetaA,
100                                 idef->iparams[i].u_b.kthetaA);
101                         break;
102                     case F_QUARTIC_ANGLES:
103                         sprintf(buf, "Q_th=%.1f_%.2f_%.2f", idef->iparams[i].qangle.theta,
104                                 idef->iparams[i].qangle.c[0], idef->iparams[i].qangle.c[1]);
105                         break;
106                     case F_TABANGLES:
107                         sprintf(buf, "Table=%d_%.2f", idef->iparams[i].tab.table,
108                                 idef->iparams[i].tab.kA);
109                         break;
110                     case F_PDIHS:
111                         sprintf(buf, "Phi=%.1f_%d_%.2f", idef->iparams[i].pdihs.phiA,
112                                 idef->iparams[i].pdihs.mult, idef->iparams[i].pdihs.cpA);
113                         break;
114                     case F_IDIHS:
115                         sprintf(buf, "Xi=%.1f_%.2f", idef->iparams[i].harmonic.rA,
116                                 idef->iparams[i].harmonic.krA);
117                         break;
118                     case F_RBDIHS:
119                         sprintf(buf, "RB-A1=%.2f", idef->iparams[i].rbdihs.rbcA[1]);
120                         break;
121                     case F_RESTRANGLES:
122                         sprintf(buf, "Theta=%.1f_%.2f", idef->iparams[i].harmonic.rA,
123                                 idef->iparams[i].harmonic.krA);
124                         break;
125                     case F_RESTRDIHS:
126                         sprintf(buf, "Theta=%.1f_%.2f", idef->iparams[i].harmonic.rA,
127                                 idef->iparams[i].harmonic.krA);
128                         break;
129                     case F_CBTDIHS:
130                         sprintf(buf, "CBT-A1=%.2f", idef->iparams[i].cbtdihs.cbtcA[1]);
131                         break;
132
133                     default:
134                         gmx_fatal(FARGS, "Unsupported function type '%s' selected",
135                                   interaction_function[ftype].longname);
136                 }
137                 grpnames[ind] = gmx_strdup(buf);
138                 ind++;
139             }
140         }
141     }
142 }
143
144 static void fill_ang(int               nft,
145                      const int*        ft,
146                      int               fac,
147                      int               nr[],
148                      int*              index[],
149                      const int         ft_ind[],
150                      const t_topology* top,
151                      gmx_bool          bNoH,
152                      real              hq)
153 {
154     int           f, ftype, i, j, indg, nr_fac;
155     gmx_bool      bUse;
156     const t_idef* idef;
157     t_atom*       atom;
158     t_iatom*      ia;
159
160
161     idef = &top->idef;
162     atom = top->atoms.atom;
163
164     for (f = 0; f < nft; f++)
165     {
166         ftype = ft[f];
167         ia    = idef->il[ftype].iatoms;
168         for (i = 0; (i < idef->il[ftype].nr);)
169         {
170             indg = ft_ind[ia[0]];
171             if (indg == -1)
172             {
173                 gmx_incons("Routine fill_ang");
174             }
175             bUse = TRUE;
176             if (bNoH)
177             {
178                 for (j = 0; j < fac; j++)
179                 {
180                     if (atom[ia[1 + j]].m < 1.5)
181                     {
182                         bUse = FALSE;
183                     }
184                 }
185             }
186             if (hq != 0.0F)
187             {
188                 for (j = 0; j < fac; j++)
189                 {
190                     if (atom[ia[1 + j]].m < 1.5 && std::abs(atom[ia[1 + j]].q) < hq)
191                     {
192                         bUse = FALSE;
193                     }
194                 }
195             }
196             if (bUse)
197             {
198                 if (nr[indg] % 1000 == 0)
199                 {
200                     srenew(index[indg], fac * (nr[indg] + 1000));
201                 }
202                 nr_fac = fac * nr[indg];
203                 for (j = 0; (j < fac); j++)
204                 {
205                     index[indg][nr_fac + j] = ia[j + 1];
206                 }
207                 nr[indg]++;
208             }
209             ia += interaction_function[ftype].nratoms + 1;
210             i += interaction_function[ftype].nratoms + 1;
211         }
212     }
213 }
214
215 static int* select_ftype(const char* opt, int* nft, int* mult)
216 {
217     int *ft = nullptr, ftype;
218
219     if (opt[0] == 'a')
220     {
221         *mult = 3;
222         for (ftype = 0; ftype < F_NRE; ftype++)
223         {
224             if ((interaction_function[ftype].flags & IF_ATYPE) || ftype == F_TABANGLES)
225             {
226                 (*nft)++;
227                 srenew(ft, *nft);
228                 ft[*nft - 1] = ftype;
229             }
230         }
231     }
232     else
233     {
234         *mult = 4;
235         *nft  = 1;
236         snew(ft, *nft);
237         switch (opt[0])
238         {
239             case 'd': ft[0] = F_PDIHS; break;
240             case 'i': ft[0] = F_IDIHS; break;
241             case 'r': ft[0] = F_RBDIHS; break;
242             default: 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 ([REF].tpx[ref]) for the",
254         "definitions of the angles, dihedrals etc."
255     };
256     static const char* opt[] = { nullptr, "angle", "dihedral", "improper", "ryckaert-bellemans",
257                                  nullptr };
258     static gmx_bool    bH    = TRUE;
259     static real        hq    = -1;
260     t_pargs            pa[]  = { { "-type", FALSE, etENUM, { opt }, "Type of angle" },
261                      { "-hyd", FALSE, etBOOL, { &bH }, "Include angles with atoms with mass < 1.5" },
262                      { "-hq",
263                        FALSE,
264                        etREAL,
265                        { &hq },
266                        "Ignore angles with atoms with mass < 1.5 and magnitude of their charge "
267                        "less than this value" } };
268
269     gmx_output_env_t* oenv;
270     FILE*             out;
271     t_topology*       top;
272     int               i, j, ntype;
273     int               nft = 0, *ft, mult = 0;
274     int**             index;
275     int*              ft_ind;
276     int*              nr;
277     char**            grpnames;
278     t_filenm fnm[] = { { efTPR, nullptr, nullptr, ffREAD }, { efNDX, nullptr, "angle", ffWRITE } };
279 #define NFILE asize(fnm)
280
281     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
282     {
283         return 0;
284     }
285
286     GMX_RELEASE_ASSERT(opt[0] != nullptr, "Options inconsistency; opt[0] is NULL");
287
288     ft = select_ftype(opt[0], &nft, &mult);
289
290     top = read_top(ftp2fn(efTPR, NFILE, fnm), nullptr);
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 }