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.
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/commandline/viewit.h"
47 #include "gromacs/compat/make_unique.h"
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/fileio/g96io.h"
50 #include "gromacs/fileio/gmxfio.h"
51 #include "gromacs/fileio/groio.h"
52 #include "gromacs/fileio/pdbio.h"
53 #include "gromacs/fileio/tngio.h"
54 #include "gromacs/fileio/tpxio.h"
55 #include "gromacs/fileio/trrio.h"
56 #include "gromacs/fileio/trxio.h"
57 #include "gromacs/fileio/xtcio.h"
58 #include "gromacs/fileio/xvgr.h"
59 #include "gromacs/gmxana/gmx_ana.h"
60 #include "gromacs/math/do_fit.h"
61 #include "gromacs/math/functions.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdtypes/md_enums.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/pbcutil/rmpbc.h"
66 #include "gromacs/topology/index.h"
67 #include "gromacs/topology/topology.h"
68 #include "gromacs/trajectory/trajectoryframe.h"
69 #include "gromacs/utility/arrayref.h"
70 #include "gromacs/utility/arraysize.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/futil.h"
73 #include "gromacs/utility/smalloc.h"
76 euSel
, euRect
, euTric
, euCompact
, euNR
80 static void calc_pbc_cluster(int ecenter
, int nrefat
, t_topology
*top
, int ePBC
,
81 rvec x
[], const int index
[], matrix box
)
83 int m
, i
, j
, j0
, j1
, jj
, ai
, aj
;
86 rvec dx
, xtest
, box_center
;
87 int nmol
, imol_center
;
89 gmx_bool
*bMol
, *bTmp
;
90 rvec
*m_com
, *m_shift
;
97 calc_box_center(ecenter
, box
, box_center
);
99 /* Initiate the pbc structure */
100 std::memset(&pbc
, 0, sizeof(pbc
));
101 set_pbc(&pbc
, ePBC
, box
);
103 /* Convert atom index to molecular */
105 molind
= top
->mols
.index
;
111 snew(bTmp
, top
->atoms
.nr
);
113 for (i
= 0; (i
< nrefat
); i
++)
115 /* Mark all molecules in the index */
118 /* Binary search assuming the molecules are sorted */
123 if (ai
< molind
[j0
+1])
127 else if (ai
>= molind
[j1
])
134 if (ai
< molind
[jj
+1])
146 /* Double check whether all atoms in all molecules that are marked are part
147 * of the cluster. Simultaneously compute the center of geometry.
149 min_dist2
= 10*gmx::square(trace(box
));
152 for (i
= 0; i
< nmol
; i
++)
154 for (j
= molind
[i
]; j
< molind
[i
+1]; j
++)
156 if (bMol
[i
] && !bTmp
[j
])
158 gmx_fatal(FARGS
, "Molecule %d marked for clustering but not atom %d in it - check your index!", i
+1, j
+1);
160 else if (!bMol
[i
] && bTmp
[j
])
162 gmx_fatal(FARGS
, "Atom %d marked for clustering but not molecule %d - this is an internal error...", j
+1, i
+1);
166 /* Make molecule whole, move 2nd and higher atom to same periodicity as 1st atom in molecule */
169 pbc_dx(&pbc
, x
[j
], x
[j
-1], dx
);
170 rvec_add(x
[j
-1], dx
, x
[j
]);
172 /* Compute center of geometry of molecule - m_com[i] was zeroed when we did snew() on it! */
173 rvec_inc(m_com
[i
], x
[j
]);
178 /* Normalize center of geometry */
179 fac
= 1.0/(molind
[i
+1]-molind
[i
]);
180 for (m
= 0; (m
< DIM
); m
++)
184 /* Determine which molecule is closest to the center of the box */
185 pbc_dx(&pbc
, box_center
, m_com
[i
], dx
);
186 tmp_r2
= iprod(dx
, dx
);
188 if (tmp_r2
< min_dist2
)
193 cluster
[ncluster
++] = i
;
200 fprintf(stderr
, "No molecules selected in the cluster\n");
203 else if (imol_center
== -1)
205 fprintf(stderr
, "No central molecules could be found\n");
210 added
[nadded
++] = imol_center
;
211 bMol
[imol_center
] = FALSE
;
213 while (nadded
< ncluster
)
215 /* Find min distance between cluster molecules and those remaining to be added */
216 min_dist2
= 10*gmx::square(trace(box
));
219 /* Loop over added mols */
220 for (i
= 0; i
< nadded
; i
++)
223 /* Loop over all mols */
224 for (j
= 0; j
< ncluster
; j
++)
227 /* check those remaining to be added */
230 pbc_dx(&pbc
, m_com
[aj
], m_com
[ai
], dx
);
231 tmp_r2
= iprod(dx
, dx
);
232 if (tmp_r2
< min_dist2
)
242 /* Add the best molecule */
243 added
[nadded
++] = jmin
;
245 /* Calculate the shift from the ai molecule */
246 pbc_dx(&pbc
, m_com
[jmin
], m_com
[imin
], dx
);
247 rvec_add(m_com
[imin
], dx
, xtest
);
248 rvec_sub(xtest
, m_com
[jmin
], m_shift
[jmin
]);
249 rvec_inc(m_com
[jmin
], m_shift
[jmin
]);
251 for (j
= molind
[jmin
]; j
< molind
[jmin
+1]; j
++)
253 rvec_inc(x
[j
], m_shift
[jmin
]);
255 fprintf(stdout
, "\rClustering iteration %d of %d...", nadded
, ncluster
);
265 fprintf(stdout
, "\n");
268 static void put_molecule_com_in_box(int unitcell_enum
, int ecenter
,
270 int natoms
, t_atom atom
[],
271 int ePBC
, matrix box
, rvec x
[])
275 rvec com
, shift
, box_center
;
280 calc_box_center(ecenter
, box
, box_center
);
281 set_pbc(&pbc
, ePBC
, box
);
284 gmx_fatal(FARGS
, "There are no molecule descriptions. I need a .tpr file for this pbc option.");
286 for (i
= 0; (i
< mols
->nr
); i
++)
291 for (j
= mols
->index
[i
]; (j
< mols
->index
[i
+1] && j
< natoms
); j
++)
294 for (d
= 0; d
< DIM
; d
++)
300 /* calculate final COM */
301 svmul(1.0/mtot
, com
, com
);
303 /* check if COM is outside box */
305 copy_rvec(com
, newCom
);
306 auto newComArrayRef
= gmx::arrayRefFromArray(&newCom
, 1);
307 switch (unitcell_enum
)
310 put_atoms_in_box(ePBC
, box
, newComArrayRef
);
313 put_atoms_in_triclinic_unitcell(ecenter
, box
, newComArrayRef
);
316 put_atoms_in_compact_unitcell(ePBC
, ecenter
, box
, newComArrayRef
);
319 rvec_sub(newCom
, com
, shift
);
320 if (norm2(shift
) > 0)
324 fprintf(debug
, "\nShifting position of molecule %d "
325 "by %8.3f %8.3f %8.3f\n", i
+1,
326 shift
[XX
], shift
[YY
], shift
[ZZ
]);
328 for (j
= mols
->index
[i
]; (j
< mols
->index
[i
+1] && j
< natoms
); j
++)
330 rvec_inc(x
[j
], shift
);
336 static void put_residue_com_in_box(int unitcell_enum
, int ecenter
,
337 int natoms
, t_atom atom
[],
338 int ePBC
, matrix box
, rvec x
[])
340 int i
, j
, res_start
, res_end
;
344 rvec box_center
, com
, shift
;
345 static const int NOTSET
= -12347;
346 calc_box_center(ecenter
, box
, box_center
);
352 for (i
= 0; i
< natoms
+1; i
++)
354 if (i
== natoms
|| (presnr
!= atom
[i
].resind
&& presnr
!= NOTSET
))
356 /* calculate final COM */
358 svmul(1.0/mtot
, com
, com
);
360 /* check if COM is outside box */
362 copy_rvec(com
, newCom
);
363 auto newComArrayRef
= gmx::arrayRefFromArray(&newCom
, 1);
364 switch (unitcell_enum
)
367 put_atoms_in_box(ePBC
, box
, newComArrayRef
);
370 put_atoms_in_triclinic_unitcell(ecenter
, box
, newComArrayRef
);
373 put_atoms_in_compact_unitcell(ePBC
, ecenter
, box
, newComArrayRef
);
376 rvec_sub(newCom
, com
, shift
);
381 fprintf(debug
, "\nShifting position of residue %d (atoms %d-%d) "
382 "by %g,%g,%g\n", atom
[res_start
].resind
+1,
383 res_start
+1, res_end
+1, shift
[XX
], shift
[YY
], shift
[ZZ
]);
385 for (j
= res_start
; j
< res_end
; j
++)
387 rvec_inc(x
[j
], shift
);
393 /* remember start of new residue */
400 for (d
= 0; d
< DIM
; d
++)
406 presnr
= atom
[i
].resind
;
411 static void center_x(int ecenter
, rvec x
[], matrix box
, int n
, int nc
, const int ci
[])
414 rvec cmin
, cmax
, box_center
, dx
;
418 copy_rvec(x
[ci
[0]], cmin
);
419 copy_rvec(x
[ci
[0]], cmax
);
420 for (i
= 0; i
< nc
; i
++)
423 for (m
= 0; m
< DIM
; m
++)
425 if (x
[ai
][m
] < cmin
[m
])
429 else if (x
[ai
][m
] > cmax
[m
])
435 calc_box_center(ecenter
, box
, box_center
);
436 for (m
= 0; m
< DIM
; m
++)
438 dx
[m
] = box_center
[m
]-(cmin
[m
]+cmax
[m
])*0.5;
441 for (i
= 0; i
< n
; i
++)
448 static void mk_filenm(char *base
, const char *ext
, int ndigit
, int file_nr
,
454 std::strcpy(out_file
, base
);
465 std::strncat(out_file
, "00000000000", ndigit
-nd
);
467 sprintf(nbuf
, "%d.", file_nr
);
468 std::strcat(out_file
, nbuf
);
469 std::strcat(out_file
, ext
);
472 static void check_trr(const char *fn
)
474 if (fn2ftp(fn
) != efTRR
)
476 gmx_fatal(FARGS
, "%s is not a trajectory file, exiting\n", fn
);
480 static void do_trunc(const char *fn
, real t0
)
493 gmx_fatal(FARGS
, "You forgot to set the truncation time");
496 /* Check whether this is a .trr file */
499 in
= gmx_trr_open(fn
, "r");
500 fp
= gmx_fio_getfp(in
);
503 fprintf(stderr
, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn
);
509 fpos
= gmx_fio_ftell(in
);
511 while (!bStop
&& gmx_trr_read_frame_header(in
, &sh
, &bOK
))
513 gmx_trr_read_frame_data(in
, &sh
, nullptr, nullptr, nullptr, nullptr);
514 fpos
= gmx_ftell(fp
);
518 gmx_fseek(fp
, fpos
, SEEK_SET
);
524 fprintf(stderr
, "Do you REALLY want to truncate this trajectory (%s) at:\n"
525 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
526 fn
, j
, t
, (long int)fpos
);
527 if (1 != scanf("%s", yesno
))
529 gmx_fatal(FARGS
, "Error reading user input");
531 if (std::strcmp(yesno
, "YES") == 0)
533 fprintf(stderr
, "Once again, I'm gonna DO this...\n");
535 if (0 != gmx_truncate(fn
, fpos
))
537 gmx_fatal(FARGS
, "Error truncating file %s", fn
);
542 fprintf(stderr
, "Ok, I'll forget about it\n");
547 fprintf(stderr
, "Already at end of file (t=%g)...\n", t
);
553 /*! \brief Read a full molecular topology if useful and available.
555 * If the input trajectory file is not in TNG format, and the output
556 * file is in TNG format, then we want to try to read a full topology
557 * (if available), so that we can write molecule information to the
558 * output file. The full topology provides better molecule information
559 * than is available from the normal t_topology data used by GROMACS
562 * Also, the t_topology is only read under (different) particular
563 * conditions. If both apply, then a .tpr file might be read
564 * twice. Trying to fix this redundancy while trjconv is still an
565 * all-purpose tool does not seem worthwhile.
567 * Because of the way gmx_prepare_tng_writing is implemented, the case
568 * where the input TNG file has no molecule information will never
569 * lead to an output TNG file having molecule information. Since
570 * molecule information will generally be present if the input TNG
571 * file was written by a GROMACS tool, this seems like reasonable
573 static std::unique_ptr
<gmx_mtop_t
>
574 read_mtop_for_tng(const char *tps_file
,
575 const char *input_file
,
576 const char *output_file
)
578 std::unique_ptr
<gmx_mtop_t
> mtop
;
580 if (fn2bTPX(tps_file
) &&
581 efTNG
!= fn2ftp(input_file
) &&
582 efTNG
== fn2ftp(output_file
))
584 int temp_natoms
= -1;
585 mtop
= gmx::compat::make_unique
<gmx_mtop_t
>();
586 read_tpx(tps_file
, nullptr, nullptr, &temp_natoms
,
587 nullptr, nullptr, mtop
.get());
593 int gmx_trjconv(int argc
, char *argv
[])
595 const char *desc
[] = {
596 "[THISMODULE] can convert trajectory files in many ways:",
598 "* from one format to another",
599 "* select a subset of atoms",
600 "* change the periodicity representation",
601 "* keep multimeric molecules together",
602 "* center atoms in the box",
603 "* fit atoms to reference structure",
604 "* reduce the number of frames",
605 "* change the timestamps of the frames ([TT]-t0[tt] and [TT]-timestep[tt])",
606 "* cut the trajectory in small subtrajectories according",
607 " to information in an index file. This allows subsequent analysis of",
608 " the subtrajectories that could, for example, be the result of a",
609 " cluster analysis. Use option [TT]-sub[tt].",
610 " This assumes that the entries in the index file are frame numbers and",
611 " dumps each group in the index file to a separate trajectory file.",
612 "* select frames within a certain range of a quantity given",
613 " in an [REF].xvg[ref] file.",
615 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
618 "The following formats are supported for input and output:",
619 "[REF].xtc[ref], [REF].trr[ref], [REF].gro[ref], [TT].g96[tt]",
620 "and [REF].pdb[ref].",
621 "The file formats are detected from the file extension.",
622 "The precision of the [REF].xtc[ref] output is taken from the",
623 "input file for [REF].xtc[ref], [REF].gro[ref] and [REF].pdb[ref],",
624 "and from the [TT]-ndec[tt] option for other input formats. The precision",
625 "is always taken from [TT]-ndec[tt], when this option is set.",
626 "All other formats have fixed precision. [REF].trr[ref]",
627 "output can be single or double precision, depending on the precision",
628 "of the [THISMODULE] binary.",
629 "Note that velocities are only supported in",
630 "[REF].trr[ref], [REF].gro[ref] and [TT].g96[tt] files.[PAR]",
632 "Option [TT]-sep[tt] can be used to write every frame to a separate",
633 "[TT].gro, .g96[tt] or [REF].pdb[ref] file. By default, all frames all written to one file.",
634 "[REF].pdb[ref] files with all frames concatenated can be viewed with",
635 "[TT]rasmol -nmrpdb[tt].[PAR]",
637 "It is possible to select part of your trajectory and write it out",
638 "to a new trajectory file in order to save disk space, e.g. for leaving",
639 "out the water from a trajectory of a protein in water.",
640 "[BB]ALWAYS[bb] put the original trajectory on tape!",
641 "We recommend to use the portable [REF].xtc[ref] format for your analysis",
642 "to save disk space and to have portable files.[PAR]",
644 "There are two options for fitting the trajectory to a reference",
645 "either for essential dynamics analysis, etc.",
646 "The first option is just plain fitting to a reference structure",
647 "in the structure file. The second option is a progressive fit",
648 "in which the first timeframe is fitted to the reference structure ",
649 "in the structure file to obtain and each subsequent timeframe is ",
650 "fitted to the previously fitted structure. This way a continuous",
651 "trajectory is generated, which might not be the case when using the",
652 "regular fit method, e.g. when your protein undergoes large",
653 "conformational transitions.[PAR]",
655 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
658 " * [TT]mol[tt] puts the center of mass of molecules in the box,",
659 " and requires a run input file to be supplied with [TT]-s[tt].",
660 " * [TT]res[tt] puts the center of mass of residues in the box.",
661 " * [TT]atom[tt] puts all the atoms in the box.",
662 " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
663 " them back. This has the effect that all molecules",
664 " will remain whole (provided they were whole in the initial",
665 " conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
666 " molecules may diffuse out of the box. The starting configuration",
667 " for this procedure is taken from the structure file, if one is",
668 " supplied, otherwise it is the first frame.",
669 " * [TT]cluster[tt] clusters all the atoms in the selected index",
670 " such that they are all closest to the center of mass of the cluster,",
671 " which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
672 " results if you in fact have a cluster. Luckily that can be checked",
673 " afterwards using a trajectory viewer. Note also that if your molecules",
674 " are broken this will not work either.",
675 " * [TT]whole[tt] only makes broken molecules whole.",
678 "Option [TT]-ur[tt] sets the unit cell representation for options",
679 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
680 "All three options give different results for triclinic boxes and",
681 "identical results for rectangular boxes.",
682 "[TT]rect[tt] is the ordinary brick shape.",
683 "[TT]tric[tt] is the triclinic unit cell.",
684 "[TT]compact[tt] puts all atoms at the closest distance from the center",
685 "of the box. This can be useful for visualizing e.g. truncated octahedra",
686 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
687 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
688 "is set differently.[PAR]",
690 "Option [TT]-center[tt] centers the system in the box. The user can",
691 "select the group which is used to determine the geometrical center.",
692 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
693 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
694 "[TT]tric[tt]: half of the sum of the box vectors,",
695 "[TT]rect[tt]: half of the box diagonal,",
696 "[TT]zero[tt]: zero.",
697 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
698 "want all molecules in the box after the centering.[PAR]",
700 "Option [TT]-box[tt] sets the size of the new box. This option only works",
701 "for leading dimensions and is thus generally only useful for rectangular boxes.",
702 "If you want to modify only some of the dimensions, e.g. when reading from",
703 "a trajectory, you can use -1 for those dimensions that should stay the same",
705 "It is not always possible to use combinations of [TT]-pbc[tt],",
706 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
707 "you want in one call to [THISMODULE]. Consider using multiple",
708 "calls, and check out the GROMACS website for suggestions.[PAR]",
710 "With [TT]-dt[tt], it is possible to reduce the number of ",
711 "frames in the output. This option relies on the accuracy of the times",
712 "in your input trajectory, so if these are inaccurate use the",
713 "[TT]-timestep[tt] option to modify the time (this can be done",
714 "simultaneously). For making smooth movies, the program [gmx-filter]",
715 "can reduce the number of frames while using low-pass frequency",
716 "filtering, this reduces aliasing of high frequency motions.[PAR]",
718 "Using [TT]-trunc[tt] [THISMODULE] can truncate [REF].trr[ref] in place, i.e.",
719 "without copying the file. This is useful when a run has crashed",
720 "during disk I/O (i.e. full disk), or when two contiguous",
721 "trajectories must be concatenated without having double frames.[PAR]",
723 "Option [TT]-dump[tt] can be used to extract a frame at or near",
724 "one specific time from your trajectory, but only works reliably",
725 "if the time interval between frames is uniform.[PAR]",
727 "Option [TT]-drop[tt] reads an [REF].xvg[ref] file with times and values.",
728 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
729 "frames with a value below and above the value of the respective options",
730 "will not be written."
746 const char *pbc_opt
[epNR
+ 1] =
748 nullptr, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
753 const char *unitcell_opt
[euNR
+1] =
754 { nullptr, "rect", "tric", "compact", nullptr };
758 ecSel
, ecTric
, ecRect
, ecZero
, ecNR
760 const char *center_opt
[ecNR
+1] =
761 { nullptr, "tric", "rect", "zero", nullptr };
767 efSel
, efNone
, efFit
, efFitXY
, efReset
, efResetXY
, efPFit
, efNR
769 const char *fit
[efNR
+ 1] =
771 nullptr, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
772 "progressive", nullptr
775 static gmx_bool bSeparate
= FALSE
, bVels
= TRUE
, bForce
= FALSE
, bCONECT
= FALSE
;
776 static gmx_bool bCenter
= FALSE
;
777 static int skip_nr
= 1, ndec
= 3, nzero
= 0;
778 static real tzero
= 0, delta_t
= 0, timestep
= 0, ttrunc
= -1, tdump
= -1, split_t
= 0;
779 static rvec newbox
= {0, 0, 0}, shift
= {0, 0, 0}, trans
= {0, 0, 0};
780 static char *exec_command
= nullptr;
781 static real dropunder
= 0, dropover
= 0;
782 static gmx_bool bRound
= FALSE
;
787 { "-skip", FALSE
, etINT
,
788 { &skip_nr
}, "Only write every nr-th frame" },
789 { "-dt", FALSE
, etTIME
,
791 "Only write frame when t MOD dt = first time (%t)" },
792 { "-round", FALSE
, etBOOL
,
793 { &bRound
}, "Round measurements to nearest picosecond"},
794 { "-dump", FALSE
, etTIME
,
795 { &tdump
}, "Dump frame nearest specified time (%t)" },
796 { "-t0", FALSE
, etTIME
,
798 "Starting time (%t) (default: don't change)" },
799 { "-timestep", FALSE
, etTIME
,
801 "Change time step between input frames (%t)" },
802 { "-pbc", FALSE
, etENUM
,
804 "PBC treatment (see help text for full description)" },
805 { "-ur", FALSE
, etENUM
,
806 { unitcell_opt
}, "Unit-cell representation" },
807 { "-center", FALSE
, etBOOL
,
808 { &bCenter
}, "Center atoms in box" },
809 { "-boxcenter", FALSE
, etENUM
,
810 { center_opt
}, "Center for -pbc and -center" },
811 { "-box", FALSE
, etRVEC
,
813 "Size for new cubic box (default: read from input)" },
814 { "-trans", FALSE
, etRVEC
,
816 "All coordinates will be translated by trans. This "
817 "can advantageously be combined with -pbc mol -ur "
819 { "-shift", FALSE
, etRVEC
,
821 "All coordinates will be shifted by framenr*shift" },
822 { "-fit", FALSE
, etENUM
,
824 "Fit molecule to ref structure in the structure file" },
825 { "-ndec", FALSE
, etINT
,
827 "Precision for .xtc and .gro writing in number of "
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
;
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
, bSetPrec
, 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], 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
,
931 top_file
= ftp2fn(efTPS
, NFILE
, fnm
);
933 /* Check command line */
934 in_file
= opt2fn("-f", NFILE
, fnm
);
938 do_trunc(in_file
, ttrunc
);
942 /* mark active cmdline options */
943 bSetBox
= opt2parg_bSet("-box", NPA
, pa
);
944 bSetTime
= opt2parg_bSet("-t0", NPA
, pa
);
945 bSetPrec
= opt2parg_bSet("-ndec", NPA
, pa
);
946 bSetUR
= opt2parg_bSet("-ur", NPA
, pa
);
947 bExec
= opt2parg_bSet("-exec", NPA
, pa
);
948 bTimeStep
= opt2parg_bSet("-timestep", NPA
, pa
);
949 bTDump
= opt2parg_bSet("-dump", NPA
, pa
);
950 bDropUnder
= opt2parg_bSet("-dropunder", NPA
, pa
);
951 bDropOver
= opt2parg_bSet("-dropover", NPA
, pa
);
952 bTrans
= opt2parg_bSet("-trans", NPA
, pa
);
953 bSplit
= (split_t
!= 0);
955 /* parse enum options */
956 fit_enum
= nenum(fit
);
957 bFit
= (fit_enum
== efFit
|| fit_enum
== efFitXY
);
958 bReset
= (fit_enum
== efReset
|| fit_enum
== efResetXY
);
959 bPFit
= fit_enum
== efPFit
;
960 pbc_enum
= nenum(pbc_opt
);
961 bPBCWhole
= pbc_enum
== epWhole
;
962 bPBCcomRes
= pbc_enum
== epComRes
;
963 bPBCcomMol
= pbc_enum
== epComMol
;
964 bPBCcomAtom
= pbc_enum
== epComAtom
;
965 bNoJump
= pbc_enum
== epNojump
;
966 bCluster
= pbc_enum
== epCluster
;
967 bPBC
= pbc_enum
!= epNone
;
968 unitcell_enum
= nenum(unitcell_opt
);
969 ecenter
= nenum(center_opt
) - ecTric
;
971 /* set and check option dependencies */
974 bFit
= TRUE
; /* for pfit, fit *must* be set */
978 bReset
= TRUE
; /* for fit, reset *must* be set */
983 nfitdim
= (fit_enum
== efFitXY
|| fit_enum
== efResetXY
) ? 2 : 3;
985 bRmPBC
= bFit
|| bPBCWhole
|| bPBCcomRes
|| bPBCcomMol
;
989 if (!(bPBCcomRes
|| bPBCcomMol
|| bPBCcomAtom
))
992 "WARNING: Option for unitcell representation (-ur %s)\n"
993 " only has effect in combination with -pbc %s, %s or %s.\n"
994 " Ingoring unitcell representation.\n\n",
995 unitcell_opt
[0], pbc_opt
[2], pbc_opt
[3], pbc_opt
[4]);
1000 gmx_fatal(FARGS
, "PBC condition treatment does not work together with rotational fit.\n"
1001 "Please do the PBC condition treatment first and then run trjconv in a second step\n"
1002 "for the rotational fit.\n"
1003 "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1007 /* ndec is in nr of decimal places, prec is a multiplication factor: */
1009 for (i
= 0; i
< ndec
; i
++)
1014 bIndex
= ftp2bSet(efNDX
, NFILE
, fnm
);
1017 /* Determine output type */
1018 out_file
= opt2fn("-o", NFILE
, fnm
);
1019 int ftp
= fn2ftp(out_file
);
1020 fprintf(stderr
, "Will write %s: %s\n", ftp2ext(ftp
), ftp2desc(ftp
));
1021 bNeedPrec
= (ftp
== efXTC
|| ftp
== efGRO
);
1022 int ftpin
= fn2ftp(in_file
);
1025 /* check if velocities are possible in input and output files */
1026 bVels
= (ftp
== efTRR
|| ftp
== efGRO
||
1027 ftp
== efG96
|| ftp
== efTNG
)
1028 && (ftpin
== efTRR
|| ftpin
== efGRO
||
1029 ftpin
== efG96
|| ftpin
== efTNG
|| ftpin
== efCPT
);
1031 if (bSeparate
|| bSplit
)
1033 outf_ext
= std::strrchr(out_file
, '.');
1034 if (outf_ext
== nullptr)
1036 gmx_fatal(FARGS
, "Output file name '%s' does not contain a '.'", out_file
);
1038 outf_base
= gmx_strdup(out_file
);
1039 outf_base
[outf_ext
- out_file
] = '\0';
1042 bSubTraj
= opt2bSet("-sub", NFILE
, fnm
);
1045 if ((ftp
!= efXTC
) && (ftp
!= efTRR
))
1047 /* It seems likely that other trajectory file types
1048 * could work here. */
1049 gmx_fatal(FARGS
, "Can only use the sub option with output file types "
1052 clust
= cluster_index(nullptr, opt2fn("-sub", NFILE
, fnm
));
1054 /* Check for number of files disabled, as FOPEN_MAX is not the correct
1055 * number to check for. In my linux box it is only 16.
1057 if (/* DISABLES CODE */ (0) && (clust
->clust
->nr
> FOPEN_MAX
-4))
1059 gmx_fatal(FARGS
, "Can not open enough (%d) files to write all the"
1060 " trajectories.\ntry splitting the index file in %d parts.\n"
1062 clust
->clust
->nr
, 1+clust
->clust
->nr
/FOPEN_MAX
, FOPEN_MAX
);
1064 gmx_warning("The -sub option could require as many open output files as there are\n"
1065 "index groups in the file (%d). If you get I/O errors opening new files,\n"
1066 "try reducing the number of index groups in the file, and perhaps\n"
1067 "using trjconv -sub several times on different chunks of your index file.\n",
1070 snew(clust_status
, clust
->clust
->nr
);
1071 snew(clust_status_id
, clust
->clust
->nr
);
1072 snew(nfwritten
, clust
->clust
->nr
);
1073 for (i
= 0; (i
< clust
->clust
->nr
); i
++)
1075 clust_status
[i
] = nullptr;
1076 clust_status_id
[i
] = -1;
1078 bSeparate
= bSplit
= FALSE
;
1083 gmx_fatal(FARGS
, "Argument for -skip (%d) needs to be greater or equal to 1.", skip_nr
);
1086 std::unique_ptr
<gmx_mtop_t
> mtop
= read_mtop_for_tng(top_file
, in_file
, out_file
);
1088 /* Determine whether to read a topology */
1089 bTPS
= (ftp2bSet(efTPS
, NFILE
, fnm
) ||
1090 bRmPBC
|| bReset
|| bPBCcomMol
|| bCluster
||
1091 (ftp
== efGRO
) || (ftp
== efPDB
) || bCONECT
);
1093 /* Determine if when can read index groups */
1094 bIndex
= (bIndex
|| bTPS
);
1098 read_tps_conf(top_file
, &top
, &ePBC
, &xp
, nullptr, top_box
,
1099 bReset
|| bPBCcomRes
);
1100 std::strncpy(top_title
, *top
.name
, 255);
1101 top_title
[255] = '\0';
1104 if (0 == top
.mols
.nr
&& (bCluster
|| bPBCcomMol
))
1106 gmx_fatal(FARGS
, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt
[pbc_enum
]);
1109 /* top_title is only used for gro and pdb,
1110 * the header in such a file is top_title, followed by
1111 * t= ... and/or step= ...
1112 * to prevent double t= or step=, remove it from top_title.
1113 * From GROMACS-2018 we only write t/step when the frame actually
1114 * has a valid time/step, so we need to check for both separately.
1116 if ((charpt
= std::strstr(top_title
, " t= ")))
1120 if ((charpt
= std::strstr(top_title
, " step= ")))
1127 gc
= gmx_conect_generate(&top
);
1131 gpbc
= gmx_rmpbc_init(&top
.idef
, ePBC
, top
.atoms
.nr
);
1135 /* get frame number index */
1137 if (opt2bSet("-fr", NFILE
, fnm
))
1139 printf("Select groups of frame number indices:\n");
1140 rd_index(opt2fn("-fr", NFILE
, fnm
), 1, &nrfri
, (int **)&frindex
, &frname
);
1143 for (i
= 0; i
< nrfri
; i
++)
1145 fprintf(debug
, "frindex[%4d]=%4d\n", i
, frindex
[i
]);
1150 /* get index groups etc. */
1153 printf("Select group for %s fit\n",
1154 bFit
? "least squares" : "translational");
1155 get_index(atoms
, ftp2fn_null(efNDX
, NFILE
, fnm
),
1156 1, &ifit
, &ind_fit
, &gn_fit
);
1162 gmx_fatal(FARGS
, "Need at least 2 atoms to fit!\n");
1166 fprintf(stderr
, "WARNING: fitting with only 2 atoms is not unique\n");
1172 printf("Select group for clustering\n");
1173 get_index(atoms
, ftp2fn_null(efNDX
, NFILE
, fnm
),
1174 1, &ifit
, &ind_fit
, &gn_fit
);
1181 printf("Select group for centering\n");
1182 get_index(atoms
, ftp2fn_null(efNDX
, NFILE
, fnm
),
1183 1, &ncent
, &cindex
, &grpnm
);
1185 printf("Select group for output\n");
1186 get_index(atoms
, ftp2fn_null(efNDX
, NFILE
, fnm
),
1187 1, &nout
, &index
, &grpnm
);
1191 /* no index file, so read natoms from TRX */
1192 if (!read_first_frame(oenv
, &trxin
, in_file
, &fr
, TRX_DONT_SKIP
))
1194 gmx_fatal(FARGS
, "Could not read a frame from %s", in_file
);
1199 snew(index
, natoms
);
1200 for (i
= 0; i
< natoms
; i
++)
1214 snew(w_rls
, atoms
->nr
);
1215 for (i
= 0; (i
< ifit
); i
++)
1217 w_rls
[ind_fit
[i
]] = atoms
->atom
[ind_fit
[i
]].m
;
1220 /* Restore reference structure and set to origin,
1221 store original location (to put structure back) */
1224 gmx_rmpbc(gpbc
, top
.atoms
.nr
, top_box
, xp
);
1226 copy_rvec(xp
[index
[0]], x_shift
);
1227 reset_x_ndim(nfitdim
, ifit
, ind_fit
, atoms
->nr
, nullptr, xp
, w_rls
);
1228 rvec_dec(x_shift
, xp
[index
[0]]);
1232 clear_rvec(x_shift
);
1235 if (bDropUnder
|| bDropOver
)
1237 /* Read the .xvg file with the drop values */
1238 fprintf(stderr
, "\nReading drop file ...");
1239 ndrop
= read_xvg(opt2fn("-drop", NFILE
, fnm
), &dropval
, &ncol
);
1240 fprintf(stderr
, " %d time points\n", ndrop
);
1241 if (ndrop
== 0 || ncol
< 2)
1243 gmx_fatal(FARGS
, "Found no data points in %s",
1244 opt2fn("-drop", NFILE
, fnm
));
1250 /* Make atoms struct for output in GRO or PDB files */
1251 if ((ftp
== efGRO
) || ((ftp
== efG96
) && bTPS
) || (ftp
== efPDB
))
1253 /* get memory for stuff to go in .pdb file, and initialize
1254 * the pdbinfo structure part if the input has it.
1256 init_t_atoms(&useatoms
, atoms
->nr
, atoms
->havePdbInfo
);
1257 sfree(useatoms
.resinfo
);
1258 useatoms
.resinfo
= atoms
->resinfo
;
1259 for (i
= 0; (i
< nout
); i
++)
1261 useatoms
.atomname
[i
] = atoms
->atomname
[index
[i
]];
1262 useatoms
.atom
[i
] = atoms
->atom
[index
[i
]];
1263 if (atoms
->havePdbInfo
)
1265 useatoms
.pdbinfo
[i
] = atoms
->pdbinfo
[index
[i
]];
1267 useatoms
.nres
= std::max(useatoms
.nres
, useatoms
.atom
[i
].resind
+1);
1271 /* select what to read */
1282 flags
= flags
| TRX_READ_V
;
1286 flags
= flags
| TRX_READ_F
;
1289 /* open trx file for reading */
1290 bHaveFirstFrame
= read_first_frame(oenv
, &trxin
, in_file
, &fr
, flags
);
1293 fprintf(stderr
, "\nPrecision of %s is %g (nm)\n", in_file
, 1/fr
.prec
);
1297 if (bSetPrec
|| !fr
.bPrec
)
1299 fprintf(stderr
, "\nSetting output precision to %g (nm)\n", 1/prec
);
1303 fprintf(stderr
, "Using output precision of %g (nm)\n", 1/prec
);
1307 if (bHaveFirstFrame
)
1311 // Determine timestep (assuming constant spacing for now) if we
1312 // need to dump frames based on time. This is required so we do not
1313 // skip the first frame if that was the one that should have been dumped
1314 double firstFrameTime
= fr
.time
;
1315 if (read_next_frame(oenv
, trxin
, &fr
))
1317 dt
= fr
.time
- firstFrameTime
;
1321 fprintf(stderr
, "Warning: Frame times are not incrementing - will dump first frame.\n");
1324 // Now close and reopen so we are at first frame again
1327 // Reopen at first frame (We already know it exists if we got here)
1328 read_first_frame(oenv
, &trxin
, in_file
, &fr
, flags
);
1331 set_trxframe_ePBC(&fr
, ePBC
);
1336 tshift
= tzero
-fr
.time
;
1346 /* check if index is meaningful */
1347 for (i
= 0; i
< nout
; i
++)
1349 if (index
[i
] >= natoms
)
1352 "Index[%d] %d is larger than the number of atoms in the\n"
1353 "trajectory file (%d). There is a mismatch in the contents\n"
1354 "of your -f, -s and/or -n files.", i
, index
[i
]+1, natoms
);
1356 bCopy
= bCopy
|| (i
!= index
[i
]);
1360 /* open output for writing */
1361 std::strcpy(filemode
, "w");
1365 trxout
= trjtools_gmx_prepare_tng_writing(out_file
,
1377 if (!bSplit
&& !bSubTraj
)
1379 trxout
= open_trx(out_file
, filemode
);
1385 if (( !bSeparate
&& !bSplit
) && !bSubTraj
)
1387 out
= gmx_ffopen(out_file
, filemode
);
1391 gmx_incons("Illegal output file format");
1407 /* Start the big loop over frames */
1413 /* Main loop over frames */
1424 /*if (frame >= clust->clust->nra)
1425 gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1426 if (frame
> clust
->maxframe
)
1432 my_clust
= clust
->inv_clust
[frame
];
1434 if ((my_clust
< 0) || (my_clust
>= clust
->clust
->nr
) ||
1443 /* generate new box */
1444 if (fr
.bBox
== FALSE
)
1448 for (m
= 0; m
< DIM
; m
++)
1452 fr
.box
[m
][m
] = newbox
[m
];
1456 if (fr
.bBox
== FALSE
)
1458 gmx_fatal(FARGS
, "Cannot preserve a box that does not exist.\n");
1466 for (i
= 0; i
< natoms
; i
++)
1468 rvec_inc(fr
.x
[i
], trans
);
1474 // If we could not read two frames or times are not incrementing
1475 // we have almost no idea what to do,
1476 // but dump the first frame so output is not broken.
1477 if (dt
<= 0 || !bDTset
)
1483 // Dump the frame if we are less than half a frame time
1484 // below it. This will also ensure we at least dump a
1485 // somewhat reasonable frame if the spacing is unequal
1486 // and we have overrun the frame time. Once we dump one
1487 // frame based on time we quit, so it does not matter
1488 // that this might be true for all subsequent frames too.
1489 bDumpFrame
= (fr
.time
> tdump
-0.5*dt
);
1497 /* determine if an atom jumped across the box and reset it if so */
1498 if (bNoJump
&& (bTPS
|| frame
!= 0))
1500 for (d
= 0; d
< DIM
; d
++)
1502 hbox
[d
] = 0.5*fr
.box
[d
][d
];
1504 for (i
= 0; i
< natoms
; i
++)
1508 rvec_dec(fr
.x
[i
], x_shift
);
1510 for (m
= DIM
-1; m
>= 0; m
--)
1514 while (fr
.x
[i
][m
]-xp
[i
][m
] <= -hbox
[m
])
1516 for (d
= 0; d
<= m
; d
++)
1518 fr
.x
[i
][d
] += fr
.box
[m
][d
];
1521 while (fr
.x
[i
][m
]-xp
[i
][m
] > hbox
[m
])
1523 for (d
= 0; d
<= m
; d
++)
1525 fr
.x
[i
][d
] -= fr
.box
[m
][d
];
1534 calc_pbc_cluster(ecenter
, ifit
, &top
, ePBC
, fr
.x
, ind_fit
, fr
.box
);
1539 /* Now modify the coords according to the flags,
1540 for normal fit, this is only done for output frames */
1543 gmx_rmpbc_trxfr(gpbc
, &fr
);
1546 reset_x_ndim(nfitdim
, ifit
, ind_fit
, natoms
, nullptr, fr
.x
, w_rls
);
1547 do_fit(natoms
, w_rls
, xp
, fr
.x
);
1550 /* store this set of coordinates for future use */
1551 if (bPFit
|| bNoJump
)
1557 for (i
= 0; (i
< natoms
); i
++)
1559 copy_rvec(fr
.x
[i
], xp
[i
]);
1560 rvec_inc(fr
.x
[i
], x_shift
);
1566 /* see if we have a frame from the frame index group */
1567 for (i
= 0; i
< nrfri
&& !bDumpFrame
; i
++)
1569 bDumpFrame
= frame
== frindex
[i
];
1572 if (debug
&& bDumpFrame
)
1574 fprintf(debug
, "dumping %d\n", frame
);
1578 ( ( !bTDump
&& !frindex
&& frame
% skip_nr
== 0 ) || bDumpFrame
);
1580 if (bWriteFrame
&& (bDropUnder
|| bDropOver
))
1582 while (dropval
[0][drop1
] < fr
.time
&& drop1
+1 < ndrop
)
1587 if (std::abs(dropval
[0][drop0
] - fr
.time
)
1588 < std::abs(dropval
[0][drop1
] - fr
.time
))
1596 if ((bDropUnder
&& dropval
[1][dropuse
] < dropunder
) ||
1597 (bDropOver
&& dropval
[1][dropuse
] > dropover
))
1599 bWriteFrame
= FALSE
;
1605 /* We should avoid modifying the input frame,
1606 * but since here we don't have the output frame yet,
1607 * we introduce a temporary output frame time variable.
1611 frout_time
= fr
.time
;
1616 frout_time
= tzero
+ frame
*timestep
;
1620 frout_time
+= tshift
;
1625 fprintf(stderr
, "\nDumping frame at t= %g %s\n",
1626 output_env_conv_time(oenv
, frout_time
), output_env_get_time_unit(oenv
).c_str());
1629 /* check for writing at each delta_t */
1630 bDoIt
= (delta_t
== 0);
1635 bDoIt
= bRmod(frout_time
, tzero
, delta_t
);
1639 /* round() is not C89 compatible, so we do this: */
1640 bDoIt
= bRmod(std::floor(frout_time
+0.5), std::floor(tzero
+0.5),
1641 std::floor(delta_t
+0.5));
1645 if (bDoIt
|| bTDump
)
1647 /* print sometimes */
1648 if ( ((outframe
% SKIP
) == 0) || (outframe
< SKIP
) )
1650 fprintf(stderr
, " -> frame %6d time %8.3f \r",
1651 outframe
, output_env_conv_time(oenv
, frout_time
));
1657 /* Now modify the coords according to the flags,
1658 for PFit we did this already! */
1662 gmx_rmpbc_trxfr(gpbc
, &fr
);
1667 reset_x_ndim(nfitdim
, ifit
, ind_fit
, natoms
, nullptr, fr
.x
, w_rls
);
1670 do_fit_ndim(nfitdim
, natoms
, w_rls
, xp
, fr
.x
);
1674 for (i
= 0; i
< natoms
; i
++)
1676 rvec_inc(fr
.x
[i
], x_shift
);
1683 center_x(ecenter
, fr
.x
, fr
.box
, natoms
, ncent
, cindex
);
1687 auto positionsArrayRef
= gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec
*>(fr
.x
), natoms
);
1690 switch (unitcell_enum
)
1693 put_atoms_in_box(ePBC
, fr
.box
, positionsArrayRef
);
1696 put_atoms_in_triclinic_unitcell(ecenter
, fr
.box
, positionsArrayRef
);
1699 put_atoms_in_compact_unitcell(ePBC
, ecenter
, fr
.box
,
1706 put_residue_com_in_box(unitcell_enum
, ecenter
,
1707 natoms
, atoms
->atom
, ePBC
, fr
.box
, fr
.x
);
1711 put_molecule_com_in_box(unitcell_enum
, ecenter
,
1713 natoms
, atoms
->atom
, ePBC
, fr
.box
, fr
.x
);
1715 /* Copy the input trxframe struct to the output trxframe struct */
1717 frout
.time
= frout_time
;
1718 frout
.bV
= (frout
.bV
&& bVels
);
1719 frout
.bF
= (frout
.bF
&& bForce
);
1720 frout
.natoms
= nout
;
1721 if (bNeedPrec
&& (bSetPrec
|| !fr
.bPrec
))
1737 for (i
= 0; i
< nout
; i
++)
1739 copy_rvec(fr
.x
[index
[i
]], frout
.x
[i
]);
1742 copy_rvec(fr
.v
[index
[i
]], frout
.v
[i
]);
1746 copy_rvec(fr
.f
[index
[i
]], frout
.f
[i
]);
1751 if (opt2parg_bSet("-shift", NPA
, pa
))
1753 for (i
= 0; i
< nout
; i
++)
1755 for (d
= 0; d
< DIM
; d
++)
1757 frout
.x
[i
][d
] += outframe
*shift
[d
];
1764 bSplitHere
= bSplit
&& bRmod(frout
.time
, tzero
, split_t
);
1768 /* round() is not C89 compatible, so we do this: */
1769 bSplitHere
= bSplit
&& bRmod(std::floor(frout
.time
+0.5),
1770 std::floor(tzero
+0.5),
1771 std::floor(split_t
+0.5));
1773 if (bSeparate
|| bSplitHere
)
1775 mk_filenm(outf_base
, ftp2ext(ftp
), nzero
, file_nr
, out_file2
);
1781 write_tng_frame(trxout
, &frout
);
1782 // TODO when trjconv behaves better: work how to read and write lambda
1792 trxout
= open_trx(out_file2
, filemode
);
1799 if (clust_status_id
[my_clust
] == -1)
1801 sprintf(buf
, "%s.%s", clust
->grpname
[my_clust
], ftp2ext(ftp
));
1802 clust_status
[my_clust
] = open_trx(buf
, "w");
1803 clust_status_id
[my_clust
] = 1;
1806 else if (clust_status_id
[my_clust
] == -2)
1808 gmx_fatal(FARGS
, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1809 clust
->grpname
[my_clust
], ntrxopen
, frame
,
1812 write_trxframe(clust_status
[my_clust
], &frout
, gc
);
1813 nfwritten
[my_clust
]++;
1814 if (nfwritten
[my_clust
] ==
1815 (clust
->clust
->index
[my_clust
+1]-
1816 clust
->clust
->index
[my_clust
]))
1818 close_trx(clust_status
[my_clust
]);
1819 clust_status
[my_clust
] = nullptr;
1820 clust_status_id
[my_clust
] = -2;
1824 gmx_fatal(FARGS
, "Less than zero open .xtc files!");
1831 write_trxframe(trxout
, &frout
, gc
);
1837 // Only add a generator statement if title is empty,
1838 // to avoid multiple generated-by statements from various programs
1839 if (std::strlen(top_title
) == 0)
1841 sprintf(top_title
, "Generated by trjconv");
1845 sprintf(timestr
, " t= %9.5f", frout
.time
);
1849 std::strcpy(timestr
, "");
1853 sprintf(stepstr
, " step= %" PRId64
, frout
.step
);
1857 std::strcpy(stepstr
, "");
1859 snprintf(title
, 256, "%s%s%s", top_title
, timestr
, stepstr
);
1860 if (bSeparate
|| bSplitHere
)
1862 out
= gmx_ffopen(out_file2
, "w");
1867 write_hconf_p(out
, title
, &useatoms
,
1868 frout
.x
, frout
.bV
? frout
.v
: nullptr, frout
.box
);
1871 fprintf(out
, "REMARK GENERATED BY TRJCONV\n");
1872 /* if reading from pdb, we want to keep the original
1873 model numbering else we write the output frame
1874 number plus one, because model 0 is not allowed in pdb */
1875 if (ftpin
== efPDB
&& fr
.bStep
&& fr
.step
> model_nr
)
1883 write_pdbfile(out
, title
, &useatoms
, frout
.x
,
1884 frout
.ePBC
, frout
.box
, ' ', model_nr
, gc
, TRUE
);
1887 const char *outputTitle
= "";
1888 if (bSeparate
|| bTDump
)
1890 outputTitle
= title
;
1893 frout
.bAtoms
= TRUE
;
1895 frout
.atoms
= &useatoms
;
1896 frout
.bStep
= FALSE
;
1897 frout
.bTime
= FALSE
;
1903 outputTitle
= title
;
1905 frout
.bAtoms
= FALSE
;
1909 write_g96_conf(out
, outputTitle
, &frout
, -1, nullptr);
1911 if (bSeparate
|| bSplitHere
)
1918 gmx_fatal(FARGS
, "DHE, ftp=%d\n", ftp
);
1920 if (bSeparate
|| bSplitHere
)
1925 /* execute command */
1929 sprintf(c
, "%s %d", exec_command
, file_nr
-1);
1930 /*fprintf(stderr,"Executing '%s'\n",c);*/
1933 gmx_fatal(FARGS
, "Error executing command: %s", c
);
1940 bHaveNextFrame
= read_next_frame(oenv
, trxin
, &fr
);
1942 while (!(bTDump
&& bDumpFrame
) && bHaveNextFrame
);
1945 if (!bHaveFirstFrame
|| (bTDump
&& !bDumpFrame
))
1947 fprintf(stderr
, "\nWARNING no output, "
1948 "last frame read at t=%g\n", fr
.time
);
1950 fprintf(stderr
, "\n");
1957 gmx_rmpbc_done(gpbc
);
1964 else if (out
!= nullptr)
1970 for (i
= 0; (i
< clust
->clust
->nr
); i
++)
1972 if (clust_status_id
[i
] >= 0)
1974 close_trx(clust_status
[i
]);
1993 do_view(oenv
, out_file
, nullptr);
1995 output_env_done(oenv
);