Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_genpr.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include "sysstuff.h"
43 #include "gromacs/commandline/pargs.h"
44 #include <string.h>
45 #include "smalloc.h"
46 #include "typedefs.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/fileio/futil.h"
49 #include "macros.h"
50 #include "vec.h"
51 #include "index.h"
52 #include "gmx_fatal.h"
53 #include "gmx_ana.h"
54
55 int gmx_genpr(int argc, char *argv[])
56 {
57     const char        *desc[] = {
58         "[THISMODULE] produces an include file for a topology containing",
59         "a list of atom numbers and three force constants for the",
60         "[IT]x[it]-, [IT]y[it]-, and [IT]z[it]-direction. A single isotropic force constant may",
61         "be given on the command line instead of three components.[PAR]",
62         "WARNING: position restraints only work for the one molecule at a time.",
63         "Position restraints are interactions within molecules, therefore",
64         "they should be included within the correct [TT][ moleculetype ][tt]",
65         "block in the topology. Since the atom numbers in every moleculetype",
66         "in the topology start at 1 and the numbers in the input file for",
67         "[THISMODULE] number consecutively from 1, [THISMODULE] will only",
68         "produce a useful file for the first molecule.[PAR]",
69         "The [TT]-of[tt] option produces an index file that can be used for",
70         "freezing atoms. In this case, the input file must be a [TT].pdb[tt] file.[PAR]",
71         "With the [TT]-disre[tt] option, half a matrix of distance restraints",
72         "is generated instead of position restraints. With this matrix, that",
73         "one typically would apply to C[GRK]alpha[grk] atoms in a protein, one can",
74         "maintain the overall conformation of a protein without tieing it to",
75         "a specific position (as with position restraints)."
76     };
77     static rvec        fc           = {1000.0, 1000.0, 1000.0};
78     static real        freeze_level = 0.0;
79     static real        disre_dist   = 0.1;
80     static real        disre_frac   = 0.0;
81     static real        disre_up2    = 1.0;
82     static gmx_bool    bDisre       = FALSE;
83     static gmx_bool    bConstr      = FALSE;
84     static real        cutoff       = -1.0;
85
86     t_pargs            pa[] = {
87         { "-fc", FALSE, etRVEC, {fc},
88           "Force constants (kJ/mol nm^2)" },
89         { "-freeze", FALSE, etREAL, {&freeze_level},
90           "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" },
91         { "-disre", FALSE, etBOOL, {&bDisre},
92           "Generate a distance restraint matrix for all the atoms in index" },
93         { "-disre_dist", FALSE, etREAL, {&disre_dist},
94           "Distance range around the actual distance for generating distance restraints" },
95         { "-disre_frac", FALSE, etREAL, {&disre_frac},
96           "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." },
97         { "-disre_up2", FALSE, etREAL, {&disre_up2},
98           "Distance between upper bound for distance restraints, and the distance at which the force becomes constant (see manual)" },
99         { "-cutoff", FALSE, etREAL, {&cutoff},
100           "Only generate distance restraints for atoms pairs within cutoff (nm)" },
101         { "-constr", FALSE, etBOOL, {&bConstr},
102           "Generate a constraint matrix rather than distance restraints. Constraints of type 2 will be generated that do generate exclusions." }
103     };
104 #define npargs asize(pa)
105
106     output_env_t     oenv;
107     t_atoms         *atoms = NULL;
108     int              i, j, k;
109     FILE            *out;
110     int              igrp;
111     real             d, dd, lo, hi;
112     atom_id         *ind_grp;
113     const char      *xfn, *nfn;
114     char            *gn_grp;
115     char             title[STRLEN];
116     matrix           box;
117     gmx_bool         bFreeze;
118     rvec             dx, *x = NULL, *v = NULL;
119
120     t_filenm         fnm[] = {
121         { efSTX, "-f",  NULL,    ffREAD },
122         { efNDX, "-n",  NULL,    ffOPTRD },
123         { efITP, "-o",  "posre", ffWRITE },
124         { efNDX, "-of", "freeze",    ffOPTWR }
125     };
126 #define NFILE asize(fnm)
127
128     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, npargs, pa,
129                            asize(desc), desc, 0, NULL, &oenv))
130     {
131         return 0;
132     }
133
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);
138
139     if (( nfn == NULL ) && ( xfn == NULL))
140     {
141         gmx_fatal(FARGS, "no index file and no structure file supplied");
142     }
143
144     if ((disre_frac < 0) || (disre_frac >= 1))
145     {
146         gmx_fatal(FARGS, "disre_frac should be between 0 and 1");
147     }
148     if (disre_dist < 0)
149     {
150         gmx_fatal(FARGS, "disre_dist should be >= 0");
151     }
152
153     if (xfn != NULL)
154     {
155         snew(atoms, 1);
156         get_stx_coordnum(xfn, &(atoms->nr));
157         init_t_atoms(atoms, atoms->nr, TRUE);
158         snew(x, atoms->nr);
159         snew(v, atoms->nr);
160         fprintf(stderr, "\nReading structure file\n");
161         read_stx_conf(xfn, title, atoms, x, v, NULL, box);
162     }
163
164     if (bFreeze)
165     {
166         if (!atoms || !atoms->pdbinfo)
167         {
168             gmx_fatal(FARGS, "No B-factors in input file %s, use a pdb file next time.",
169                       xfn);
170         }
171
172         out = opt2FILE("-of", NFILE, fnm, "w");
173         fprintf(out, "[ freeze ]\n");
174         for (i = 0; (i < atoms->nr); i++)
175         {
176             if (atoms->pdbinfo[i].bfac <= freeze_level)
177             {
178                 fprintf(out, "%d\n", i+1);
179             }
180         }
181         ffclose(out);
182     }
183     else if ((bDisre || bConstr) && x)
184     {
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);
188
189         out = ftp2FILE(efITP, NFILE, fnm, "w");
190         if (bConstr)
191         {
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");
195         }
196         else
197         {
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");
202         }
203         for (i = k = 0; i < igrp; i++)
204         {
205             for (j = i+1; j < igrp; j++, k++)
206             {
207                 rvec_sub(x[ind_grp[i]], x[ind_grp[j]], dx);
208                 d = norm(dx);
209                 if (bConstr)
210                 {
211                     fprintf(out, "%5d %5d %1d %10g\n", ind_grp[i]+1, ind_grp[j]+1, 2, d);
212                 }
213                 else
214                 {
215                     if (cutoff < 0 || d < cutoff)
216                     {
217                         if (disre_frac > 0)
218                         {
219                             dd = min(disre_dist, disre_frac*d);
220                         }
221                         else
222                         {
223                             dd = disre_dist;
224                         }
225                         lo = max(0, d-dd);
226                         hi = d+dd;
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,
229                                 lo, hi, hi+1, 1.0);
230                     }
231                 }
232             }
233         }
234         ffclose(out);
235     }
236     else
237     {
238         printf("Select group to position restrain\n");
239         get_index(atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
240
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++)
246         {
247             fprintf(out, "%4d %4d %10g %10g %10g\n",
248                     ind_grp[i]+1, 1, fc[XX], fc[YY], fc[ZZ]);
249         }
250         ffclose(out);
251     }
252     if (xfn)
253     {
254         sfree(x);
255         sfree(v);
256     }
257
258     return 0;
259 }