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.
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/fileio/gmxfio.h"
46 #include "gromacs/fileio/pdbio.h"
47 #include "gromacs/fileio/tngio_for_tools.h"
48 #include "gromacs/fileio/tpxio.h"
49 #include "gromacs/fileio/trnio.h"
50 #include "gromacs/fileio/trxio.h"
51 #include "gromacs/fileio/xtcio.h"
52 #include "gromacs/fileio/xvgr.h"
53 #include "gromacs/gmxana/gmx_ana.h"
54 #include "gromacs/legacyheaders/copyrite.h"
55 #include "gromacs/legacyheaders/macros.h"
56 #include "gromacs/legacyheaders/names.h"
57 #include "gromacs/legacyheaders/typedefs.h"
58 #include "gromacs/legacyheaders/viewit.h"
59 #include "gromacs/math/do_fit.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/pbcutil/pbc.h"
62 #include "gromacs/pbcutil/rmpbc.h"
63 #include "gromacs/topology/index.h"
64 #include "gromacs/topology/topology.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/futil.h"
67 #include "gromacs/utility/smalloc.h"
70 euSel, euRect, euTric, euCompact, euNR
74 static void calc_pbc_cluster(int ecenter, int nrefat, t_topology *top, int ePBC,
75 rvec x[], atom_id index[], matrix box)
77 int m, i, j, j0, j1, jj, ai, aj;
80 rvec dx, xtest, box_center;
81 int nmol, imol_center;
83 gmx_bool *bMol, *bTmp;
84 rvec *m_com, *m_shift;
92 calc_box_center(ecenter, box, box_center);
94 /* Initiate the pbc structure */
95 memset(&pbc, 0, sizeof(pbc));
96 set_pbc(&pbc, ePBC, box);
98 /* Convert atom index to molecular */
100 molind = top->mols.index;
106 snew(bTmp, top->atoms.nr);
108 for (i = 0; (i < nrefat); i++)
110 /* Mark all molecules in the index */
113 /* Binary search assuming the molecules are sorted */
118 if (ai < molind[j0+1])
122 else if (ai >= molind[j1])
129 if (ai < molind[jj+1])
141 /* Double check whether all atoms in all molecules that are marked are part
142 * of the cluster. Simultaneously compute the center of geometry.
144 min_dist2 = 10*sqr(trace(box));
147 for (i = 0; i < nmol; i++)
149 for (j = molind[i]; j < molind[i+1]; j++)
151 if (bMol[i] && !bTmp[j])
153 gmx_fatal(FARGS, "Molecule %d marked for clustering but not atom %d in it - check your index!", i+1, j+1);
155 else if (!bMol[i] && bTmp[j])
157 gmx_fatal(FARGS, "Atom %d marked for clustering but not molecule %d - this is an internal error...", j+1, i+1);
161 /* Make molecule whole, move 2nd and higher atom to same periodicity as 1st atom in molecule */
164 pbc_dx(&pbc, x[j], x[j-1], dx);
165 rvec_add(x[j-1], dx, x[j]);
167 /* Compute center of geometry of molecule - m_com[i] was zeroed when we did snew() on it! */
168 rvec_inc(m_com[i], x[j]);
173 /* Normalize center of geometry */
174 fac = 1.0/(molind[i+1]-molind[i]);
175 for (m = 0; (m < DIM); m++)
179 /* Determine which molecule is closest to the center of the box */
180 pbc_dx(&pbc, box_center, m_com[i], dx);
181 tmp_r2 = iprod(dx, dx);
183 if (tmp_r2 < min_dist2)
188 cluster[ncluster++] = i;
195 fprintf(stderr, "No molecules selected in the cluster\n");
198 else if (imol_center == -1)
200 fprintf(stderr, "No central molecules could be found\n");
205 added[nadded++] = imol_center;
206 bMol[imol_center] = FALSE;
208 while (nadded < ncluster)
210 /* Find min distance between cluster molecules and those remaining to be added */
211 min_dist2 = 10*sqr(trace(box));
214 /* Loop over added mols */
215 for (i = 0; i < nadded; i++)
218 /* Loop over all mols */
219 for (j = 0; j < ncluster; j++)
222 /* check those remaining to be added */
225 pbc_dx(&pbc, m_com[aj], m_com[ai], dx);
226 tmp_r2 = iprod(dx, dx);
227 if (tmp_r2 < min_dist2)
237 /* Add the best molecule */
238 added[nadded++] = jmin;
240 /* Calculate the shift from the ai molecule */
241 pbc_dx(&pbc, m_com[jmin], m_com[imin], dx);
242 rvec_add(m_com[imin], dx, xtest);
243 rvec_sub(xtest, m_com[jmin], m_shift[jmin]);
244 rvec_inc(m_com[jmin], m_shift[jmin]);
246 for (j = molind[jmin]; j < molind[jmin+1]; j++)
248 rvec_inc(x[j], m_shift[jmin]);
250 fprintf(stdout, "\rClustering iteration %d of %d...", nadded, ncluster);
260 fprintf(stdout, "\n");
263 static void put_molecule_com_in_box(int unitcell_enum, int ecenter,
265 int natoms, t_atom atom[],
266 int ePBC, matrix box, rvec x[])
270 rvec com, new_com, shift, dx, box_center;
275 calc_box_center(ecenter, box, box_center);
276 set_pbc(&pbc, ePBC, box);
279 gmx_fatal(FARGS, "There are no molecule descriptions. I need a .tpr file for this pbc option.");
281 for (i = 0; (i < mols->nr); i++)
286 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
289 for (d = 0; d < DIM; d++)
295 /* calculate final COM */
296 svmul(1.0/mtot, com, com);
298 /* check if COM is outside box */
299 copy_rvec(com, new_com);
300 switch (unitcell_enum)
303 put_atoms_in_box(ePBC, box, 1, &new_com);
306 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
309 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
312 rvec_sub(new_com, com, shift);
313 if (norm2(shift) > 0)
317 fprintf(debug, "\nShifting position of molecule %d "
318 "by %8.3f %8.3f %8.3f\n", i+1,
319 shift[XX], shift[YY], shift[ZZ]);
321 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
323 rvec_inc(x[j], shift);
329 static void put_residue_com_in_box(int unitcell_enum, int ecenter,
330 int natoms, t_atom atom[],
331 int ePBC, matrix box, rvec x[])
333 atom_id i, j, res_start, res_end, res_nat;
337 rvec box_center, com, new_com, shift;
339 calc_box_center(ecenter, box, box_center);
345 for (i = 0; i < natoms+1; i++)
347 if (i == natoms || (presnr != atom[i].resind && presnr != NOTSET))
349 /* calculate final COM */
351 res_nat = res_end - res_start;
352 svmul(1.0/mtot, com, com);
354 /* check if COM is outside box */
355 copy_rvec(com, new_com);
356 switch (unitcell_enum)
359 put_atoms_in_box(ePBC, box, 1, &new_com);
362 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
365 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
368 rvec_sub(new_com, com, shift);
373 fprintf(debug, "\nShifting position of residue %d (atoms %d-%d) "
374 "by %g,%g,%g\n", atom[res_start].resind+1,
375 res_start+1, res_end+1, shift[XX], shift[YY], shift[ZZ]);
377 for (j = res_start; j < res_end; j++)
379 rvec_inc(x[j], shift);
385 /* remember start of new residue */
392 for (d = 0; d < DIM; d++)
398 presnr = atom[i].resind;
403 static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, atom_id ci[])
406 rvec cmin, cmax, box_center, dx;
410 copy_rvec(x[ci[0]], cmin);
411 copy_rvec(x[ci[0]], cmax);
412 for (i = 0; i < nc; i++)
415 for (m = 0; m < DIM; m++)
417 if (x[ai][m] < cmin[m])
421 else if (x[ai][m] > cmax[m])
427 calc_box_center(ecenter, box, box_center);
428 for (m = 0; m < DIM; m++)
430 dx[m] = box_center[m]-(cmin[m]+cmax[m])*0.5;
433 for (i = 0; i < n; i++)
440 static void mk_filenm(char *base, const char *ext, int ndigit, int file_nr,
446 strcpy(out_file, base);
457 strncat(out_file, "00000000000", ndigit-nd);
459 sprintf(nbuf, "%d.", file_nr);
460 strcat(out_file, nbuf);
461 strcat(out_file, ext);
464 void check_trn(const char *fn)
466 if (fn2ftp(fn) != efTRR)
468 gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
472 void do_trunc(const char *fn, real t0)
485 gmx_fatal(FARGS, "You forgot to set the truncation time");
488 /* Check whether this is a .trr file */
491 in = open_trn(fn, "r");
492 fp = gmx_fio_getfp(in);
495 fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
501 fpos = gmx_fio_ftell(in);
503 while (!bStop && fread_trnheader(in, &sh, &bOK))
505 fread_htrn(in, &sh, NULL, NULL, NULL, NULL);
506 fpos = gmx_ftell(fp);
510 gmx_fseek(fp, fpos, SEEK_SET);
516 fprintf(stderr, "Do you REALLY want to truncate this trajectory (%s) at:\n"
517 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
518 fn, j, t, (long int)fpos);
519 if (1 != scanf("%s", yesno))
521 gmx_fatal(FARGS, "Error reading user input");
523 if (strcmp(yesno, "YES") == 0)
525 fprintf(stderr, "Once again, I'm gonna DO this...\n");
527 if (0 != gmx_truncate(fn, fpos))
529 gmx_fatal(FARGS, "Error truncating file %s", fn);
534 fprintf(stderr, "Ok, I'll forget about it\n");
539 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
545 /*! \brief Read a full molecular topology if useful and available.
547 * If the input trajectory file is not in TNG format, and the output
548 * file is in TNG format, then we want to try to read a full topology
549 * (if available), so that we can write molecule information to the
550 * output file. The full topology provides better molecule information
551 * than is available from the normal t_topology data used by GROMACS
554 * Also, the t_topology is only read under (different) particular
555 * conditions. If both apply, then a .tpr file might be read
556 * twice. Trying to fix this redundancy while trjconv is still an
557 * all-purpose tool does not seem worthwhile.
559 * Because of the way gmx_prepare_tng_writing is implemented, the case
560 * where the input TNG file has no molecule information will never
561 * lead to an output TNG file having molecule information. Since
562 * molecule information will generally be present if the input TNG
563 * file was written by a GROMACS tool, this seems like reasonable
565 static gmx_mtop_t *read_mtop_for_tng(const char *tps_file,
566 const char *input_file,
567 const char *output_file)
569 gmx_mtop_t *mtop = NULL;
571 if (fn2bTPX(tps_file) &&
572 efTNG != fn2ftp(input_file) &&
573 efTNG == fn2ftp(output_file))
575 int temp_natoms = -1;
577 read_tpx(tps_file, NULL, NULL, &temp_natoms,
578 NULL, NULL, NULL, mtop);
584 int gmx_trjconv(int argc, char *argv[])
586 const char *desc[] = {
587 "[THISMODULE] can convert trajectory files in many ways:[BR]",
588 "* from one format to another[BR]",
589 "* select a subset of atoms[BR]",
590 "* change the periodicity representation[BR]",
591 "* keep multimeric molecules together[BR]",
592 "* center atoms in the box[BR]",
593 "* fit atoms to reference structure[BR]",
594 "* reduce the number of frames[BR]",
595 "* change the timestamps of the frames ",
596 "([TT]-t0[tt] and [TT]-timestep[tt])[BR]",
597 "* cut the trajectory in small subtrajectories according",
598 "to information in an index file. This allows subsequent analysis of",
599 "the subtrajectories that could, for example, be the result of a",
600 "cluster analysis. Use option [TT]-sub[tt].",
601 "This assumes that the entries in the index file are frame numbers and",
602 "dumps each group in the index file to a separate trajectory file.[BR]",
603 "* select frames within a certain range of a quantity given",
604 "in an [TT].xvg[tt] file.[PAR]",
606 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
609 "The following formats are supported for input and output:",
610 "[TT].xtc[tt], [TT].trr[tt], [TT].gro[tt], [TT].g96[tt]",
612 "The file formats are detected from the file extension.",
613 "The precision of [TT].xtc[tt] and [TT].gro[tt] output is taken from the",
614 "input file for [TT].xtc[tt], [TT].gro[tt] and [TT].pdb[tt],",
615 "and from the [TT]-ndec[tt] option for other input formats. The precision",
616 "is always taken from [TT]-ndec[tt], when this option is set.",
617 "All other formats have fixed precision. [TT].trr[tt]",
618 "output can be single or double precision, depending on the precision",
619 "of the [THISMODULE] binary.",
620 "Note that velocities are only supported in",
621 "[TT].trr[tt], [TT].gro[tt] and [TT].g96[tt] files.[PAR]",
623 "Option [TT]-sep[tt] can be used to write every frame to a separate",
624 "[TT].gro, .g96[tt] or [TT].pdb[tt] file. By default, all frames all written to one file.",
625 "[TT].pdb[tt] files with all frames concatenated can be viewed with",
626 "[TT]rasmol -nmrpdb[tt].[PAR]",
628 "It is possible to select part of your trajectory and write it out",
629 "to a new trajectory file in order to save disk space, e.g. for leaving",
630 "out the water from a trajectory of a protein in water.",
631 "[BB]ALWAYS[bb] put the original trajectory on tape!",
632 "We recommend to use the portable [TT].xtc[tt] format for your analysis",
633 "to save disk space and to have portable files.[PAR]",
635 "There are two options for fitting the trajectory to a reference",
636 "either for essential dynamics analysis, etc.",
637 "The first option is just plain fitting to a reference structure",
638 "in the structure file. The second option is a progressive fit",
639 "in which the first timeframe is fitted to the reference structure ",
640 "in the structure file to obtain and each subsequent timeframe is ",
641 "fitted to the previously fitted structure. This way a continuous",
642 "trajectory is generated, which might not be the case when using the",
643 "regular fit method, e.g. when your protein undergoes large",
644 "conformational transitions.[PAR]",
646 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
648 "[TT]* mol[tt] puts the center of mass of molecules in the box,",
649 "and requires a run input file to be supplied with [TT]-s[tt].[BR]",
650 "[TT]* res[tt] puts the center of mass of residues in the box.[BR]",
651 "[TT]* atom[tt] puts all the atoms in the box.[BR]",
652 "[TT]* nojump[tt] checks if atoms jump across the box and then puts",
653 "them back. This has the effect that all molecules",
654 "will remain whole (provided they were whole in the initial",
655 "conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
656 "molecules may diffuse out of the box. The starting configuration",
657 "for this procedure is taken from the structure file, if one is",
658 "supplied, otherwise it is the first frame.[BR]",
659 "[TT]* cluster[tt] clusters all the atoms in the selected index",
660 "such that they are all closest to the center of mass of the cluster,",
661 "which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
662 "results if you in fact have a cluster. Luckily that can be checked",
663 "afterwards using a trajectory viewer. Note also that if your molecules",
664 "are broken this will not work either.[BR]",
665 "The separate option [TT]-clustercenter[tt] can be used to specify an",
666 "approximate center for the cluster. This is useful e.g. if you have",
667 "two big vesicles, and you want to maintain their relative positions.[BR]",
668 "[TT]* whole[tt] only makes broken molecules whole.[PAR]",
670 "Option [TT]-ur[tt] sets the unit cell representation for options",
671 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
672 "All three options give different results for triclinic boxes and",
673 "identical results for rectangular boxes.",
674 "[TT]rect[tt] is the ordinary brick shape.",
675 "[TT]tric[tt] is the triclinic unit cell.",
676 "[TT]compact[tt] puts all atoms at the closest distance from the center",
677 "of the box. This can be useful for visualizing e.g. truncated octahedra",
678 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
679 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
680 "is set differently.[PAR]",
682 "Option [TT]-center[tt] centers the system in the box. The user can",
683 "select the group which is used to determine the geometrical center.",
684 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
685 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
686 "[TT]tric[tt]: half of the sum of the box vectors,",
687 "[TT]rect[tt]: half of the box diagonal,",
688 "[TT]zero[tt]: zero.",
689 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
690 "want all molecules in the box after the centering.[PAR]",
692 "Option [TT]-box[tt] sets the size of the new box. If you want to"
693 "modify only some of the dimensions, e.g. when reading from a trajectory,"
694 "you can use -1 for those dimensions that should stay the same"
696 "It is not always possible to use combinations of [TT]-pbc[tt],",
697 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
698 "you want in one call to [THISMODULE]. Consider using multiple",
699 "calls, and check out the GROMACS website for suggestions.[PAR]",
701 "With [TT]-dt[tt], it is possible to reduce the number of ",
702 "frames in the output. This option relies on the accuracy of the times",
703 "in your input trajectory, so if these are inaccurate use the",
704 "[TT]-timestep[tt] option to modify the time (this can be done",
705 "simultaneously). For making smooth movies, the program [gmx-filter]",
706 "can reduce the number of frames while using low-pass frequency",
707 "filtering, this reduces aliasing of high frequency motions.[PAR]",
709 "Using [TT]-trunc[tt] [THISMODULE] can truncate [TT].trr[tt] in place, i.e.",
710 "without copying the file. This is useful when a run has crashed",
711 "during disk I/O (i.e. full disk), or when two contiguous",
712 "trajectories must be concatenated without having double frames.[PAR]",
714 "Option [TT]-dump[tt] can be used to extract a frame at or near",
715 "one specific time from your trajectory.[PAR]",
717 "Option [TT]-drop[tt] reads an [TT].xvg[tt] file with times and values.",
718 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
719 "frames with a value below and above the value of the respective options",
720 "will not be written."
736 const char *pbc_opt[epNR + 1] =
738 NULL, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
743 const char *unitcell_opt[euNR+1] =
744 { NULL, "rect", "tric", "compact", NULL };
748 ecSel, ecTric, ecRect, ecZero, ecNR
750 const char *center_opt[ecNR+1] =
751 { NULL, "tric", "rect", "zero", NULL };
757 efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
759 const char *fit[efNR + 1] =
761 NULL, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
765 static gmx_bool bSeparate = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
766 static gmx_bool bCenter = FALSE;
767 static int skip_nr = 1, ndec = 3, nzero = 0;
768 static real tzero = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
769 static rvec newbox = {0, 0, 0}, shift = {0, 0, 0}, trans = {0, 0, 0};
770 static char *exec_command = NULL;
771 static real dropunder = 0, dropover = 0;
772 static gmx_bool bRound = FALSE;
777 { "-skip", FALSE, etINT,
778 { &skip_nr }, "Only write every nr-th frame" },
779 { "-dt", FALSE, etTIME,
781 "Only write frame when t MOD dt = first time (%t)" },
782 { "-round", FALSE, etBOOL,
783 { &bRound }, "Round measurements to nearest picosecond"},
784 { "-dump", FALSE, etTIME,
785 { &tdump }, "Dump frame nearest specified time (%t)" },
786 { "-t0", FALSE, etTIME,
788 "Starting time (%t) (default: don't change)" },
789 { "-timestep", FALSE, etTIME,
791 "Change time step between input frames (%t)" },
792 { "-pbc", FALSE, etENUM,
794 "PBC treatment (see help text for full description)" },
795 { "-ur", FALSE, etENUM,
796 { unitcell_opt }, "Unit-cell representation" },
797 { "-center", FALSE, etBOOL,
798 { &bCenter }, "Center atoms in box" },
799 { "-boxcenter", FALSE, etENUM,
800 { center_opt }, "Center for -pbc and -center" },
801 { "-box", FALSE, etRVEC,
803 "Size for new cubic box (default: read from input)" },
804 { "-trans", FALSE, etRVEC,
806 "All coordinates will be translated by trans. This "
807 "can advantageously be combined with -pbc mol -ur "
809 { "-shift", FALSE, etRVEC,
811 "All coordinates will be shifted by framenr*shift" },
812 { "-fit", FALSE, etENUM,
814 "Fit molecule to ref structure in the structure file" },
815 { "-ndec", FALSE, etINT,
817 "Precision for .xtc and .gro writing in number of "
819 { "-vel", FALSE, etBOOL,
820 { &bVels }, "Read and write velocities if possible" },
821 { "-force", FALSE, etBOOL,
822 { &bForce }, "Read and write forces if possible" },
823 { "-trunc", FALSE, etTIME,
825 "Truncate input trajectory file after this time (%t)" },
826 { "-exec", FALSE, etSTR,
828 "Execute command for every output frame with the "
829 "frame number as argument" },
830 { "-split", FALSE, etTIME,
832 "Start writing new file when t MOD split = first "
834 { "-sep", FALSE, etBOOL,
836 "Write each frame to a separate .gro, .g96 or .pdb "
838 { "-nzero", FALSE, etINT,
840 "If the -sep flag is set, use these many digits "
841 "for the file numbers and prepend zeros as needed" },
842 { "-dropunder", FALSE, etREAL,
843 { &dropunder }, "Drop all frames below this value" },
844 { "-dropover", FALSE, etREAL,
845 { &dropover }, "Drop all frames above this value" },
846 { "-conect", FALSE, etBOOL,
848 "Add conect records when writing [TT].pdb[tt] files. Useful "
849 "for visualization of non-standard molecules, e.g. "
850 "coarse grained ones" }
852 #define NPA asize(pa)
855 t_trxstatus *trxout = NULL;
857 int ftp, ftpin = 0, file_nr;
858 t_trxframe fr, frout;
860 rvec *xmem = NULL, *vmem = NULL, *fmem = NULL;
861 rvec *xp = NULL, x_shift, hbox, box_center, dx;
862 real xtcpr, lambda, *w_rls = NULL;
863 int m, i, d, frame, outframe, natoms, nout, ncent, nre, newstep = 0, model_nr;
866 gmx_mtop_t *mtop = NULL;
867 gmx_conect gc = NULL;
869 t_atoms *atoms = NULL, useatoms;
871 atom_id *index, *cindex;
875 int ifit, irms, my_clust = -1;
876 atom_id *ind_fit, *ind_rms;
877 char *gn_fit, *gn_rms;
878 t_cluster_ndx *clust = NULL;
879 t_trxstatus **clust_status = NULL;
880 int *clust_status_id = NULL;
882 int *nfwritten = NULL;
883 int ndrop = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
885 real tshift = 0, t0 = -1, dt = 0.001, prec;
886 gmx_bool bFit, bFitXY, bPFit, bReset;
888 gmx_rmpbc_t gpbc = NULL;
889 gmx_bool bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
890 gmx_bool bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
891 gmx_bool bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetPrec, bNeedPrec;
892 gmx_bool bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
893 gmx_bool bSubTraj = FALSE, bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
894 gmx_bool bWriteFrame, bSplitHere;
895 const char *top_file, *in_file, *out_file = NULL;
896 char out_file2[256], *charpt;
897 char *outf_base = NULL;
898 const char *outf_ext = NULL;
899 char top_title[256], title[256], command[256], filemode[5];
901 gmx_bool bWarnCompact = FALSE;
906 { efTRX, "-f", NULL, ffREAD },
907 { efTRO, "-o", NULL, ffWRITE },
908 { efTPS, NULL, NULL, ffOPTRD },
909 { efNDX, NULL, NULL, ffOPTRD },
910 { efNDX, "-fr", "frames", ffOPTRD },
911 { efNDX, "-sub", "cluster", ffOPTRD },
912 { efXVG, "-drop", "drop", ffOPTRD }
914 #define NFILE asize(fnm)
916 if (!parse_common_args(&argc, argv,
917 PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW |
919 NFILE, fnm, NPA, pa, asize(desc), desc,
925 top_file = ftp2fn(efTPS, NFILE, fnm);
928 /* Check command line */
929 in_file = opt2fn("-f", NFILE, fnm);
933 do_trunc(in_file, ttrunc);
937 /* mark active cmdline options */
938 bSetBox = opt2parg_bSet("-box", NPA, pa);
939 bSetTime = opt2parg_bSet("-t0", NPA, pa);
940 bSetPrec = opt2parg_bSet("-ndec", NPA, pa);
941 bSetUR = opt2parg_bSet("-ur", NPA, pa);
942 bExec = opt2parg_bSet("-exec", NPA, pa);
943 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
944 bTDump = opt2parg_bSet("-dump", NPA, pa);
945 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
946 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
947 bTrans = opt2parg_bSet("-trans", NPA, pa);
948 bSplit = (split_t != 0);
950 /* parse enum options */
951 fit_enum = nenum(fit);
952 bFit = (fit_enum == efFit || fit_enum == efFitXY);
953 bFitXY = fit_enum == efFitXY;
954 bReset = (fit_enum == efReset || fit_enum == efResetXY);
955 bPFit = fit_enum == efPFit;
956 pbc_enum = nenum(pbc_opt);
957 bPBCWhole = pbc_enum == epWhole;
958 bPBCcomRes = pbc_enum == epComRes;
959 bPBCcomMol = pbc_enum == epComMol;
960 bPBCcomAtom = pbc_enum == epComAtom;
961 bNoJump = pbc_enum == epNojump;
962 bCluster = pbc_enum == epCluster;
963 bPBC = pbc_enum != epNone;
964 unitcell_enum = nenum(unitcell_opt);
965 ecenter = nenum(center_opt) - ecTric;
967 /* set and check option dependencies */
970 bFit = TRUE; /* for pfit, fit *must* be set */
974 bReset = TRUE; /* for fit, reset *must* be set */
979 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
981 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
985 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
988 "WARNING: Option for unitcell representation (-ur %s)\n"
989 " only has effect in combination with -pbc %s, %s or %s.\n"
990 " Ingoring unitcell representation.\n\n",
991 unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
997 gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
998 "Please do the PBC condition treatment first and then run trjconv in a second step\n"
999 "for the rotational fit.\n"
1000 "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1004 /* ndec is in nr of decimal places, prec is a multiplication factor: */
1006 for (i = 0; i < ndec; i++)
1011 bIndex = ftp2bSet(efNDX, NFILE, fnm);
1014 /* Determine output type */
1015 out_file = opt2fn("-o", NFILE, fnm);
1016 ftp = fn2ftp(out_file);
1017 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
1018 bNeedPrec = (ftp == efXTC || ftp == efGRO);
1021 /* check if velocities are possible in input and output files */
1022 ftpin = fn2ftp(in_file);
1023 bVels = (ftp == efTRR || ftp == efGRO ||
1024 ftp == efG96 || ftp == efTNG)
1025 && (ftpin == efTRR || ftpin == efGRO ||
1026 ftpin == efG96 || ftpin == efTNG || ftpin == efCPT);
1028 if (bSeparate || bSplit)
1030 outf_ext = strrchr(out_file, '.');
1031 if (outf_ext == NULL)
1033 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1035 outf_base = gmx_strdup(out_file);
1036 outf_base[outf_ext - out_file] = '\0';
1039 bSubTraj = opt2bSet("-sub", NFILE, fnm);
1042 if ((ftp != efXTC) && (ftp != efTRR))
1044 /* It seems likely that other trajectory file types
1045 * could work here. */
1046 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1049 clust = cluster_index(NULL, opt2fn("-sub", NFILE, fnm));
1051 /* Check for number of files disabled, as FOPEN_MAX is not the correct
1052 * number to check for. In my linux box it is only 16.
1054 if (0 && (clust->clust->nr > FOPEN_MAX-4))
1056 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1057 " trajectories.\ntry splitting the index file in %d parts.\n"
1059 clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1061 gmx_warning("The -sub option could require as many open output files as there are\n"
1062 "index groups in the file (%d). If you get I/O errors opening new files,\n"
1063 "try reducing the number of index groups in the file, and perhaps\n"
1064 "using trjconv -sub several times on different chunks of your index file.\n",
1067 snew(clust_status, clust->clust->nr);
1068 snew(clust_status_id, clust->clust->nr);
1069 snew(nfwritten, clust->clust->nr);
1070 for (i = 0; (i < clust->clust->nr); i++)
1072 clust_status[i] = NULL;
1073 clust_status_id[i] = -1;
1075 bSeparate = bSplit = FALSE;
1082 mtop = read_mtop_for_tng(top_file, in_file, out_file);
1084 /* Determine whether to read a topology */
1085 bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1086 bRmPBC || bReset || bPBCcomMol || bCluster ||
1087 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1089 /* Determine if when can read index groups */
1090 bIndex = (bIndex || bTPS);
1094 read_tps_conf(top_file, top_title, &top, &ePBC, &xp, NULL, top_box,
1095 bReset || bPBCcomRes);
1098 if (0 == top.mols.nr && (bCluster || bPBCcomMol))
1100 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1103 /* top_title is only used for gro and pdb,
1104 * the header in such a file is top_title t= ...
1105 * to prevent a double t=, remove it from top_title
1107 if ((charpt = strstr(top_title, " t= ")))
1114 gc = gmx_conect_generate(&top);
1118 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
1122 /* get frame number index */
1124 if (opt2bSet("-fr", NFILE, fnm))
1126 printf("Select groups of frame number indices:\n");
1127 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, (atom_id **)&frindex, &frname);
1130 for (i = 0; i < nrfri; i++)
1132 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1137 /* get index groups etc. */
1140 printf("Select group for %s fit\n",
1141 bFit ? "least squares" : "translational");
1142 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1143 1, &ifit, &ind_fit, &gn_fit);
1149 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1153 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1159 printf("Select group for clustering\n");
1160 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1161 1, &ifit, &ind_fit, &gn_fit);
1168 printf("Select group for centering\n");
1169 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1170 1, &ncent, &cindex, &grpnm);
1172 printf("Select group for output\n");
1173 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1174 1, &nout, &index, &grpnm);
1178 /* no index file, so read natoms from TRX */
1179 if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
1181 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1186 snew(index, natoms);
1187 for (i = 0; i < natoms; i++)
1201 snew(w_rls, atoms->nr);
1202 for (i = 0; (i < ifit); i++)
1204 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1207 /* Restore reference structure and set to origin,
1208 store original location (to put structure back) */
1211 gmx_rmpbc(gpbc, top.atoms.nr, top_box, xp);
1213 copy_rvec(xp[index[0]], x_shift);
1214 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, NULL, xp, w_rls);
1215 rvec_dec(x_shift, xp[index[0]]);
1219 clear_rvec(x_shift);
1222 if (bDropUnder || bDropOver)
1224 /* Read the .xvg file with the drop values */
1225 fprintf(stderr, "\nReading drop file ...");
1226 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1227 fprintf(stderr, " %d time points\n", ndrop);
1228 if (ndrop == 0 || ncol < 2)
1230 gmx_fatal(FARGS, "Found no data points in %s",
1231 opt2fn("-drop", NFILE, fnm));
1237 /* Make atoms struct for output in GRO or PDB files */
1238 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1240 /* get memory for stuff to go in .pdb file, and initialize
1241 * the pdbinfo structure part if the input has it.
1243 init_t_atoms(&useatoms, atoms->nr, (atoms->pdbinfo != NULL));
1244 sfree(useatoms.resinfo);
1245 useatoms.resinfo = atoms->resinfo;
1246 for (i = 0; (i < nout); i++)
1248 useatoms.atomname[i] = atoms->atomname[index[i]];
1249 useatoms.atom[i] = atoms->atom[index[i]];
1250 if (atoms->pdbinfo != NULL)
1252 useatoms.pdbinfo[i] = atoms->pdbinfo[index[i]];
1254 useatoms.nres = max(useatoms.nres, useatoms.atom[i].resind+1);
1258 /* select what to read */
1269 flags = flags | TRX_READ_V;
1273 flags = flags | TRX_READ_F;
1276 /* open trx file for reading */
1277 bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
1280 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1284 if (bSetPrec || !fr.bPrec)
1286 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1290 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1294 if (bHaveFirstFrame)
1296 set_trxframe_ePBC(&fr, ePBC);
1302 tshift = tzero-fr.time;
1312 /* check if index is meaningful */
1313 for (i = 0; i < nout; i++)
1315 if (index[i] >= natoms)
1318 "Index[%d] %d is larger than the number of atoms in the\n"
1319 "trajectory file (%d). There is a mismatch in the contents\n"
1320 "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1322 bCopy = bCopy || (i != index[i]);
1326 /* open output for writing */
1327 strcpy(filemode, "w");
1331 trjtools_gmx_prepare_tng_writing(out_file,
1344 if (!bSplit && !bSubTraj)
1346 trxout = open_trx(out_file, filemode);
1352 if (( !bSeparate && !bSplit ) && !bSubTraj)
1354 out = gmx_ffopen(out_file, filemode);
1358 gmx_incons("Illegal output file format");
1374 /* Start the big loop over frames */
1381 /* Main loop over frames */
1392 /*if (frame >= clust->clust->nra)
1393 gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1394 if (frame > clust->maxframe)
1400 my_clust = clust->inv_clust[frame];
1402 if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1403 (my_clust == NO_ATID))
1411 /* generate new box */
1413 for (m = 0; m < DIM; m++)
1417 fr.box[m][m] = newbox[m];
1424 for (i = 0; i < natoms; i++)
1426 rvec_inc(fr.x[i], trans);
1432 /* determine timestep */
1445 /* This is not very elegant, as one can not dump a frame after
1446 * a timestep with is more than twice as small as the first one. */
1447 bDumpFrame = (fr.time > tdump-0.5*dt) && (fr.time <= tdump+0.5*dt);
1454 /* determine if an atom jumped across the box and reset it if so */
1455 if (bNoJump && (bTPS || frame != 0))
1457 for (d = 0; d < DIM; d++)
1459 hbox[d] = 0.5*fr.box[d][d];
1461 for (i = 0; i < natoms; i++)
1465 rvec_dec(fr.x[i], x_shift);
1467 for (m = DIM-1; m >= 0; m--)
1471 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1473 for (d = 0; d <= m; d++)
1475 fr.x[i][d] += fr.box[m][d];
1478 while (fr.x[i][m]-xp[i][m] > hbox[m])
1480 for (d = 0; d <= m; d++)
1482 fr.x[i][d] -= fr.box[m][d];
1491 calc_pbc_cluster(ecenter, ifit, &top, ePBC, fr.x, ind_fit, fr.box);
1496 /* Now modify the coords according to the flags,
1497 for normal fit, this is only done for output frames */
1500 gmx_rmpbc_trxfr(gpbc, &fr);
1503 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1504 do_fit(natoms, w_rls, xp, fr.x);
1507 /* store this set of coordinates for future use */
1508 if (bPFit || bNoJump)
1514 for (i = 0; (i < natoms); i++)
1516 copy_rvec(fr.x[i], xp[i]);
1517 rvec_inc(fr.x[i], x_shift);
1523 /* see if we have a frame from the frame index group */
1524 for (i = 0; i < nrfri && !bDumpFrame; i++)
1526 bDumpFrame = frame == frindex[i];
1529 if (debug && bDumpFrame)
1531 fprintf(debug, "dumping %d\n", frame);
1535 ( ( !bTDump && !frindex && frame % skip_nr == 0 ) || bDumpFrame );
1537 if (bWriteFrame && (bDropUnder || bDropOver))
1539 while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1544 if (fabs(dropval[0][drop0] - fr.time)
1545 < fabs(dropval[0][drop1] - fr.time))
1553 if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1554 (bDropOver && dropval[1][dropuse] > dropover))
1556 bWriteFrame = FALSE;
1562 /* We should avoid modifying the input frame,
1563 * but since here we don't have the output frame yet,
1564 * we introduce a temporary output frame time variable.
1568 frout_time = fr.time;
1573 frout_time = tzero + frame*timestep;
1578 frout_time += tshift;
1583 fprintf(stderr, "\nDumping frame at t= %g %s\n",
1584 output_env_conv_time(oenv, frout_time), output_env_get_time_unit(oenv));
1587 /* check for writing at each delta_t */
1588 bDoIt = (delta_t == 0);
1593 bDoIt = bRmod(frout_time, tzero, delta_t);
1597 /* round() is not C89 compatible, so we do this: */
1598 bDoIt = bRmod(floor(frout_time+0.5), floor(tzero+0.5),
1599 floor(delta_t+0.5));
1603 if (bDoIt || bTDump)
1605 /* print sometimes */
1606 if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1608 fprintf(stderr, " -> frame %6d time %8.3f \r",
1609 outframe, output_env_conv_time(oenv, frout_time));
1614 /* Now modify the coords according to the flags,
1615 for PFit we did this already! */
1619 gmx_rmpbc_trxfr(gpbc, &fr);
1624 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1627 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1631 for (i = 0; i < natoms; i++)
1633 rvec_inc(fr.x[i], x_shift);
1640 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1646 switch (unitcell_enum)
1649 put_atoms_in_box(ePBC, fr.box, natoms, fr.x);
1652 put_atoms_in_triclinic_unitcell(ecenter, fr.box, natoms, fr.x);
1655 warn = put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1657 if (warn && !bWarnCompact)
1659 fprintf(stderr, "\n%s\n", warn);
1660 bWarnCompact = TRUE;
1667 put_residue_com_in_box(unitcell_enum, ecenter,
1668 natoms, atoms->atom, ePBC, fr.box, fr.x);
1672 put_molecule_com_in_box(unitcell_enum, ecenter,
1674 natoms, atoms->atom, ePBC, fr.box, fr.x);
1676 /* Copy the input trxframe struct to the output trxframe struct */
1678 frout.time = frout_time;
1679 frout.bV = (frout.bV && bVels);
1680 frout.bF = (frout.bF && bForce);
1681 frout.natoms = nout;
1682 if (bNeedPrec && (bSetPrec || !fr.bPrec))
1698 for (i = 0; i < nout; i++)
1700 copy_rvec(fr.x[index[i]], frout.x[i]);
1703 copy_rvec(fr.v[index[i]], frout.v[i]);
1707 copy_rvec(fr.f[index[i]], frout.f[i]);
1712 if (opt2parg_bSet("-shift", NPA, pa))
1714 for (i = 0; i < nout; i++)
1716 for (d = 0; d < DIM; d++)
1718 frout.x[i][d] += outframe*shift[d];
1725 bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1729 /* round() is not C89 compatible, so we do this: */
1730 bSplitHere = bSplit && bRmod(floor(frout.time+0.5),
1732 floor(split_t+0.5));
1734 if (bSeparate || bSplitHere)
1736 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1742 write_tng_frame(trxout, &frout);
1743 // TODO when trjconv behaves better: work how to read and write lambda
1753 trxout = open_trx(out_file2, filemode);
1760 if (clust_status_id[my_clust] == -1)
1762 sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1763 clust_status[my_clust] = open_trx(buf, "w");
1764 clust_status_id[my_clust] = 1;
1767 else if (clust_status_id[my_clust] == -2)
1769 gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1770 clust->grpname[my_clust], ntrxopen, frame,
1773 write_trxframe(clust_status[my_clust], &frout, gc);
1774 nfwritten[my_clust]++;
1775 if (nfwritten[my_clust] ==
1776 (clust->clust->index[my_clust+1]-
1777 clust->clust->index[my_clust]))
1779 close_trx(clust_status[my_clust]);
1780 clust_status[my_clust] = NULL;
1781 clust_status_id[my_clust] = -2;
1785 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1792 write_trxframe(trxout, &frout, gc);
1798 sprintf(title, "Generated by trjconv : %s t= %9.5f",
1799 top_title, frout.time);
1800 if (bSeparate || bSplitHere)
1802 out = gmx_ffopen(out_file2, "w");
1807 write_hconf_p(out, title, &useatoms, prec2ndec(frout.prec),
1808 frout.x, frout.bV ? frout.v : NULL, frout.box);
1811 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1812 sprintf(title, "%s t= %9.5f", top_title, frout.time);
1813 /* if reading from pdb, we want to keep the original
1814 model numbering else we write the output frame
1815 number plus one, because model 0 is not allowed in pdb */
1816 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1824 write_pdbfile(out, title, &useatoms, frout.x,
1825 frout.ePBC, frout.box, ' ', model_nr, gc, TRUE);
1828 frout.title = title;
1829 if (bSeparate || bTDump)
1831 frout.bTitle = TRUE;
1834 frout.bAtoms = TRUE;
1836 frout.atoms = &useatoms;
1837 frout.bStep = FALSE;
1838 frout.bTime = FALSE;
1842 frout.bTitle = (outframe == 0);
1843 frout.bAtoms = FALSE;
1847 write_g96_conf(out, &frout, -1, NULL);
1849 if (bSeparate || bSplitHere)
1856 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1858 if (bSeparate || bSplitHere)
1863 /* execute command */
1867 sprintf(c, "%s %d", exec_command, file_nr-1);
1868 /*fprintf(stderr,"Executing '%s'\n",c);*/
1871 gmx_fatal(FARGS, "Error executing command: %s", c);
1878 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1880 while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1883 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1885 fprintf(stderr, "\nWARNING no output, "
1886 "last frame read at t=%g\n", fr.time);
1888 fprintf(stderr, "\n");
1895 gmx_rmpbc_done(gpbc);
1902 else if (out != NULL)
1908 for (i = 0; (i < clust->clust->nr); i++)
1910 if (clust_status_id[i] >= 0)
1912 close_trx(clust_status[i]);
1920 do_view(oenv, out_file, NULL);