Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_morph.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 "gromacs/commandline/pargs.h"
40 #include "gromacs/fileio/confio.h"
41 #include "gromacs/fileio/trxio.h"
42 #include "gromacs/fileio/xvgr.h"
43 #include "gromacs/gmxana/gmx_ana.h"
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/legacyheaders/viewit.h"
46 #include "gromacs/math/do_fit.h"
47 #include "gromacs/topology/atoms.h"
48 #include "gromacs/topology/index.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/smalloc.h"
51
52 static real dointerp(int n, rvec x1[], rvec x2[], rvec xx[],
53                      int I, int N, real first, real last)
54 {
55     int    i, j;
56     double fac, fac0, fac1;
57
58     fac  = first + (I*(last-first))/(N-1);
59     fac0 = 1-fac;
60     fac1 = fac;
61     for (i = 0; (i < n); i++)
62     {
63         for (j = 0; (j < DIM); j++)
64         {
65             xx[i][j] = fac0*x1[i][j] + fac1*x2[i][j];
66         }
67     }
68
69     return fac;
70 }
71
72 int gmx_morph(int argc, char *argv[])
73 {
74     const char      *desc[] = {
75         "[THISMODULE] does a linear interpolation of conformations in order to",
76         "create intermediates. Of course these are completely unphysical, but",
77         "that you may try to justify yourself. Output is in the form of a ",
78         "generic trajectory. The number of intermediates can be controlled with",
79         "the [TT]-ninterm[tt] flag. The first and last flag correspond to the way of",
80         "interpolating: 0 corresponds to input structure 1 while",
81         "1 corresponds to input structure 2.",
82         "If you specify [TT]-first[tt] < 0 or [TT]-last[tt] > 1 extrapolation will be",
83         "on the path from input structure x[SUB]1[sub] to x[SUB]2[sub]. In general, the coordinates",
84         "of the intermediate x(i) out of N total intermediates correspond to:[PAR]",
85         "x(i) = x[SUB]1[sub] + (first+(i/(N-1))*(last-first))*(x[SUB]2[sub]-x[SUB]1[sub])[PAR]",
86         "Finally the RMSD with respect to both input structures can be computed",
87         "if explicitly selected ([TT]-or[tt] option). In that case, an index file may be",
88         "read to select the group from which the RMS is computed."
89     };
90     t_filenm         fnm[] = {
91         { efSTX, "-f1", "conf1",  ffREAD },
92         { efSTX, "-f2", "conf2",  ffREAD },
93         { efTRX, "-o",  "interm", ffWRITE },
94         { efXVG, "-or", "rms-interm", ffOPTWR },
95         { efNDX, "-n",  "index",  ffOPTRD }
96     };
97 #define NFILE asize(fnm)
98     static  int      ninterm = 11;
99     static  real     first   = 0.0;
100     static  real     last    = 1.0;
101     static  gmx_bool bFit    = TRUE;
102     t_pargs          pa []   = {
103         { "-ninterm", FALSE, etINT,  {&ninterm},
104           "Number of intermediates" },
105         { "-first",   FALSE, etREAL, {&first},
106           "Corresponds to first generated structure (0 is input x[SUB]1[sub], see above)" },
107         { "-last",    FALSE, etREAL, {&last},
108           "Corresponds to last generated structure (1 is input x[SUB]2[sub], see above)" },
109         { "-fit",     FALSE, etBOOL, {&bFit},
110           "Do a least squares fit of the second to the first structure before interpolating" }
111     };
112     const char      *leg[] = { "Ref = 1\\Sst\\N conf", "Ref = 2\\Snd\\N conf" };
113     FILE            *fp    = NULL;
114     int              i, isize, is_lsq, nat1, nat2;
115     t_trxstatus     *status;
116     atom_id         *index, *index_lsq, *index_all, *dummy;
117     t_atoms          atoms;
118     rvec            *x1, *x2, *xx, *v;
119     matrix           box;
120     real             rms1, rms2, fac, *mass;
121     char             title[STRLEN], *grpname;
122     gmx_bool         bRMS;
123     output_env_t     oenv;
124
125     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW,
126                            NFILE, fnm, asize(pa), pa, asize(desc), desc,
127                            0, NULL, &oenv))
128     {
129         return 0;
130     }
131     get_stx_coordnum (opt2fn("-f1", NFILE, fnm), &nat1);
132     get_stx_coordnum (opt2fn("-f2", NFILE, fnm), &nat2);
133     if (nat1 != nat2)
134     {
135         gmx_fatal(FARGS, "Number of atoms in first structure is %d, in second %d",
136                   nat1, nat2);
137     }
138
139     init_t_atoms(&atoms, nat1, TRUE);
140     snew(x1, nat1);
141     snew(x2, nat1);
142     snew(xx, nat1);
143     snew(v, nat1);
144
145     read_stx_conf(opt2fn("-f1", NFILE, fnm), title, &atoms, x1, v, NULL, box);
146     read_stx_conf(opt2fn("-f2", NFILE, fnm), title, &atoms, x2, v, NULL, box);
147
148     snew(mass, nat1);
149     snew(index_all, nat1);
150     for (i = 0; (i < nat1); i++)
151     {
152         mass[i]      = 1;
153         index_all[i] = i;
154     }
155     if (bFit)
156     {
157         printf("Select group for LSQ superposition:\n");
158         get_index(&atoms, opt2fn_null("-n", NFILE, fnm), 1, &is_lsq, &index_lsq,
159                   &grpname);
160         reset_x(is_lsq, index_lsq, nat1, index_all, x1, mass);
161         reset_x(is_lsq, index_lsq, nat1, index_all, x2, mass);
162         do_fit(nat1, mass, x1, x2);
163     }
164
165     bRMS = opt2bSet("-or", NFILE, fnm);
166     if (bRMS)
167     {
168         fp = xvgropen(opt2fn("-or", NFILE, fnm), "RMSD", "Conf", "(nm)", oenv);
169         xvgr_legend(fp, asize(leg), leg, oenv);
170         printf("Select group for RMSD calculation:\n");
171         get_index(&atoms, opt2fn_null("-n", NFILE, fnm), 1, &isize, &index, &grpname);
172         printf("You selected group %s, containing %d atoms\n", grpname, isize);
173         rms1 = rmsdev_ind(isize, index, mass, x1, x2);
174         fprintf(stderr, "RMSD between input conformations is %g nm\n", rms1);
175     }
176
177     snew(dummy, nat1);
178     for (i = 0; (i < nat1); i++)
179     {
180         dummy[i] = i;
181     }
182     status = open_trx(ftp2fn(efTRX, NFILE, fnm), "w");
183
184     for (i = 0; (i < ninterm); i++)
185     {
186         fac = dointerp(nat1, x1, x2, xx, i, ninterm, first, last);
187         write_trx(status, nat1, dummy, &atoms, i, fac, box, xx, NULL, NULL);
188         if (bRMS)
189         {
190             rms1 = rmsdev_ind(isize, index, mass, x1, xx);
191             rms2 = rmsdev_ind(isize, index, mass, x2, xx);
192             fprintf(fp, "%10g  %10g  %10g\n", fac, rms1, rms2);
193         }
194     }
195
196     close_trx(status);
197
198     if (bRMS)
199     {
200         gmx_ffclose(fp);
201         do_view(oenv, opt2fn("-or", NFILE, fnm), "-nxy");
202     }
203
204     return 0;
205 }