Refactor preprocessing atom types.
[gromacs.git] / src / gromacs / gmxpreprocess / genrestr.cpp
bloba70cbb5c1c1b4b161996583a57422c318bf664ab
1 /*
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,2019, 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.
37 #include "gmxpre.h"
39 #include "genrestr.h"
41 #include <cmath>
42 #include <cstring>
44 #include <algorithm>
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/topology/index.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/utility/arraysize.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
57 int gmx_genrestr(int argc, char *argv[])
59 const char *desc[] = {
60 "[THISMODULE] 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 based on",
63 "the contents of the [TT]-f[tt] file. A single isotropic force constant may",
64 "be given on the command line instead of three components.[PAR]",
65 "WARNING: Position restraints are interactions within molecules, therefore",
66 "they must be included within the correct [TT][ moleculetype ][tt]",
67 "block in the topology. The atom indices within the",
68 "[TT][ position_restraints ][tt] block must be within the range of the",
69 "atom indices for that molecule type. Since the atom numbers in every",
70 "moleculetype in the topology start at 1 and the numbers in the input file",
71 "for [THISMODULE] number consecutively from 1, [THISMODULE] will only",
72 "produce a useful file for the first molecule. You may wish to",
73 "edit the resulting index file to remove the lines for later atoms,",
74 "or construct a suitable index group to provide",
75 "as input to [THISMODULE].[PAR]",
76 "The [TT]-of[tt] option produces an index file that can be used for",
77 "freezing atoms. In this case, the input file must be a [REF].pdb[ref] file.[PAR]",
78 "With the [TT]-disre[tt] option, half a matrix of distance restraints",
79 "is generated instead of position restraints. With this matrix, that",
80 "one typically would apply to C[GRK]alpha[grk] atoms in a protein, one can",
81 "maintain the overall conformation of a protein without tieing it to",
82 "a specific position (as with position restraints)."
84 static rvec fc = {1000.0, 1000.0, 1000.0};
85 static real freeze_level = 0.0;
86 static real disre_dist = 0.1;
87 static real disre_frac = 0.0;
88 static real disre_up2 = 1.0;
89 static gmx_bool bDisre = FALSE;
90 static gmx_bool bConstr = FALSE;
91 static real cutoff = -1.0;
93 t_pargs pa[] = {
94 { "-fc", FALSE, etRVEC, {fc},
95 "Force constants (kJ/mol nm^2)" },
96 { "-freeze", FALSE, etREAL, {&freeze_level},
97 "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" },
98 { "-disre", FALSE, etBOOL, {&bDisre},
99 "Generate a distance restraint matrix for all the atoms in index" },
100 { "-disre_dist", FALSE, etREAL, {&disre_dist},
101 "Distance range around the actual distance for generating distance restraints" },
102 { "-disre_frac", FALSE, etREAL, {&disre_frac},
103 "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." },
104 { "-disre_up2", FALSE, etREAL, {&disre_up2},
105 "Distance between upper bound for distance restraints, and the distance at which the force becomes constant (see manual)" },
106 { "-cutoff", FALSE, etREAL, {&cutoff},
107 "Only generate distance restraints for atoms pairs within cutoff (nm)" },
108 { "-constr", FALSE, etBOOL, {&bConstr},
109 "Generate a constraint matrix rather than distance restraints. Constraints of type 2 will be generated that do generate exclusions." }
111 #define npargs asize(pa)
113 gmx_output_env_t *oenv;
114 t_atoms *atoms = nullptr;
115 int i, j, k;
116 FILE *out;
117 int igrp;
118 real d, dd, lo, hi;
119 int *ind_grp;
120 const char *xfn, *nfn;
121 char *gn_grp;
122 matrix box;
123 gmx_bool bFreeze;
124 rvec dx, *x = nullptr, *v = nullptr;
126 t_filenm fnm[] = {
127 { efSTX, "-f", nullptr, ffREAD },
128 { efNDX, "-n", nullptr, ffOPTRD },
129 { efITP, "-o", "posre", ffWRITE },
130 { efNDX, "-of", "freeze", ffOPTWR }
132 #define NFILE asize(fnm)
134 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, npargs, pa,
135 asize(desc), desc, 0, nullptr, &oenv))
137 return 0;
140 bFreeze = opt2bSet("-of", NFILE, fnm) || opt2parg_bSet("-freeze", asize(pa), pa);
141 bDisre = bDisre || opt2parg_bSet("-disre_dist", npargs, pa);
142 xfn = opt2fn_null("-f", NFILE, fnm);
143 nfn = opt2fn_null("-n", NFILE, fnm);
145 if (( nfn == nullptr ) && ( xfn == nullptr))
147 gmx_fatal(FARGS, "no index file and no structure file supplied");
150 if ((disre_frac < 0) || (disre_frac >= 1))
152 gmx_fatal(FARGS, "disre_frac should be between 0 and 1");
154 if (disre_dist < 0)
156 gmx_fatal(FARGS, "disre_dist should be >= 0");
159 const char *title = "";
160 if (xfn != nullptr)
162 fprintf(stderr, "\nReading structure file\n");
163 t_topology *top = nullptr;
164 snew(top, 1);
165 read_tps_conf(xfn, top, nullptr, &x, &v, box, FALSE);
166 title = *top->name;
167 atoms = &top->atoms;
168 if (atoms->pdbinfo == nullptr)
170 snew(atoms->pdbinfo, atoms->nr);
174 if (bFreeze)
176 if (!atoms || !atoms->pdbinfo)
178 gmx_fatal(FARGS, "No B-factors in input file %s, use a pdb file next time.",
179 xfn);
182 out = opt2FILE("-of", NFILE, fnm, "w");
183 fprintf(out, "[ freeze ]\n");
184 for (i = 0; (i < atoms->nr); i++)
186 if (atoms->pdbinfo[i].bfac <= freeze_level)
188 fprintf(out, "%d\n", i+1);
191 gmx_ffclose(out);
193 else if ((bDisre || bConstr) && x)
195 printf("Select group to generate %s matrix from\n",
196 bConstr ? "constraint" : "distance restraint");
197 get_index(atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
199 out = ftp2FILE(efITP, NFILE, fnm, "w");
200 if (bConstr)
202 fprintf(out, "; constraints for %s of %s\n\n", gn_grp, title);
203 fprintf(out, "[ constraints ]\n");
204 fprintf(out, ";%4s %5s %1s %10s\n", "i", "j", "tp", "dist");
206 else
208 fprintf(out, "; distance restraints for %s of %s\n\n", gn_grp, title);
209 fprintf(out, "[ distance_restraints ]\n");
210 fprintf(out, ";%4s %5s %1s %5s %10s %10s %10s %10s %10s\n", "i", "j", "?",
211 "label", "funct", "lo", "up1", "up2", "weight");
213 for (i = k = 0; i < igrp; i++)
215 for (j = i+1; j < igrp; j++, k++)
217 rvec_sub(x[ind_grp[i]], x[ind_grp[j]], dx);
218 d = norm(dx);
219 if (bConstr)
221 fprintf(out, "%5d %5d %1d %10g\n", ind_grp[i]+1, ind_grp[j]+1, 2, d);
223 else
225 if (cutoff < 0 || d < cutoff)
227 if (disre_frac > 0)
229 dd = std::min(disre_dist, disre_frac*d);
231 else
233 dd = disre_dist;
235 lo = std::max(0.0_real, d-dd);
236 hi = d+dd;
237 fprintf(out, "%5d %5d %1d %5d %10d %10g %10g %10g %10g\n",
238 ind_grp[i]+1, ind_grp[j]+1, 1, k, 1,
239 lo, hi, hi+disre_up2, 1.0);
244 gmx_ffclose(out);
246 else
248 printf("Select group to position restrain\n");
249 get_index(atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
251 out = ftp2FILE(efITP, NFILE, fnm, "w");
252 fprintf(out, "; position restraints for %s of %s\n\n", gn_grp, title);
253 fprintf(out, "[ position_restraints ]\n");
254 fprintf(out, ";%3s %5s %9s %10s %10s\n", "i", "funct", "fcx", "fcy", "fcz");
255 for (i = 0; i < igrp; i++)
257 fprintf(out, "%4d %4d %10g %10g %10g\n",
258 ind_grp[i]+1, 1, fc[XX], fc[YY], fc[ZZ]);
260 gmx_ffclose(out);
262 if (xfn)
264 sfree(x);
265 sfree(v);
268 return 0;