3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
66 static int ocomp(const void *a, const void *b)
83 int gmx_trjorder(int argc, char *argv[])
85 const char *desc[] = {
86 "[TT]trjorder[tt] orders molecules according to the smallest distance",
87 "to atoms in a reference group",
88 "or on z-coordinate (with option [TT]-z[tt]).",
89 "With distance ordering, it will ask for a group of reference",
90 "atoms and a group of molecules. For each frame of the trajectory",
91 "the selected molecules will be reordered according to the shortest",
92 "distance between atom number [TT]-da[tt] in the molecule and all the",
93 "atoms in the reference group. The center of mass of the molecules can",
94 "be used instead of a reference atom by setting [TT]-da[tt] to 0.",
95 "All atoms in the trajectory are written",
96 "to the output trajectory.[PAR]",
97 "[TT]trjorder[tt] can be useful for e.g. analyzing the n waters closest to a",
99 "In that case the reference group would be the protein and the group",
100 "of molecules would consist of all the water atoms. When an index group",
101 "of the first n waters is made, the ordered trajectory can be used",
102 "with any Gromacs program to analyze the n closest waters.",
104 "If the output file is a [TT].pdb[tt] file, the distance to the reference target",
105 "will be stored in the B-factor field in order to color with e.g. Rasmol.",
107 "With option [TT]-nshell[tt] the number of molecules within a shell",
108 "of radius [TT]-r[tt] around the reference group are printed."
110 static int na = 3, ref_a = 1;
111 static real rcut = 0;
112 static gmx_bool bCOM = FALSE, bZ = FALSE;
114 { "-na", FALSE, etINT, {&na},
115 "Number of atoms in a molecule" },
116 { "-da", FALSE, etINT, {&ref_a},
117 "Atom used for the distance calculation, 0 is COM" },
118 { "-com", FALSE, etBOOL, {&bCOM},
119 "Use the distance to the center of mass of the reference group" },
120 { "-r", FALSE, etREAL, {&rcut},
121 "Cutoff used for the distance calculation when computing the number of molecules in a shell around e.g. a protein" },
122 { "-z", FALSE, etBOOL, {&bZ},
123 "Order molecules on z-coordinate" }
128 gmx_bool bNShell, bPDBout;
131 rvec *x, *xsol, xcom, dx;
135 real t, totmass, mass, rcut2 = 0, n2;
136 int natoms, nwat, ncut;
137 char **grpname, title[256];
138 int i, j, d, *isize, isize_ref = 0, isize_sol;
139 atom_id sa, sr, *swi, **index, *ind_ref = NULL, *ind_sol;
142 { efTRX, "-f", NULL, ffREAD },
143 { efTPS, NULL, NULL, ffREAD },
144 { efNDX, NULL, NULL, ffOPTRD },
145 { efTRO, "-o", "ordered", ffOPTWR },
146 { efXVG, "-nshell", "nshell", ffOPTWR }
148 #define NFILE asize(fnm)
150 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_BE_NICE,
151 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
156 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x, NULL, box, TRUE);
159 /* get index groups */
160 printf("Select %sa group of molecules to be ordered:\n",
161 bZ ? "" : "a group of reference atoms and ");
165 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), bZ ? 1 : 2,
166 isize, index, grpname);
170 isize_ref = isize[0];
171 isize_sol = isize[1];
177 isize_sol = isize[0];
181 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
182 if (natoms > top.atoms.nr)
184 gmx_fatal(FARGS, "Number of atoms in the run input file is larger than in the trjactory");
186 for (i = 0; (i < 2); i++)
188 for (j = 0; (j < isize[i]); j++)
190 if (index[i][j] > natoms)
192 gmx_fatal(FARGS, "An atom number in group %s is larger than the number of atoms in the trajectory");
197 if ((isize_sol % na) != 0)
199 gmx_fatal(FARGS, "Number of atoms in the molecule group (%d) is not a multiple of na (%d)",
206 gmx_fatal(FARGS, "The reference atom can not be larger than the number of atoms in a molecule");
212 for (i = 0; (i < natoms); i++)
219 bNShell = ((opt2bSet("-nshell", NFILE, fnm)) ||
220 (opt2parg_bSet("-r", asize(pa), pa)));
225 fp = xvgropen(opt2fn("-nshell", NFILE, fnm), "Number of molecules",
226 "Time (ps)", "N", oenv);
227 printf("Will compute the number of molecules within a radius of %g\n",
230 if (!bNShell || opt2bSet("-o", NFILE, fnm))
232 bPDBout = (fn2ftp(opt2fn("-o", NFILE, fnm)) == efPDB);
233 if (bPDBout && !top.atoms.pdbinfo)
235 fprintf(stderr, "Creating pdbfino records\n");
236 snew(top.atoms.pdbinfo, top.atoms.nr);
238 out = open_trx(opt2fn("-o", NFILE, fnm), "w");
240 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
243 gmx_rmpbc(gpbc, natoms, box, x);
244 set_pbc(&pbc, ePBC, box);
248 /* Calculate the COM of all solvent molecules */
249 for (i = 0; i < nwat; i++)
253 for (j = 0; j < na; j++)
255 sa = ind_sol[i*na+j];
256 mass = top.atoms.atom[sa].m;
258 for (d = 0; d < DIM; d++)
260 xsol[i][d] += mass*x[sa][d];
263 svmul(1/totmass, xsol[i], xsol[i]);
268 /* Copy the reference atom of all solvent molecules */
269 for (i = 0; i < nwat; i++)
271 copy_rvec(x[ind_sol[i*na+ref_a]], xsol[i]);
277 for (i = 0; (i < nwat); i++)
281 order[i].d2 = xsol[i][ZZ];
288 for (i = 0; i < isize_ref; i++)
290 mass = top.atoms.atom[ind_ref[i]].m;
292 for (j = 0; j < DIM; j++)
294 xcom[j] += mass*x[ind_ref[i]][j];
297 svmul(1/totmass, xcom, xcom);
298 for (i = 0; (i < nwat); i++)
301 pbc_dx(&pbc, xcom, xsol[i], dx);
303 order[i].d2 = norm2(dx);
308 /* Set distance to first atom */
309 for (i = 0; (i < nwat); i++)
312 pbc_dx(&pbc, x[ind_ref[0]], xsol[i], dx);
314 order[i].d2 = norm2(dx);
316 for (j = 1; (j < isize_ref); j++)
319 for (i = 0; (i < nwat); i++)
322 pbc_dx(&pbc, x[sr], xsol[i], dx);
324 if (n2 < order[i].d2)
335 for (i = 0; (i < nwat); i++)
337 if (order[i].d2 <= rcut2)
342 fprintf(fp, "%10.3f %8d\n", t, ncut);
346 qsort(order, nwat, sizeof(*order), ocomp);
347 for (i = 0; (i < nwat); i++)
349 for (j = 0; (j < na); j++)
351 swi[ind_sol[na*i]+j] = order[i].i+j;
355 /* Store the distance as the B-factor */
358 for (i = 0; (i < nwat); i++)
360 for (j = 0; (j < na); j++)
362 top.atoms.pdbinfo[order[i].i+j].bfac = sqrt(order[i].d2);
366 write_trx(out, natoms, swi, &top.atoms, 0, t, box, x, NULL, NULL);
369 while (read_next_x(oenv, status, &t, x, box));
379 gmx_rmpbc_done(gpbc);