2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
40 #include "gromacs/fileio/confio.h"
42 #include "gromacs/topology/index.h"
44 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/xvgr.h"
48 #include "gromacs/math/do_fit.h"
49 #include "gromacs/topology/atoms.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/smalloc.h"
53 static real dointerp(int n, rvec x1[], rvec x2[], rvec xx[],
54 int I, int N, real first, real last)
57 double fac, fac0, fac1;
59 fac = first + (I*(last-first))/(N-1);
62 for (i = 0; (i < n); i++)
64 for (j = 0; (j < DIM); j++)
66 xx[i][j] = fac0*x1[i][j] + fac1*x2[i][j];
73 int gmx_morph(int argc, char *argv[])
75 const char *desc[] = {
76 "[THISMODULE] does a linear interpolation of conformations in order to",
77 "create intermediates. Of course these are completely unphysical, but",
78 "that you may try to justify yourself. Output is in the form of a ",
79 "generic trajectory. The number of intermediates can be controlled with",
80 "the [TT]-ninterm[tt] flag. The first and last flag correspond to the way of",
81 "interpolating: 0 corresponds to input structure 1 while",
82 "1 corresponds to input structure 2.",
83 "If you specify [TT]-first[tt] < 0 or [TT]-last[tt] > 1 extrapolation will be",
84 "on the path from input structure x[SUB]1[sub] to x[SUB]2[sub]. In general, the coordinates",
85 "of the intermediate x(i) out of N total intermediates correspond to:[PAR]",
86 "x(i) = x[SUB]1[sub] + (first+(i/(N-1))*(last-first))*(x[SUB]2[sub]-x[SUB]1[sub])[PAR]",
87 "Finally the RMSD with respect to both input structures can be computed",
88 "if explicitly selected ([TT]-or[tt] option). In that case, an index file may be",
89 "read to select the group from which the RMS is computed."
92 { efSTX, "-f1", "conf1", ffREAD },
93 { efSTX, "-f2", "conf2", ffREAD },
94 { efTRX, "-o", "interm", ffWRITE },
95 { efXVG, "-or", "rms-interm", ffOPTWR },
96 { efNDX, "-n", "index", ffOPTRD }
98 #define NFILE asize(fnm)
99 static int ninterm = 11;
100 static real first = 0.0;
101 static real last = 1.0;
102 static gmx_bool bFit = TRUE;
104 { "-ninterm", FALSE, etINT, {&ninterm},
105 "Number of intermediates" },
106 { "-first", FALSE, etREAL, {&first},
107 "Corresponds to first generated structure (0 is input x[SUB]1[sub], see above)" },
108 { "-last", FALSE, etREAL, {&last},
109 "Corresponds to last generated structure (1 is input x[SUB]2[sub], see above)" },
110 { "-fit", FALSE, etBOOL, {&bFit},
111 "Do a least squares fit of the second to the first structure before interpolating" }
113 const char *leg[] = { "Ref = 1\\Sst\\N conf", "Ref = 2\\Snd\\N conf" };
115 int i, isize, is_lsq, nat1, nat2;
117 atom_id *index, *index_lsq, *index_all, *dummy;
119 rvec *x1, *x2, *xx, *v;
121 real rms1, rms2, fac, *mass;
122 char title[STRLEN], *grpname;
126 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW,
127 NFILE, fnm, asize(pa), pa, asize(desc), desc,
132 get_stx_coordnum (opt2fn("-f1", NFILE, fnm), &nat1);
133 get_stx_coordnum (opt2fn("-f2", NFILE, fnm), &nat2);
136 gmx_fatal(FARGS, "Number of atoms in first structure is %d, in second %d",
140 init_t_atoms(&atoms, nat1, TRUE);
146 read_stx_conf(opt2fn("-f1", NFILE, fnm), title, &atoms, x1, v, NULL, box);
147 read_stx_conf(opt2fn("-f2", NFILE, fnm), title, &atoms, x2, v, NULL, box);
150 snew(index_all, nat1);
151 for (i = 0; (i < nat1); i++)
158 printf("Select group for LSQ superposition:\n");
159 get_index(&atoms, opt2fn_null("-n", NFILE, fnm), 1, &is_lsq, &index_lsq,
161 reset_x(is_lsq, index_lsq, nat1, index_all, x1, mass);
162 reset_x(is_lsq, index_lsq, nat1, index_all, x2, mass);
163 do_fit(nat1, mass, x1, x2);
166 bRMS = opt2bSet("-or", NFILE, fnm);
169 fp = xvgropen(opt2fn("-or", NFILE, fnm), "RMSD", "Conf", "(nm)", oenv);
170 xvgr_legend(fp, asize(leg), leg, oenv);
171 printf("Select group for RMSD calculation:\n");
172 get_index(&atoms, opt2fn_null("-n", NFILE, fnm), 1, &isize, &index, &grpname);
173 printf("You selected group %s, containing %d atoms\n", grpname, isize);
174 rms1 = rmsdev_ind(isize, index, mass, x1, x2);
175 fprintf(stderr, "RMSD between input conformations is %g nm\n", rms1);
179 for (i = 0; (i < nat1); i++)
183 status = open_trx(ftp2fn(efTRX, NFILE, fnm), "w");
185 for (i = 0; (i < ninterm); i++)
187 fac = dointerp(nat1, x1, x2, xx, i, ninterm, first, last);
188 write_trx(status, nat1, dummy, &atoms, i, fac, box, xx, NULL, NULL);
191 rms1 = rmsdev_ind(isize, index, mass, x1, xx);
192 rms2 = rmsdev_ind(isize, index, mass, x2, xx);
193 fprintf(fp, "%10g %10g %10g\n", fac, rms1, rms2);
202 do_view(oenv, opt2fn("-or", NFILE, fnm), "-nxy");