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) 2013, 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.
47 #include "gromacs/fileio/confio.h"
49 #include "gromacs/fileio/tpxio.h"
50 #include "gromacs/fileio/trxio.h"
53 #include "hackblock.h"
55 int gmx_protonate(int argc, char *argv[])
57 const char *desc[] = {
58 "[THISMODULE] reads (a) conformation(s) and adds all missing",
59 "hydrogens as defined in [TT]gmx2.ff/aminoacids.hdb[tt]. If only [TT]-s[tt] is",
60 "specified, this conformation will be protonated, if also [TT]-f[tt]",
61 "is specified, the conformation(s) will be read from this file, ",
62 "which can be either a single conformation or a trajectory.",
64 "If a [TT].pdb[tt] file is supplied, residue names might not correspond to",
65 "to the GROMACS naming conventions, in which case these residues will",
66 "probably not be properly protonated.",
68 "If an index file is specified, please note that the atom numbers",
69 "should correspond to the [BB]protonated[bb] state."
77 t_atoms *atoms, *iatoms;
84 int nidx, natoms, natoms_out;
87 gmx_bool bReadMultiple;
90 const char *bugs[] = {
91 "For the moment, only .pdb files are accepted to the -s flag"
95 { efTPS, NULL, NULL, ffREAD },
96 { efTRX, "-f", NULL, ffOPTRD },
97 { efNDX, NULL, NULL, ffOPTRD },
98 { efTRO, "-o", "protonated", ffWRITE }
100 #define NFILE asize(fnm)
102 if (!parse_common_args(&argc, argv, PCA_CAN_TIME,
103 NFILE, fnm, 0, NULL, asize(desc), desc, asize(bugs), bugs, &oenv))
108 infile = opt2fn("-s", NFILE, fnm);
109 read_tps_conf(infile, title, &top, &ePBC, &x, NULL, box, FALSE);
110 atoms = &(top.atoms);
111 printf("Select group to process:\n");
112 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidx, &index, &grpnm);
113 bReadMultiple = opt2bSet("-f", NFILE, fnm);
116 infile = opt2fn("-f", NFILE, fnm);
117 if (!read_first_frame(oenv, &status, infile, &fr, TRX_NEED_X ) )
119 gmx_fatal(FARGS, "cannot read coordinate file %s", infile);
125 clear_trxframe(&fr, TRUE);
126 fr.natoms = atoms->nr;
132 copy_mat(box, fr.box);
139 gmx_fatal(FARGS, "no atoms in coordinate file %s", infile);
142 if (natoms > atoms->nr)
144 gmx_fatal(FARGS, "topology with %d atoms does not match "
145 "coordinates with %d atoms", atoms->nr, natoms);
148 for (i = 0; i < nidx; i++)
150 if (index[i] > natoms)
152 gmx_fatal(FARGS, "An atom number in group %s is larger than the number of "
153 "atoms (%d) in the coordinate file %s", grpnm, natoms, infile);
157 /* get indexed copy of atoms */
159 init_t_atoms(iatoms, nidx, FALSE);
160 snew(iatoms->atom, iatoms->nr);
162 for (i = 0; i < nidx; i++)
164 iatoms->atom[i] = atoms->atom[index[i]];
165 iatoms->atomname[i] = atoms->atomname[index[i]];
166 if (i > 0 && (atoms->atom[index[i]].resind != atoms->atom[index[i-1]].resind) )
170 iatoms->atom[i].resind = resind;
171 iatoms->resinfo[resind] = atoms->resinfo[atoms->atom[index[i]].resind];
172 /* allocate some space for the rtp name and copy from name */
173 snew(iatoms->resinfo[resind].rtp, 1);
174 *iatoms->resinfo[resind].rtp = gmx_strdup(*atoms->resinfo[resind].name);
176 iatoms->nres = max(iatoms->nres, iatoms->atom[i].resind+1);
179 init_t_protonate(&protdata);
181 out = open_trx(opt2fn("-o", NFILE, fnm), "w");
188 fprintf(debug, "FRAME %d (%d %g)\n", frame, fr.step, fr.time);
190 /* get indexed copy of x */
191 for (i = 0; i < nidx; i++)
193 copy_rvec(fr.x[index[i]], ix[i]);
196 natoms_out = protonate(&iatoms, &ix, &protdata);
198 /* setup output frame */
200 frout.natoms = natoms_out;
202 frout.atoms = iatoms;
208 write_trxframe(out, &frout, NULL);
211 while (bReadMultiple && read_next_frame(oenv, status, &fr) );