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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
64 static void get_refx(output_env_t oenv,const char *trxfn,int nfitdim,int skip,
66 gmx_bool bMW,t_topology *top,int ePBC,rvec *x_ref)
68 int natoms,nfr_all,nfr,i,j,a,r,c,min_fr;
71 double tot_mass,msd,*srmsd,min_srmsd,srmsd_tot;
76 gmx_rmpbc_t gpbc=NULL;
83 natoms = read_first_x(oenv,&status,trxfn,&ti[nfr],&x,box);
89 if (index[a] >= natoms)
91 gmx_fatal(FARGS,"Atom index (%d) is larger than the number of atoms in the trajecory (%d)",index[a]+1,natoms);
93 w_rls[a] = (bMW ? top->atoms.atom[index[a]].m : 1.0);
96 gpbc = gmx_rmpbc_init(&top->idef,ePBC,natoms,box);
100 if (nfr_all % skip == 0)
102 gmx_rmpbc(gpbc,natoms,box,x);
106 copy_rvec(x[index[i]],xi[nfr][i]);
108 reset_x(gnx,NULL,gnx,NULL,xi[nfr],w_rls);
118 while(read_next_x(oenv,status,&ti[nfr],natoms,x,box));
122 gmx_rmpbc_done(gpbc);
127 printf("\rProcessing frame %d of %d",i,nfr);
128 for(j=i+1; j<nfr; j++)
130 calc_fit_R(nfitdim,gnx,w_rls,xi[i],xi[j],R);
140 xf += R[r][c]*xi[j][a][c];
142 msd += w_rls[a]*sqr(xi[i][a][r] - xf);
146 srmsd[i] += sqrt(msd);
147 srmsd[j] += sqrt(msd);
154 min_srmsd = GMX_REAL_MAX;
160 srmsd[i] /= (nfr - 1);
161 if (srmsd[i] < min_srmsd)
163 min_srmsd = srmsd[i];
167 srmsd_tot += srmsd[i];
171 printf("Average RMSD between all structures: %.3f\n",srmsd_tot/nfr);
172 printf("Structure with lowest RMSD to all others: time %g, av. RMSD %.3f\n",
177 copy_rvec(xi[min_fr][a],x_ref[index[a]]);
183 int gmx_rotmat(int argc,char *argv[])
185 const char *desc[] = {
186 "[TT]g_rotmat[tt] plots the rotation matrix required for least squares fitting",
187 "a conformation onto the reference conformation provided with",
188 "[TT]-s[tt]. Translation is removed before fitting.",
189 "The output are the three vectors that give the new directions",
190 "of the x, y and z directions of the reference conformation,",
191 "for example: (zx,zy,zz) is the orientation of the reference",
192 "z-axis in the trajectory frame.",
194 "This tool is useful for, for instance,",
195 "determining the orientation of a molecule",
196 "at an interface, possibly on a trajectory produced with",
197 "[TT]trjconv -fit rotxy+transxy[tt] to remove the rotation",
198 "in the [IT]x-y[it] plane.",
200 "Option [TT]-ref[tt] determines a reference structure for fitting,",
201 "instead of using the structure from [TT]-s[tt]. The structure with",
202 "the lowest sum of RMSD's to all other structures is used.",
203 "Since the computational cost of this procedure grows with",
204 "the square of the number of frames, the [TT]-skip[tt] option",
205 "can be useful. A full fit or only a fit in the [IT]x-y[it] plane can",
208 "Option [TT]-fitxy[tt] fits in the [IT]x-y[it] plane before determining",
209 "the rotation matrix."
211 const char *reffit[] =
212 { NULL, "none", "xyz", "xy", NULL };
214 static gmx_bool bFitXY=FALSE,bMW=TRUE;
216 { "-ref", FALSE, etENUM, {reffit},
217 "Determine the optimal reference structure" },
218 { "-skip", FALSE, etINT, {&skip},
219 "Use every nr-th frame for [TT]-ref[tt]" },
220 { "-fitxy", FALSE, etBOOL, {&bFitXY},
221 "Fit the x/y rotation before determining the rotation" },
222 { "-mw", FALSE, etBOOL, {&bMW},
223 "Use mass weighted fitting" }
233 char *grpname,title[256];
235 gmx_rmpbc_t gpbc=NULL;
239 const char *leg[] = { "xx", "xy", "xz", "yx", "yy", "yz", "zx", "zy", "zz" };
240 #define NLEG asize(leg)
242 { efTRX, "-f", NULL, ffREAD },
243 { efTPS, NULL, NULL, ffREAD },
244 { efNDX, NULL, NULL, ffOPTRD },
245 { efXVG, NULL, "rotmat", ffWRITE }
247 #define NFILE asize(fnm)
249 CopyRight(stderr,argv[0]);
251 parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
252 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
254 read_tps_conf(ftp2fn(efTPS,NFILE,fnm),title,&top,&ePBC,&x_ref,NULL,box,bMW);
256 gpbc = gmx_rmpbc_init(&top.idef,ePBC,top.atoms.nr,box);
258 gmx_rmpbc(gpbc,top.atoms.nr,box,x_ref);
260 get_index(&top.atoms,ftp2fn_null(efNDX,NFILE,fnm),1,&gnx,&index,&grpname);
262 if (reffit[0][0] != 'n')
264 get_refx(oenv,ftp2fn(efTRX,NFILE,fnm),reffit[0][2]=='z' ? 3 : 2,skip,
265 gnx,index,bMW,&top,ePBC,x_ref);
268 natoms = read_first_x(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
273 if (index[i] >= natoms)
275 gmx_fatal(FARGS,"Atom index (%d) is larger than the number of atoms in the trajecory (%d)",index[i]+1,natoms);
277 w_rls[index[i]] = (bMW ? top.atoms.atom[index[i]].m : 1.0);
280 if (reffit[0][0] == 'n')
282 reset_x(gnx,index,natoms,NULL,x_ref,w_rls);
285 out = xvgropen(ftp2fn(efXVG,NFILE,fnm),
286 "Fit matrix","Time (ps)","",oenv);
287 xvgr_legend(out,NLEG,leg,oenv);
291 gmx_rmpbc(gpbc,natoms,box,x);
293 reset_x(gnx,index,natoms,NULL,x,w_rls);
297 do_fit_ndim(2,natoms,w_rls,x_ref,x);
300 calc_fit_R(DIM,natoms,w_rls,x_ref,x,R);
303 "%7g %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f\n",
305 R[XX][XX],R[XX][YY],R[XX][ZZ],
306 R[YY][XX],R[YY][YY],R[YY][ZZ],
307 R[ZZ][XX],R[ZZ][YY],R[ZZ][ZZ]);
309 while(read_next_x(oenv,status,&t,natoms,x,box));
311 gmx_rmpbc_done(gpbc);
317 do_view(oenv,ftp2fn(efXVG,NFILE,fnm),"-nxy");