3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
50 #include "gmx_fatal.h"
53 int gmx_genpr(int argc, char *argv[])
55 const char *desc[] = {
56 "[TT]genrestr[tt] produces an include file for a topology containing",
57 "a list of atom numbers and three force constants for the",
58 "[IT]x[it]-, [IT]y[it]-, and [IT]z[it]-direction. A single isotropic force constant may",
59 "be given on the command line instead of three components.[PAR]",
60 "WARNING: position restraints only work for the one molecule at a time.",
61 "Position restraints are interactions within molecules, therefore",
62 "they should be included within the correct [TT][ moleculetype ][tt]",
63 "block in the topology. Since the atom numbers in every moleculetype",
64 "in the topology start at 1 and the numbers in the input file for",
65 "[TT]genrestr[tt] number consecutively from 1, [TT]genrestr[tt] will only",
66 "produce a useful file for the first molecule.[PAR]",
67 "The [TT]-of[tt] option produces an index file that can be used for",
68 "freezing atoms. In this case, the input file must be a [TT].pdb[tt] file.[PAR]",
69 "With the [TT]-disre[tt] option, half a matrix of distance restraints",
70 "is generated instead of position restraints. With this matrix, that",
71 "one typically would apply to C[GRK]alpha[grk] atoms in a protein, one can",
72 "maintain the overall conformation of a protein without tieing it to",
73 "a specific position (as with position restraints)."
75 static rvec fc = {1000.0, 1000.0, 1000.0};
76 static real freeze_level = 0.0;
77 static real disre_dist = 0.1;
78 static real disre_frac = 0.0;
79 static real disre_up2 = 1.0;
80 static gmx_bool bDisre = FALSE;
81 static gmx_bool bConstr = FALSE;
82 static real cutoff = -1.0;
85 { "-fc", FALSE, etRVEC, {fc},
86 "Force constants (kJ/mol nm^2)" },
87 { "-freeze", FALSE, etREAL, {&freeze_level},
88 "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" },
89 { "-disre", FALSE, etBOOL, {&bDisre},
90 "Generate a distance restraint matrix for all the atoms in index" },
91 { "-disre_dist", FALSE, etREAL, {&disre_dist},
92 "Distance range around the actual distance for generating distance restraints" },
93 { "-disre_frac", FALSE, etREAL, {&disre_frac},
94 "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." },
95 { "-disre_up2", FALSE, etREAL, {&disre_up2},
96 "Distance between upper bound for distance restraints, and the distance at which the force becomes constant (see manual)" },
97 { "-cutoff", FALSE, etREAL, {&cutoff},
98 "Only generate distance restraints for atoms pairs within cutoff (nm)" },
99 { "-constr", FALSE, etBOOL, {&bConstr},
100 "Generate a constraint matrix rather than distance restraints. Constraints of type 2 will be generated that do generate exclusions." }
102 #define npargs asize(pa)
105 t_atoms *atoms = NULL;
111 const char *xfn, *nfn;
116 rvec dx, *x = NULL, *v = NULL;
119 { efSTX, "-f", NULL, ffREAD },
120 { efNDX, "-n", NULL, ffOPTRD },
121 { efITP, "-o", "posre", ffWRITE },
122 { efNDX, "-of", "freeze", ffOPTWR }
124 #define NFILE asize(fnm)
126 parse_common_args(&argc, argv, 0, NFILE, fnm, npargs, pa,
127 asize(desc), desc, 0, NULL, &oenv);
129 bFreeze = opt2bSet("-of", NFILE, fnm) || opt2parg_bSet("-freeze", asize(pa), pa);
130 bDisre = bDisre || opt2parg_bSet("-disre_dist", npargs, pa);
131 xfn = opt2fn_null("-f", NFILE, fnm);
132 nfn = opt2fn_null("-n", NFILE, fnm);
134 if (( nfn == NULL ) && ( xfn == NULL))
136 gmx_fatal(FARGS, "no index file and no structure file supplied");
139 if ((disre_frac < 0) || (disre_frac >= 1))
141 gmx_fatal(FARGS, "disre_frac should be between 0 and 1");
145 gmx_fatal(FARGS, "disre_dist should be >= 0");
151 get_stx_coordnum(xfn, &(atoms->nr));
152 init_t_atoms(atoms, atoms->nr, TRUE);
155 fprintf(stderr, "\nReading structure file\n");
156 read_stx_conf(xfn, title, atoms, x, v, NULL, box);
161 if (!atoms || !atoms->pdbinfo)
163 gmx_fatal(FARGS, "No B-factors in input file %s, use a pdb file next time.",
167 out = opt2FILE("-of", NFILE, fnm, "w");
168 fprintf(out, "[ freeze ]\n");
169 for (i = 0; (i < atoms->nr); i++)
171 if (atoms->pdbinfo[i].bfac <= freeze_level)
173 fprintf(out, "%d\n", i+1);
178 else if ((bDisre || bConstr) && x)
180 printf("Select group to generate %s matrix from\n",
181 bConstr ? "constraint" : "distance restraint");
182 get_index(atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
184 out = ftp2FILE(efITP, NFILE, fnm, "w");
187 fprintf(out, "; constraints for %s of %s\n\n", gn_grp, title);
188 fprintf(out, "[ constraints ]\n");
189 fprintf(out, ";%4s %5s %1s %10s\n", "i", "j", "tp", "dist");
193 fprintf(out, "; distance restraints for %s of %s\n\n", gn_grp, title);
194 fprintf(out, "[ distance_restraints ]\n");
195 fprintf(out, ";%4s %5s %1s %5s %10s %10s %10s %10s %10s\n", "i", "j", "?",
196 "label", "funct", "lo", "up1", "up2", "weight");
198 for (i = k = 0; i < igrp; i++)
200 for (j = i+1; j < igrp; j++, k++)
202 rvec_sub(x[ind_grp[i]], x[ind_grp[j]], dx);
206 fprintf(out, "%5d %5d %1d %10g\n", ind_grp[i]+1, ind_grp[j]+1, 2, d);
210 if (cutoff < 0 || d < cutoff)
214 dd = min(disre_dist, disre_frac*d);
222 fprintf(out, "%5d %5d %1d %5d %10d %10g %10g %10g %10g\n",
223 ind_grp[i]+1, ind_grp[j]+1, 1, k, 1,
233 printf("Select group to position restrain\n");
234 get_index(atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
236 out = ftp2FILE(efITP, NFILE, fnm, "w");
237 fprintf(out, "; position restraints for %s of %s\n\n", gn_grp, title);
238 fprintf(out, "[ position_restraints ]\n");
239 fprintf(out, ";%3s %5s %9s %10s %10s\n", "i", "funct", "fcx", "fcy", "fcz");
240 for (i = 0; i < igrp; i++)
242 fprintf(out, "%4d %4d %10g %10g %10g\n",
243 ind_grp[i]+1, 1, fc[XX], fc[YY], fc[ZZ]);