Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_genpr.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <math.h>
40 #include "sysstuff.h"
41 #include "gromacs/commandline/pargs.h"
42 #include <string.h>
43 #include "smalloc.h"
44 #include "typedefs.h"
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/futil.h"
47 #include "macros.h"
48 #include "vec.h"
49 #include "index.h"
50 #include "gmx_fatal.h"
51 #include "gmx_ana.h"
52
53 int gmx_genpr(int argc, char *argv[])
54 {
55     const char        *desc[] = {
56         "[THISMODULE] produces an include file for a topology containing",
57         "a list of atom numbers and three force constants for the",
58         "[IT]x[it]-, [IT]y[it]-, and [IT]z[it]-direction. A single isotropic force constant may",
59         "be given on the command line instead of three components.[PAR]",
60         "WARNING: position restraints only work for the one molecule at a time.",
61         "Position restraints are interactions within molecules, therefore",
62         "they should be included within the correct [TT][ moleculetype ][tt]",
63         "block in the topology. Since the atom numbers in every moleculetype",
64         "in the topology start at 1 and the numbers in the input file for",
65         "[THISMODULE] number consecutively from 1, [THISMODULE] will only",
66         "produce a useful file for the first molecule.[PAR]",
67         "The [TT]-of[tt] option produces an index file that can be used for",
68         "freezing atoms. In this case, the input file must be a [TT].pdb[tt] file.[PAR]",
69         "With the [TT]-disre[tt] option, half a matrix of distance restraints",
70         "is generated instead of position restraints. With this matrix, that",
71         "one typically would apply to C[GRK]alpha[grk] atoms in a protein, one can",
72         "maintain the overall conformation of a protein without tieing it to",
73         "a specific position (as with position restraints)."
74     };
75     static rvec        fc           = {1000.0, 1000.0, 1000.0};
76     static real        freeze_level = 0.0;
77     static real        disre_dist   = 0.1;
78     static real        disre_frac   = 0.0;
79     static real        disre_up2    = 1.0;
80     static gmx_bool    bDisre       = FALSE;
81     static gmx_bool    bConstr      = FALSE;
82     static real        cutoff       = -1.0;
83
84     t_pargs            pa[] = {
85         { "-fc", FALSE, etRVEC, {fc},
86           "Force constants (kJ/mol nm^2)" },
87         { "-freeze", FALSE, etREAL, {&freeze_level},
88           "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" },
89         { "-disre", FALSE, etBOOL, {&bDisre},
90           "Generate a distance restraint matrix for all the atoms in index" },
91         { "-disre_dist", FALSE, etREAL, {&disre_dist},
92           "Distance range around the actual distance for generating distance restraints" },
93         { "-disre_frac", FALSE, etREAL, {&disre_frac},
94           "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." },
95         { "-disre_up2", FALSE, etREAL, {&disre_up2},
96           "Distance between upper bound for distance restraints, and the distance at which the force becomes constant (see manual)" },
97         { "-cutoff", FALSE, etREAL, {&cutoff},
98           "Only generate distance restraints for atoms pairs within cutoff (nm)" },
99         { "-constr", FALSE, etBOOL, {&bConstr},
100           "Generate a constraint matrix rather than distance restraints. Constraints of type 2 will be generated that do generate exclusions." }
101     };
102 #define npargs asize(pa)
103
104     output_env_t     oenv;
105     t_atoms         *atoms = NULL;
106     int              i, j, k;
107     FILE            *out;
108     int              igrp;
109     real             d, dd, lo, hi;
110     atom_id         *ind_grp;
111     const char      *xfn, *nfn;
112     char            *gn_grp;
113     char             title[STRLEN];
114     matrix           box;
115     gmx_bool         bFreeze;
116     rvec             dx, *x = NULL, *v = NULL;
117
118     t_filenm         fnm[] = {
119         { efSTX, "-f",  NULL,    ffREAD },
120         { efNDX, "-n",  NULL,    ffOPTRD },
121         { efITP, "-o",  "posre", ffWRITE },
122         { efNDX, "-of", "freeze",    ffOPTWR }
123     };
124 #define NFILE asize(fnm)
125
126     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, npargs, pa,
127                            asize(desc), desc, 0, NULL, &oenv))
128     {
129         return 0;
130     }
131
132     bFreeze = opt2bSet("-of", NFILE, fnm) || opt2parg_bSet("-freeze", asize(pa), pa);
133     bDisre  = bDisre || opt2parg_bSet("-disre_dist", npargs, pa);
134     xfn     = opt2fn_null("-f", NFILE, fnm);
135     nfn     = opt2fn_null("-n", NFILE, fnm);
136
137     if (( nfn == NULL ) && ( xfn == NULL))
138     {
139         gmx_fatal(FARGS, "no index file and no structure file supplied");
140     }
141
142     if ((disre_frac < 0) || (disre_frac >= 1))
143     {
144         gmx_fatal(FARGS, "disre_frac should be between 0 and 1");
145     }
146     if (disre_dist < 0)
147     {
148         gmx_fatal(FARGS, "disre_dist should be >= 0");
149     }
150
151     if (xfn != NULL)
152     {
153         snew(atoms, 1);
154         get_stx_coordnum(xfn, &(atoms->nr));
155         init_t_atoms(atoms, atoms->nr, TRUE);
156         snew(x, atoms->nr);
157         snew(v, atoms->nr);
158         fprintf(stderr, "\nReading structure file\n");
159         read_stx_conf(xfn, title, atoms, x, v, NULL, box);
160     }
161
162     if (bFreeze)
163     {
164         if (!atoms || !atoms->pdbinfo)
165         {
166             gmx_fatal(FARGS, "No B-factors in input file %s, use a pdb file next time.",
167                       xfn);
168         }
169
170         out = opt2FILE("-of", NFILE, fnm, "w");
171         fprintf(out, "[ freeze ]\n");
172         for (i = 0; (i < atoms->nr); i++)
173         {
174             if (atoms->pdbinfo[i].bfac <= freeze_level)
175             {
176                 fprintf(out, "%d\n", i+1);
177             }
178         }
179         ffclose(out);
180     }
181     else if ((bDisre || bConstr) && x)
182     {
183         printf("Select group to generate %s matrix from\n",
184                bConstr ? "constraint" : "distance restraint");
185         get_index(atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
186
187         out = ftp2FILE(efITP, NFILE, fnm, "w");
188         if (bConstr)
189         {
190             fprintf(out, "; constraints for %s of %s\n\n", gn_grp, title);
191             fprintf(out, "[ constraints ]\n");
192             fprintf(out, ";%4s %5s %1s %10s\n", "i", "j", "tp", "dist");
193         }
194         else
195         {
196             fprintf(out, "; distance restraints for %s of %s\n\n", gn_grp, title);
197             fprintf(out, "[ distance_restraints ]\n");
198             fprintf(out, ";%4s %5s %1s %5s %10s %10s %10s %10s %10s\n", "i", "j", "?",
199                     "label", "funct", "lo", "up1", "up2", "weight");
200         }
201         for (i = k = 0; i < igrp; i++)
202         {
203             for (j = i+1; j < igrp; j++, k++)
204             {
205                 rvec_sub(x[ind_grp[i]], x[ind_grp[j]], dx);
206                 d = norm(dx);
207                 if (bConstr)
208                 {
209                     fprintf(out, "%5d %5d %1d %10g\n", ind_grp[i]+1, ind_grp[j]+1, 2, d);
210                 }
211                 else
212                 {
213                     if (cutoff < 0 || d < cutoff)
214                     {
215                         if (disre_frac > 0)
216                         {
217                             dd = min(disre_dist, disre_frac*d);
218                         }
219                         else
220                         {
221                             dd = disre_dist;
222                         }
223                         lo = max(0, d-dd);
224                         hi = d+dd;
225                         fprintf(out, "%5d %5d %1d %5d %10d %10g %10g %10g %10g\n",
226                                 ind_grp[i]+1, ind_grp[j]+1, 1, k, 1,
227                                 lo, hi, hi+1, 1.0);
228                     }
229                 }
230             }
231         }
232         ffclose(out);
233     }
234     else
235     {
236         printf("Select group to position restrain\n");
237         get_index(atoms, nfn, 1, &igrp, &ind_grp, &gn_grp);
238
239         out = ftp2FILE(efITP, NFILE, fnm, "w");
240         fprintf(out, "; position restraints for %s of %s\n\n", gn_grp, title);
241         fprintf(out, "[ position_restraints ]\n");
242         fprintf(out, ";%3s %5s %9s %10s %10s\n", "i", "funct", "fcx", "fcy", "fcz");
243         for (i = 0; i < igrp; i++)
244         {
245             fprintf(out, "%4d %4d %10g %10g %10g\n",
246                     ind_grp[i]+1, 1, fc[XX], fc[YY], fc[ZZ]);
247         }
248         ffclose(out);
249     }
250     if (xfn)
251     {
252         sfree(x);
253         sfree(v);
254     }
255
256     return 0;
257 }