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,2019, 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.
48 #include "gromacs/commandline/pargs.h"
49 #include "gromacs/commandline/viewit.h"
50 #include "gromacs/fileio/confio.h"
51 #include "gromacs/fileio/g96io.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/groio.h"
54 #include "gromacs/fileio/pdbio.h"
55 #include "gromacs/fileio/tngio.h"
56 #include "gromacs/fileio/tpxio.h"
57 #include "gromacs/fileio/trrio.h"
58 #include "gromacs/fileio/trxio.h"
59 #include "gromacs/fileio/xtcio.h"
60 #include "gromacs/fileio/xvgr.h"
61 #include "gromacs/math/do_fit.h"
62 #include "gromacs/math/functions.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/pbcutil/rmpbc.h"
67 #include "gromacs/topology/index.h"
68 #include "gromacs/topology/topology.h"
69 #include "gromacs/trajectory/trajectoryframe.h"
70 #include "gromacs/utility/arrayref.h"
71 #include "gromacs/utility/arraysize.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/futil.h"
74 #include "gromacs/utility/smalloc.h"
77 euSel, euRect, euTric, euCompact, euNR
81 static void calc_pbc_cluster(int ecenter, int nrefat, t_topology *top, int ePBC,
82 rvec x[], const int index[], matrix box)
84 int m, i, j, j0, j1, jj, ai, aj;
87 rvec dx, xtest, box_center;
88 int nmol, imol_center;
90 gmx_bool *bMol, *bTmp;
91 rvec *m_com, *m_shift;
98 calc_box_center(ecenter, box, box_center);
100 /* Initiate the pbc structure */
101 std::memset(&pbc, 0, sizeof(pbc));
102 set_pbc(&pbc, ePBC, box);
104 /* Convert atom index to molecular */
106 molind = top->mols.index;
112 snew(bTmp, top->atoms.nr);
114 for (i = 0; (i < nrefat); i++)
116 /* Mark all molecules in the index */
119 /* Binary search assuming the molecules are sorted */
124 if (ai < molind[j0+1])
128 else if (ai >= molind[j1])
135 if (ai < molind[jj+1])
147 /* Double check whether all atoms in all molecules that are marked are part
148 * of the cluster. Simultaneously compute the center of geometry.
150 min_dist2 = 10*gmx::square(trace(box));
153 for (i = 0; i < nmol; i++)
155 for (j = molind[i]; j < molind[i+1]; j++)
157 if (bMol[i] && !bTmp[j])
159 gmx_fatal(FARGS, "Molecule %d marked for clustering but not atom %d in it - check your index!", i+1, j+1);
161 else if (!bMol[i] && bTmp[j])
163 gmx_fatal(FARGS, "Atom %d marked for clustering but not molecule %d - this is an internal error...", j+1, i+1);
167 /* Make molecule whole, move 2nd and higher atom to same periodicity as 1st atom in molecule */
170 pbc_dx(&pbc, x[j], x[j-1], dx);
171 rvec_add(x[j-1], dx, x[j]);
173 /* Compute center of geometry of molecule - m_com[i] was zeroed when we did snew() on it! */
174 rvec_inc(m_com[i], x[j]);
179 /* Normalize center of geometry */
180 fac = 1.0/(molind[i+1]-molind[i]);
181 for (m = 0; (m < DIM); m++)
185 /* Determine which molecule is closest to the center of the box */
186 pbc_dx(&pbc, box_center, m_com[i], dx);
187 tmp_r2 = iprod(dx, dx);
189 if (tmp_r2 < min_dist2)
194 cluster[ncluster++] = i;
201 fprintf(stderr, "No molecules selected in the cluster\n");
204 else if (imol_center == -1)
206 fprintf(stderr, "No central molecules could be found\n");
211 added[nadded++] = imol_center;
212 bMol[imol_center] = FALSE;
214 while (nadded < ncluster)
216 /* Find min distance between cluster molecules and those remaining to be added */
217 min_dist2 = 10*gmx::square(trace(box));
220 /* Loop over added mols */
221 for (i = 0; i < nadded; i++)
224 /* Loop over all mols */
225 for (j = 0; j < ncluster; j++)
228 /* check those remaining to be added */
231 pbc_dx(&pbc, m_com[aj], m_com[ai], dx);
232 tmp_r2 = iprod(dx, dx);
233 if (tmp_r2 < min_dist2)
243 /* Add the best molecule */
244 added[nadded++] = jmin;
246 /* Calculate the shift from the ai molecule */
247 pbc_dx(&pbc, m_com[jmin], m_com[imin], dx);
248 rvec_add(m_com[imin], dx, xtest);
249 rvec_sub(xtest, m_com[jmin], m_shift[jmin]);
250 rvec_inc(m_com[jmin], m_shift[jmin]);
252 for (j = molind[jmin]; j < molind[jmin+1]; j++)
254 rvec_inc(x[j], m_shift[jmin]);
256 fprintf(stdout, "\rClustering iteration %d of %d...", nadded, ncluster);
266 fprintf(stdout, "\n");
269 static void put_molecule_com_in_box(int unitcell_enum, int ecenter,
271 int natoms, t_atom atom[],
272 int ePBC, matrix box, rvec x[])
276 rvec com, shift, box_center;
281 calc_box_center(ecenter, box, box_center);
282 set_pbc(&pbc, ePBC, box);
285 gmx_fatal(FARGS, "There are no molecule descriptions. I need a .tpr file for this pbc option.");
287 for (i = 0; (i < mols->nr); i++)
292 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
295 for (d = 0; d < DIM; d++)
301 /* calculate final COM */
302 svmul(1.0/mtot, com, com);
304 /* check if COM is outside box */
306 copy_rvec(com, newCom);
307 auto newComArrayRef = gmx::arrayRefFromArray(&newCom, 1);
308 switch (unitcell_enum)
311 put_atoms_in_box(ePBC, box, newComArrayRef);
314 put_atoms_in_triclinic_unitcell(ecenter, box, newComArrayRef);
317 put_atoms_in_compact_unitcell(ePBC, ecenter, box, newComArrayRef);
320 rvec_sub(newCom, com, shift);
321 if (norm2(shift) > 0)
325 fprintf(debug, "\nShifting position of molecule %d "
326 "by %8.3f %8.3f %8.3f\n", i+1,
327 shift[XX], shift[YY], shift[ZZ]);
329 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
331 rvec_inc(x[j], shift);
337 static void put_residue_com_in_box(int unitcell_enum, int ecenter,
338 int natoms, t_atom atom[],
339 int ePBC, matrix box, rvec x[])
341 int i, j, res_start, res_end;
345 rvec box_center, com, shift;
346 static const int NOTSET = -12347;
347 calc_box_center(ecenter, box, box_center);
353 for (i = 0; i < natoms+1; i++)
355 if (i == natoms || (presnr != atom[i].resind && presnr != NOTSET))
357 /* calculate final COM */
359 svmul(1.0/mtot, com, com);
361 /* check if COM is outside box */
363 copy_rvec(com, newCom);
364 auto newComArrayRef = gmx::arrayRefFromArray(&newCom, 1);
365 switch (unitcell_enum)
368 put_atoms_in_box(ePBC, box, newComArrayRef);
371 put_atoms_in_triclinic_unitcell(ecenter, box, newComArrayRef);
374 put_atoms_in_compact_unitcell(ePBC, ecenter, box, newComArrayRef);
377 rvec_sub(newCom, com, shift);
378 if (norm2(shift) != 0.0f)
382 fprintf(debug, "\nShifting position of residue %d (atoms %d-%d) "
383 "by %g,%g,%g\n", atom[res_start].resind+1,
384 res_start+1, res_end+1, shift[XX], shift[YY], shift[ZZ]);
386 for (j = res_start; j < res_end; j++)
388 rvec_inc(x[j], shift);
394 /* remember start of new residue */
401 for (d = 0; d < DIM; d++)
407 presnr = atom[i].resind;
412 static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, const int ci[])
415 rvec cmin, cmax, box_center, dx;
419 copy_rvec(x[ci[0]], cmin);
420 copy_rvec(x[ci[0]], cmax);
421 for (i = 0; i < nc; i++)
424 for (m = 0; m < DIM; m++)
426 if (x[ai][m] < cmin[m])
430 else if (x[ai][m] > cmax[m])
436 calc_box_center(ecenter, box, box_center);
437 for (m = 0; m < DIM; m++)
439 dx[m] = box_center[m]-(cmin[m]+cmax[m])*0.5;
442 for (i = 0; i < n; i++)
449 static void mk_filenm(char *base, const char *ext, int ndigit, int file_nr,
455 std::strcpy(out_file, base);
466 std::strncat(out_file, "00000000000", ndigit-nd);
468 sprintf(nbuf, "%d.", file_nr);
469 std::strcat(out_file, nbuf);
470 std::strcat(out_file, ext);
473 static void check_trr(const char *fn)
475 if (fn2ftp(fn) != efTRR)
477 gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
481 static void do_trunc(const char *fn, real t0)
494 gmx_fatal(FARGS, "You forgot to set the truncation time");
497 /* Check whether this is a .trr file */
500 in = gmx_trr_open(fn, "r");
501 fp = gmx_fio_getfp(in);
504 fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
510 fpos = gmx_fio_ftell(in);
512 while (!bStop && gmx_trr_read_frame_header(in, &sh, &bOK))
514 gmx_trr_read_frame_data(in, &sh, nullptr, nullptr, nullptr, nullptr);
515 fpos = gmx_ftell(fp);
519 gmx_fseek(fp, fpos, SEEK_SET);
525 fprintf(stderr, "Do you REALLY want to truncate this trajectory (%s) at:\n"
526 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
527 fn, j, t, static_cast<long int>(fpos));
528 if (1 != scanf("%s", yesno))
530 gmx_fatal(FARGS, "Error reading user input");
532 if (std::strcmp(yesno, "YES") == 0)
534 fprintf(stderr, "Once again, I'm gonna DO this...\n");
536 if (0 != gmx_truncate(fn, fpos))
538 gmx_fatal(FARGS, "Error truncating file %s", fn);
543 fprintf(stderr, "Ok, I'll forget about it\n");
548 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
554 /*! \brief Read a full molecular topology if useful and available.
556 * If the input trajectory file is not in TNG format, and the output
557 * file is in TNG format, then we want to try to read a full topology
558 * (if available), so that we can write molecule information to the
559 * output file. The full topology provides better molecule information
560 * than is available from the normal t_topology data used by GROMACS
563 * Also, the t_topology is only read under (different) particular
564 * conditions. If both apply, then a .tpr file might be read
565 * twice. Trying to fix this redundancy while trjconv is still an
566 * all-purpose tool does not seem worthwhile.
568 * Because of the way gmx_prepare_tng_writing is implemented, the case
569 * where the input TNG file has no molecule information will never
570 * lead to an output TNG file having molecule information. Since
571 * molecule information will generally be present if the input TNG
572 * file was written by a GROMACS tool, this seems like reasonable
574 static std::unique_ptr<gmx_mtop_t>
575 read_mtop_for_tng(const char *tps_file,
576 const char *input_file,
577 const char *output_file)
579 std::unique_ptr<gmx_mtop_t> mtop;
581 if (fn2bTPX(tps_file) &&
582 efTNG != fn2ftp(input_file) &&
583 efTNG == fn2ftp(output_file))
585 int temp_natoms = -1;
586 mtop = std::make_unique<gmx_mtop_t>();
587 read_tpx(tps_file, nullptr, nullptr, &temp_natoms,
588 nullptr, nullptr, mtop.get());
594 int gmx_trjconv(int argc, char *argv[])
596 const char *desc[] = {
597 "[THISMODULE] can convert trajectory files in many ways:",
599 "* from one format to another",
600 "* select a subset of atoms",
601 "* change the periodicity representation",
602 "* keep multimeric molecules together",
603 "* center atoms in the box",
604 "* fit atoms to reference structure",
605 "* reduce the number of frames",
606 "* change the timestamps of the frames ([TT]-t0[tt] and [TT]-timestep[tt])",
607 "* cut the trajectory in small subtrajectories according",
608 " to information in an index file. This allows subsequent analysis of",
609 " the subtrajectories that could, for example, be the result of a",
610 " cluster analysis. Use option [TT]-sub[tt].",
611 " This assumes that the entries in the index file are frame numbers and",
612 " dumps each group in the index file to a separate trajectory file.",
613 "* select frames within a certain range of a quantity given",
614 " in an [REF].xvg[ref] file.",
616 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
619 "The following formats are supported for input and output:",
620 "[REF].xtc[ref], [REF].trr[ref], [REF].gro[ref], [TT].g96[tt]",
621 "and [REF].pdb[ref].",
622 "The file formats are detected from the file extension.",
623 "The precision of the [REF].xtc[ref] output is taken from the",
624 "input file for [REF].xtc[ref], [REF].gro[ref] and [REF].pdb[ref],",
625 "and from the [TT]-ndec[tt] option for other input formats. The precision",
626 "is always taken from [TT]-ndec[tt], when this option is set.",
627 "All other formats have fixed precision. [REF].trr[ref]",
628 "output can be single or double precision, depending on the precision",
629 "of the [THISMODULE] binary.",
630 "Note that velocities are only supported in",
631 "[REF].trr[ref], [REF].gro[ref] and [TT].g96[tt] files.[PAR]",
633 "Option [TT]-sep[tt] can be used to write every frame to a separate",
634 "[TT].gro, .g96[tt] or [REF].pdb[ref] file. By default, all frames all written to one file.",
635 "[REF].pdb[ref] files with all frames concatenated can be viewed with",
636 "[TT]rasmol -nmrpdb[tt].[PAR]",
638 "It is possible to select part of your trajectory and write it out",
639 "to a new trajectory file in order to save disk space, e.g. for leaving",
640 "out the water from a trajectory of a protein in water.",
641 "[BB]ALWAYS[bb] put the original trajectory on tape!",
642 "We recommend to use the portable [REF].xtc[ref] format for your analysis",
643 "to save disk space and to have portable files.[PAR]",
645 "There are two options for fitting the trajectory to a reference",
646 "either for essential dynamics analysis, etc.",
647 "The first option is just plain fitting to a reference structure",
648 "in the structure file. The second option is a progressive fit",
649 "in which the first timeframe is fitted to the reference structure ",
650 "in the structure file to obtain and each subsequent timeframe is ",
651 "fitted to the previously fitted structure. This way a continuous",
652 "trajectory is generated, which might not be the case when using the",
653 "regular fit method, e.g. when your protein undergoes large",
654 "conformational transitions.[PAR]",
656 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
659 " * [TT]mol[tt] puts the center of mass of molecules in the box,",
660 " and requires a run input file to be supplied with [TT]-s[tt].",
661 " * [TT]res[tt] puts the center of mass of residues in the box.",
662 " * [TT]atom[tt] puts all the atoms in the box.",
663 " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
664 " them back. This has the effect that all molecules",
665 " will remain whole (provided they were whole in the initial",
666 " conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
667 " molecules may diffuse out of the box. The starting configuration",
668 " for this procedure is taken from the structure file, if one is",
669 " supplied, otherwise it is the first frame.",
670 " * [TT]cluster[tt] clusters all the atoms in the selected index",
671 " such that they are all closest to the center of mass of the cluster,",
672 " which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
673 " results if you in fact have a cluster. Luckily that can be checked",
674 " afterwards using a trajectory viewer. Note also that if your molecules",
675 " are broken this will not work either.",
676 " * [TT]whole[tt] only makes broken molecules whole.",
679 "Option [TT]-ur[tt] sets the unit cell representation for options",
680 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
681 "All three options give different results for triclinic boxes and",
682 "identical results for rectangular boxes.",
683 "[TT]rect[tt] is the ordinary brick shape.",
684 "[TT]tric[tt] is the triclinic unit cell.",
685 "[TT]compact[tt] puts all atoms at the closest distance from the center",
686 "of the box. This can be useful for visualizing e.g. truncated octahedra",
687 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
688 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
689 "is set differently.[PAR]",
691 "Option [TT]-center[tt] centers the system in the box. The user can",
692 "select the group which is used to determine the geometrical center.",
693 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
694 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
695 "[TT]tric[tt]: half of the sum of the box vectors,",
696 "[TT]rect[tt]: half of the box diagonal,",
697 "[TT]zero[tt]: zero.",
698 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
699 "want all molecules in the box after the centering.[PAR]",
701 "Option [TT]-box[tt] sets the size of the new box. This option only works",
702 "for leading dimensions and is thus generally only useful for rectangular boxes.",
703 "If you want to modify only some of the dimensions, e.g. when reading from",
704 "a trajectory, you can use -1 for those dimensions that should stay the same",
706 "It is not always possible to use combinations of [TT]-pbc[tt],",
707 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
708 "you want in one call to [THISMODULE]. Consider using multiple",
709 "calls, and check out the GROMACS website for suggestions.[PAR]",
711 "With [TT]-dt[tt], it is possible to reduce the number of ",
712 "frames in the output. This option relies on the accuracy of the times",
713 "in your input trajectory, so if these are inaccurate use the",
714 "[TT]-timestep[tt] option to modify the time (this can be done",
715 "simultaneously). For making smooth movies, the program [gmx-filter]",
716 "can reduce the number of frames while using low-pass frequency",
717 "filtering, this reduces aliasing of high frequency motions.[PAR]",
719 "Using [TT]-trunc[tt] [THISMODULE] can truncate [REF].trr[ref] in place, i.e.",
720 "without copying the file. This is useful when a run has crashed",
721 "during disk I/O (i.e. full disk), or when two contiguous",
722 "trajectories must be concatenated without having double frames.[PAR]",
724 "Option [TT]-dump[tt] can be used to extract a frame at or near",
725 "one specific time from your trajectory, but only works reliably",
726 "if the time interval between frames is uniform.[PAR]",
728 "Option [TT]-drop[tt] reads an [REF].xvg[ref] file with times and values.",
729 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
730 "frames with a value below and above the value of the respective options",
731 "will not be written."
747 const char *pbc_opt[epNR + 1] =
749 nullptr, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
754 const char *unitcell_opt[euNR+1] =
755 { nullptr, "rect", "tric", "compact", nullptr };
759 ecSel, ecTric, ecRect, ecZero, ecNR
761 const char *center_opt[ecNR+1] =
762 { nullptr, "tric", "rect", "zero", nullptr };
768 efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
770 const char *fit[efNR + 1] =
772 nullptr, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
773 "progressive", nullptr
776 static gmx_bool bSeparate = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
777 static gmx_bool bCenter = FALSE;
778 static int skip_nr = 1, ndec = 3, nzero = 0;
779 static real tzero = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
780 static rvec newbox = {0, 0, 0}, shift = {0, 0, 0}, trans = {0, 0, 0};
781 static char *exec_command = nullptr;
782 static real dropunder = 0, dropover = 0;
783 static gmx_bool bRound = FALSE;
788 { "-skip", FALSE, etINT,
789 { &skip_nr }, "Only write every nr-th frame" },
790 { "-dt", FALSE, etTIME,
792 "Only write frame when t MOD dt = first time (%t)" },
793 { "-round", FALSE, etBOOL,
794 { &bRound }, "Round measurements to nearest picosecond"},
795 { "-dump", FALSE, etTIME,
796 { &tdump }, "Dump frame nearest specified time (%t)" },
797 { "-t0", FALSE, etTIME,
799 "Starting time (%t) (default: don't change)" },
800 { "-timestep", FALSE, etTIME,
802 "Change time step between input frames (%t)" },
803 { "-pbc", FALSE, etENUM,
805 "PBC treatment (see help text for full description)" },
806 { "-ur", FALSE, etENUM,
807 { unitcell_opt }, "Unit-cell representation" },
808 { "-center", FALSE, etBOOL,
809 { &bCenter }, "Center atoms in box" },
810 { "-boxcenter", FALSE, etENUM,
811 { center_opt }, "Center for -pbc and -center" },
812 { "-box", FALSE, etRVEC,
814 "Size for new cubic box (default: read from input)" },
815 { "-trans", FALSE, etRVEC,
817 "All coordinates will be translated by trans. This "
818 "can advantageously be combined with -pbc mol -ur "
820 { "-shift", FALSE, etRVEC,
822 "All coordinates will be shifted by framenr*shift" },
823 { "-fit", FALSE, etENUM,
825 "Fit molecule to ref structure in the structure file" },
826 { "-ndec", FALSE, etINT,
828 "Number of decimal places to write to .xtc output" },
829 { "-vel", FALSE, etBOOL,
830 { &bVels }, "Read and write velocities if possible" },
831 { "-force", FALSE, etBOOL,
832 { &bForce }, "Read and write forces if possible" },
833 { "-trunc", FALSE, etTIME,
835 "Truncate input trajectory file after this time (%t)" },
836 { "-exec", FALSE, etSTR,
838 "Execute command for every output frame with the "
839 "frame number as argument" },
840 { "-split", FALSE, etTIME,
842 "Start writing new file when t MOD split = first "
844 { "-sep", FALSE, etBOOL,
846 "Write each frame to a separate .gro, .g96 or .pdb "
848 { "-nzero", FALSE, etINT,
850 "If the -sep flag is set, use these many digits "
851 "for the file numbers and prepend zeros as needed" },
852 { "-dropunder", FALSE, etREAL,
853 { &dropunder }, "Drop all frames below this value" },
854 { "-dropover", FALSE, etREAL,
855 { &dropover }, "Drop all frames above this value" },
856 { "-conect", FALSE, etBOOL,
858 "Add conect records when writing [REF].pdb[ref] files. Useful "
859 "for visualization of non-standard molecules, e.g. "
860 "coarse grained ones" }
862 #define NPA asize(pa)
865 t_trxstatus *trxout = nullptr;
868 t_trxframe fr, frout;
870 rvec *xmem = nullptr, *vmem = nullptr, *fmem = nullptr;
871 rvec *xp = nullptr, x_shift, hbox;
872 real *w_rls = nullptr;
873 int m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
875 t_topology *top = nullptr;
876 gmx_conect gc = nullptr;
878 t_atoms *atoms = nullptr, useatoms;
880 int *index = nullptr, *cindex = nullptr;
881 char *grpnm = nullptr;
884 int ifit, my_clust = -1;
887 t_cluster_ndx *clust = nullptr;
888 t_trxstatus **clust_status = nullptr;
889 int *clust_status_id = nullptr;
891 int *nfwritten = nullptr;
892 int ndrop = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
894 real tshift = 0, dt = -1, prec;
895 gmx_bool bFit, bPFit, bReset;
897 gmx_rmpbc_t gpbc = nullptr;
898 gmx_bool bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
899 gmx_bool bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
900 gmx_bool bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetXtcPrec, bNeedPrec;
901 gmx_bool bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
902 gmx_bool bSubTraj = FALSE, bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
903 gmx_bool bWriteFrame, bSplitHere;
904 const char *top_file, *in_file, *out_file = nullptr;
905 char out_file2[256], *charpt;
906 char *outf_base = nullptr;
907 const char *outf_ext = nullptr;
908 char top_title[256], timestr[32], stepstr[32], filemode[5];
909 gmx_output_env_t *oenv;
912 { efTRX, "-f", nullptr, ffREAD },
913 { efTRO, "-o", nullptr, ffWRITE },
914 { efTPS, nullptr, nullptr, ffOPTRD },
915 { efNDX, nullptr, nullptr, ffOPTRD },
916 { efNDX, "-fr", "frames", ffOPTRD },
917 { efNDX, "-sub", "cluster", ffOPTRD },
918 { efXVG, "-drop", "drop", ffOPTRD }
920 #define NFILE asize(fnm)
922 if (!parse_common_args(&argc, argv,
923 PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW |
925 NFILE, fnm, NPA, pa, asize(desc), desc,
930 fprintf(stdout, "Note that major changes are planned in future for "
931 "trjconv, to improve usability and utility.\n");
933 top_file = ftp2fn(efTPS, NFILE, fnm);
935 /* Check command line */
936 in_file = opt2fn("-f", NFILE, fnm);
940 do_trunc(in_file, ttrunc);
944 /* mark active cmdline options */
945 bSetBox = opt2parg_bSet("-box", NPA, pa);
946 bSetTime = opt2parg_bSet("-t0", NPA, pa);
947 bSetXtcPrec = opt2parg_bSet("-ndec", NPA, pa);
948 bSetUR = opt2parg_bSet("-ur", NPA, pa);
949 bExec = opt2parg_bSet("-exec", NPA, pa);
950 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
951 bTDump = opt2parg_bSet("-dump", NPA, pa);
952 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
953 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
954 bTrans = opt2parg_bSet("-trans", NPA, pa);
955 bSplit = (split_t != 0);
957 /* parse enum options */
958 fit_enum = nenum(fit);
959 bFit = (fit_enum == efFit || fit_enum == efFitXY);
960 bReset = (fit_enum == efReset || fit_enum == efResetXY);
961 bPFit = fit_enum == efPFit;
962 pbc_enum = nenum(pbc_opt);
963 bPBCWhole = pbc_enum == epWhole;
964 bPBCcomRes = pbc_enum == epComRes;
965 bPBCcomMol = pbc_enum == epComMol;
966 bPBCcomAtom = pbc_enum == epComAtom;
967 bNoJump = pbc_enum == epNojump;
968 bCluster = pbc_enum == epCluster;
969 bPBC = pbc_enum != epNone;
970 unitcell_enum = nenum(unitcell_opt);
971 ecenter = nenum(center_opt) - ecTric;
973 /* set and check option dependencies */
976 bFit = TRUE; /* for pfit, fit *must* be set */
980 bReset = TRUE; /* for fit, reset *must* be set */
985 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
987 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
991 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
994 "WARNING: Option for unitcell representation (-ur %s)\n"
995 " only has effect in combination with -pbc %s, %s or %s.\n"
996 " Ingoring unitcell representation.\n\n",
997 unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
1002 gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
1003 "Please do the PBC condition treatment first and then run trjconv in a second step\n"
1004 "for the rotational fit.\n"
1005 "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1009 /* ndec for XTC writing is in nr of decimal places, prec is a multiplication factor: */
1011 for (i = 0; i < ndec; i++)
1016 bIndex = ftp2bSet(efNDX, NFILE, fnm);
1019 /* Determine output type */
1020 out_file = opt2fn("-o", NFILE, fnm);
1021 int ftp = fn2ftp(out_file);
1022 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
1023 bNeedPrec = (ftp == efXTC);
1024 int ftpin = fn2ftp(in_file);
1027 /* check if velocities are possible in input and output files */
1028 bVels = (ftp == efTRR || ftp == efGRO ||
1029 ftp == efG96 || ftp == efTNG)
1030 && (ftpin == efTRR || ftpin == efGRO ||
1031 ftpin == efG96 || ftpin == efTNG || ftpin == efCPT);
1033 if (bSeparate || bSplit)
1035 outf_ext = std::strrchr(out_file, '.');
1036 if (outf_ext == nullptr)
1038 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1040 outf_base = gmx_strdup(out_file);
1041 outf_base[outf_ext - out_file] = '\0';
1044 bSubTraj = opt2bSet("-sub", NFILE, fnm);
1047 if ((ftp != efXTC) && (ftp != efTRR))
1049 /* It seems likely that other trajectory file types
1050 * could work here. */
1051 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1054 clust = cluster_index(nullptr, opt2fn("-sub", NFILE, fnm));
1056 /* Check for number of files disabled, as FOPEN_MAX is not the correct
1057 * number to check for. In my linux box it is only 16.
1059 if (/* DISABLES CODE */ (false))
1060 //if (clust->clust->nr > FOPEN_MAX-4)
1062 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1063 " trajectories.\ntry splitting the index file in %d parts.\n"
1065 clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1067 gmx_warning("The -sub option could require as many open output files as there are\n"
1068 "index groups in the file (%d). If you get I/O errors opening new files,\n"
1069 "try reducing the number of index groups in the file, and perhaps\n"
1070 "using trjconv -sub several times on different chunks of your index file.\n",
1073 snew(clust_status, clust->clust->nr);
1074 snew(clust_status_id, clust->clust->nr);
1075 snew(nfwritten, clust->clust->nr);
1076 for (i = 0; (i < clust->clust->nr); i++)
1078 clust_status[i] = nullptr;
1079 clust_status_id[i] = -1;
1081 bSeparate = bSplit = FALSE;
1086 gmx_fatal(FARGS, "Argument for -skip (%d) needs to be greater or equal to 1.", skip_nr);
1089 std::unique_ptr<gmx_mtop_t> mtop = read_mtop_for_tng(top_file, in_file, out_file);
1091 /* Determine whether to read a topology */
1092 bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1093 bRmPBC || bReset || bPBCcomMol || bCluster ||
1094 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1096 /* Determine if when can read index groups */
1097 bIndex = (bIndex || bTPS);
1102 read_tps_conf(top_file, top, &ePBC, &xp, nullptr, top_box,
1103 bReset || bPBCcomRes);
1104 std::strncpy(top_title, *top->name, 255);
1105 top_title[255] = '\0';
1106 atoms = &top->atoms;
1108 if (0 == top->mols.nr && (bCluster || bPBCcomMol))
1110 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1113 /* top_title is only used for gro and pdb,
1114 * the header in such a file is top_title, followed by
1115 * t= ... and/or step= ...
1116 * to prevent double t= or step=, remove it from top_title.
1117 * From GROMACS-2018 we only write t/step when the frame actually
1118 * has a valid time/step, so we need to check for both separately.
1120 if ((charpt = std::strstr(top_title, " t= ")))
1124 if ((charpt = std::strstr(top_title, " step= ")))
1131 gc = gmx_conect_generate(top);
1135 gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
1139 /* get frame number index */
1141 if (opt2bSet("-fr", NFILE, fnm))
1143 printf("Select groups of frame number indices:\n");
1144 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, &frindex, &frname);
1147 for (i = 0; i < nrfri; i++)
1149 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1154 /* get index groups etc. */
1157 printf("Select group for %s fit\n",
1158 bFit ? "least squares" : "translational");
1159 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1160 1, &ifit, &ind_fit, &gn_fit);
1166 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1170 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1176 printf("Select group for clustering\n");
1177 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1178 1, &ifit, &ind_fit, &gn_fit);
1185 printf("Select group for centering\n");
1186 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1187 1, &ncent, &cindex, &grpnm);
1189 printf("Select group for output\n");
1190 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1191 1, &nout, &index, &grpnm);
1195 /* no index file, so read natoms from TRX */
1196 if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
1198 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1203 snew(index, natoms);
1204 for (i = 0; i < natoms; i++)
1218 snew(w_rls, atoms->nr);
1219 for (i = 0; (i < ifit); i++)
1221 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1224 /* Restore reference structure and set to origin,
1225 store original location (to put structure back) */
1228 gmx_rmpbc(gpbc, top->atoms.nr, top_box, xp);
1230 copy_rvec(xp[index[0]], x_shift);
1231 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, nullptr, xp, w_rls);
1232 rvec_dec(x_shift, xp[index[0]]);
1236 clear_rvec(x_shift);
1239 if (bDropUnder || bDropOver)
1241 /* Read the .xvg file with the drop values */
1242 fprintf(stderr, "\nReading drop file ...");
1243 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1244 fprintf(stderr, " %d time points\n", ndrop);
1245 if (ndrop == 0 || ncol < 2)
1247 gmx_fatal(FARGS, "Found no data points in %s",
1248 opt2fn("-drop", NFILE, fnm));
1254 /* Make atoms struct for output in GRO or PDB files */
1255 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1257 /* get memory for stuff to go in .pdb file, and initialize
1258 * the pdbinfo structure part if the input has it.
1260 init_t_atoms(&useatoms, atoms->nr, atoms->havePdbInfo);
1261 sfree(useatoms.resinfo);
1262 useatoms.resinfo = atoms->resinfo;
1263 for (i = 0; (i < nout); i++)
1265 useatoms.atomname[i] = atoms->atomname[index[i]];
1266 useatoms.atom[i] = atoms->atom[index[i]];
1267 if (atoms->havePdbInfo)
1269 useatoms.pdbinfo[i] = atoms->pdbinfo[index[i]];
1271 useatoms.nres = std::max(useatoms.nres, useatoms.atom[i].resind+1);
1275 /* select what to read */
1286 flags = flags | TRX_READ_V;
1290 flags = flags | TRX_READ_F;
1293 /* open trx file for reading */
1294 bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
1297 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1301 if (bSetXtcPrec || !fr.bPrec)
1303 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1307 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1311 if (bHaveFirstFrame)
1315 // Determine timestep (assuming constant spacing for now) if we
1316 // need to dump frames based on time. This is required so we do not
1317 // skip the first frame if that was the one that should have been dumped
1318 double firstFrameTime = fr.time;
1319 if (read_next_frame(oenv, trxin, &fr))
1321 dt = fr.time - firstFrameTime;
1325 fprintf(stderr, "Warning: Frame times are not incrementing - will dump first frame.\n");
1328 // Now close and reopen so we are at first frame again
1331 // Reopen at first frame (We already know it exists if we got here)
1332 read_first_frame(oenv, &trxin, in_file, &fr, flags);
1335 set_trxframe_ePBC(&fr, ePBC);
1340 tshift = tzero-fr.time;
1350 /* check if index is meaningful */
1351 for (i = 0; i < nout; i++)
1353 if (index[i] >= natoms)
1356 "Index[%d] %d is larger than the number of atoms in the\n"
1357 "trajectory file (%d). There is a mismatch in the contents\n"
1358 "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1360 bCopy = bCopy || (i != index[i]);
1364 /* open output for writing */
1365 std::strcpy(filemode, "w");
1369 trxout = trjtools_gmx_prepare_tng_writing(out_file,
1375 gmx::arrayRefFromArray(index, nout),
1381 if (!bSplit && !bSubTraj)
1383 trxout = open_trx(out_file, filemode);
1389 if (( !bSeparate && !bSplit ) && !bSubTraj)
1391 out = gmx_ffopen(out_file, filemode);
1395 gmx_incons("Illegal output file format");
1411 /* Start the big loop over frames */
1417 /* Main loop over frames */
1428 /*if (frame >= clust->clust->nra)
1429 gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1430 if (frame > clust->maxframe)
1436 my_clust = clust->inv_clust[frame];
1438 if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1447 /* generate new box */
1452 for (m = 0; m < DIM; m++)
1456 fr.box[m][m] = newbox[m];
1462 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1470 for (i = 0; i < natoms; i++)
1472 rvec_inc(fr.x[i], trans);
1478 // If we could not read two frames or times are not incrementing
1479 // we have almost no idea what to do,
1480 // but dump the first frame so output is not broken.
1481 if (dt <= 0 || !bDTset)
1487 // Dump the frame if we are less than half a frame time
1488 // below it. This will also ensure we at least dump a
1489 // somewhat reasonable frame if the spacing is unequal
1490 // and we have overrun the frame time. Once we dump one
1491 // frame based on time we quit, so it does not matter
1492 // that this might be true for all subsequent frames too.
1493 bDumpFrame = (fr.time > tdump-0.5*dt);
1501 /* determine if an atom jumped across the box and reset it if so */
1502 if (bNoJump && (bTPS || frame != 0))
1504 for (d = 0; d < DIM; d++)
1506 hbox[d] = 0.5*fr.box[d][d];
1508 for (i = 0; i < natoms; i++)
1512 rvec_dec(fr.x[i], x_shift);
1514 for (m = DIM-1; m >= 0; m--)
1518 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1520 for (d = 0; d <= m; d++)
1522 fr.x[i][d] += fr.box[m][d];
1525 while (fr.x[i][m]-xp[i][m] > hbox[m])
1527 for (d = 0; d <= m; d++)
1529 fr.x[i][d] -= fr.box[m][d];
1538 calc_pbc_cluster(ecenter, ifit, top, ePBC, fr.x, ind_fit, fr.box);
1543 /* Now modify the coords according to the flags,
1544 for normal fit, this is only done for output frames */
1547 gmx_rmpbc_trxfr(gpbc, &fr);
1550 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1551 do_fit(natoms, w_rls, xp, fr.x);
1554 /* store this set of coordinates for future use */
1555 if (bPFit || bNoJump)
1561 for (i = 0; (i < natoms); i++)
1563 copy_rvec(fr.x[i], xp[i]);
1564 rvec_inc(fr.x[i], x_shift);
1570 /* see if we have a frame from the frame index group */
1571 for (i = 0; i < nrfri && !bDumpFrame; i++)
1573 bDumpFrame = frame == frindex[i];
1576 if (debug && bDumpFrame)
1578 fprintf(debug, "dumping %d\n", frame);
1582 ( ( !bTDump && (frindex == nullptr) && frame % skip_nr == 0 ) || bDumpFrame );
1584 if (bWriteFrame && (bDropUnder || bDropOver))
1586 while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1591 if (std::abs(dropval[0][drop0] - fr.time)
1592 < std::abs(dropval[0][drop1] - fr.time))
1600 if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1601 (bDropOver && dropval[1][dropuse] > dropover))
1603 bWriteFrame = FALSE;
1609 /* We should avoid modifying the input frame,
1610 * but since here we don't have the output frame yet,
1611 * we introduce a temporary output frame time variable.
1615 frout_time = fr.time;
1620 frout_time = tzero + frame*timestep;
1624 frout_time += tshift;
1629 fprintf(stderr, "\nDumping frame at t= %g %s\n",
1630 output_env_conv_time(oenv, frout_time), output_env_get_time_unit(oenv).c_str());
1633 /* check for writing at each delta_t */
1634 bDoIt = (delta_t == 0);
1639 bDoIt = bRmod(frout_time, tzero, delta_t);
1643 /* round() is not C89 compatible, so we do this: */
1644 bDoIt = bRmod(std::floor(frout_time+0.5), std::floor(tzero+0.5),
1645 std::floor(delta_t+0.5));
1649 if (bDoIt || bTDump)
1651 /* print sometimes */
1652 if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1654 fprintf(stderr, " -> frame %6d time %8.3f \r",
1655 outframe, output_env_conv_time(oenv, frout_time));
1661 /* Now modify the coords according to the flags,
1662 for PFit we did this already! */
1666 gmx_rmpbc_trxfr(gpbc, &fr);
1671 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1674 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1678 for (i = 0; i < natoms; i++)
1680 rvec_inc(fr.x[i], x_shift);
1687 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1691 auto positionsArrayRef = gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec *>(fr.x), natoms);
1694 switch (unitcell_enum)
1697 put_atoms_in_box(ePBC, fr.box, positionsArrayRef);
1700 put_atoms_in_triclinic_unitcell(ecenter, fr.box, positionsArrayRef);
1703 put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1710 put_residue_com_in_box(unitcell_enum, ecenter,
1711 natoms, atoms->atom, ePBC, fr.box, fr.x);
1715 put_molecule_com_in_box(unitcell_enum, ecenter,
1717 natoms, atoms->atom, ePBC, fr.box, fr.x);
1719 /* Copy the input trxframe struct to the output trxframe struct */
1721 frout.time = frout_time;
1722 frout.bV = (frout.bV && bVels);
1723 frout.bF = (frout.bF && bForce);
1724 frout.natoms = nout;
1725 if (bNeedPrec && (bSetXtcPrec || !fr.bPrec))
1741 for (i = 0; i < nout; i++)
1743 copy_rvec(fr.x[index[i]], frout.x[i]);
1746 copy_rvec(fr.v[index[i]], frout.v[i]);
1750 copy_rvec(fr.f[index[i]], frout.f[i]);
1755 if (opt2parg_bSet("-shift", NPA, pa))
1757 for (i = 0; i < nout; i++)
1759 for (d = 0; d < DIM; d++)
1761 frout.x[i][d] += outframe*shift[d];
1768 bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1772 /* round() is not C89 compatible, so we do this: */
1773 bSplitHere = bSplit && bRmod(std::floor(frout.time+0.5),
1774 std::floor(tzero+0.5),
1775 std::floor(split_t+0.5));
1777 if (bSeparate || bSplitHere)
1779 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1786 write_tng_frame(trxout, &frout);
1787 // TODO when trjconv behaves better: work how to read and write lambda
1797 trxout = open_trx(out_file2, filemode);
1804 if (clust_status_id[my_clust] == -1)
1806 sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1807 clust_status[my_clust] = open_trx(buf, "w");
1808 clust_status_id[my_clust] = 1;
1811 else if (clust_status_id[my_clust] == -2)
1813 gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1814 clust->grpname[my_clust], ntrxopen, frame,
1817 write_trxframe(clust_status[my_clust], &frout, gc);
1818 nfwritten[my_clust]++;
1819 if (nfwritten[my_clust] ==
1820 (clust->clust->index[my_clust+1]-
1821 clust->clust->index[my_clust]))
1823 close_trx(clust_status[my_clust]);
1824 clust_status[my_clust] = nullptr;
1825 clust_status_id[my_clust] = -2;
1829 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1836 write_trxframe(trxout, &frout, gc);
1842 // Only add a generator statement if title is empty,
1843 // to avoid multiple generated-by statements from various programs
1844 if (std::strlen(top_title) == 0)
1846 sprintf(top_title, "Generated by trjconv");
1850 sprintf(timestr, " t= %9.5f", frout.time);
1854 std::strcpy(timestr, "");
1858 sprintf(stepstr, " step= %" PRId64, frout.step);
1862 std::strcpy(stepstr, "");
1864 title = gmx::formatString("%s%s%s", top_title, timestr, stepstr);
1865 if (bSeparate || bSplitHere)
1867 out = gmx_ffopen(out_file2, "w");
1872 write_hconf_p(out, title.c_str(), &useatoms,
1873 frout.x, frout.bV ? frout.v : nullptr, frout.box);
1876 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1877 /* if reading from pdb, we want to keep the original
1878 model numbering else we write the output frame
1879 number plus one, because model 0 is not allowed in pdb */
1880 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1888 write_pdbfile(out, title.c_str(), &useatoms, frout.x,
1889 frout.ePBC, frout.box, ' ', model_nr, gc);
1892 const char *outputTitle = "";
1893 if (bSeparate || bTDump)
1895 outputTitle = title.c_str();
1898 frout.bAtoms = TRUE;
1900 frout.atoms = &useatoms;
1901 frout.bStep = FALSE;
1902 frout.bTime = FALSE;
1908 outputTitle = title.c_str();
1910 frout.bAtoms = FALSE;
1914 write_g96_conf(out, outputTitle, &frout, -1, nullptr);
1916 if (bSeparate || bSplitHere)
1923 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1925 if (bSeparate || bSplitHere)
1930 /* execute command */
1934 sprintf(c, "%s %d", exec_command, file_nr-1);
1935 /*fprintf(stderr,"Executing '%s'\n",c);*/
1938 gmx_fatal(FARGS, "Error executing command: %s", c);
1945 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1947 while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1950 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1952 fprintf(stderr, "\nWARNING no output, "
1953 "last frame read at t=%g\n", fr.time);
1955 fprintf(stderr, "\n");
1962 gmx_rmpbc_done(gpbc);
1969 else if (out != nullptr)
1975 for (i = 0; (i < clust->clust->nr); i++)
1977 if (clust_status_id[i] >= 0)
1979 close_trx(clust_status[i]);
1999 do_view(oenv, out_file, nullptr);
2001 output_env_done(oenv);