Remove unused code detected by PGI compiler
[gromacs.git] / src / gromacs / gmxana / gmx_trjorder.cpp
blobc9767a9ece3a6ee2b4a039606bdfc508744ad117
1 /*
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,2015, 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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdlib>
41 #include <cstring>
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/fileio/trx.h"
46 #include "gromacs/fileio/trxio.h"
47 #include "gromacs/fileio/txtdump.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/gmxana/gmx_ana.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/pbcutil/pbc.h"
52 #include "gromacs/pbcutil/rmpbc.h"
53 #include "gromacs/topology/index.h"
54 #include "gromacs/topology/topology.h"
55 #include "gromacs/utility/arraysize.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/smalloc.h"
60 typedef struct {
61 int i;
62 real d2;
63 } t_order;
65 t_order *order;
67 static int ocomp(const void *a, const void *b)
69 t_order *oa, *ob;
71 oa = (t_order *)a;
72 ob = (t_order *)b;
74 if (oa->d2 < ob->d2)
76 return -1;
78 else
80 return 1;
84 int gmx_trjorder(int argc, char *argv[])
86 const char *desc[] = {
87 "[THISMODULE] orders molecules according to the smallest distance",
88 "to atoms in a reference group",
89 "or on z-coordinate (with option [TT]-z[tt]).",
90 "With distance ordering, it will ask for a group of reference",
91 "atoms and a group of molecules. For each frame of the trajectory",
92 "the selected molecules will be reordered according to the shortest",
93 "distance between atom number [TT]-da[tt] in the molecule and all the",
94 "atoms in the reference group. The center of mass of the molecules can",
95 "be used instead of a reference atom by setting [TT]-da[tt] to 0.",
96 "All atoms in the trajectory are written",
97 "to the output trajectory.[PAR]",
98 "[THISMODULE] can be useful for e.g. analyzing the n waters closest to a",
99 "protein.",
100 "In that case the reference group would be the protein and the group",
101 "of molecules would consist of all the water atoms. When an index group",
102 "of the first n waters is made, the ordered trajectory can be used",
103 "with any GROMACS program to analyze the n closest waters.",
104 "[PAR]",
105 "If the output file is a [REF].pdb[ref] file, the distance to the reference target",
106 "will be stored in the B-factor field in order to color with e.g. Rasmol.",
107 "[PAR]",
108 "With option [TT]-nshell[tt] the number of molecules within a shell",
109 "of radius [TT]-r[tt] around the reference group are printed."
111 static int na = 3, ref_a = 1;
112 static real rcut = 0;
113 static gmx_bool bCOM = FALSE, bZ = FALSE;
114 t_pargs pa[] = {
115 { "-na", FALSE, etINT, {&na},
116 "Number of atoms in a molecule" },
117 { "-da", FALSE, etINT, {&ref_a},
118 "Atom used for the distance calculation, 0 is COM" },
119 { "-com", FALSE, etBOOL, {&bCOM},
120 "Use the distance to the center of mass of the reference group" },
121 { "-r", FALSE, etREAL, {&rcut},
122 "Cutoff used for the distance calculation when computing the number of molecules in a shell around e.g. a protein" },
123 { "-z", FALSE, etBOOL, {&bZ},
124 "Order molecules on z-coordinate" }
126 FILE *fp;
127 t_trxstatus *out;
128 t_trxstatus *status;
129 gmx_bool bNShell, bPDBout;
130 t_topology top;
131 int ePBC;
132 rvec *x, *xsol, xcom, dx;
133 matrix box;
134 t_pbc pbc;
135 gmx_rmpbc_t gpbc;
136 real t, totmass, mass, rcut2 = 0, n2;
137 int natoms, nwat, ncut;
138 char **grpname;
139 int i, j, d, *isize, isize_ref = 0, isize_sol;
140 int sa, sr, *swi, **index, *ind_ref = NULL, *ind_sol;
141 gmx_output_env_t *oenv;
142 t_filenm fnm[] = {
143 { efTRX, "-f", NULL, ffREAD },
144 { efTPS, NULL, NULL, ffREAD },
145 { efNDX, NULL, NULL, ffOPTRD },
146 { efTRO, "-o", "ordered", ffOPTWR },
147 { efXVG, "-nshell", "nshell", ffOPTWR }
149 #define NFILE asize(fnm)
151 if (!parse_common_args(&argc, argv, PCA_CAN_TIME,
152 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
154 return 0;
157 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, &x, NULL, box, TRUE);
158 sfree(x);
160 /* get index groups */
161 printf("Select %sa group of molecules to be ordered:\n",
162 bZ ? "" : "a group of reference atoms and ");
163 snew(grpname, 2);
164 snew(index, 2);
165 snew(isize, 2);
166 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), bZ ? 1 : 2,
167 isize, index, grpname);
169 if (!bZ)
171 isize_ref = isize[0];
172 isize_sol = isize[1];
173 ind_ref = index[0];
174 ind_sol = index[1];
176 else
178 isize_sol = isize[0];
179 ind_sol = index[0];
182 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
183 if (natoms > top.atoms.nr)
185 gmx_fatal(FARGS, "Number of atoms in the run input file is larger than in the trjactory");
187 for (i = 0; (i < 2); i++)
189 for (j = 0; (j < isize[i]); j++)
191 if (index[i][j] > natoms)
193 gmx_fatal(FARGS, "An atom number in group %s is larger than the number of atoms in the trajectory");
198 if ((isize_sol % na) != 0)
200 gmx_fatal(FARGS, "Number of atoms in the molecule group (%d) is not a multiple of na (%d)",
201 isize[1], na);
204 nwat = isize_sol/na;
205 if (ref_a > na)
207 gmx_fatal(FARGS, "The reference atom can not be larger than the number of atoms in a molecule");
209 ref_a--;
210 snew(xsol, nwat);
211 snew(order, nwat);
212 snew(swi, natoms);
213 for (i = 0; (i < natoms); i++)
215 swi[i] = i;
218 out = NULL;
219 fp = NULL;
220 bNShell = ((opt2bSet("-nshell", NFILE, fnm)) ||
221 (opt2parg_bSet("-r", asize(pa), pa)));
222 bPDBout = FALSE;
223 if (bNShell)
225 rcut2 = rcut*rcut;
226 fp = xvgropen(opt2fn("-nshell", NFILE, fnm), "Number of molecules",
227 "Time (ps)", "N", oenv);
228 printf("Will compute the number of molecules within a radius of %g\n",
229 rcut);
231 if (!bNShell || opt2bSet("-o", NFILE, fnm))
233 bPDBout = (fn2ftp(opt2fn("-o", NFILE, fnm)) == efPDB);
234 if (bPDBout && !top.atoms.pdbinfo)
236 fprintf(stderr, "Creating pdbfino records\n");
237 snew(top.atoms.pdbinfo, top.atoms.nr);
239 out = open_trx(opt2fn("-o", NFILE, fnm), "w");
241 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
244 gmx_rmpbc(gpbc, natoms, box, x);
245 set_pbc(&pbc, ePBC, box);
247 if (ref_a == -1)
249 /* Calculate the COM of all solvent molecules */
250 for (i = 0; i < nwat; i++)
252 totmass = 0;
253 clear_rvec(xsol[i]);
254 for (j = 0; j < na; j++)
256 sa = ind_sol[i*na+j];
257 mass = top.atoms.atom[sa].m;
258 totmass += mass;
259 for (d = 0; d < DIM; d++)
261 xsol[i][d] += mass*x[sa][d];
264 svmul(1.0/totmass, xsol[i], xsol[i]);
267 else
269 /* Copy the reference atom of all solvent molecules */
270 for (i = 0; i < nwat; i++)
272 copy_rvec(x[ind_sol[i*na+ref_a]], xsol[i]);
276 if (bZ)
278 for (i = 0; (i < nwat); i++)
280 sa = ind_sol[na*i];
281 order[i].i = sa;
282 order[i].d2 = xsol[i][ZZ];
285 else if (bCOM)
287 totmass = 0;
288 clear_rvec(xcom);
289 for (i = 0; i < isize_ref; i++)
291 mass = top.atoms.atom[ind_ref[i]].m;
292 totmass += mass;
293 for (j = 0; j < DIM; j++)
295 xcom[j] += mass*x[ind_ref[i]][j];
298 svmul(1/totmass, xcom, xcom);
299 for (i = 0; (i < nwat); i++)
301 sa = ind_sol[na*i];
302 pbc_dx(&pbc, xcom, xsol[i], dx);
303 order[i].i = sa;
304 order[i].d2 = norm2(dx);
307 else
309 /* Set distance to first atom */
310 for (i = 0; (i < nwat); i++)
312 sa = ind_sol[na*i];
313 pbc_dx(&pbc, x[ind_ref[0]], xsol[i], dx);
314 order[i].i = sa;
315 order[i].d2 = norm2(dx);
317 for (j = 1; (j < isize_ref); j++)
319 sr = ind_ref[j];
320 for (i = 0; (i < nwat); i++)
322 pbc_dx(&pbc, x[sr], xsol[i], dx);
323 n2 = norm2(dx);
324 if (n2 < order[i].d2)
326 order[i].d2 = n2;
332 if (bNShell)
334 ncut = 0;
335 for (i = 0; (i < nwat); i++)
337 if (order[i].d2 <= rcut2)
339 ncut++;
342 fprintf(fp, "%10.3f %8d\n", t, ncut);
344 if (out)
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 */
356 if (bPDBout)
358 for (i = 0; (i < nwat); i++)
360 for (j = 0; (j < na); j++)
362 top.atoms.pdbinfo[order[i].i+j].bfac = std::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));
370 close_trj(status);
371 if (out)
373 close_trx(out);
375 if (fp)
377 xvgrclose(fp);
379 gmx_rmpbc_done(gpbc);
381 return 0;