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