Remove unused function generate_excls and make clean_excls static
[gromacs.git] / src / gromacs / gmxpreprocess / genrestr.cpp
blob3d462299cfdc7e4458868d3eb896e913abc1ea41
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/mtop_util.h"
51 #include "gromacs/topology/topology.h"
52 #include "gromacs/utility/arraysize.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/utility/smalloc.h"
58 int gmx_genrestr(int argc, char *argv[])
60 const char *desc[] = {
61 "[THISMODULE] produces an #include file for a topology containing",
62 "a list of atom numbers and three force constants for the",
63 "[IT]x[it]-, [IT]y[it]-, and [IT]z[it]-direction based on",
64 "the contents of the [TT]-f[tt] file. A single isotropic force constant may",
65 "be given on the command line instead of three components.[PAR]",
66 "WARNING: Position restraints are interactions within molecules, therefore",
67 "they must be included within the correct [TT][ moleculetype ][tt]",
68 "block in the topology. The atom indices within the",
69 "[TT][ position_restraints ][tt] block must be within the range of the",
70 "atom indices for that molecule type. Since the atom numbers in every",
71 "moleculetype in the topology start at 1 and the numbers in the input file",
72 "for [THISMODULE] number consecutively from 1, [THISMODULE] will only",
73 "produce a useful file for the first molecule. You may wish to",
74 "edit the resulting index file to remove the lines for later atoms,",
75 "or construct a suitable index group to provide",
76 "as input to [THISMODULE].[PAR]",
77 "The [TT]-of[tt] option produces an index file that can be used for",
78 "freezing atoms. In this case, the input file must be a [REF].pdb[ref] file.[PAR]",
79 "With the [TT]-disre[tt] option, half a matrix of distance restraints",
80 "is generated instead of position restraints. With this matrix, that",
81 "one typically would apply to C[GRK]alpha[grk] atoms in a protein, one can",
82 "maintain the overall conformation of a protein without tieing it to",
83 "a specific position (as with position restraints)."
85 static rvec fc = {1000.0, 1000.0, 1000.0};
86 static real freeze_level = 0.0;
87 static real disre_dist = 0.1;
88 static real disre_frac = 0.0;
89 static real disre_up2 = 1.0;
90 static gmx_bool bDisre = FALSE;
91 static gmx_bool bConstr = FALSE;
92 static real cutoff = -1.0;
94 t_pargs pa[] = {
95 { "-fc", FALSE, etRVEC, {fc},
96 "Force constants (kJ/mol nm^2)" },
97 { "-freeze", FALSE, etREAL, {&freeze_level},
98 "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" },
99 { "-disre", FALSE, etBOOL, {&bDisre},
100 "Generate a distance restraint matrix for all the atoms in index" },
101 { "-disre_dist", FALSE, etREAL, {&disre_dist},
102 "Distance range around the actual distance for generating distance restraints" },
103 { "-disre_frac", FALSE, etREAL, {&disre_frac},
104 "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." },
105 { "-disre_up2", FALSE, etREAL, {&disre_up2},
106 "Distance between upper bound for distance restraints, and the distance at which the force becomes constant (see manual)" },
107 { "-cutoff", FALSE, etREAL, {&cutoff},
108 "Only generate distance restraints for atoms pairs within cutoff (nm)" },
109 { "-constr", FALSE, etBOOL, {&bConstr},
110 "Generate a constraint matrix rather than distance restraints. Constraints of type 2 will be generated that do generate exclusions." }
112 #define npargs asize(pa)
114 gmx_output_env_t *oenv;
115 t_atoms atoms;
116 int i, j, k;
117 FILE *out;
118 int igrp;
119 real d, dd, lo, hi;
120 int *ind_grp;
121 const char *xfn, *nfn;
122 char *gn_grp;
123 matrix box;
124 gmx_bool bFreeze;
125 rvec dx, *x = nullptr, *v = nullptr;
127 t_filenm fnm[] = {
128 { efSTX, "-f", nullptr, ffREAD },
129 { efNDX, "-n", nullptr, ffOPTRD },
130 { efITP, "-o", "posre", ffWRITE },
131 { efNDX, "-of", "freeze", ffOPTWR }
133 #define NFILE asize(fnm)
135 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, npargs, pa,
136 asize(desc), desc, 0, nullptr, &oenv))
138 return 0;
141 bFreeze = opt2bSet("-of", NFILE, fnm) || opt2parg_bSet("-freeze", asize(pa), pa);
142 bDisre = bDisre || opt2parg_bSet("-disre_dist", npargs, pa);
143 xfn = opt2fn_null("-f", NFILE, fnm);
144 nfn = opt2fn_null("-n", NFILE, fnm);
146 if (( nfn == nullptr ) && ( xfn == nullptr))
148 gmx_fatal(FARGS, "no index file and no structure file supplied");
151 if ((disre_frac < 0) || (disre_frac >= 1))
153 gmx_fatal(FARGS, "disre_frac should be between 0 and 1");
155 if (disre_dist < 0)
157 gmx_fatal(FARGS, "disre_dist should be >= 0");
160 const char *title = "";
161 bool haveTopology = false;
162 if (xfn != nullptr)
164 fprintf(stderr, "\nReading structure file\n");
165 gmx_mtop_t mtop;
166 readConfAndTopology(xfn, &haveTopology, &mtop, nullptr, &x, &v, box);
167 title = *mtop.name;
168 atoms = gmx_mtop_global_atoms(&mtop);
169 if (atoms.pdbinfo == nullptr)
171 snew(atoms.pdbinfo, atoms.nr);
173 haveTopology = true;
176 if (bFreeze)
178 if (!haveTopology || !atoms.pdbinfo)
180 gmx_fatal(FARGS, "No B-factors in input file %s, use a pdb file next time.",
181 xfn);
184 out = opt2FILE("-of", NFILE, fnm, "w");
185 fprintf(out, "[ freeze ]\n");
186 for (i = 0; (i < atoms.nr); i++)
188 if (atoms.pdbinfo[i].bfac <= freeze_level)
190 fprintf(out, "%d\n", i+1);
193 gmx_ffclose(out);
195 else if ((bDisre || bConstr) && x)
197 printf("Select group to generate %s matrix from\n",
198 bConstr ? "constraint" : "distance restraint");
199 get_index(&atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
201 out = ftp2FILE(efITP, NFILE, fnm, "w");
202 if (bConstr)
204 fprintf(out, "; constraints for %s of %s\n\n", gn_grp, title);
205 fprintf(out, "[ constraints ]\n");
206 fprintf(out, ";%4s %5s %1s %10s\n", "i", "j", "tp", "dist");
208 else
210 fprintf(out, "; distance restraints for %s of %s\n\n", gn_grp, title);
211 fprintf(out, "[ distance_restraints ]\n");
212 fprintf(out, ";%4s %5s %1s %5s %10s %10s %10s %10s %10s\n", "i", "j", "?",
213 "label", "funct", "lo", "up1", "up2", "weight");
215 for (i = k = 0; i < igrp; i++)
217 for (j = i+1; j < igrp; j++, k++)
219 rvec_sub(x[ind_grp[i]], x[ind_grp[j]], dx);
220 d = norm(dx);
221 if (bConstr)
223 fprintf(out, "%5d %5d %1d %10g\n", ind_grp[i]+1, ind_grp[j]+1, 2, d);
225 else
227 if (cutoff < 0 || d < cutoff)
229 if (disre_frac > 0)
231 dd = std::min(disre_dist, disre_frac*d);
233 else
235 dd = disre_dist;
237 lo = std::max(0.0_real, d-dd);
238 hi = d+dd;
239 fprintf(out, "%5d %5d %1d %5d %10d %10g %10g %10g %10g\n",
240 ind_grp[i]+1, ind_grp[j]+1, 1, k, 1,
241 lo, hi, hi+disre_up2, 1.0);
246 gmx_ffclose(out);
248 else
250 printf("Select group to position restrain\n");
251 get_index(&atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
253 out = ftp2FILE(efITP, NFILE, fnm, "w");
254 fprintf(out, "; position restraints for %s of %s\n\n", gn_grp, title);
255 fprintf(out, "[ position_restraints ]\n");
256 fprintf(out, ";%3s %5s %9s %10s %10s\n", "i", "funct", "fcx", "fcy", "fcz");
257 for (i = 0; i < igrp; i++)
259 fprintf(out, "%4d %4d %10g %10g %10g\n",
260 ind_grp[i]+1, 1, fc[XX], fc[YY], fc[ZZ]);
262 gmx_ffclose(out);
264 if (xfn)
266 sfree(x);
267 sfree(v);
268 done_atom(&atoms);
271 return 0;