Merge branch release-4-6 into release-5-0
[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,2014, 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 "gromacs/utility/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 based on",
61         "the contents of the [TT]-f[tt] file. A single isotropic force constant may",
62         "be given on the command line instead of three components.[PAR]",
63         "WARNING: Position restraints are interactions within molecules, therefore",
64         "they must be included within the correct [TT][ moleculetype ][tt]",
65         "block in the topology. The atom indices within the",
66         "[TT][ position_restraints ][tt] block must be within the range of the",
67         "atom indices for that molecule type. Since the atom numbers in every",
68         "moleculetype in the topology start at 1 and the numbers in the input file",
69         "for [THISMODULE] number consecutively from 1, [THISMODULE] will only",
70         "produce a useful file for the first molecule. You may wish to",
71         "edit the resulting index file to remove the lines for later atoms,",
72         "or construct a suitable index group to provide",
73         "as input to [THISMODULE].[PAR]",
74         "The [TT]-of[tt] option produces an index file that can be used for",
75         "freezing atoms. In this case, the input file must be a [TT].pdb[tt] file.[PAR]",
76         "With the [TT]-disre[tt] option, half a matrix of distance restraints",
77         "is generated instead of position restraints. With this matrix, that",
78         "one typically would apply to C[GRK]alpha[grk] atoms in a protein, one can",
79         "maintain the overall conformation of a protein without tieing it to",
80         "a specific position (as with position restraints)."
81     };
82     static rvec        fc           = {1000.0, 1000.0, 1000.0};
83     static real        freeze_level = 0.0;
84     static real        disre_dist   = 0.1;
85     static real        disre_frac   = 0.0;
86     static real        disre_up2    = 1.0;
87     static gmx_bool    bDisre       = FALSE;
88     static gmx_bool    bConstr      = FALSE;
89     static real        cutoff       = -1.0;
90
91     t_pargs            pa[] = {
92         { "-fc", FALSE, etRVEC, {fc},
93           "Force constants (kJ/mol nm^2)" },
94         { "-freeze", FALSE, etREAL, {&freeze_level},
95           "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" },
96         { "-disre", FALSE, etBOOL, {&bDisre},
97           "Generate a distance restraint matrix for all the atoms in index" },
98         { "-disre_dist", FALSE, etREAL, {&disre_dist},
99           "Distance range around the actual distance for generating distance restraints" },
100         { "-disre_frac", FALSE, etREAL, {&disre_frac},
101           "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." },
102         { "-disre_up2", FALSE, etREAL, {&disre_up2},
103           "Distance between upper bound for distance restraints, and the distance at which the force becomes constant (see manual)" },
104         { "-cutoff", FALSE, etREAL, {&cutoff},
105           "Only generate distance restraints for atoms pairs within cutoff (nm)" },
106         { "-constr", FALSE, etBOOL, {&bConstr},
107           "Generate a constraint matrix rather than distance restraints. Constraints of type 2 will be generated that do generate exclusions." }
108     };
109 #define npargs asize(pa)
110
111     output_env_t     oenv;
112     t_atoms         *atoms = NULL;
113     int              i, j, k;
114     FILE            *out;
115     int              igrp;
116     real             d, dd, lo, hi;
117     atom_id         *ind_grp;
118     const char      *xfn, *nfn;
119     char            *gn_grp;
120     char             title[STRLEN];
121     matrix           box;
122     gmx_bool         bFreeze;
123     rvec             dx, *x = NULL, *v = NULL;
124
125     t_filenm         fnm[] = {
126         { efSTX, "-f",  NULL,    ffREAD },
127         { efNDX, "-n",  NULL,    ffOPTRD },
128         { efITP, "-o",  "posre", ffWRITE },
129         { efNDX, "-of", "freeze",    ffOPTWR }
130     };
131 #define NFILE asize(fnm)
132
133     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, npargs, pa,
134                            asize(desc), desc, 0, NULL, &oenv))
135     {
136         return 0;
137     }
138
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);
143
144     if (( nfn == NULL ) && ( xfn == NULL))
145     {
146         gmx_fatal(FARGS, "no index file and no structure file supplied");
147     }
148
149     if ((disre_frac < 0) || (disre_frac >= 1))
150     {
151         gmx_fatal(FARGS, "disre_frac should be between 0 and 1");
152     }
153     if (disre_dist < 0)
154     {
155         gmx_fatal(FARGS, "disre_dist should be >= 0");
156     }
157
158     if (xfn != NULL)
159     {
160         snew(atoms, 1);
161         get_stx_coordnum(xfn, &(atoms->nr));
162         init_t_atoms(atoms, atoms->nr, TRUE);
163         snew(x, atoms->nr);
164         snew(v, atoms->nr);
165         fprintf(stderr, "\nReading structure file\n");
166         read_stx_conf(xfn, title, atoms, x, v, NULL, box);
167     }
168
169     if (bFreeze)
170     {
171         if (!atoms || !atoms->pdbinfo)
172         {
173             gmx_fatal(FARGS, "No B-factors in input file %s, use a pdb file next time.",
174                       xfn);
175         }
176
177         out = opt2FILE("-of", NFILE, fnm, "w");
178         fprintf(out, "[ freeze ]\n");
179         for (i = 0; (i < atoms->nr); i++)
180         {
181             if (atoms->pdbinfo[i].bfac <= freeze_level)
182             {
183                 fprintf(out, "%d\n", i+1);
184             }
185         }
186         gmx_ffclose(out);
187     }
188     else if ((bDisre || bConstr) && x)
189     {
190         printf("Select group to generate %s matrix from\n",
191                bConstr ? "constraint" : "distance restraint");
192         get_index(atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
193
194         out = ftp2FILE(efITP, NFILE, fnm, "w");
195         if (bConstr)
196         {
197             fprintf(out, "; constraints for %s of %s\n\n", gn_grp, title);
198             fprintf(out, "[ constraints ]\n");
199             fprintf(out, ";%4s %5s %1s %10s\n", "i", "j", "tp", "dist");
200         }
201         else
202         {
203             fprintf(out, "; distance restraints for %s of %s\n\n", gn_grp, title);
204             fprintf(out, "[ distance_restraints ]\n");
205             fprintf(out, ";%4s %5s %1s %5s %10s %10s %10s %10s %10s\n", "i", "j", "?",
206                     "label", "funct", "lo", "up1", "up2", "weight");
207         }
208         for (i = k = 0; i < igrp; i++)
209         {
210             for (j = i+1; j < igrp; j++, k++)
211             {
212                 rvec_sub(x[ind_grp[i]], x[ind_grp[j]], dx);
213                 d = norm(dx);
214                 if (bConstr)
215                 {
216                     fprintf(out, "%5d %5d %1d %10g\n", ind_grp[i]+1, ind_grp[j]+1, 2, d);
217                 }
218                 else
219                 {
220                     if (cutoff < 0 || d < cutoff)
221                     {
222                         if (disre_frac > 0)
223                         {
224                             dd = min(disre_dist, disre_frac*d);
225                         }
226                         else
227                         {
228                             dd = disre_dist;
229                         }
230                         lo = max(0, d-dd);
231                         hi = d+dd;
232                         fprintf(out, "%5d %5d %1d %5d %10d %10g %10g %10g %10g\n",
233                                 ind_grp[i]+1, ind_grp[j]+1, 1, k, 1,
234                                 lo, hi, hi+disre_up2, 1.0);
235                     }
236                 }
237             }
238         }
239         gmx_ffclose(out);
240     }
241     else
242     {
243         printf("Select group to position restrain\n");
244         get_index(atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
245
246         out = ftp2FILE(efITP, NFILE, fnm, "w");
247         fprintf(out, "; position restraints for %s of %s\n\n", gn_grp, title);
248         fprintf(out, "[ position_restraints ]\n");
249         fprintf(out, ";%3s %5s %9s %10s %10s\n", "i", "funct", "fcx", "fcy", "fcz");
250         for (i = 0; i < igrp; i++)
251         {
252             fprintf(out, "%4d %4d %10g %10g %10g\n",
253                     ind_grp[i]+1, 1, fc[XX], fc[YY], fc[ZZ]);
254         }
255         gmx_ffclose(out);
256     }
257     if (xfn)
258     {
259         sfree(x);
260         sfree(v);
261     }
262
263     return 0;
264 }