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, 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/fileio/confio.h"
47 #include "gromacs/fileio/gmxfio.h"
48 #include "gromacs/fileio/pdbio.h"
49 #include "gromacs/fileio/tngio.h"
50 #include "gromacs/fileio/tngio_for_tools.h"
51 #include "gromacs/fileio/trxio.h"
52 #include "gromacs/fileio/xtcio.h"
53 #include "gromacs/fileio/xvgr.h"
54 #include "gromacs/gmxana/gmx_ana.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdtypes/md_enums.h"
57 #include "gromacs/topology/index.h"
58 #include "gromacs/trajectory/trajectoryframe.h"
59 #include "gromacs/utility/arraysize.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
64 #define TIME_EXPLICIT 0
65 #define TIME_CONTINUE 1
70 #define FLAGS (TRX_READ_X | TRX_READ_V | TRX_READ_F)
72 static void scan_trj_files(char **fnms
, int nfiles
, real
*readtime
,
73 real
*timestep
, int imax
,
74 const gmx_output_env_t
*oenv
)
76 /* Check start time of all files */
82 for (i
= 0; i
< nfiles
; i
++)
84 ok
= read_first_frame(oenv
, &status
, fnms
[i
], &fr
, FLAGS
);
88 gmx_fatal(FARGS
, "\nCouldn't read frame from file." );
92 readtime
[i
] = fr
.time
;
97 fprintf(stderr
, "\nWARNING: Couldn't find a time in the frame.\n");
108 if (natoms
!= fr
.natoms
)
110 gmx_fatal(FARGS
, "\nDifferent numbers of atoms (%d/%d) in files",
116 if (fr
.natoms
<= imax
)
118 gmx_fatal(FARGS
, "\nNot enough atoms (%d) for index group (%d)",
123 ok
= read_next_frame(oenv
, status
, &fr
);
126 timestep
[i
] = fr
.time
- readtime
[i
];
147 fprintf(stderr
, "\n");
151 static void sort_files(char **fnms
, real
*settime
, int nfile
)
157 for (i
= 0; i
< nfile
; i
++)
160 for (j
= i
+ 1; j
< nfile
; j
++)
162 if (settime
[j
] < settime
[minidx
])
169 timeswap
= settime
[i
];
170 settime
[i
] = settime
[minidx
];
171 settime
[minidx
] = timeswap
;
173 fnms
[i
] = fnms
[minidx
];
174 fnms
[minidx
] = chptr
;
179 static void edit_files(char **fnms
, int nfiles
, real
*readtime
, real
*timestep
,
180 real
*settime
, int *cont_type
, gmx_bool bSetTime
,
181 gmx_bool bSort
, const gmx_output_env_t
*oenv
)
185 char inputstring
[STRLEN
], *chptr
;
189 fprintf(stderr
, "\n\nEnter the new start time (%s) for each file.\n"
190 "There are two special options, both disable sorting:\n\n"
191 "c (continue) - The start time is taken from the end\n"
192 "of the previous file. Use it when your continuation run\n"
193 "restarts with t=0.\n\n"
194 "l (last) - The time in this file will be changed the\n"
195 "same amount as in the previous. Use it when the time in the\n"
196 "new run continues from the end of the previous one,\n"
197 "since this takes possible overlap into account.\n\n",
198 output_env_get_time_unit(oenv
));
202 " File Current start (%s) New start (%s)\n"
203 "---------------------------------------------------------\n",
204 output_env_get_time_unit(oenv
), output_env_get_time_unit(oenv
));
206 for (i
= 0; i
< nfiles
; i
++)
208 fprintf(stderr
, "%25s %10.3f %s ", fnms
[i
],
209 output_env_conv_time(oenv
, readtime
[i
]), output_env_get_time_unit(oenv
));
213 if (NULL
== fgets(inputstring
, STRLEN
- 1, stdin
))
215 gmx_fatal(FARGS
, "Error reading user input" );
218 inputstring
[std::strlen(inputstring
)-1] = 0;
220 if (inputstring
[0] == 'c' || inputstring
[0] == 'C')
222 cont_type
[i
] = TIME_CONTINUE
;
225 settime
[i
] = FLT_MAX
;
227 else if (inputstring
[0] == 'l' ||
228 inputstring
[0] == 'L')
230 cont_type
[i
] = TIME_LAST
;
233 settime
[i
] = FLT_MAX
;
237 settime
[i
] = strtod(inputstring
, &chptr
)*
238 output_env_get_time_invfactor(oenv
);
239 if (chptr
== inputstring
)
241 fprintf(stderr
, "'%s' not recognized as a floating point number, 'c' or 'l'. "
242 "Try again: ", inputstring
);
246 cont_type
[i
] = TIME_EXPLICIT
;
253 if (cont_type
[0] != TIME_EXPLICIT
)
255 cont_type
[0] = TIME_EXPLICIT
;
261 for (i
= 0; i
< nfiles
; i
++)
263 settime
[i
] = readtime
[i
];
268 fprintf(stderr
, "Sorting disabled.\n");
272 sort_files(fnms
, settime
, nfiles
);
274 /* Write out the new order and start times */
275 fprintf(stderr
, "\nSummary of files and start times used:\n\n"
276 " File Start time Time step\n"
277 "---------------------------------------------------------\n");
278 for (i
= 0; i
< nfiles
; i
++)
280 switch (cont_type
[i
])
283 fprintf(stderr
, "%25s %10.3f %s %10.3f %s",
285 output_env_conv_time(oenv
, settime
[i
]), output_env_get_time_unit(oenv
),
286 output_env_conv_time(oenv
, timestep
[i
]), output_env_get_time_unit(oenv
));
288 cont_type
[i
-1] == TIME_EXPLICIT
&& settime
[i
] == settime
[i
-1])
290 fprintf(stderr
, " WARNING: same Start time as previous");
292 fprintf(stderr
, "\n");
295 fprintf(stderr
, "%25s Continue from last file\n", fnms
[i
]);
298 fprintf(stderr
, "%25s Change by same amount as last file\n",
303 fprintf(stderr
, "\n");
305 settime
[nfiles
] = FLT_MAX
;
306 cont_type
[nfiles
] = TIME_EXPLICIT
;
307 readtime
[nfiles
] = FLT_MAX
;
310 static void do_demux(int nset
, char *fnms
[], char *fnms_out
[], int nval
,
311 real
**value
, real
*time
, real dt_remd
, int isize
,
312 int index
[], real dt
, const gmx_output_env_t
*oenv
)
314 int i
, j
, k
, natoms
, nnn
;
315 t_trxstatus
**fp_in
, **fp_out
;
316 gmx_bool bCont
, *bSet
;
317 real t
, first_time
= 0;
325 for (i
= 0; (i
< nset
); i
++)
327 nnn
= read_first_frame(oenv
, &(fp_in
[i
]), fnms
[i
], &(trx
[i
]),
332 first_time
= trx
[i
].time
;
334 else if (natoms
!= nnn
)
336 gmx_fatal(FARGS
, "Trajectory file %s has %d atoms while previous trajs had %d atoms", fnms
[i
], nnn
, natoms
);
342 else if (t
!= trx
[i
].time
)
344 gmx_fatal(FARGS
, "Trajectory file %s has time %f while previous trajs had time %f", fnms
[i
], trx
[i
].time
, t
);
349 for (i
= 0; (i
< nset
); i
++)
351 fp_out
[i
] = open_trx(fnms_out
[i
], "w");
354 if (std::round(time
[k
] - t
) != 0)
356 gmx_fatal(FARGS
, "First time in demuxing table does not match trajectories");
360 while ((k
+1 < nval
) && ((trx
[0].time
- time
[k
+1]) > dt_remd
*0.1))
366 fprintf(debug
, "trx[0].time = %g, time[k] = %g\n", trx
[0].time
, time
[k
]);
368 for (i
= 0; (i
< nset
); i
++)
372 for (i
= 0; (i
< nset
); i
++)
374 j
= std::round(value
[i
][k
]);
375 range_check(j
, 0, nset
);
378 gmx_fatal(FARGS
, "Demuxing the same replica %d twice at time %f",
383 if (dt
== 0 || bRmod(trx
[i
].time
, first_time
, dt
))
387 write_trxframe_indexed(fp_out
[j
], &trx
[i
], isize
, index
, NULL
);
391 write_trxframe(fp_out
[j
], &trx
[i
], NULL
);
397 for (i
= 0; (i
< nset
); i
++)
399 bCont
= bCont
&& read_next_frame(oenv
, fp_in
[i
], &trx
[i
]);
404 for (i
= 0; (i
< nset
); i
++)
407 close_trx(fp_out
[i
]);
411 int gmx_trjcat(int argc
, char *argv
[])
415 "[THISMODULE] concatenates several input trajectory files in sorted order. ",
416 "In case of double time frames the one in the later file is used. ",
417 "By specifying [TT]-settime[tt] you will be asked for the start time ",
418 "of each file. The input files are taken from the command line, ",
419 "such that a command like [TT]gmx trjcat -f *.trr -o fixed.trr[tt] should do ",
420 "the trick. Using [TT]-cat[tt], you can simply paste several files ",
421 "together without removal of frames with identical time stamps.[PAR]",
422 "One important option is inferred when the output file is amongst the",
423 "input files. In that case that particular file will be appended to",
424 "which implies you do not need to store double the amount of data.",
425 "Obviously the file to append to has to be the one with lowest starting",
426 "time since one can only append at the end of a file.[PAR]",
427 "If the [TT]-demux[tt] option is given, the N trajectories that are",
428 "read, are written in another order as specified in the [REF].xvg[ref] file.",
429 "The [REF].xvg[ref] file should contain something like::",
434 "The first number is the time, and subsequent numbers point to",
435 "trajectory indices.",
436 "The frames corresponding to the numbers present at the first line",
437 "are collected into the output trajectory. If the number of frames in",
438 "the trajectory does not match that in the [REF].xvg[ref] file then the program",
439 "tries to be smart. Beware."
441 static gmx_bool bVels
= TRUE
;
442 static gmx_bool bCat
= FALSE
;
443 static gmx_bool bSort
= TRUE
;
444 static gmx_bool bKeepLast
= FALSE
;
445 static gmx_bool bKeepLastAppend
= FALSE
;
446 static gmx_bool bOverwrite
= FALSE
;
447 static gmx_bool bSetTime
= FALSE
;
448 static gmx_bool bDeMux
;
449 static real begin
= -1;
450 static real end
= -1;
456 { "-b", FALSE
, etTIME
,
457 { &begin
}, "First time to use (%t)" },
458 { "-e", FALSE
, etTIME
,
459 { &end
}, "Last time to use (%t)" },
460 { "-dt", FALSE
, etTIME
,
461 { &dt
}, "Only write frame when t MOD dt = first time (%t)" },
462 { "-vel", FALSE
, etBOOL
,
463 { &bVels
}, "Read and write velocities if possible" },
464 { "-settime", FALSE
, etBOOL
,
465 { &bSetTime
}, "Change starting time interactively" },
466 { "-sort", FALSE
, etBOOL
,
467 { &bSort
}, "Sort trajectory files (not frames)" },
468 { "-keeplast", FALSE
, etBOOL
,
469 { &bKeepLast
}, "Keep overlapping frames at end of trajectory" },
470 { "-overwrite", FALSE
, etBOOL
,
471 { &bOverwrite
}, "Overwrite overlapping frames during appending" },
472 { "-cat", FALSE
, etBOOL
,
473 { &bCat
}, "Do not discard double time frames" }
475 #define npargs asize(pa)
476 int ftpin
, i
, frame
, frame_out
;
477 t_trxstatus
*status
, *trxout
= NULL
;
479 t_trxframe fr
, frout
;
480 char **fnms
, **fnms_out
, *out_file
;
482 gmx_bool bNewFile
, bIndex
, bWrite
;
483 int nfile_in
, nfile_out
, *cont_type
;
484 real
*readtime
, *timest
, *settime
;
485 real first_time
= 0, lasttime
, last_ok_t
= -1, timestep
;
486 gmx_bool lastTimeSet
= FALSE
;
487 real last_frame_time
, searchtime
;
489 int *index
= NULL
, imax
;
491 real
**val
= NULL
, *t
= NULL
, dt_remd
;
492 int n
, nset
, ftpout
= -1, prevEndStep
= 0, filetype
;
494 gmx_output_env_t
*oenv
;
497 { efTRX
, "-f", NULL
, ffRDMULT
},
498 { efTRO
, "-o", NULL
, ffWRMULT
},
499 { efNDX
, "-n", "index", ffOPTRD
},
500 { efXVG
, "-demux", "remd", ffOPTRD
}
503 #define NFILE asize(fnm)
505 if (!parse_common_args(&argc
, argv
, PCA_TIME_UNIT
, NFILE
, fnm
,
506 asize(pa
), pa
, asize(desc
), desc
, 0, NULL
, &oenv
))
511 bIndex
= ftp2bSet(efNDX
, NFILE
, fnm
);
512 bDeMux
= ftp2bSet(efXVG
, NFILE
, fnm
);
513 bSort
= bSort
&& !bDeMux
;
518 printf("Select group for output\n");
519 rd_index(ftp2fn(efNDX
, NFILE
, fnm
), 1, &isize
, &index
, &grpname
);
522 for (i
= 1; i
< isize
; i
++)
524 imax
= std::max(imax
, index
[i
]);
531 val
= read_xvg_time(opt2fn("-demux", NFILE
, fnm
), TRUE
,
532 opt2parg_bSet("-b", npargs
, pa
), begin
,
533 opt2parg_bSet("-e", npargs
, pa
), end
, 1, &nset
, &n
,
535 printf("Read %d sets of %d points, dt = %g\n\n", nset
, n
, dt_remd
);
538 fprintf(debug
, "Dump of replica_index.xvg\n");
539 for (i
= 0; (i
< n
); i
++)
541 fprintf(debug
, "%10g", t
[i
]);
542 for (j
= 0; (j
< nset
); j
++)
544 fprintf(debug
, " %3d", static_cast<int>(std::round(val
[j
][i
])));
546 fprintf(debug
, "\n");
551 nfile_in
= opt2fns(&fnms
, "-f", NFILE
, fnm
);
554 gmx_fatal(FARGS
, "No input files!" );
557 if (bDeMux
&& (nfile_in
!= nset
))
559 gmx_fatal(FARGS
, "You have specified %d files and %d entries in the demux table", nfile_in
, nset
);
562 ftpin
= fn2ftp(fnms
[0]);
564 for (i
= 1; i
< nfile_in
; i
++)
566 if (ftpin
!= fn2ftp(fnms
[i
]))
568 gmx_fatal(FARGS
, "All input files must be of the same format");
572 nfile_out
= opt2fns(&fnms_out
, "-o", NFILE
, fnm
);
575 gmx_fatal(FARGS
, "No output files!");
577 if ((nfile_out
> 1) && !bDeMux
)
579 gmx_fatal(FARGS
, "Don't know what to do with more than 1 output file if not demultiplexing");
581 else if (bDeMux
&& (nfile_out
!= nset
) && (nfile_out
!= 1))
583 gmx_fatal(FARGS
, "Number of output files should be 1 or %d (#input files), not %d", nset
, nfile_out
);
587 if (nfile_out
!= nset
)
589 char *buf
= gmx_strdup(fnms_out
[0]);
590 snew(fnms_out
, nset
);
591 for (i
= 0; (i
< nset
); i
++)
593 snew(fnms_out
[i
], std::strlen(buf
)+32);
594 sprintf(fnms_out
[i
], "%d_%s", i
, buf
);
598 do_demux(nfile_in
, fnms
, fnms_out
, n
, val
, t
, dt_remd
, isize
, index
, dt
, oenv
);
602 snew(readtime
, nfile_in
+1);
603 snew(timest
, nfile_in
+1);
604 scan_trj_files(fnms
, nfile_in
, readtime
, timest
, imax
, oenv
);
606 snew(settime
, nfile_in
+1);
607 snew(cont_type
, nfile_in
+1);
608 edit_files(fnms
, nfile_in
, readtime
, timest
, settime
, cont_type
, bSetTime
, bSort
,
611 /* Check whether the output file is amongst the input files
612 * This has to be done after sorting etc.
614 out_file
= fnms_out
[0];
615 ftpout
= fn2ftp(out_file
);
617 for (i
= 0; ((i
< nfile_in
) && (n_append
== -1)); i
++)
619 if (std::strcmp(fnms
[i
], out_file
) == 0)
626 fprintf(stderr
, "Will append to %s rather than creating a new file\n",
629 else if (n_append
!= -1)
631 gmx_fatal(FARGS
, "Can only append to the first file which is %s (not %s)",
635 /* Not checking input format, could be dangerous :-) */
636 /* Not checking output format, equally dangerous :-) */
640 /* the default is not to change the time at all,
641 * but this is overridden by the edit_files routine
651 gmx_fatal(FARGS
, "When writing TNG the input file format must also be TNG");
655 trjtools_gmx_prepare_tng_writing(out_file
, 'w', NULL
, &trxout
,
656 fnms
[0], isize
, NULL
, index
, grpname
);
660 trjtools_gmx_prepare_tng_writing(out_file
, 'w', NULL
, &trxout
,
661 fnms
[0], -1, NULL
, NULL
, NULL
);
666 trxout
= open_trx(out_file
, "w");
668 std::memset(&frout
, 0, sizeof(frout
));
674 if (!read_first_frame(oenv
, &status
, out_file
, &fr
, FLAGS
))
676 gmx_fatal(FARGS
, "Reading first frame from %s", out_file
);
679 stfio
= trx_get_fileio(status
);
680 if (!bKeepLast
&& !bOverwrite
)
682 fprintf(stderr
, "\n\nWARNING: Appending without -overwrite implies -keeplast "
683 "between the first two files. \n"
684 "If the trajectories have an overlap and have not been written binary \n"
685 "reproducible this will produce an incorrect trajectory!\n\n");
687 filetype
= gmx_fio_getftp(stfio
);
688 /* Fails if last frame is incomplete
689 * We can't do anything about it without overwriting
691 if (filetype
== efXTC
|| filetype
== efTNG
)
693 lasttime
= trx_get_time_of_final_frame(status
);
698 while (read_next_frame(oenv
, status
, &fr
))
705 bKeepLastAppend
= TRUE
;
707 trxout
= open_trx(out_file
, "a");
711 if (gmx_fio_getftp(stfio
) != efXTC
)
713 gmx_fatal(FARGS
, "Overwrite only supported for XTC." );
715 last_frame_time
= trx_get_time_of_final_frame(status
);
717 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
718 * or when seek time = 0 */
719 if (nfile_in
> 1 && settime
[1] < last_frame_time
+timest
[0]*0.5)
721 /* Jump to one time-frame before the start of next
723 searchtime
= settime
[1]-timest
[0]*1.25;
727 searchtime
= last_frame_time
;
729 if (xtc_seek_time(stfio
, searchtime
, fr
.natoms
, TRUE
))
731 gmx_fatal(FARGS
, "Error seeking to append position.");
733 read_next_frame(oenv
, status
, &fr
);
734 if (std::abs(searchtime
- fr
.time
) > timest
[0]*0.5)
736 gmx_fatal(FARGS
, "Error seeking: attempted to seek to %f but got %f.",
737 searchtime
, fr
.time
);
741 fpos
= gmx_fio_ftell(stfio
);
743 trxout
= open_trx(out_file
, "r+");
744 if (gmx_fio_seek(trx_get_fileio(trxout
), fpos
))
746 gmx_fatal(FARGS
, "Error seeking to append position.");
751 printf("\n Will append after %f \n", lasttime
);
755 /* Lets stitch up some files */
756 timestep
= timest
[0];
757 for (i
= n_append
+1; (i
< nfile_in
); i
++)
761 /* set the next time from the last frame in previous file */
764 /* When writing TNG the step determine which frame to write. Use an
765 * offset to be able to increase steps properly when changing files. */
768 prevEndStep
= frout
.step
;
773 if (cont_type
[i
] == TIME_CONTINUE
)
776 begin
+= 0.5*timestep
;
777 settime
[i
] = frout
.time
;
778 cont_type
[i
] = TIME_EXPLICIT
;
780 else if (cont_type
[i
] == TIME_LAST
)
783 begin
+= 0.5*timestep
;
785 /* Or, if the time in the next part should be changed by the
786 * same amount, start at half a timestep from the last time
787 * so we dont repeat frames.
789 /* I don't understand the comment above, but for all the cases
790 * I tried the code seems to work properly. B. Hess 2008-4-2.
793 /* Or, if time is set explicitly, we check for overlap/gap */
794 if (cont_type
[i
] == TIME_EXPLICIT
)
796 if ( ( i
< nfile_in
) &&
797 ( frout
.time
< settime
[i
]-1.5*timestep
) )
799 fprintf(stderr
, "WARNING: Frames around t=%f %s have a different "
800 "spacing than the rest,\n"
801 "might be a gap or overlap that couldn't be corrected "
802 "automatically.\n", output_env_conv_time(oenv
, frout
.time
),
803 output_env_get_time_unit(oenv
));
808 /* if we don't have a timestep in the current file, use the old one */
811 timestep
= timest
[i
];
813 read_first_frame(oenv
, &status
, fnms
[i
], &fr
, FLAGS
);
817 fprintf(stderr
, "\nWARNING: Couldn't find a time in the frame.\n");
820 if (cont_type
[i
] == TIME_EXPLICIT
)
822 t_corr
= settime
[i
]-fr
.time
;
824 /* t_corr is the amount we want to change the time.
825 * If the user has chosen not to change the time for
826 * this part of the trajectory t_corr remains at
827 * the value it had in the last part, changing this
828 * by the same amount.
829 * If no value was given for the first trajectory part
830 * we let the time start at zero, see the edit_files routine.
841 printf("lasttime %g\n", lasttime
);
845 /* copy the input frame to the output frame */
847 /* set the new time by adding the correct calculated above */
848 frout
.time
+= t_corr
;
851 frout
.step
+= prevEndStep
;
853 /* quit if we have reached the end of what should be written */
854 if ((end
> 0) && (frout
.time
> end
+GMX_REAL_EPS
))
860 /* determine if we should write this frame (dt is handled elsewhere) */
861 if (bCat
) /* write all frames of all files */
865 else if (bKeepLast
|| (bKeepLastAppend
&& i
== 1))
866 /* write till last frame of this traj
867 and skip first frame(s) of next traj */
869 bWrite
= ( frout
.time
> lasttime
+0.5*timestep
);
871 else /* write till first frame of next traj */
873 bWrite
= ( frout
.time
< settime
[i
+1]-0.5*timestep
);
876 if (bWrite
&& (frout
.time
>= begin
) )
881 first_time
= frout
.time
;
883 lasttime
= frout
.time
;
885 if (dt
== 0 || bRmod(frout
.time
, first_time
, dt
))
888 last_ok_t
= frout
.time
;
891 fprintf(stderr
, "\nContinue writing frames from %s t=%g %s, "
893 fnms
[i
], output_env_conv_time(oenv
, frout
.time
), output_env_get_time_unit(oenv
),
900 write_trxframe_indexed(trxout
, &frout
, isize
, index
,
905 write_trxframe(trxout
, &frout
, NULL
);
907 if ( ((frame
% 10) == 0) || (frame
< 10) )
909 fprintf(stderr
, " -> frame %6d time %8.3f %s \r",
910 frame_out
, output_env_conv_time(oenv
, frout
.time
), output_env_get_time_unit(oenv
));
915 while (read_next_frame(oenv
, status
, &fr
));
923 fprintf(stderr
, "\nLast frame written was %d, time %f %s\n",
924 frame
, output_env_conv_time(oenv
, last_ok_t
), output_env_get_time_unit(oenv
));