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 by the GROMACS development team.
7 * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
11 *
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
16 *
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
21 *
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 *
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
34 *
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
37 */
38 #include "gmxpre.h"
39
40 #include "mk_angndx.h"
41
42 #include <cmath>
43
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/topology/ifunc.h"
47 #include "gromacs/topology/topology.h"
48 #include "gromacs/utility/arraysize.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/futil.h"
52 #include "gromacs/utility/gmxassert.h"
53 #include "gromacs/utility/smalloc.h"
54
calc_ntype(int nft,const int * ft,const t_idef * idef)55 static int calc_ntype(int nft, const int* ft, const t_idef* idef)
56 {
57 int i, f, nf = 0;
58
59 for (i = 0; (i < idef->ntypes); i++)
60 {
61 for (f = 0; f < nft; f++)
62 {
63 if (idef->functype[i] == ft[f])
64 {
65 nf++;
66 }
67 }
68 }
69
70 return nf;
71 }
72
fill_ft_ind(int nft,const int * ft,const t_idef * idef,int ft_ind[],char * grpnames[])73 static void fill_ft_ind(int nft, const int* ft, const t_idef* idef, 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 // Fall through intended
124 case F_RESTRDIHS:
125 sprintf(buf, "Theta=%.1f_%.2f", idef->iparams[i].harmonic.rA,
126 idef->iparams[i].harmonic.krA);
127 break;
128 case F_CBTDIHS:
129 sprintf(buf, "CBT-A1=%.2f", idef->iparams[i].cbtdihs.cbtcA[1]);
130 break;
131
132 default:
133 gmx_fatal(FARGS, "Unsupported function type '%s' selected",
134 interaction_function[ftype].longname);
135 }
136 grpnames[ind] = gmx_strdup(buf);
137 ind++;
138 }
139 }
140 }
141 }
142
fill_ang(int nft,const int * ft,int fac,int nr[],int * index[],const int ft_ind[],const t_topology * top,gmx_bool bNoH,real hq)143 static void fill_ang(int nft,
144 const int* ft,
145 int fac,
146 int nr[],
147 int* index[],
148 const int ft_ind[],
149 const t_topology* top,
150 gmx_bool bNoH,
151 real hq)
152 {
153 int f, ftype, i, j, indg, nr_fac;
154 gmx_bool bUse;
155 const t_idef* idef;
156 t_atom* atom;
157 t_iatom* ia;
158
159
160 idef = &top->idef;
161 atom = top->atoms.atom;
162
163 for (f = 0; f < nft; f++)
164 {
165 ftype = ft[f];
166 ia = idef->il[ftype].iatoms;
167 for (i = 0; (i < idef->il[ftype].nr);)
168 {
169 indg = ft_ind[ia[0]];
170 if (indg == -1)
171 {
172 gmx_incons("Routine fill_ang");
173 }
174 bUse = TRUE;
175 if (bNoH)
176 {
177 for (j = 0; j < fac; j++)
178 {
179 if (atom[ia[1 + j]].m < 1.5)
180 {
181 bUse = FALSE;
182 }
183 }
184 }
185 if (hq != 0.0F)
186 {
187 for (j = 0; j < fac; j++)
188 {
189 if (atom[ia[1 + j]].m < 1.5 && std::abs(atom[ia[1 + j]].q) < hq)
190 {
191 bUse = FALSE;
192 }
193 }
194 }
195 if (bUse)
196 {
197 if (nr[indg] % 1000 == 0)
198 {
199 srenew(index[indg], fac * (nr[indg] + 1000));
200 }
201 nr_fac = fac * nr[indg];
202 for (j = 0; (j < fac); j++)
203 {
204 index[indg][nr_fac + j] = ia[j + 1];
205 }
206 nr[indg]++;
207 }
208 ia += interaction_function[ftype].nratoms + 1;
209 i += interaction_function[ftype].nratoms + 1;
210 }
211 }
212 }
213
select_ftype(const char * opt,int * nft,int * mult)214 static int* select_ftype(const char* opt, int* nft, int* mult)
215 {
216 int *ft = nullptr, ftype;
217
218 if (opt[0] == 'a')
219 {
220 *mult = 3;
221 for (ftype = 0; ftype < F_NRE; ftype++)
222 {
223 if ((interaction_function[ftype].flags & IF_ATYPE) || ftype == F_TABANGLES)
224 {
225 (*nft)++;
226 srenew(ft, *nft);
227 ft[*nft - 1] = ftype;
228 }
229 }
230 }
231 else
232 {
233 *mult = 4;
234 *nft = 1;
235 snew(ft, *nft);
236 switch (opt[0])
237 {
238 case 'd': ft[0] = F_PDIHS; break;
239 case 'i': ft[0] = F_IDIHS; break;
240 case 'r': ft[0] = F_RBDIHS; break;
241 default: break;
242 }
243 }
244
245 return ft;
246 }
247
gmx_mk_angndx(int argc,char * argv[])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 ([REF].tpx[ref]) for the",
253 "definitions of the angles, dihedrals etc."
254 };
255 static const char* opt[] = { nullptr, "angle", "dihedral", "improper", "ryckaert-bellemans",
256 nullptr };
257 static gmx_bool bH = TRUE;
258 static real hq = -1;
259 t_pargs pa[] = { { "-type", FALSE, etENUM, { opt }, "Type of angle" },
260 { "-hyd", FALSE, etBOOL, { &bH }, "Include angles with atoms with mass < 1.5" },
261 { "-hq",
262 FALSE,
263 etREAL,
264 { &hq },
265 "Ignore angles with atoms with mass < 1.5 and magnitude of their charge "
266 "less than this value" } };
267
268 gmx_output_env_t* oenv;
269 FILE* out;
270 t_topology* top;
271 int i, j, ntype;
272 int nft = 0, *ft, mult = 0;
273 int** index;
274 int* ft_ind;
275 int* nr;
276 char** grpnames;
277 t_filenm fnm[] = { { efTPR, nullptr, nullptr, ffREAD }, { efNDX, nullptr, "angle", ffWRITE } };
278 #define NFILE asize(fnm)
279
280 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
281 {
282 return 0;
283 }
284
285 GMX_RELEASE_ASSERT(opt[0] != nullptr, "Options inconsistency; opt[0] is NULL");
286
287 ft = select_ftype(opt[0], &nft, &mult);
288
289 top = read_top(ftp2fn(efTPR, NFILE, fnm), nullptr);
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 }
321