2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
41 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/fileio/trxio.h"
43 #include "gromacs/gmxana/gmx_ana.h"
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/legacyheaders/typedefs.h"
46 #include "gromacs/utility/cstringutil.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/futil.h"
49 #include "gromacs/utility/gmxassert.h"
50 #include "gromacs/utility/smalloc.h"
52 static int calc_ntype(int nft
, int *ft
, t_idef
*idef
)
56 for (i
= 0; (i
< idef
->ntypes
); i
++)
58 for (f
= 0; f
< nft
; f
++)
60 if (idef
->functype
[i
] == ft
[f
])
70 static void fill_ft_ind(int nft
, int *ft
, t_idef
*idef
,
71 int ft_ind
[], char *grpnames
[])
74 int i
, f
, ftype
, ind
= 0;
76 /* Loop over all the function types in the topology */
77 for (i
= 0; (i
< idef
->ntypes
); i
++)
80 /* Check all the selected function types */
81 for (f
= 0; f
< nft
; f
++)
84 if (idef
->functype
[i
] == ftype
)
90 sprintf(buf
, "Theta=%.1f_%.2f", idef
->iparams
[i
].harmonic
.rA
,
91 idef
->iparams
[i
].harmonic
.krA
);
94 sprintf(buf
, "Cos_th=%.1f_%.2f", idef
->iparams
[i
].harmonic
.rA
,
95 idef
->iparams
[i
].harmonic
.krA
);
98 sprintf(buf
, "UB_th=%.1f_%.2f2f", idef
->iparams
[i
].u_b
.thetaA
,
99 idef
->iparams
[i
].u_b
.kthetaA
);
101 case F_QUARTIC_ANGLES
:
102 sprintf(buf
, "Q_th=%.1f_%.2f_%.2f", idef
->iparams
[i
].qangle
.theta
,
103 idef
->iparams
[i
].qangle
.c
[0], idef
->iparams
[i
].qangle
.c
[1]);
106 sprintf(buf
, "Table=%d_%.2f", idef
->iparams
[i
].tab
.table
,
107 idef
->iparams
[i
].tab
.kA
);
110 sprintf(buf
, "Phi=%.1f_%d_%.2f", idef
->iparams
[i
].pdihs
.phiA
,
111 idef
->iparams
[i
].pdihs
.mult
, idef
->iparams
[i
].pdihs
.cpA
);
114 sprintf(buf
, "Xi=%.1f_%.2f", idef
->iparams
[i
].harmonic
.rA
,
115 idef
->iparams
[i
].harmonic
.krA
);
118 sprintf(buf
, "RB-A1=%.2f", idef
->iparams
[i
].rbdihs
.rbcA
[1]);
121 sprintf(buf
, "Theta=%.1f_%.2f", idef
->iparams
[i
].harmonic
.rA
,
122 idef
->iparams
[i
].harmonic
.krA
);
125 sprintf(buf
, "Theta=%.1f_%.2f", idef
->iparams
[i
].harmonic
.rA
,
126 idef
->iparams
[i
].harmonic
.krA
);
129 sprintf(buf
, "CBT-A1=%.2f", idef
->iparams
[i
].cbtdihs
.cbtcA
[1]);
133 gmx_fatal(FARGS
, "Unsupported function type '%s' selected",
134 interaction_function
[ftype
].longname
);
136 grpnames
[ind
] = gmx_strdup(buf
);
143 static void fill_ang(int nft
, int *ft
, int fac
,
144 int nr
[], int *index
[], int ft_ind
[], t_topology
*top
,
145 gmx_bool bNoH
, real hq
)
147 int f
, ftype
, i
, j
, indg
, nr_fac
;
155 atom
= top
->atoms
.atom
;
157 for (f
= 0; f
< nft
; f
++)
160 ia
= idef
->il
[ftype
].iatoms
;
161 for (i
= 0; (i
< idef
->il
[ftype
].nr
); )
163 indg
= ft_ind
[ia
[0]];
166 gmx_incons("Routine fill_ang");
171 for (j
= 0; j
< fac
; j
++)
173 if (atom
[ia
[1+j
]].m
< 1.5)
181 for (j
= 0; j
< fac
; j
++)
183 if (atom
[ia
[1+j
]].m
< 1.5 && std::abs(atom
[ia
[1+j
]].q
) < hq
)
191 if (nr
[indg
] % 1000 == 0)
193 srenew(index
[indg
], fac
*(nr
[indg
]+1000));
195 nr_fac
= fac
*nr
[indg
];
196 for (j
= 0; (j
< fac
); j
++)
198 index
[indg
][nr_fac
+j
] = ia
[j
+1];
202 ia
+= interaction_function
[ftype
].nratoms
+1;
203 i
+= interaction_function
[ftype
].nratoms
+1;
208 static int *select_ftype(const char *opt
, int *nft
, int *mult
)
210 int *ft
= NULL
, ftype
;
215 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
217 if ((interaction_function
[ftype
].flags
& IF_ATYPE
) ||
218 ftype
== F_TABANGLES
)
250 int gmx_mk_angndx(int argc
, char *argv
[])
252 static const char *desc
[] = {
253 "[THISMODULE] makes an index file for calculation of",
254 "angle distributions etc. It uses a run input file ([REF].tpx[ref]) for the",
255 "definitions of the angles, dihedrals etc."
257 static const char *opt
[] = { NULL
, "angle", "dihedral", "improper", "ryckaert-bellemans", NULL
};
258 static gmx_bool bH
= TRUE
;
261 { "-type", FALSE
, etENUM
, {opt
},
263 { "-hyd", FALSE
, etBOOL
, {&bH
},
264 "Include angles with atoms with mass < 1.5" },
265 { "-hq", FALSE
, etREAL
, {&hq
},
266 "Ignore angles with atoms with mass < 1.5 and magnitude of their charge less than this value" }
273 int nft
= 0, *ft
, mult
= 0;
279 { efTPR
, NULL
, NULL
, ffREAD
},
280 { efNDX
, NULL
, "angle", ffWRITE
}
282 #define NFILE asize(fnm)
284 if (!parse_common_args(&argc
, argv
, 0, NFILE
, fnm
, asize(pa
), pa
,
285 asize(desc
), desc
, 0, NULL
, &oenv
))
290 GMX_RELEASE_ASSERT(opt
[0] != 0, "Options inconsistency; opt[0] is NULL");
292 ft
= select_ftype(opt
[0], &nft
, &mult
);
294 top
= read_top(ftp2fn(efTPR
, NFILE
, fnm
), NULL
);
296 ntype
= calc_ntype(nft
, ft
, &(top
->idef
));
297 snew(grpnames
, ntype
);
298 snew(ft_ind
, top
->idef
.ntypes
);
299 fill_ft_ind(nft
, ft
, &top
->idef
, ft_ind
, grpnames
);
303 fill_ang(nft
, ft
, mult
, nr
, index
, ft_ind
, top
, !bH
, hq
);
305 out
= ftp2FILE(efNDX
, NFILE
, fnm
, "w");
306 for (i
= 0; (i
< ntype
); i
++)
310 fprintf(out
, "[ %s ]\n", grpnames
[i
]);
311 for (j
= 0; (j
< nr
[i
]*mult
); j
++)
313 fprintf(out
, " %5d", index
[i
][j
]+1);