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,2014,2015,2017,2018, 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.
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/gmxana/gmx_ana.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/topology/index.h"
49 #include "gromacs/topology/topology.h"
50 #include "gromacs/utility/arraysize.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/smalloc.h"
56 int gmx_genpr(int argc
, char *argv
[])
58 const char *desc
[] = {
59 "[THISMODULE] produces an #include file for a topology containing",
60 "a list of atom numbers and three force constants for the",
61 "[IT]x[it]-, [IT]y[it]-, and [IT]z[it]-direction based on",
62 "the contents of the [TT]-f[tt] file. A single isotropic force constant may",
63 "be given on the command line instead of three components.[PAR]",
64 "WARNING: Position restraints are interactions within molecules, therefore",
65 "they must be included within the correct [TT][ moleculetype ][tt]",
66 "block in the topology. The atom indices within the",
67 "[TT][ position_restraints ][tt] block must be within the range of the",
68 "atom indices for that molecule type. Since the atom numbers in every",
69 "moleculetype in the topology start at 1 and the numbers in the input file",
70 "for [THISMODULE] number consecutively from 1, [THISMODULE] will only",
71 "produce a useful file for the first molecule. You may wish to",
72 "edit the resulting index file to remove the lines for later atoms,",
73 "or construct a suitable index group to provide",
74 "as input to [THISMODULE].[PAR]",
75 "The [TT]-of[tt] option produces an index file that can be used for",
76 "freezing atoms. In this case, the input file must be a [REF].pdb[ref] file.[PAR]",
77 "With the [TT]-disre[tt] option, half a matrix of distance restraints",
78 "is generated instead of position restraints. With this matrix, that",
79 "one typically would apply to C[GRK]alpha[grk] atoms in a protein, one can",
80 "maintain the overall conformation of a protein without tieing it to",
81 "a specific position (as with position restraints)."
83 static rvec fc
= {1000.0, 1000.0, 1000.0};
84 static real freeze_level
= 0.0;
85 static real disre_dist
= 0.1;
86 static real disre_frac
= 0.0;
87 static real disre_up2
= 1.0;
88 static gmx_bool bDisre
= FALSE
;
89 static gmx_bool bConstr
= FALSE
;
90 static real cutoff
= -1.0;
93 { "-fc", FALSE
, etRVEC
, {fc
},
94 "Force constants (kJ/mol nm^2)" },
95 { "-freeze", FALSE
, etREAL
, {&freeze_level
},
96 "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" },
97 { "-disre", FALSE
, etBOOL
, {&bDisre
},
98 "Generate a distance restraint matrix for all the atoms in index" },
99 { "-disre_dist", FALSE
, etREAL
, {&disre_dist
},
100 "Distance range around the actual distance for generating distance restraints" },
101 { "-disre_frac", FALSE
, etREAL
, {&disre_frac
},
102 "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." },
103 { "-disre_up2", FALSE
, etREAL
, {&disre_up2
},
104 "Distance between upper bound for distance restraints, and the distance at which the force becomes constant (see manual)" },
105 { "-cutoff", FALSE
, etREAL
, {&cutoff
},
106 "Only generate distance restraints for atoms pairs within cutoff (nm)" },
107 { "-constr", FALSE
, etBOOL
, {&bConstr
},
108 "Generate a constraint matrix rather than distance restraints. Constraints of type 2 will be generated that do generate exclusions." }
110 #define npargs asize(pa)
112 gmx_output_env_t
*oenv
;
113 t_atoms
*atoms
= nullptr;
119 const char *xfn
, *nfn
;
123 rvec dx
, *x
= nullptr, *v
= nullptr;
126 { efSTX
, "-f", nullptr, ffREAD
},
127 { efNDX
, "-n", nullptr, ffOPTRD
},
128 { efITP
, "-o", "posre", ffWRITE
},
129 { efNDX
, "-of", "freeze", ffOPTWR
}
131 #define NFILE asize(fnm)
133 if (!parse_common_args(&argc
, argv
, 0, NFILE
, fnm
, npargs
, pa
,
134 asize(desc
), desc
, 0, nullptr, &oenv
))
139 bFreeze
= opt2bSet("-of", NFILE
, fnm
) || opt2parg_bSet("-freeze", asize(pa
), pa
);
140 bDisre
= bDisre
|| opt2parg_bSet("-disre_dist", npargs
, pa
);
141 xfn
= opt2fn_null("-f", NFILE
, fnm
);
142 nfn
= opt2fn_null("-n", NFILE
, fnm
);
144 if (( nfn
== nullptr ) && ( xfn
== nullptr))
146 gmx_fatal(FARGS
, "no index file and no structure file supplied");
149 if ((disre_frac
< 0) || (disre_frac
>= 1))
151 gmx_fatal(FARGS
, "disre_frac should be between 0 and 1");
155 gmx_fatal(FARGS
, "disre_dist should be >= 0");
158 const char *title
= "";
161 fprintf(stderr
, "\nReading structure file\n");
162 t_topology
*top
= nullptr;
164 read_tps_conf(xfn
, top
, nullptr, &x
, &v
, box
, FALSE
);
167 if (atoms
->pdbinfo
== nullptr)
169 snew(atoms
->pdbinfo
, atoms
->nr
);
175 if (!atoms
|| !atoms
->pdbinfo
)
177 gmx_fatal(FARGS
, "No B-factors in input file %s, use a pdb file next time.",
181 out
= opt2FILE("-of", NFILE
, fnm
, "w");
182 fprintf(out
, "[ freeze ]\n");
183 for (i
= 0; (i
< atoms
->nr
); i
++)
185 if (atoms
->pdbinfo
[i
].bfac
<= freeze_level
)
187 fprintf(out
, "%d\n", i
+1);
192 else if ((bDisre
|| bConstr
) && x
)
194 printf("Select group to generate %s matrix from\n",
195 bConstr
? "constraint" : "distance restraint");
196 get_index(atoms
, nfn
, 1, &igrp
, &ind_grp
, &gn_grp
);
198 out
= ftp2FILE(efITP
, NFILE
, fnm
, "w");
201 fprintf(out
, "; constraints for %s of %s\n\n", gn_grp
, title
);
202 fprintf(out
, "[ constraints ]\n");
203 fprintf(out
, ";%4s %5s %1s %10s\n", "i", "j", "tp", "dist");
207 fprintf(out
, "; distance restraints for %s of %s\n\n", gn_grp
, title
);
208 fprintf(out
, "[ distance_restraints ]\n");
209 fprintf(out
, ";%4s %5s %1s %5s %10s %10s %10s %10s %10s\n", "i", "j", "?",
210 "label", "funct", "lo", "up1", "up2", "weight");
212 for (i
= k
= 0; i
< igrp
; i
++)
214 for (j
= i
+1; j
< igrp
; j
++, k
++)
216 rvec_sub(x
[ind_grp
[i
]], x
[ind_grp
[j
]], dx
);
220 fprintf(out
, "%5d %5d %1d %10g\n", ind_grp
[i
]+1, ind_grp
[j
]+1, 2, d
);
224 if (cutoff
< 0 || d
< cutoff
)
228 dd
= std::min(disre_dist
, disre_frac
*d
);
234 lo
= std::max(0.0_real
, d
-dd
);
236 fprintf(out
, "%5d %5d %1d %5d %10d %10g %10g %10g %10g\n",
237 ind_grp
[i
]+1, ind_grp
[j
]+1, 1, k
, 1,
238 lo
, hi
, hi
+disre_up2
, 1.0);
247 printf("Select group to position restrain\n");
248 get_index(atoms
, nfn
, 1, &igrp
, &ind_grp
, &gn_grp
);
250 out
= ftp2FILE(efITP
, NFILE
, fnm
, "w");
251 fprintf(out
, "; position restraints for %s of %s\n\n", gn_grp
, title
);
252 fprintf(out
, "[ position_restraints ]\n");
253 fprintf(out
, ";%3s %5s %9s %10s %10s\n", "i", "funct", "fcx", "fcy", "fcz");
254 for (i
= 0; i
< igrp
; i
++)
256 fprintf(out
, "%4d %4d %10g %10g %10g\n",
257 ind_grp
[i
]+1, 1, fc
[XX
], fc
[YY
], fc
[ZZ
]);