Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / gmxana / gmx_traj.cpp
blob606abcadff6e716b699335ac4451fd8ee69aed7e
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,2016,2017,2018, 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 <algorithm>
44 #include <string>
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/commandline/viewit.h"
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/fileio/trxio.h"
50 #include "gromacs/fileio/xvgr.h"
51 #include "gromacs/gmxana/gmx_ana.h"
52 #include "gromacs/linearalgebra/nrjac.h"
53 #include "gromacs/math/functions.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdtypes/md_enums.h"
57 #include "gromacs/pbcutil/pbc.h"
58 #include "gromacs/pbcutil/rmpbc.h"
59 #include "gromacs/topology/index.h"
60 #include "gromacs/topology/topology.h"
61 #include "gromacs/trajectory/trajectoryframe.h"
62 #include "gromacs/utility/arraysize.h"
63 #include "gromacs/utility/cstringutil.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/futil.h"
66 #include "gromacs/utility/smalloc.h"
67 #include "gromacs/utility/stringutil.h"
69 static void low_print_data(FILE *fp, real time, rvec x[], int n, const int *index,
70 const gmx_bool bDim[], const char *sffmt)
72 int i, ii, d;
74 fprintf(fp, " %g", time);
75 for (i = 0; i < n; i++)
77 if (index != nullptr)
79 ii = index[i];
81 else
83 ii = i;
85 for (d = 0; d < DIM; d++)
87 if (bDim[d])
89 fprintf(fp, sffmt, x[ii][d]);
92 if (bDim[DIM])
94 fprintf(fp, sffmt, norm(x[ii]));
97 fprintf(fp, "\n");
100 static void average_data(rvec x[], rvec xav[], const real *mass,
101 int ngrps, const int isize[], int **index)
103 int g, i, ind, d;
104 real m;
105 rvec tmp;
106 double sum[DIM], mtot;
108 for (g = 0; g < ngrps; g++)
110 clear_dvec(sum);
111 clear_rvec(xav[g]);
112 mtot = 0;
113 for (i = 0; i < isize[g]; i++)
115 ind = index[g][i];
116 if (mass != nullptr)
118 m = mass[ind];
119 svmul(m, x[ind], tmp);
120 for (d = 0; d < DIM; d++)
122 sum[d] += tmp[d];
124 mtot += m;
126 else
128 for (d = 0; d < DIM; d++)
130 sum[d] += x[ind][d];
134 if (mass != nullptr)
136 for (d = 0; d < DIM; d++)
138 xav[g][d] = sum[d]/mtot;
141 else
143 /* mass=NULL, so these are forces: sum only (do not average) */
144 for (d = 0; d < DIM; d++)
146 xav[g][d] = sum[d];
152 static void print_data(FILE *fp, real time, rvec x[], real *mass, gmx_bool bCom,
153 int ngrps, int isize[], int **index, gmx_bool bDim[],
154 const char *sffmt)
156 static rvec *xav = nullptr;
158 if (bCom)
160 if (xav == nullptr)
162 snew(xav, ngrps);
164 average_data(x, xav, mass, ngrps, isize, index);
165 low_print_data(fp, time, xav, ngrps, nullptr, bDim, sffmt);
167 else
169 low_print_data(fp, time, x, isize[0], index[0], bDim, sffmt);
173 static void write_trx_x(t_trxstatus *status, const t_trxframe *fr, real *mass, gmx_bool bCom,
174 int ngrps, int isize[], int **index)
176 static rvec *xav = nullptr;
177 static t_atoms *atoms = nullptr;
178 t_trxframe fr_av;
179 int i;
181 if (bCom)
183 if (xav == nullptr)
185 snew(xav, ngrps);
186 snew(atoms, 1);
187 *atoms = *fr->atoms;
188 snew(atoms->atom, ngrps);
189 atoms->nr = ngrps;
190 /* Note that atom and residue names will be the ones
191 * of the first atom in each group.
193 for (i = 0; i < ngrps; i++)
195 atoms->atom[i] = fr->atoms->atom[index[i][0]];
196 atoms->atomname[i] = fr->atoms->atomname[index[i][0]];
199 average_data(fr->x, xav, mass, ngrps, isize, index);
200 fr_av = *fr;
201 fr_av.natoms = ngrps;
202 fr_av.atoms = atoms;
203 fr_av.x = xav;
204 write_trxframe(status, &fr_av, nullptr);
206 else
208 write_trxframe_indexed(status, fr, isize[0], index[0], nullptr);
212 static void make_legend(FILE *fp, int ngrps, int isize, int index[],
213 char **name, gmx_bool bCom, gmx_bool bMol, const gmx_bool bDim[],
214 const gmx_output_env_t *oenv)
216 char **leg;
217 const char *dimtxt[] = { " X", " Y", " Z", "" };
218 int n, i, j, d;
220 if (bCom)
222 n = ngrps;
224 else
226 n = isize;
229 snew(leg, 4*n);
230 j = 0;
231 for (i = 0; i < n; i++)
233 for (d = 0; d <= DIM; d++)
235 if (bDim[d])
237 snew(leg[j], STRLEN);
238 if (bMol)
240 sprintf(leg[j], "mol %d%s", index[i]+1, dimtxt[d]);
242 else if (bCom)
244 sprintf(leg[j], "%s%s", name[i], dimtxt[d]);
246 else
248 sprintf(leg[j], "atom %d%s", index[i]+1, dimtxt[d]);
250 j++;
254 xvgr_legend(fp, j, leg, oenv);
256 for (i = 0; i < j; i++)
258 sfree(leg[i]);
260 sfree(leg);
263 static real ekrot(rvec x[], rvec v[], const real mass[], int isize, const int index[])
265 real TCM[5][5], L[5][5];
266 double tm, m0, lxx, lxy, lxz, lyy, lyz, lzz, ekrot;
267 rvec a0, ocm;
268 dvec dx, b0;
269 dvec xcm, vcm, acm;
270 int i, j, m, n;
272 clear_dvec(xcm);
273 clear_dvec(vcm);
274 clear_dvec(acm);
275 tm = 0.0;
276 for (i = 0; i < isize; i++)
278 j = index[i];
279 m0 = mass[j];
280 tm += m0;
281 cprod(x[j], v[j], a0);
282 for (m = 0; (m < DIM); m++)
284 xcm[m] += m0*x[j][m]; /* c.o.m. position */
285 vcm[m] += m0*v[j][m]; /* c.o.m. velocity */
286 acm[m] += m0*a0[m]; /* rotational velocity around c.o.m. */
289 dcprod(xcm, vcm, b0);
290 for (m = 0; (m < DIM); m++)
292 xcm[m] /= tm;
293 vcm[m] /= tm;
294 acm[m] -= b0[m]/tm;
297 lxx = lxy = lxz = lyy = lyz = lzz = 0.0;
298 for (i = 0; i < isize; i++)
300 j = index[i];
301 m0 = mass[j];
302 for (m = 0; m < DIM; m++)
304 dx[m] = x[j][m] - xcm[m];
306 lxx += dx[XX]*dx[XX]*m0;
307 lxy += dx[XX]*dx[YY]*m0;
308 lxz += dx[XX]*dx[ZZ]*m0;
309 lyy += dx[YY]*dx[YY]*m0;
310 lyz += dx[YY]*dx[ZZ]*m0;
311 lzz += dx[ZZ]*dx[ZZ]*m0;
314 L[XX][XX] = lyy + lzz;
315 L[YY][XX] = -lxy;
316 L[ZZ][XX] = -lxz;
317 L[XX][YY] = -lxy;
318 L[YY][YY] = lxx + lzz;
319 L[ZZ][YY] = -lyz;
320 L[XX][ZZ] = -lxz;
321 L[YY][ZZ] = -lyz;
322 L[ZZ][ZZ] = lxx + lyy;
324 m_inv_gen(&L[0][0], DIM, &TCM[0][0]);
326 /* Compute omega (hoeksnelheid) */
327 clear_rvec(ocm);
328 ekrot = 0;
329 for (m = 0; m < DIM; m++)
331 for (n = 0; n < DIM; n++)
333 ocm[m] += TCM[m][n]*acm[n];
335 ekrot += 0.5*ocm[m]*acm[m];
338 return ekrot;
341 static real ektrans(rvec v[], const real mass[], int isize, const int index[])
343 dvec mvcom;
344 double mtot;
345 int i, j, d;
347 clear_dvec(mvcom);
348 mtot = 0;
349 for (i = 0; i < isize; i++)
351 j = index[i];
352 for (d = 0; d < DIM; d++)
354 mvcom[d] += mass[j]*v[j][d];
356 mtot += mass[j];
359 return dnorm2(mvcom)/(mtot*2);
362 static real temp(rvec v[], const real mass[], int isize, const int index[])
364 double ekin2;
365 int i, j;
367 ekin2 = 0;
368 for (i = 0; i < isize; i++)
370 j = index[i];
371 ekin2 += mass[j]*norm2(v[j]);
374 return ekin2/(3*isize*BOLTZ);
377 static void remove_jump(matrix box, int natoms, rvec xp[], rvec x[])
379 rvec hbox;
380 int d, i, m;
382 for (d = 0; d < DIM; d++)
384 hbox[d] = 0.5*box[d][d];
386 for (i = 0; i < natoms; i++)
388 for (m = DIM-1; m >= 0; m--)
390 while (x[i][m]-xp[i][m] <= -hbox[m])
392 for (d = 0; d <= m; d++)
394 x[i][d] += box[m][d];
397 while (x[i][m]-xp[i][m] > hbox[m])
399 for (d = 0; d <= m; d++)
401 x[i][d] -= box[m][d];
408 static void write_pdb_bfac(const char *fname, const char *xname,
409 const char *title, t_atoms *atoms, int ePBC, matrix box,
410 int isize, int *index, int nfr_x, rvec *x,
411 int nfr_v, rvec *sum,
412 const gmx_bool bDim[], real scale_factor,
413 const gmx_output_env_t *oenv)
415 FILE *fp;
416 real max, len2, scale;
417 int maxi;
418 int i, m, onedim;
420 if ((nfr_x == 0) || (nfr_v == 0))
422 fprintf(stderr, "No frames found for %s, will not write %s\n",
423 title, fname);
425 else
427 fprintf(stderr, "Used %d frames for %s\n", nfr_x, "coordinates");
428 fprintf(stderr, "Used %d frames for %s\n", nfr_v, title);
429 onedim = -1;
430 if (!bDim[DIM])
432 m = 0;
433 for (i = 0; i < DIM; i++)
435 if (bDim[i])
437 onedim = i;
438 m++;
441 if (m != 1)
443 onedim = -1;
446 scale = 1.0/nfr_v;
447 for (i = 0; i < isize; i++)
449 svmul(scale, sum[index[i]], sum[index[i]]);
452 fp = xvgropen(xname, title, "Atom", "Spatial component", oenv);
453 for (i = 0; i < isize; i++)
455 fprintf(fp, "%-5d %10.3f %10.3f %10.3f\n", 1+i,
456 sum[index[i]][XX], sum[index[i]][YY], sum[index[i]][ZZ]);
458 xvgrclose(fp);
459 max = 0;
460 maxi = 0;
461 for (i = 0; i < isize; i++)
463 len2 = 0;
464 for (m = 0; m < DIM; m++)
466 if (bDim[m] || bDim[DIM])
468 len2 += gmx::square(sum[index[i]][m]);
471 if (len2 > max)
473 max = len2;
474 maxi = index[i];
477 if (scale_factor != 0)
479 scale = scale_factor;
481 else
483 if (max == 0)
485 scale = 1;
487 else
489 scale = 10.0/std::sqrt(max);
493 printf("Maximum %s is %g on atom %d %s, res. %s %d\n",
494 title, std::sqrt(max), maxi+1, *(atoms->atomname[maxi]),
495 *(atoms->resinfo[atoms->atom[maxi].resind].name),
496 atoms->resinfo[atoms->atom[maxi].resind].nr);
498 if (atoms->pdbinfo == nullptr)
500 snew(atoms->pdbinfo, atoms->nr);
502 atoms->havePdbInfo = TRUE;
504 if (onedim == -1)
506 for (i = 0; i < isize; i++)
508 len2 = 0;
509 for (m = 0; m < DIM; m++)
511 if (bDim[m] || bDim[DIM])
513 len2 += gmx::square(sum[index[i]][m]);
516 atoms->pdbinfo[index[i]].bfac = std::sqrt(len2)*scale;
519 else
521 for (i = 0; i < isize; i++)
523 atoms->pdbinfo[index[i]].bfac = sum[index[i]][onedim]*scale;
526 write_sto_conf_indexed(fname, title, atoms, x, nullptr, ePBC, box, isize, index);
530 static void update_histo(int gnx, const int index[], rvec v[],
531 int *nhisto, int **histo, real binwidth)
533 int i, m, in, nnn;
534 real vn, vnmax;
536 if (*histo == nullptr)
538 vnmax = 0;
539 for (i = 0; (i < gnx); i++)
541 vn = norm(v[index[i]]);
542 vnmax = std::max(vn, vnmax);
544 vnmax *= 2;
545 *nhisto = static_cast<int>(1+(vnmax/binwidth));
546 snew(*histo, *nhisto);
548 for (i = 0; (i < gnx); i++)
550 vn = norm(v[index[i]]);
551 in = static_cast<int>(vn/binwidth);
552 if (in >= *nhisto)
554 nnn = in+100;
555 fprintf(stderr, "Extending histogram from %d to %d\n", *nhisto, nnn);
557 srenew(*histo, nnn);
558 for (m = *nhisto; (m < nnn); m++)
560 (*histo)[m] = 0;
562 *nhisto = nnn;
564 (*histo)[in]++;
568 static void print_histo(const char *fn, int nhisto, int histo[], real binwidth,
569 const gmx_output_env_t *oenv)
571 FILE *fp;
572 int i;
574 fp = xvgropen(fn, "Velocity distribution", "V (nm/ps)", "arbitrary units",
575 oenv);
576 for (i = 0; (i < nhisto); i++)
578 fprintf(fp, "%10.3e %10d\n", i*binwidth, histo[i]);
580 xvgrclose(fp);
583 int gmx_traj(int argc, char *argv[])
585 const char *desc[] = {
586 "[THISMODULE] plots coordinates, velocities, forces and/or the box.",
587 "With [TT]-com[tt] the coordinates, velocities and forces are",
588 "calculated for the center of mass of each group.",
589 "When [TT]-mol[tt] is set, the numbers in the index file are",
590 "interpreted as molecule numbers and the same procedure as with",
591 "[TT]-com[tt] is used for each molecule.[PAR]",
592 "Option [TT]-ot[tt] plots the temperature of each group,",
593 "provided velocities are present in the trajectory file.",
594 "No corrections are made for constrained degrees of freedom!",
595 "This implies [TT]-com[tt].[PAR]",
596 "Options [TT]-ekt[tt] and [TT]-ekr[tt] plot the translational and",
597 "rotational kinetic energy of each group,",
598 "provided velocities are present in the trajectory file.",
599 "This implies [TT]-com[tt].[PAR]",
600 "Options [TT]-cv[tt] and [TT]-cf[tt] write the average velocities",
601 "and average forces as temperature factors to a [REF].pdb[ref] file with",
602 "the average coordinates or the coordinates at [TT]-ctime[tt].",
603 "The temperature factors are scaled such that the maximum is 10.",
604 "The scaling can be changed with the option [TT]-scale[tt].",
605 "To get the velocities or forces of one",
606 "frame set both [TT]-b[tt] and [TT]-e[tt] to the time of",
607 "desired frame. When averaging over frames you might need to use",
608 "the [TT]-nojump[tt] option to obtain the correct average coordinates.",
609 "If you select either of these option the average force and velocity",
610 "for each atom are written to an [REF].xvg[ref] file as well",
611 "(specified with [TT]-av[tt] or [TT]-af[tt]).[PAR]",
612 "Option [TT]-vd[tt] computes a velocity distribution, i.e. the",
613 "norm of the vector is plotted. In addition in the same graph",
614 "the kinetic energy distribution is given.",
616 "See [gmx-trajectory] for plotting similar data for selections."
618 static gmx_bool bMol = FALSE, bCom = FALSE, bPBC = TRUE, bNoJump = FALSE;
619 static gmx_bool bX = TRUE, bY = TRUE, bZ = TRUE, bNorm = FALSE, bFP = FALSE;
620 static int ngroups = 1;
621 static real ctime = -1, scale = 0, binwidth = 1;
622 t_pargs pa[] = {
623 { "-com", FALSE, etBOOL, {&bCom},
624 "Plot data for the com of each group" },
625 { "-pbc", FALSE, etBOOL, {&bPBC},
626 "Make molecules whole for COM" },
627 { "-mol", FALSE, etBOOL, {&bMol},
628 "Index contains molecule numbers instead of atom numbers" },
629 { "-nojump", FALSE, etBOOL, {&bNoJump},
630 "Remove jumps of atoms across the box" },
631 { "-x", FALSE, etBOOL, {&bX},
632 "Plot X-component" },
633 { "-y", FALSE, etBOOL, {&bY},
634 "Plot Y-component" },
635 { "-z", FALSE, etBOOL, {&bZ},
636 "Plot Z-component" },
637 { "-ng", FALSE, etINT, {&ngroups},
638 "Number of groups to consider" },
639 { "-len", FALSE, etBOOL, {&bNorm},
640 "Plot vector length" },
641 { "-fp", FALSE, etBOOL, {&bFP},
642 "Full precision output" },
643 { "-bin", FALSE, etREAL, {&binwidth},
644 "Binwidth for velocity histogram (nm/ps)" },
645 { "-ctime", FALSE, etREAL, {&ctime},
646 "Use frame at this time for x in [TT]-cv[tt] and [TT]-cf[tt] instead of the average x" },
647 { "-scale", FALSE, etREAL, {&scale},
648 "Scale factor for [REF].pdb[ref] output, 0 is autoscale" }
650 FILE *outx = nullptr, *outv = nullptr, *outf = nullptr, *outb = nullptr, *outt = nullptr;
651 FILE *outekt = nullptr, *outekr = nullptr;
652 t_topology top;
653 int ePBC;
654 real *mass, time;
655 const char *indexfn;
656 t_trxframe fr;
657 int flags, nvhisto = 0, *vhisto = nullptr;
658 rvec *xtop, *xp = nullptr;
659 rvec *sumx = nullptr, *sumv = nullptr, *sumf = nullptr;
660 matrix topbox;
661 t_trxstatus *status;
662 t_trxstatus *status_out = nullptr;
663 gmx_rmpbc_t gpbc = nullptr;
664 int i, j;
665 int nr_xfr, nr_vfr, nr_ffr;
666 char **grpname;
667 int *isize0, *isize;
668 int **index0, **index;
669 int *atndx;
670 t_block *mols;
671 gmx_bool bTop, bOX, bOXT, bOV, bOF, bOB, bOT, bEKT, bEKR, bCV, bCF;
672 gmx_bool bDim[4], bDum[4], bVD;
673 char sffmt[STRLEN];
674 const char *box_leg[6] = { "XX", "YY", "ZZ", "YX", "ZX", "ZY" };
675 gmx_output_env_t *oenv;
677 t_filenm fnm[] = {
678 { efTRX, "-f", nullptr, ffREAD },
679 { efTPS, nullptr, nullptr, ffREAD },
680 { efNDX, nullptr, nullptr, ffOPTRD },
681 { efXVG, "-ox", "coord", ffOPTWR },
682 { efTRX, "-oxt", "coord", ffOPTWR },
683 { efXVG, "-ov", "veloc", ffOPTWR },
684 { efXVG, "-of", "force", ffOPTWR },
685 { efXVG, "-ob", "box", ffOPTWR },
686 { efXVG, "-ot", "temp", ffOPTWR },
687 { efXVG, "-ekt", "ektrans", ffOPTWR },
688 { efXVG, "-ekr", "ekrot", ffOPTWR },
689 { efXVG, "-vd", "veldist", ffOPTWR },
690 { efPDB, "-cv", "veloc", ffOPTWR },
691 { efPDB, "-cf", "force", ffOPTWR },
692 { efXVG, "-av", "all_veloc", ffOPTWR },
693 { efXVG, "-af", "all_force", ffOPTWR }
695 #define NFILE asize(fnm)
697 if (!parse_common_args(&argc, argv,
698 PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW,
699 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
701 return 0;
704 if (bMol)
706 fprintf(stderr, "Interpreting indexfile entries as molecules.\n"
707 "Using center of mass.\n");
710 bOX = opt2bSet("-ox", NFILE, fnm);
711 bOXT = opt2bSet("-oxt", NFILE, fnm);
712 bOV = opt2bSet("-ov", NFILE, fnm);
713 bOF = opt2bSet("-of", NFILE, fnm);
714 bOB = opt2bSet("-ob", NFILE, fnm);
715 bOT = opt2bSet("-ot", NFILE, fnm);
716 bEKT = opt2bSet("-ekt", NFILE, fnm);
717 bEKR = opt2bSet("-ekr", NFILE, fnm);
718 bCV = opt2bSet("-cv", NFILE, fnm) || opt2bSet("-av", NFILE, fnm);
719 bCF = opt2bSet("-cf", NFILE, fnm) || opt2bSet("-af", NFILE, fnm);
720 bVD = opt2bSet("-vd", NFILE, fnm) || opt2parg_bSet("-bin", asize(pa), pa);
721 if (bMol || bOT || bEKT || bEKR)
723 bCom = TRUE;
726 bDim[XX] = bX;
727 bDim[YY] = bY;
728 bDim[ZZ] = bZ;
729 bDim[DIM] = bNorm;
731 if (bFP)
733 sprintf(sffmt, "\t%s", gmx_real_fullprecision_pfmt);
735 else
737 sprintf(sffmt, "\t%%g");
739 std::string sffmt6 = gmx::formatString("%s%s%s%s%s%s", sffmt, sffmt, sffmt, sffmt, sffmt, sffmt);
741 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC,
742 &xtop, nullptr, topbox,
743 bCom && (bOX || bOXT || bOV || bOT || bEKT || bEKR));
744 sfree(xtop);
745 if ((bMol || bCV || bCF) && !bTop)
747 gmx_fatal(FARGS, "Need a run input file for option -mol, -cv or -cf");
750 if (bMol)
752 indexfn = ftp2fn(efNDX, NFILE, fnm);
754 else
756 indexfn = ftp2fn_null(efNDX, NFILE, fnm);
759 if (!(bCom && !bMol))
761 ngroups = 1;
763 snew(grpname, ngroups);
764 snew(isize0, ngroups);
765 snew(index0, ngroups);
766 get_index(&(top.atoms), indexfn, ngroups, isize0, index0, grpname);
768 if (bMol)
770 mols = &(top.mols);
771 atndx = mols->index;
772 ngroups = isize0[0];
773 snew(isize, ngroups);
774 snew(index, ngroups);
775 for (i = 0; i < ngroups; i++)
777 if (index0[0][i] < 0 || index0[0][i] >= mols->nr)
779 gmx_fatal(FARGS, "Molecule index (%d) is out of range (%d-%d)",
780 index0[0][i]+1, 1, mols->nr);
782 isize[i] = atndx[index0[0][i]+1] - atndx[index0[0][i]];
783 snew(index[i], isize[i]);
784 for (j = 0; j < isize[i]; j++)
786 index[i][j] = atndx[index0[0][i]] + j;
790 else
792 isize = isize0;
793 index = index0;
795 if (bCom)
797 snew(mass, top.atoms.nr);
798 for (i = 0; i < top.atoms.nr; i++)
800 mass[i] = top.atoms.atom[i].m;
803 else
805 mass = nullptr;
808 flags = 0;
809 std::string label(output_env_get_xvgr_tlabel(oenv));
810 if (bOX)
812 flags = flags | TRX_READ_X;
813 outx = xvgropen(opt2fn("-ox", NFILE, fnm),
814 bCom ? "Center of mass" : "Coordinate",
815 label, "Coordinate (nm)", oenv);
816 make_legend(outx, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
818 if (bOXT)
820 flags = flags | TRX_READ_X;
821 status_out = open_trx(opt2fn("-oxt", NFILE, fnm), "w");
823 if (bOV)
825 flags = flags | TRX_READ_V;
826 outv = xvgropen(opt2fn("-ov", NFILE, fnm),
827 bCom ? "Center of mass velocity" : "Velocity",
828 label, "Velocity (nm/ps)", oenv);
829 make_legend(outv, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
831 if (bOF)
833 flags = flags | TRX_READ_F;
834 outf = xvgropen(opt2fn("-of", NFILE, fnm), "Force",
835 label, "Force (kJ mol\\S-1\\N nm\\S-1\\N)",
836 oenv);
837 make_legend(outf, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
839 if (bOB)
841 outb = xvgropen(opt2fn("-ob", NFILE, fnm), "Box vector elements",
842 label, "(nm)", oenv);
844 xvgr_legend(outb, 6, box_leg, oenv);
846 if (bOT)
848 bDum[XX] = FALSE;
849 bDum[YY] = FALSE;
850 bDum[ZZ] = FALSE;
851 bDum[DIM] = TRUE;
852 flags = flags | TRX_READ_V;
853 outt = xvgropen(opt2fn("-ot", NFILE, fnm), "Temperature",
854 label, "(K)", oenv);
855 make_legend(outt, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
857 if (bEKT)
859 bDum[XX] = FALSE;
860 bDum[YY] = FALSE;
861 bDum[ZZ] = FALSE;
862 bDum[DIM] = TRUE;
863 flags = flags | TRX_READ_V;
864 outekt = xvgropen(opt2fn("-ekt", NFILE, fnm), "Center of mass translation",
865 label, "Energy (kJ mol\\S-1\\N)", oenv);
866 make_legend(outekt, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
868 if (bEKR)
870 bDum[XX] = FALSE;
871 bDum[YY] = FALSE;
872 bDum[ZZ] = FALSE;
873 bDum[DIM] = TRUE;
874 flags = flags | TRX_READ_X | TRX_READ_V;
875 outekr = xvgropen(opt2fn("-ekr", NFILE, fnm), "Center of mass rotation",
876 label, "Energy (kJ mol\\S-1\\N)", oenv);
877 make_legend(outekr, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
879 if (bVD)
881 flags = flags | TRX_READ_V;
883 if (bCV)
885 flags = flags | TRX_READ_X | TRX_READ_V;
887 if (bCF)
889 flags = flags | TRX_READ_X | TRX_READ_F;
891 if ((flags == 0) && !bOB)
893 fprintf(stderr, "Please select one or more output file options\n");
894 exit(0);
897 read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
900 if ((bOV || bOF) && fn2ftp(ftp2fn(efTRX, NFILE, fnm)) == efXTC)
902 gmx_fatal(FARGS, "Cannot extract velocities or forces since your input XTC file does not contain them.");
905 if (bCV || bCF)
907 snew(sumx, fr.natoms);
909 if (bCV)
911 snew(sumv, fr.natoms);
913 if (bCF)
915 snew(sumf, fr.natoms);
917 nr_xfr = 0;
918 nr_vfr = 0;
919 nr_ffr = 0;
921 if (bCom && bPBC)
923 gpbc = gmx_rmpbc_init(&top.idef, ePBC, fr.natoms);
928 time = output_env_conv_time(oenv, fr.time);
930 if (fr.bX && bNoJump && fr.bBox)
932 if (xp)
934 remove_jump(fr.box, fr.natoms, xp, fr.x);
936 else
938 snew(xp, fr.natoms);
940 for (i = 0; i < fr.natoms; i++)
942 copy_rvec(fr.x[i], xp[i]);
946 if (fr.bX && bCom && bPBC)
948 gmx_rmpbc_trxfr(gpbc, &fr);
951 if (bVD && fr.bV)
953 update_histo(isize[0], index[0], fr.v, &nvhisto, &vhisto, binwidth);
956 if (bOX && fr.bX)
958 print_data(outx, time, fr.x, mass, bCom, ngroups, isize, index, bDim, sffmt);
960 if (bOXT && fr.bX)
962 t_trxframe frout = fr;
963 if (!frout.bAtoms)
965 frout.atoms = &top.atoms;
966 frout.bAtoms = TRUE;
968 frout.bV = FALSE;
969 frout.bF = FALSE;
970 write_trx_x(status_out, &frout, mass, bCom, ngroups, isize, index);
972 if (bOV && fr.bV)
974 print_data(outv, time, fr.v, mass, bCom, ngroups, isize, index, bDim, sffmt);
976 if (bOF && fr.bF)
978 print_data(outf, time, fr.f, nullptr, bCom, ngroups, isize, index, bDim, sffmt);
980 if (bOB && fr.bBox)
982 fprintf(outb, "\t%g", fr.time);
983 fprintf(outb, sffmt6.c_str(),
984 fr.box[XX][XX], fr.box[YY][YY], fr.box[ZZ][ZZ],
985 fr.box[YY][XX], fr.box[ZZ][XX], fr.box[ZZ][YY]);
986 fprintf(outb, "\n");
988 if (bOT && fr.bV)
990 fprintf(outt, " %g", time);
991 for (i = 0; i < ngroups; i++)
993 fprintf(outt, sffmt, temp(fr.v, mass, isize[i], index[i]));
995 fprintf(outt, "\n");
997 if (bEKT && fr.bV)
999 fprintf(outekt, " %g", time);
1000 for (i = 0; i < ngroups; i++)
1002 fprintf(outekt, sffmt, ektrans(fr.v, mass, isize[i], index[i]));
1004 fprintf(outekt, "\n");
1006 if (bEKR && fr.bX && fr.bV)
1008 fprintf(outekr, " %g", time);
1009 for (i = 0; i < ngroups; i++)
1011 fprintf(outekr, sffmt, ekrot(fr.x, fr.v, mass, isize[i], index[i]));
1013 fprintf(outekr, "\n");
1015 if ((bCV || bCF) && fr.bX &&
1016 (ctime < 0 || (fr.time >= ctime*0.999999 &&
1017 fr.time <= ctime*1.000001)))
1019 for (i = 0; i < fr.natoms; i++)
1021 rvec_inc(sumx[i], fr.x[i]);
1023 nr_xfr++;
1025 if (bCV && fr.bV)
1027 for (i = 0; i < fr.natoms; i++)
1029 rvec_inc(sumv[i], fr.v[i]);
1031 nr_vfr++;
1033 if (bCF && fr.bF)
1035 for (i = 0; i < fr.natoms; i++)
1037 rvec_inc(sumf[i], fr.f[i]);
1039 nr_ffr++;
1043 while (read_next_frame(oenv, status, &fr));
1045 if (gpbc != nullptr)
1047 gmx_rmpbc_done(gpbc);
1050 /* clean up a bit */
1051 close_trx(status);
1053 if (bOX)
1055 xvgrclose(outx);
1057 if (bOXT)
1059 close_trx(status_out);
1061 if (bOV)
1063 xvgrclose(outv);
1065 if (bOF)
1067 xvgrclose(outf);
1069 if (bOB)
1071 xvgrclose(outb);
1073 if (bOT)
1075 xvgrclose(outt);
1077 if (bEKT)
1079 xvgrclose(outekt);
1081 if (bEKR)
1083 xvgrclose(outekr);
1086 if (bVD)
1088 print_histo(opt2fn("-vd", NFILE, fnm), nvhisto, vhisto, binwidth, oenv);
1091 if (bCV || bCF)
1093 if (nr_xfr > 1)
1095 if (ePBC != epbcNONE && !bNoJump)
1097 fprintf(stderr, "\nWARNING: More than one frame was used for option -cv or -cf\n"
1098 "If atoms jump across the box you should use the -nojump or -ctime option\n\n");
1100 for (i = 0; i < isize[0]; i++)
1102 svmul(1.0/nr_xfr, sumx[index[0][i]], sumx[index[0][i]]);
1105 else if (nr_xfr == 0)
1107 fprintf(stderr, "\nWARNING: No coordinate frames found for option -cv or -cf\n\n");
1110 if (bCV)
1112 write_pdb_bfac(opt2fn("-cv", NFILE, fnm),
1113 opt2fn("-av", NFILE, fnm), "average velocity", &(top.atoms),
1114 ePBC, topbox, isize[0], index[0], nr_xfr, sumx,
1115 nr_vfr, sumv, bDim, scale, oenv);
1117 if (bCF)
1119 write_pdb_bfac(opt2fn("-cf", NFILE, fnm),
1120 opt2fn("-af", NFILE, fnm), "average force", &(top.atoms),
1121 ePBC, topbox, isize[0], index[0], nr_xfr, sumx,
1122 nr_ffr, sumf, bDim, scale, oenv);
1125 /* view it */
1126 view_all(oenv, NFILE, fnm);
1128 done_top(&top);
1129 // Free index and isize only if they are distinct from index0 and isize0
1130 if (bMol)
1132 for (int i = 0; i < ngroups; i++)
1134 sfree(index[i]);
1136 sfree(index);
1137 sfree(isize);
1139 for (int i = 0; i < ngroups; i++)
1141 sfree(index0[i]);
1142 sfree(grpname[i]);
1144 sfree(index0);
1145 sfree(isize0);
1146 sfree(grpname);
1147 done_frame(&fr);
1148 output_env_done(oenv);
1150 return 0;