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, 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/fileio/confio.h"
48 #include "gromacs/fileio/g96io.h"
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/fileio/groio.h"
51 #include "gromacs/fileio/pdbio.h"
52 #include "gromacs/fileio/tngio_for_tools.h"
53 #include "gromacs/fileio/tpxio.h"
54 #include "gromacs/fileio/trrio.h"
55 #include "gromacs/fileio/trxio.h"
56 #include "gromacs/fileio/xtcio.h"
57 #include "gromacs/fileio/xvgr.h"
58 #include "gromacs/gmxana/gmx_ana.h"
59 #include "gromacs/math/do_fit.h"
60 #include "gromacs/math/functions.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/pbcutil/rmpbc.h"
65 #include "gromacs/topology/index.h"
66 #include "gromacs/topology/topology.h"
67 #include "gromacs/trajectory/trajectoryframe.h"
68 #include "gromacs/utility/arraysize.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/futil.h"
71 #include "gromacs/utility/smalloc.h"
74 euSel
, euRect
, euTric
, euCompact
, euNR
78 static void calc_pbc_cluster(int ecenter
, int nrefat
, t_topology
*top
, int ePBC
,
79 rvec x
[], int index
[], matrix box
)
81 int m
, i
, j
, j0
, j1
, jj
, ai
, aj
;
84 rvec dx
, xtest
, box_center
;
85 int nmol
, imol_center
;
87 gmx_bool
*bMol
, *bTmp
;
88 rvec
*m_com
, *m_shift
;
95 calc_box_center(ecenter
, box
, box_center
);
97 /* Initiate the pbc structure */
98 std::memset(&pbc
, 0, sizeof(pbc
));
99 set_pbc(&pbc
, ePBC
, box
);
101 /* Convert atom index to molecular */
103 molind
= top
->mols
.index
;
109 snew(bTmp
, top
->atoms
.nr
);
111 for (i
= 0; (i
< nrefat
); i
++)
113 /* Mark all molecules in the index */
116 /* Binary search assuming the molecules are sorted */
121 if (ai
< molind
[j0
+1])
125 else if (ai
>= molind
[j1
])
132 if (ai
< molind
[jj
+1])
144 /* Double check whether all atoms in all molecules that are marked are part
145 * of the cluster. Simultaneously compute the center of geometry.
147 min_dist2
= 10*gmx::square(trace(box
));
150 for (i
= 0; i
< nmol
; i
++)
152 for (j
= molind
[i
]; j
< molind
[i
+1]; j
++)
154 if (bMol
[i
] && !bTmp
[j
])
156 gmx_fatal(FARGS
, "Molecule %d marked for clustering but not atom %d in it - check your index!", i
+1, j
+1);
158 else if (!bMol
[i
] && bTmp
[j
])
160 gmx_fatal(FARGS
, "Atom %d marked for clustering but not molecule %d - this is an internal error...", j
+1, i
+1);
164 /* Make molecule whole, move 2nd and higher atom to same periodicity as 1st atom in molecule */
167 pbc_dx(&pbc
, x
[j
], x
[j
-1], dx
);
168 rvec_add(x
[j
-1], dx
, x
[j
]);
170 /* Compute center of geometry of molecule - m_com[i] was zeroed when we did snew() on it! */
171 rvec_inc(m_com
[i
], x
[j
]);
176 /* Normalize center of geometry */
177 fac
= 1.0/(molind
[i
+1]-molind
[i
]);
178 for (m
= 0; (m
< DIM
); m
++)
182 /* Determine which molecule is closest to the center of the box */
183 pbc_dx(&pbc
, box_center
, m_com
[i
], dx
);
184 tmp_r2
= iprod(dx
, dx
);
186 if (tmp_r2
< min_dist2
)
191 cluster
[ncluster
++] = i
;
198 fprintf(stderr
, "No molecules selected in the cluster\n");
201 else if (imol_center
== -1)
203 fprintf(stderr
, "No central molecules could be found\n");
208 added
[nadded
++] = imol_center
;
209 bMol
[imol_center
] = FALSE
;
211 while (nadded
< ncluster
)
213 /* Find min distance between cluster molecules and those remaining to be added */
214 min_dist2
= 10*gmx::square(trace(box
));
217 /* Loop over added mols */
218 for (i
= 0; i
< nadded
; i
++)
221 /* Loop over all mols */
222 for (j
= 0; j
< ncluster
; j
++)
225 /* check those remaining to be added */
228 pbc_dx(&pbc
, m_com
[aj
], m_com
[ai
], dx
);
229 tmp_r2
= iprod(dx
, dx
);
230 if (tmp_r2
< min_dist2
)
240 /* Add the best molecule */
241 added
[nadded
++] = jmin
;
243 /* Calculate the shift from the ai molecule */
244 pbc_dx(&pbc
, m_com
[jmin
], m_com
[imin
], dx
);
245 rvec_add(m_com
[imin
], dx
, xtest
);
246 rvec_sub(xtest
, m_com
[jmin
], m_shift
[jmin
]);
247 rvec_inc(m_com
[jmin
], m_shift
[jmin
]);
249 for (j
= molind
[jmin
]; j
< molind
[jmin
+1]; j
++)
251 rvec_inc(x
[j
], m_shift
[jmin
]);
253 fprintf(stdout
, "\rClustering iteration %d of %d...", nadded
, ncluster
);
263 fprintf(stdout
, "\n");
266 static void put_molecule_com_in_box(int unitcell_enum
, int ecenter
,
268 int natoms
, t_atom atom
[],
269 int ePBC
, matrix box
, rvec x
[])
273 rvec com
, new_com
, shift
, box_center
;
278 calc_box_center(ecenter
, box
, box_center
);
279 set_pbc(&pbc
, ePBC
, box
);
282 gmx_fatal(FARGS
, "There are no molecule descriptions. I need a .tpr file for this pbc option.");
284 for (i
= 0; (i
< mols
->nr
); i
++)
289 for (j
= mols
->index
[i
]; (j
< mols
->index
[i
+1] && j
< natoms
); j
++)
292 for (d
= 0; d
< DIM
; d
++)
298 /* calculate final COM */
299 svmul(1.0/mtot
, com
, com
);
301 /* check if COM is outside box */
302 copy_rvec(com
, new_com
);
303 switch (unitcell_enum
)
306 put_atoms_in_box(ePBC
, box
, 1, &new_com
);
309 put_atoms_in_triclinic_unitcell(ecenter
, box
, 1, &new_com
);
312 put_atoms_in_compact_unitcell(ePBC
, ecenter
, box
, 1, &new_com
);
315 rvec_sub(new_com
, com
, shift
);
316 if (norm2(shift
) > 0)
320 fprintf(debug
, "\nShifting position of molecule %d "
321 "by %8.3f %8.3f %8.3f\n", i
+1,
322 shift
[XX
], shift
[YY
], shift
[ZZ
]);
324 for (j
= mols
->index
[i
]; (j
< mols
->index
[i
+1] && j
< natoms
); j
++)
326 rvec_inc(x
[j
], shift
);
332 static void put_residue_com_in_box(int unitcell_enum
, int ecenter
,
333 int natoms
, t_atom atom
[],
334 int ePBC
, matrix box
, rvec x
[])
336 int i
, j
, res_start
, res_end
;
340 rvec box_center
, com
, new_com
, shift
;
341 static const int NOTSET
= -12347;
342 calc_box_center(ecenter
, box
, box_center
);
348 for (i
= 0; i
< natoms
+1; i
++)
350 if (i
== natoms
|| (presnr
!= atom
[i
].resind
&& presnr
!= NOTSET
))
352 /* calculate final COM */
354 svmul(1.0/mtot
, com
, com
);
356 /* check if COM is outside box */
357 copy_rvec(com
, new_com
);
358 switch (unitcell_enum
)
361 put_atoms_in_box(ePBC
, box
, 1, &new_com
);
364 put_atoms_in_triclinic_unitcell(ecenter
, box
, 1, &new_com
);
367 put_atoms_in_compact_unitcell(ePBC
, ecenter
, box
, 1, &new_com
);
370 rvec_sub(new_com
, com
, shift
);
375 fprintf(debug
, "\nShifting position of residue %d (atoms %d-%d) "
376 "by %g,%g,%g\n", atom
[res_start
].resind
+1,
377 res_start
+1, res_end
+1, shift
[XX
], shift
[YY
], shift
[ZZ
]);
379 for (j
= res_start
; j
< res_end
; j
++)
381 rvec_inc(x
[j
], shift
);
387 /* remember start of new residue */
394 for (d
= 0; d
< DIM
; d
++)
400 presnr
= atom
[i
].resind
;
405 static void center_x(int ecenter
, rvec x
[], matrix box
, int n
, int nc
, int ci
[])
408 rvec cmin
, cmax
, box_center
, dx
;
412 copy_rvec(x
[ci
[0]], cmin
);
413 copy_rvec(x
[ci
[0]], cmax
);
414 for (i
= 0; i
< nc
; i
++)
417 for (m
= 0; m
< DIM
; m
++)
419 if (x
[ai
][m
] < cmin
[m
])
423 else if (x
[ai
][m
] > cmax
[m
])
429 calc_box_center(ecenter
, box
, box_center
);
430 for (m
= 0; m
< DIM
; m
++)
432 dx
[m
] = box_center
[m
]-(cmin
[m
]+cmax
[m
])*0.5;
435 for (i
= 0; i
< n
; i
++)
442 static void mk_filenm(char *base
, const char *ext
, int ndigit
, int file_nr
,
448 std::strcpy(out_file
, base
);
459 std::strncat(out_file
, "00000000000", ndigit
-nd
);
461 sprintf(nbuf
, "%d.", file_nr
);
462 std::strcat(out_file
, nbuf
);
463 std::strcat(out_file
, ext
);
466 void check_trr(const char *fn
)
468 if (fn2ftp(fn
) != efTRR
)
470 gmx_fatal(FARGS
, "%s is not a trajectory file, exiting\n", fn
);
474 void do_trunc(const char *fn
, real t0
)
487 gmx_fatal(FARGS
, "You forgot to set the truncation time");
490 /* Check whether this is a .trr file */
493 in
= gmx_trr_open(fn
, "r");
494 fp
= gmx_fio_getfp(in
);
497 fprintf(stderr
, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn
);
503 fpos
= gmx_fio_ftell(in
);
505 while (!bStop
&& gmx_trr_read_frame_header(in
, &sh
, &bOK
))
507 gmx_trr_read_frame_data(in
, &sh
, NULL
, NULL
, NULL
, NULL
);
508 fpos
= gmx_ftell(fp
);
512 gmx_fseek(fp
, fpos
, SEEK_SET
);
518 fprintf(stderr
, "Do you REALLY want to truncate this trajectory (%s) at:\n"
519 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
520 fn
, j
, t
, (long int)fpos
);
521 if (1 != scanf("%s", yesno
))
523 gmx_fatal(FARGS
, "Error reading user input");
525 if (std::strcmp(yesno
, "YES") == 0)
527 fprintf(stderr
, "Once again, I'm gonna DO this...\n");
529 if (0 != gmx_truncate(fn
, fpos
))
531 gmx_fatal(FARGS
, "Error truncating file %s", fn
);
536 fprintf(stderr
, "Ok, I'll forget about it\n");
541 fprintf(stderr
, "Already at end of file (t=%g)...\n", t
);
547 /*! \brief Read a full molecular topology if useful and available.
549 * If the input trajectory file is not in TNG format, and the output
550 * file is in TNG format, then we want to try to read a full topology
551 * (if available), so that we can write molecule information to the
552 * output file. The full topology provides better molecule information
553 * than is available from the normal t_topology data used by GROMACS
556 * Also, the t_topology is only read under (different) particular
557 * conditions. If both apply, then a .tpr file might be read
558 * twice. Trying to fix this redundancy while trjconv is still an
559 * all-purpose tool does not seem worthwhile.
561 * Because of the way gmx_prepare_tng_writing is implemented, the case
562 * where the input TNG file has no molecule information will never
563 * lead to an output TNG file having molecule information. Since
564 * molecule information will generally be present if the input TNG
565 * file was written by a GROMACS tool, this seems like reasonable
567 static gmx_mtop_t
*read_mtop_for_tng(const char *tps_file
,
568 const char *input_file
,
569 const char *output_file
)
571 gmx_mtop_t
*mtop
= NULL
;
573 if (fn2bTPX(tps_file
) &&
574 efTNG
!= fn2ftp(input_file
) &&
575 efTNG
== fn2ftp(output_file
))
577 int temp_natoms
= -1;
579 read_tpx(tps_file
, NULL
, NULL
, &temp_natoms
,
586 int gmx_trjconv(int argc
, char *argv
[])
588 const char *desc
[] = {
589 "[THISMODULE] can convert trajectory files in many ways:",
591 "* from one format to another",
592 "* select a subset of atoms",
593 "* change the periodicity representation",
594 "* keep multimeric molecules together",
595 "* center atoms in the box",
596 "* fit atoms to reference structure",
597 "* reduce the number of frames",
598 "* change the timestamps of the frames ([TT]-t0[tt] and [TT]-timestep[tt])",
599 "* cut the trajectory in small subtrajectories according",
600 " to information in an index file. This allows subsequent analysis of",
601 " the subtrajectories that could, for example, be the result of a",
602 " cluster analysis. Use option [TT]-sub[tt].",
603 " This assumes that the entries in the index file are frame numbers and",
604 " dumps each group in the index file to a separate trajectory file.",
605 "* select frames within a certain range of a quantity given",
606 " in an [REF].xvg[ref] file.",
608 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
611 "The following formats are supported for input and output:",
612 "[REF].xtc[ref], [REF].trr[ref], [REF].gro[ref], [TT].g96[tt]",
613 "and [REF].pdb[ref].",
614 "The file formats are detected from the file extension.",
615 "The precision of [REF].xtc[ref] and [REF].gro[ref] output is taken from the",
616 "input file for [REF].xtc[ref], [REF].gro[ref] and [REF].pdb[ref],",
617 "and from the [TT]-ndec[tt] option for other input formats. The precision",
618 "is always taken from [TT]-ndec[tt], when this option is set.",
619 "All other formats have fixed precision. [REF].trr[ref]",
620 "output can be single or double precision, depending on the precision",
621 "of the [THISMODULE] binary.",
622 "Note that velocities are only supported in",
623 "[REF].trr[ref], [REF].gro[ref] and [TT].g96[tt] files.[PAR]",
625 "Option [TT]-sep[tt] can be used to write every frame to a separate",
626 "[TT].gro, .g96[tt] or [REF].pdb[ref] file. By default, all frames all written to one file.",
627 "[REF].pdb[ref] files with all frames concatenated can be viewed with",
628 "[TT]rasmol -nmrpdb[tt].[PAR]",
630 "It is possible to select part of your trajectory and write it out",
631 "to a new trajectory file in order to save disk space, e.g. for leaving",
632 "out the water from a trajectory of a protein in water.",
633 "[BB]ALWAYS[bb] put the original trajectory on tape!",
634 "We recommend to use the portable [REF].xtc[ref] format for your analysis",
635 "to save disk space and to have portable files.[PAR]",
637 "There are two options for fitting the trajectory to a reference",
638 "either for essential dynamics analysis, etc.",
639 "The first option is just plain fitting to a reference structure",
640 "in the structure file. The second option is a progressive fit",
641 "in which the first timeframe is fitted to the reference structure ",
642 "in the structure file to obtain and each subsequent timeframe is ",
643 "fitted to the previously fitted structure. This way a continuous",
644 "trajectory is generated, which might not be the case when using the",
645 "regular fit method, e.g. when your protein undergoes large",
646 "conformational transitions.[PAR]",
648 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
651 " * [TT]mol[tt] puts the center of mass of molecules in the box,",
652 " and requires a run input file to be supplied with [TT]-s[tt].",
653 " * [TT]res[tt] puts the center of mass of residues in the box.",
654 " * [TT]atom[tt] puts all the atoms in the box.",
655 " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
656 " them back. This has the effect that all molecules",
657 " will remain whole (provided they were whole in the initial",
658 " conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
659 " molecules may diffuse out of the box. The starting configuration",
660 " for this procedure is taken from the structure file, if one is",
661 " supplied, otherwise it is the first frame.",
662 " * [TT]cluster[tt] clusters all the atoms in the selected index",
663 " such that they are all closest to the center of mass of the cluster,",
664 " which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
665 " results if you in fact have a cluster. Luckily that can be checked",
666 " afterwards using a trajectory viewer. Note also that if your molecules",
667 " are broken this will not work either.",
668 " * [TT]whole[tt] only makes broken molecules whole.",
671 "Option [TT]-ur[tt] sets the unit cell representation for options",
672 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
673 "All three options give different results for triclinic boxes and",
674 "identical results for rectangular boxes.",
675 "[TT]rect[tt] is the ordinary brick shape.",
676 "[TT]tric[tt] is the triclinic unit cell.",
677 "[TT]compact[tt] puts all atoms at the closest distance from the center",
678 "of the box. This can be useful for visualizing e.g. truncated octahedra",
679 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
680 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
681 "is set differently.[PAR]",
683 "Option [TT]-center[tt] centers the system in the box. The user can",
684 "select the group which is used to determine the geometrical center.",
685 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
686 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
687 "[TT]tric[tt]: half of the sum of the box vectors,",
688 "[TT]rect[tt]: half of the box diagonal,",
689 "[TT]zero[tt]: zero.",
690 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
691 "want all molecules in the box after the centering.[PAR]",
693 "Option [TT]-box[tt] sets the size of the new box. This option only works",
694 "for leading dimensions and is thus generally only useful for rectangular boxes.",
695 "If you want to modify only some of the dimensions, e.g. when reading from",
696 "a trajectory, you can use -1 for those dimensions that should stay the same",
698 "It is not always possible to use combinations of [TT]-pbc[tt],",
699 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
700 "you want in one call to [THISMODULE]. Consider using multiple",
701 "calls, and check out the GROMACS website for suggestions.[PAR]",
703 "With [TT]-dt[tt], it is possible to reduce the number of ",
704 "frames in the output. This option relies on the accuracy of the times",
705 "in your input trajectory, so if these are inaccurate use the",
706 "[TT]-timestep[tt] option to modify the time (this can be done",
707 "simultaneously). For making smooth movies, the program [gmx-filter]",
708 "can reduce the number of frames while using low-pass frequency",
709 "filtering, this reduces aliasing of high frequency motions.[PAR]",
711 "Using [TT]-trunc[tt] [THISMODULE] can truncate [REF].trr[ref] in place, i.e.",
712 "without copying the file. This is useful when a run has crashed",
713 "during disk I/O (i.e. full disk), or when two contiguous",
714 "trajectories must be concatenated without having double frames.[PAR]",
716 "Option [TT]-dump[tt] can be used to extract a frame at or near",
717 "one specific time from your trajectory, but only works reliably",
718 "if the time interval between frames is uniform.[PAR]",
720 "Option [TT]-drop[tt] reads an [REF].xvg[ref] file with times and values.",
721 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
722 "frames with a value below and above the value of the respective options",
723 "will not be written."
739 const char *pbc_opt
[epNR
+ 1] =
741 NULL
, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
746 const char *unitcell_opt
[euNR
+1] =
747 { NULL
, "rect", "tric", "compact", NULL
};
751 ecSel
, ecTric
, ecRect
, ecZero
, ecNR
753 const char *center_opt
[ecNR
+1] =
754 { NULL
, "tric", "rect", "zero", NULL
};
760 efSel
, efNone
, efFit
, efFitXY
, efReset
, efResetXY
, efPFit
, efNR
762 const char *fit
[efNR
+ 1] =
764 NULL
, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
768 static gmx_bool bSeparate
= FALSE
, bVels
= TRUE
, bForce
= FALSE
, bCONECT
= FALSE
;
769 static gmx_bool bCenter
= FALSE
;
770 static int skip_nr
= 1, ndec
= 3, nzero
= 0;
771 static real tzero
= 0, delta_t
= 0, timestep
= 0, ttrunc
= -1, tdump
= -1, split_t
= 0;
772 static rvec newbox
= {0, 0, 0}, shift
= {0, 0, 0}, trans
= {0, 0, 0};
773 static char *exec_command
= NULL
;
774 static real dropunder
= 0, dropover
= 0;
775 static gmx_bool bRound
= FALSE
;
780 { "-skip", FALSE
, etINT
,
781 { &skip_nr
}, "Only write every nr-th frame" },
782 { "-dt", FALSE
, etTIME
,
784 "Only write frame when t MOD dt = first time (%t)" },
785 { "-round", FALSE
, etBOOL
,
786 { &bRound
}, "Round measurements to nearest picosecond"},
787 { "-dump", FALSE
, etTIME
,
788 { &tdump
}, "Dump frame nearest specified time (%t)" },
789 { "-t0", FALSE
, etTIME
,
791 "Starting time (%t) (default: don't change)" },
792 { "-timestep", FALSE
, etTIME
,
794 "Change time step between input frames (%t)" },
795 { "-pbc", FALSE
, etENUM
,
797 "PBC treatment (see help text for full description)" },
798 { "-ur", FALSE
, etENUM
,
799 { unitcell_opt
}, "Unit-cell representation" },
800 { "-center", FALSE
, etBOOL
,
801 { &bCenter
}, "Center atoms in box" },
802 { "-boxcenter", FALSE
, etENUM
,
803 { center_opt
}, "Center for -pbc and -center" },
804 { "-box", FALSE
, etRVEC
,
806 "Size for new cubic box (default: read from input)" },
807 { "-trans", FALSE
, etRVEC
,
809 "All coordinates will be translated by trans. This "
810 "can advantageously be combined with -pbc mol -ur "
812 { "-shift", FALSE
, etRVEC
,
814 "All coordinates will be shifted by framenr*shift" },
815 { "-fit", FALSE
, etENUM
,
817 "Fit molecule to ref structure in the structure file" },
818 { "-ndec", FALSE
, etINT
,
820 "Precision for .xtc and .gro writing in number of "
822 { "-vel", FALSE
, etBOOL
,
823 { &bVels
}, "Read and write velocities if possible" },
824 { "-force", FALSE
, etBOOL
,
825 { &bForce
}, "Read and write forces if possible" },
826 { "-trunc", FALSE
, etTIME
,
828 "Truncate input trajectory file after this time (%t)" },
829 { "-exec", FALSE
, etSTR
,
831 "Execute command for every output frame with the "
832 "frame number as argument" },
833 { "-split", FALSE
, etTIME
,
835 "Start writing new file when t MOD split = first "
837 { "-sep", FALSE
, etBOOL
,
839 "Write each frame to a separate .gro, .g96 or .pdb "
841 { "-nzero", FALSE
, etINT
,
843 "If the -sep flag is set, use these many digits "
844 "for the file numbers and prepend zeros as needed" },
845 { "-dropunder", FALSE
, etREAL
,
846 { &dropunder
}, "Drop all frames below this value" },
847 { "-dropover", FALSE
, etREAL
,
848 { &dropover
}, "Drop all frames above this value" },
849 { "-conect", FALSE
, etBOOL
,
851 "Add conect records when writing [REF].pdb[ref] files. Useful "
852 "for visualization of non-standard molecules, e.g. "
853 "coarse grained ones" }
855 #define NPA asize(pa)
858 t_trxstatus
*trxout
= NULL
;
861 t_trxframe fr
, frout
;
863 rvec
*xmem
= NULL
, *vmem
= NULL
, *fmem
= NULL
;
864 rvec
*xp
= NULL
, x_shift
, hbox
;
866 int m
, i
, d
, frame
, outframe
, natoms
, nout
, ncent
, newstep
= 0, model_nr
;
869 gmx_mtop_t
*mtop
= NULL
;
870 gmx_conect gc
= NULL
;
872 t_atoms
*atoms
= NULL
, useatoms
;
878 int ifit
, my_clust
= -1;
881 t_cluster_ndx
*clust
= NULL
;
882 t_trxstatus
**clust_status
= NULL
;
883 int *clust_status_id
= NULL
;
885 int *nfwritten
= NULL
;
886 int ndrop
= 0, ncol
, drop0
= 0, drop1
= 0, dropuse
= 0;
888 real tshift
= 0, t0
= -1, dt
= 0.001, prec
;
889 gmx_bool bFit
, bPFit
, bReset
;
891 gmx_rmpbc_t gpbc
= NULL
;
892 gmx_bool bRmPBC
, bPBCWhole
, bPBCcomRes
, bPBCcomMol
, bPBCcomAtom
, bPBC
, bNoJump
, bCluster
;
893 gmx_bool bCopy
, bDoIt
, bIndex
, bTDump
, bSetTime
, bTPS
= FALSE
, bDTset
= FALSE
;
894 gmx_bool bExec
, bTimeStep
= FALSE
, bDumpFrame
= FALSE
, bSetPrec
, bNeedPrec
;
895 gmx_bool bHaveFirstFrame
, bHaveNextFrame
, bSetBox
, bSetUR
, bSplit
= FALSE
;
896 gmx_bool bSubTraj
= FALSE
, bDropUnder
= FALSE
, bDropOver
= FALSE
, bTrans
= FALSE
;
897 gmx_bool bWriteFrame
, bSplitHere
;
898 const char *top_file
, *in_file
, *out_file
= NULL
;
899 char out_file2
[256], *charpt
;
900 char *outf_base
= NULL
;
901 const char *outf_ext
= NULL
;
902 char top_title
[256], title
[256], filemode
[5];
903 gmx_output_env_t
*oenv
;
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 bReset
= (fit_enum
== efReset
|| fit_enum
== efResetXY
);
954 bPFit
= fit_enum
== efPFit
;
955 pbc_enum
= nenum(pbc_opt
);
956 bPBCWhole
= pbc_enum
== epWhole
;
957 bPBCcomRes
= pbc_enum
== epComRes
;
958 bPBCcomMol
= pbc_enum
== epComMol
;
959 bPBCcomAtom
= pbc_enum
== epComAtom
;
960 bNoJump
= pbc_enum
== epNojump
;
961 bCluster
= pbc_enum
== epCluster
;
962 bPBC
= pbc_enum
!= epNone
;
963 unitcell_enum
= nenum(unitcell_opt
);
964 ecenter
= nenum(center_opt
) - ecTric
;
966 /* set and check option dependencies */
969 bFit
= TRUE
; /* for pfit, fit *must* be set */
973 bReset
= TRUE
; /* for fit, reset *must* be set */
978 nfitdim
= (fit_enum
== efFitXY
|| fit_enum
== efResetXY
) ? 2 : 3;
980 bRmPBC
= bFit
|| bPBCWhole
|| bPBCcomRes
|| bPBCcomMol
;
984 if (!(bPBCcomRes
|| bPBCcomMol
|| bPBCcomAtom
))
987 "WARNING: Option for unitcell representation (-ur %s)\n"
988 " only has effect in combination with -pbc %s, %s or %s.\n"
989 " Ingoring unitcell representation.\n\n",
990 unitcell_opt
[0], pbc_opt
[2], pbc_opt
[3], pbc_opt
[4]);
995 gmx_fatal(FARGS
, "PBC condition treatment does not work together with rotational fit.\n"
996 "Please do the PBC condition treatment first and then run trjconv in a second step\n"
997 "for the rotational fit.\n"
998 "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1002 /* ndec is in nr of decimal places, prec is a multiplication factor: */
1004 for (i
= 0; i
< ndec
; i
++)
1009 bIndex
= ftp2bSet(efNDX
, NFILE
, fnm
);
1012 /* Determine output type */
1013 out_file
= opt2fn("-o", NFILE
, fnm
);
1014 int ftp
= fn2ftp(out_file
);
1015 fprintf(stderr
, "Will write %s: %s\n", ftp2ext(ftp
), ftp2desc(ftp
));
1016 bNeedPrec
= (ftp
== efXTC
|| ftp
== efGRO
);
1017 int ftpin
= fn2ftp(in_file
);
1020 /* check if velocities are possible in input and output files */
1021 bVels
= (ftp
== efTRR
|| ftp
== efGRO
||
1022 ftp
== efG96
|| ftp
== efTNG
)
1023 && (ftpin
== efTRR
|| ftpin
== efGRO
||
1024 ftpin
== efG96
|| ftpin
== efTNG
|| ftpin
== efCPT
);
1026 if (bSeparate
|| bSplit
)
1028 outf_ext
= std::strrchr(out_file
, '.');
1029 if (outf_ext
== NULL
)
1031 gmx_fatal(FARGS
, "Output file name '%s' does not contain a '.'", out_file
);
1033 outf_base
= gmx_strdup(out_file
);
1034 outf_base
[outf_ext
- out_file
] = '\0';
1037 bSubTraj
= opt2bSet("-sub", NFILE
, fnm
);
1040 if ((ftp
!= efXTC
) && (ftp
!= efTRR
))
1042 /* It seems likely that other trajectory file types
1043 * could work here. */
1044 gmx_fatal(FARGS
, "Can only use the sub option with output file types "
1047 clust
= cluster_index(NULL
, opt2fn("-sub", NFILE
, fnm
));
1049 /* Check for number of files disabled, as FOPEN_MAX is not the correct
1050 * number to check for. In my linux box it is only 16.
1052 if (0 && (clust
->clust
->nr
> FOPEN_MAX
-4))
1054 gmx_fatal(FARGS
, "Can not open enough (%d) files to write all the"
1055 " trajectories.\ntry splitting the index file in %d parts.\n"
1057 clust
->clust
->nr
, 1+clust
->clust
->nr
/FOPEN_MAX
, FOPEN_MAX
);
1059 gmx_warning("The -sub option could require as many open output files as there are\n"
1060 "index groups in the file (%d). If you get I/O errors opening new files,\n"
1061 "try reducing the number of index groups in the file, and perhaps\n"
1062 "using trjconv -sub several times on different chunks of your index file.\n",
1065 snew(clust_status
, clust
->clust
->nr
);
1066 snew(clust_status_id
, clust
->clust
->nr
);
1067 snew(nfwritten
, clust
->clust
->nr
);
1068 for (i
= 0; (i
< clust
->clust
->nr
); i
++)
1070 clust_status
[i
] = NULL
;
1071 clust_status_id
[i
] = -1;
1073 bSeparate
= bSplit
= FALSE
;
1080 mtop
= read_mtop_for_tng(top_file
, in_file
, out_file
);
1082 /* Determine whether to read a topology */
1083 bTPS
= (ftp2bSet(efTPS
, NFILE
, fnm
) ||
1084 bRmPBC
|| bReset
|| bPBCcomMol
|| bCluster
||
1085 (ftp
== efGRO
) || (ftp
== efPDB
) || bCONECT
);
1087 /* Determine if when can read index groups */
1088 bIndex
= (bIndex
|| bTPS
);
1092 read_tps_conf(top_file
, &top
, &ePBC
, &xp
, NULL
, top_box
,
1093 bReset
|| bPBCcomRes
);
1094 std::strncpy(top_title
, *top
.name
, 255);
1095 top_title
[255] = '\0';
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
= std::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
, (int **)&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
= std::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 std::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
) ||
1411 /* generate new box */
1412 if (fr
.bBox
== FALSE
)
1416 for (m
= 0; m
< DIM
; m
++)
1420 fr
.box
[m
][m
] = newbox
[m
];
1424 if (fr
.bBox
== FALSE
)
1426 gmx_fatal(FARGS
, "Cannot preserve a box that does not exist.\n");
1434 for (i
= 0; i
< natoms
; i
++)
1436 rvec_inc(fr
.x
[i
], trans
);
1442 /* determine timestep */
1455 /* This is not very elegant, as one can not dump a frame after
1456 * a timestep with is more than twice as small as the first one. */
1457 bDumpFrame
= (fr
.time
> tdump
-0.5*dt
) && (fr
.time
<= tdump
+0.5*dt
);
1464 /* determine if an atom jumped across the box and reset it if so */
1465 if (bNoJump
&& (bTPS
|| frame
!= 0))
1467 for (d
= 0; d
< DIM
; d
++)
1469 hbox
[d
] = 0.5*fr
.box
[d
][d
];
1471 for (i
= 0; i
< natoms
; i
++)
1475 rvec_dec(fr
.x
[i
], x_shift
);
1477 for (m
= DIM
-1; m
>= 0; m
--)
1481 while (fr
.x
[i
][m
]-xp
[i
][m
] <= -hbox
[m
])
1483 for (d
= 0; d
<= m
; d
++)
1485 fr
.x
[i
][d
] += fr
.box
[m
][d
];
1488 while (fr
.x
[i
][m
]-xp
[i
][m
] > hbox
[m
])
1490 for (d
= 0; d
<= m
; d
++)
1492 fr
.x
[i
][d
] -= fr
.box
[m
][d
];
1501 calc_pbc_cluster(ecenter
, ifit
, &top
, ePBC
, fr
.x
, ind_fit
, fr
.box
);
1506 /* Now modify the coords according to the flags,
1507 for normal fit, this is only done for output frames */
1510 gmx_rmpbc_trxfr(gpbc
, &fr
);
1513 reset_x_ndim(nfitdim
, ifit
, ind_fit
, natoms
, NULL
, fr
.x
, w_rls
);
1514 do_fit(natoms
, w_rls
, xp
, fr
.x
);
1517 /* store this set of coordinates for future use */
1518 if (bPFit
|| bNoJump
)
1524 for (i
= 0; (i
< natoms
); i
++)
1526 copy_rvec(fr
.x
[i
], xp
[i
]);
1527 rvec_inc(fr
.x
[i
], x_shift
);
1533 /* see if we have a frame from the frame index group */
1534 for (i
= 0; i
< nrfri
&& !bDumpFrame
; i
++)
1536 bDumpFrame
= frame
== frindex
[i
];
1539 if (debug
&& bDumpFrame
)
1541 fprintf(debug
, "dumping %d\n", frame
);
1545 ( ( !bTDump
&& !frindex
&& frame
% skip_nr
== 0 ) || bDumpFrame
);
1547 if (bWriteFrame
&& (bDropUnder
|| bDropOver
))
1549 while (dropval
[0][drop1
] < fr
.time
&& drop1
+1 < ndrop
)
1554 if (std::abs(dropval
[0][drop0
] - fr
.time
)
1555 < std::abs(dropval
[0][drop1
] - fr
.time
))
1563 if ((bDropUnder
&& dropval
[1][dropuse
] < dropunder
) ||
1564 (bDropOver
&& dropval
[1][dropuse
] > dropover
))
1566 bWriteFrame
= FALSE
;
1572 /* We should avoid modifying the input frame,
1573 * but since here we don't have the output frame yet,
1574 * we introduce a temporary output frame time variable.
1578 frout_time
= fr
.time
;
1583 frout_time
= tzero
+ frame
*timestep
;
1588 frout_time
+= tshift
;
1593 fprintf(stderr
, "\nDumping frame at t= %g %s\n",
1594 output_env_conv_time(oenv
, frout_time
), output_env_get_time_unit(oenv
));
1597 /* check for writing at each delta_t */
1598 bDoIt
= (delta_t
== 0);
1603 bDoIt
= bRmod(frout_time
, tzero
, delta_t
);
1607 /* round() is not C89 compatible, so we do this: */
1608 bDoIt
= bRmod(std::floor(frout_time
+0.5), std::floor(tzero
+0.5),
1609 std::floor(delta_t
+0.5));
1613 if (bDoIt
|| bTDump
)
1615 /* print sometimes */
1616 if ( ((outframe
% SKIP
) == 0) || (outframe
< SKIP
) )
1618 fprintf(stderr
, " -> frame %6d time %8.3f \r",
1619 outframe
, output_env_conv_time(oenv
, frout_time
));
1625 /* Now modify the coords according to the flags,
1626 for PFit we did this already! */
1630 gmx_rmpbc_trxfr(gpbc
, &fr
);
1635 reset_x_ndim(nfitdim
, ifit
, ind_fit
, natoms
, NULL
, fr
.x
, w_rls
);
1638 do_fit_ndim(nfitdim
, natoms
, w_rls
, xp
, fr
.x
);
1642 for (i
= 0; i
< natoms
; i
++)
1644 rvec_inc(fr
.x
[i
], x_shift
);
1651 center_x(ecenter
, fr
.x
, fr
.box
, natoms
, ncent
, cindex
);
1657 switch (unitcell_enum
)
1660 put_atoms_in_box(ePBC
, fr
.box
, natoms
, fr
.x
);
1663 put_atoms_in_triclinic_unitcell(ecenter
, fr
.box
, natoms
, fr
.x
);
1666 put_atoms_in_compact_unitcell(ePBC
, ecenter
, fr
.box
,
1673 put_residue_com_in_box(unitcell_enum
, ecenter
,
1674 natoms
, atoms
->atom
, ePBC
, fr
.box
, fr
.x
);
1678 put_molecule_com_in_box(unitcell_enum
, ecenter
,
1680 natoms
, atoms
->atom
, ePBC
, fr
.box
, fr
.x
);
1682 /* Copy the input trxframe struct to the output trxframe struct */
1684 frout
.time
= frout_time
;
1685 frout
.bV
= (frout
.bV
&& bVels
);
1686 frout
.bF
= (frout
.bF
&& bForce
);
1687 frout
.natoms
= nout
;
1688 if (bNeedPrec
&& (bSetPrec
|| !fr
.bPrec
))
1704 for (i
= 0; i
< nout
; i
++)
1706 copy_rvec(fr
.x
[index
[i
]], frout
.x
[i
]);
1709 copy_rvec(fr
.v
[index
[i
]], frout
.v
[i
]);
1713 copy_rvec(fr
.f
[index
[i
]], frout
.f
[i
]);
1718 if (opt2parg_bSet("-shift", NPA
, pa
))
1720 for (i
= 0; i
< nout
; i
++)
1722 for (d
= 0; d
< DIM
; d
++)
1724 frout
.x
[i
][d
] += outframe
*shift
[d
];
1731 bSplitHere
= bSplit
&& bRmod(frout
.time
, tzero
, split_t
);
1735 /* round() is not C89 compatible, so we do this: */
1736 bSplitHere
= bSplit
&& bRmod(std::floor(frout
.time
+0.5),
1737 std::floor(tzero
+0.5),
1738 std::floor(split_t
+0.5));
1740 if (bSeparate
|| bSplitHere
)
1742 mk_filenm(outf_base
, ftp2ext(ftp
), nzero
, file_nr
, out_file2
);
1748 write_tng_frame(trxout
, &frout
);
1749 // TODO when trjconv behaves better: work how to read and write lambda
1759 trxout
= open_trx(out_file2
, filemode
);
1766 if (clust_status_id
[my_clust
] == -1)
1768 sprintf(buf
, "%s.%s", clust
->grpname
[my_clust
], ftp2ext(ftp
));
1769 clust_status
[my_clust
] = open_trx(buf
, "w");
1770 clust_status_id
[my_clust
] = 1;
1773 else if (clust_status_id
[my_clust
] == -2)
1775 gmx_fatal(FARGS
, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1776 clust
->grpname
[my_clust
], ntrxopen
, frame
,
1779 write_trxframe(clust_status
[my_clust
], &frout
, gc
);
1780 nfwritten
[my_clust
]++;
1781 if (nfwritten
[my_clust
] ==
1782 (clust
->clust
->index
[my_clust
+1]-
1783 clust
->clust
->index
[my_clust
]))
1785 close_trx(clust_status
[my_clust
]);
1786 clust_status
[my_clust
] = NULL
;
1787 clust_status_id
[my_clust
] = -2;
1791 gmx_fatal(FARGS
, "Less than zero open .xtc files!");
1798 write_trxframe(trxout
, &frout
, gc
);
1804 sprintf(title
, "Generated by trjconv : %s t= %9.5f",
1805 top_title
, frout
.time
);
1806 if (bSeparate
|| bSplitHere
)
1808 out
= gmx_ffopen(out_file2
, "w");
1813 write_hconf_p(out
, title
, &useatoms
,
1814 frout
.x
, frout
.bV
? frout
.v
: NULL
, frout
.box
);
1817 fprintf(out
, "REMARK GENERATED BY TRJCONV\n");
1818 sprintf(title
, "%s t= %9.5f", top_title
, frout
.time
);
1819 /* if reading from pdb, we want to keep the original
1820 model numbering else we write the output frame
1821 number plus one, because model 0 is not allowed in pdb */
1822 if (ftpin
== efPDB
&& fr
.bStep
&& fr
.step
> model_nr
)
1830 write_pdbfile(out
, title
, &useatoms
, frout
.x
,
1831 frout
.ePBC
, frout
.box
, ' ', model_nr
, gc
, TRUE
);
1834 frout
.title
= title
;
1835 if (bSeparate
|| bTDump
)
1837 frout
.bTitle
= TRUE
;
1840 frout
.bAtoms
= TRUE
;
1842 frout
.atoms
= &useatoms
;
1843 frout
.bStep
= FALSE
;
1844 frout
.bTime
= FALSE
;
1848 frout
.bTitle
= (outframe
== 0);
1849 frout
.bAtoms
= FALSE
;
1853 write_g96_conf(out
, &frout
, -1, NULL
);
1855 if (bSeparate
|| bSplitHere
)
1862 gmx_fatal(FARGS
, "DHE, ftp=%d\n", ftp
);
1864 if (bSeparate
|| bSplitHere
)
1869 /* execute command */
1873 sprintf(c
, "%s %d", exec_command
, file_nr
-1);
1874 /*fprintf(stderr,"Executing '%s'\n",c);*/
1877 gmx_fatal(FARGS
, "Error executing command: %s", c
);
1884 bHaveNextFrame
= read_next_frame(oenv
, trxin
, &fr
);
1886 while (!(bTDump
&& bDumpFrame
) && bHaveNextFrame
);
1889 if (!bHaveFirstFrame
|| (bTDump
&& !bDumpFrame
))
1891 fprintf(stderr
, "\nWARNING no output, "
1892 "last frame read at t=%g\n", fr
.time
);
1894 fprintf(stderr
, "\n");
1901 gmx_rmpbc_done(gpbc
);
1908 else if (out
!= NULL
)
1914 for (i
= 0; (i
< clust
->clust
->nr
); i
++)
1916 if (clust_status_id
[i
] >= 0)
1918 close_trx(clust_status
[i
]);
1926 do_view(oenv
, out_file
, NULL
);