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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
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.
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.
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.
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.
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.
54 #include "gmx_fatal.h"
57 int gmx_genpr(int argc, char *argv[])
59 const char *desc[] = {
60 "[TT]genrestr[tt] produces an include file for a topology containing",
61 "a list of atom numbers and three force constants for the",
62 "[IT]x[it]-, [IT]y[it]-, and [IT]z[it]-direction. A single isotropic force constant may",
63 "be given on the command line instead of three components.[PAR]",
64 "WARNING: position restraints only work for the one molecule at a time.",
65 "Position restraints are interactions within molecules, therefore",
66 "they should be included within the correct [TT][ moleculetype ][tt]",
67 "block in the topology. Since the atom numbers in every moleculetype",
68 "in the topology start at 1 and the numbers in the input file for",
69 "[TT]genrestr[tt] number consecutively from 1, [TT]genrestr[tt] will only",
70 "produce a useful file for the first molecule.[PAR]",
71 "The [TT]-of[tt] option produces an index file that can be used for",
72 "freezing atoms. In this case, the input file must be a [TT].pdb[tt] file.[PAR]",
73 "With the [TT]-disre[tt] option, half a matrix of distance restraints",
74 "is generated instead of position restraints. With this matrix, that",
75 "one typically would apply to C[GRK]alpha[grk] atoms in a protein, one can",
76 "maintain the overall conformation of a protein without tieing it to",
77 "a specific position (as with position restraints)."
79 static rvec fc = {1000.0, 1000.0, 1000.0};
80 static real freeze_level = 0.0;
81 static real disre_dist = 0.1;
82 static real disre_frac = 0.0;
83 static real disre_up2 = 1.0;
84 static gmx_bool bDisre = FALSE;
85 static gmx_bool bConstr = FALSE;
86 static real cutoff = -1.0;
89 { "-fc", FALSE, etRVEC, {fc},
90 "Force constants (kJ/mol nm^2)" },
91 { "-freeze", FALSE, etREAL, {&freeze_level},
92 "If the [TT]-of[tt] option or this one is given an index file will be written containing atom numbers of all atoms that have a B-factor less than the level given here" },
93 { "-disre", FALSE, etBOOL, {&bDisre},
94 "Generate a distance restraint matrix for all the atoms in index" },
95 { "-disre_dist", FALSE, etREAL, {&disre_dist},
96 "Distance range around the actual distance for generating distance restraints" },
97 { "-disre_frac", FALSE, etREAL, {&disre_frac},
98 "Fraction of distance to be used as interval rather than a fixed distance. If the fraction of the distance that you specify here is less than the distance given in the previous option, that one is used instead." },
99 { "-disre_up2", FALSE, etREAL, {&disre_up2},
100 "Distance between upper bound for distance restraints, and the distance at which the force becomes constant (see manual)" },
101 { "-cutoff", FALSE, etREAL, {&cutoff},
102 "Only generate distance restraints for atoms pairs within cutoff (nm)" },
103 { "-constr", FALSE, etBOOL, {&bConstr},
104 "Generate a constraint matrix rather than distance restraints. Constraints of type 2 will be generated that do generate exclusions." }
106 #define npargs asize(pa)
109 t_atoms *atoms = NULL;
115 const char *xfn, *nfn;
120 rvec dx, *x = NULL, *v = NULL;
123 { efSTX, "-f", NULL, ffREAD },
124 { efNDX, "-n", NULL, ffOPTRD },
125 { efITP, "-o", "posre", ffWRITE },
126 { efNDX, "-of", "freeze", ffOPTWR }
128 #define NFILE asize(fnm)
130 CopyRight(stderr, argv[0]);
131 parse_common_args(&argc, argv, 0, NFILE, fnm, npargs, pa,
132 asize(desc), desc, 0, NULL, &oenv);
134 bFreeze = opt2bSet("-of", NFILE, fnm) || opt2parg_bSet("-freeze", asize(pa), pa);
135 bDisre = bDisre || opt2parg_bSet("-disre_dist", npargs, pa);
136 xfn = opt2fn_null("-f", NFILE, fnm);
137 nfn = opt2fn_null("-n", NFILE, fnm);
139 if (( nfn == NULL ) && ( xfn == NULL))
141 gmx_fatal(FARGS, "no index file and no structure file supplied");
144 if ((disre_frac < 0) || (disre_frac >= 1))
146 gmx_fatal(FARGS, "disre_frac should be between 0 and 1");
150 gmx_fatal(FARGS, "disre_dist should be >= 0");
156 get_stx_coordnum(xfn, &(atoms->nr));
157 init_t_atoms(atoms, atoms->nr, TRUE);
160 fprintf(stderr, "\nReading structure file\n");
161 read_stx_conf(xfn, title, atoms, x, v, NULL, box);
166 if (!atoms || !atoms->pdbinfo)
168 gmx_fatal(FARGS, "No B-factors in input file %s, use a pdb file next time.",
172 out = opt2FILE("-of", NFILE, fnm, "w");
173 fprintf(out, "[ freeze ]\n");
174 for (i = 0; (i < atoms->nr); i++)
176 if (atoms->pdbinfo[i].bfac <= freeze_level)
178 fprintf(out, "%d\n", i+1);
183 else if ((bDisre || bConstr) && x)
185 printf("Select group to generate %s matrix from\n",
186 bConstr ? "constraint" : "distance restraint");
187 get_index(atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
189 out = ftp2FILE(efITP, NFILE, fnm, "w");
192 fprintf(out, "; constraints for %s of %s\n\n", gn_grp, title);
193 fprintf(out, "[ constraints ]\n");
194 fprintf(out, ";%4s %5s %1s %10s\n", "i", "j", "tp", "dist");
198 fprintf(out, "; distance restraints for %s of %s\n\n", gn_grp, title);
199 fprintf(out, "[ distance_restraints ]\n");
200 fprintf(out, ";%4s %5s %1s %5s %10s %10s %10s %10s %10s\n", "i", "j", "?",
201 "label", "funct", "lo", "up1", "up2", "weight");
203 for (i = k = 0; i < igrp; i++)
205 for (j = i+1; j < igrp; j++, k++)
207 rvec_sub(x[ind_grp[i]], x[ind_grp[j]], dx);
211 fprintf(out, "%5d %5d %1d %10g\n", ind_grp[i]+1, ind_grp[j]+1, 2, d);
215 if (cutoff < 0 || d < cutoff)
219 dd = min(disre_dist, disre_frac*d);
227 fprintf(out, "%5d %5d %1d %5d %10d %10g %10g %10g %10g\n",
228 ind_grp[i]+1, ind_grp[j]+1, 1, k, 1,
238 printf("Select group to position restrain\n");
239 get_index(atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
241 out = ftp2FILE(efITP, NFILE, fnm, "w");
242 fprintf(out, "; position restraints for %s of %s\n\n", gn_grp, title);
243 fprintf(out, "[ position_restraints ]\n");
244 fprintf(out, ";%3s %5s %9s %10s %10s\n", "i", "funct", "fcx", "fcy", "fcz");
245 for (i = 0; i < igrp; i++)
247 fprintf(out, "%4d %4d %10g %10g %10g\n",
248 ind_grp[i]+1, 1, fc[XX], fc[YY], fc[ZZ]);