2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
46 #include "gromacs/utility/smalloc.h"
48 #include "gromacs/fileio/gmxfio.h"
49 #include "gromacs/fileio/tpxio.h"
50 #include "gromacs/fileio/trxio.h"
51 #include "gromacs/fileio/trnio.h"
52 #include "gromacs/fileio/tngio.h"
53 #include "gromacs/fileio/tngio_for_tools.h"
54 #include "gromacs/commandline/pargs.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/fileio/pdbio.h"
57 #include "gromacs/fileio/confio.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/fileio/xtcio.h"
62 #include "gromacs/fileio/xvgr.h"
65 #include "gromacs/utility/fatalerror.h"
67 #define TIME_EXPLICIT 0
68 #define TIME_CONTINUE 1
73 #define FLAGS (TRX_READ_X | TRX_READ_V | TRX_READ_F)
75 static void scan_trj_files(char **fnms
, int nfiles
, real
*readtime
,
76 real
*timestep
, atom_id imax
,
77 const output_env_t oenv
)
79 /* Check start time of all files */
86 for (i
= 0; i
< nfiles
; i
++)
88 ok
= read_first_frame(oenv
, &status
, fnms
[i
], &fr
, FLAGS
);
92 gmx_fatal(FARGS
, "\nCouldn't read frame from file." );
96 readtime
[i
] = fr
.time
;
101 fprintf(stderr
, "\nWARNING: Couldn't find a time in the frame.\n");
112 if (natoms
!= fr
.natoms
)
114 gmx_fatal(FARGS
, "\nDifferent numbers of atoms (%d/%d) in files",
120 if (fr
.natoms
<= imax
)
122 gmx_fatal(FARGS
, "\nNot enough atoms (%d) for index group (%d)",
127 ok
= read_next_frame(oenv
, status
, &fr
);
130 timestep
[i
] = fr
.time
- readtime
[i
];
151 fprintf(stderr
, "\n");
155 static void sort_files(char **fnms
, real
*settime
, int nfile
)
161 for (i
= 0; i
< nfile
; i
++)
164 for (j
= i
+ 1; j
< nfile
; j
++)
166 if (settime
[j
] < settime
[minidx
])
173 timeswap
= settime
[i
];
174 settime
[i
] = settime
[minidx
];
175 settime
[minidx
] = timeswap
;
177 fnms
[i
] = fnms
[minidx
];
178 fnms
[minidx
] = chptr
;
183 static void edit_files(char **fnms
, int nfiles
, real
*readtime
, real
*timestep
,
184 real
*settime
, int *cont_type
, gmx_bool bSetTime
,
185 gmx_bool bSort
, const output_env_t oenv
)
189 char inputstring
[STRLEN
], *chptr
;
193 fprintf(stderr
, "\n\nEnter the new start time (%s) for each file.\n"
194 "There are two special options, both disable sorting:\n\n"
195 "c (continue) - The start time is taken from the end\n"
196 "of the previous file. Use it when your continuation run\n"
197 "restarts with t=0.\n\n"
198 "l (last) - The time in this file will be changed the\n"
199 "same amount as in the previous. Use it when the time in the\n"
200 "new run continues from the end of the previous one,\n"
201 "since this takes possible overlap into account.\n\n",
202 output_env_get_time_unit(oenv
));
206 " File Current start (%s) New start (%s)\n"
207 "---------------------------------------------------------\n",
208 output_env_get_time_unit(oenv
), output_env_get_time_unit(oenv
));
210 for (i
= 0; i
< nfiles
; i
++)
212 fprintf(stderr
, "%25s %10.3f %s ", fnms
[i
],
213 output_env_conv_time(oenv
, readtime
[i
]), output_env_get_time_unit(oenv
));
217 if (NULL
== fgets(inputstring
, STRLEN
- 1, stdin
))
219 gmx_fatal(FARGS
, "Error reading user input" );
222 inputstring
[strlen(inputstring
)-1] = 0;
224 if (inputstring
[0] == 'c' || inputstring
[0] == 'C')
226 cont_type
[i
] = TIME_CONTINUE
;
229 settime
[i
] = FLT_MAX
;
231 else if (inputstring
[0] == 'l' ||
232 inputstring
[0] == 'L')
234 cont_type
[i
] = TIME_LAST
;
237 settime
[i
] = FLT_MAX
;
241 settime
[i
] = strtod(inputstring
, &chptr
)*
242 output_env_get_time_invfactor(oenv
);
243 if (chptr
== inputstring
)
245 fprintf(stderr
, "'%s' not recognized as a floating point number, 'c' or 'l'. "
246 "Try again: ", inputstring
);
250 cont_type
[i
] = TIME_EXPLICIT
;
257 if (cont_type
[0] != TIME_EXPLICIT
)
259 cont_type
[0] = TIME_EXPLICIT
;
265 for (i
= 0; i
< nfiles
; i
++)
267 settime
[i
] = readtime
[i
];
272 fprintf(stderr
, "Sorting disabled.\n");
276 sort_files(fnms
, settime
, nfiles
);
278 /* Write out the new order and start times */
279 fprintf(stderr
, "\nSummary of files and start times used:\n\n"
280 " File Start time Time step\n"
281 "---------------------------------------------------------\n");
282 for (i
= 0; i
< nfiles
; i
++)
284 switch (cont_type
[i
])
287 fprintf(stderr
, "%25s %10.3f %s %10.3f %s",
289 output_env_conv_time(oenv
, settime
[i
]), output_env_get_time_unit(oenv
),
290 output_env_conv_time(oenv
, timestep
[i
]), output_env_get_time_unit(oenv
));
292 cont_type
[i
-1] == TIME_EXPLICIT
&& settime
[i
] == settime
[i
-1])
294 fprintf(stderr
, " WARNING: same Start time as previous");
296 fprintf(stderr
, "\n");
299 fprintf(stderr
, "%25s Continue from last file\n", fnms
[i
]);
302 fprintf(stderr
, "%25s Change by same amount as last file\n",
307 fprintf(stderr
, "\n");
309 settime
[nfiles
] = FLT_MAX
;
310 cont_type
[nfiles
] = TIME_EXPLICIT
;
311 readtime
[nfiles
] = FLT_MAX
;
314 static void do_demux(int nset
, char *fnms
[], char *fnms_out
[], int nval
,
315 real
**value
, real
*time
, real dt_remd
, int isize
,
316 atom_id index
[], real dt
, const output_env_t oenv
)
318 int i
, j
, k
, natoms
, nnn
;
319 t_trxstatus
**fp_in
, **fp_out
;
320 gmx_bool bCont
, *bSet
;
321 real t
, first_time
= 0;
329 for (i
= 0; (i
< nset
); i
++)
331 nnn
= read_first_frame(oenv
, &(fp_in
[i
]), fnms
[i
], &(trx
[i
]),
336 first_time
= trx
[i
].time
;
338 else if (natoms
!= nnn
)
340 gmx_fatal(FARGS
, "Trajectory file %s has %d atoms while previous trajs had %d atoms", fnms
[i
], nnn
, natoms
);
346 else if (t
!= trx
[i
].time
)
348 gmx_fatal(FARGS
, "Trajectory file %s has time %f while previous trajs had time %f", fnms
[i
], trx
[i
].time
, t
);
353 for (i
= 0; (i
< nset
); i
++)
355 fp_out
[i
] = open_trx(fnms_out
[i
], "w");
358 if (gmx_nint(time
[k
] - t
) != 0)
360 gmx_fatal(FARGS
, "First time in demuxing table does not match trajectories");
364 while ((k
+1 < nval
) && ((trx
[0].time
- time
[k
+1]) > dt_remd
*0.1))
370 fprintf(debug
, "trx[0].time = %g, time[k] = %g\n", trx
[0].time
, time
[k
]);
372 for (i
= 0; (i
< nset
); i
++)
376 for (i
= 0; (i
< nset
); i
++)
378 j
= gmx_nint(value
[i
][k
]);
379 range_check(j
, 0, nset
);
382 gmx_fatal(FARGS
, "Demuxing the same replica %d twice at time %f",
387 if (dt
== 0 || bRmod(trx
[i
].time
, first_time
, dt
))
391 write_trxframe_indexed(fp_out
[j
], &trx
[i
], isize
, index
, NULL
);
395 write_trxframe(fp_out
[j
], &trx
[i
], NULL
);
401 for (i
= 0; (i
< nset
); i
++)
403 bCont
= bCont
&& read_next_frame(oenv
, fp_in
[i
], &trx
[i
]);
408 for (i
= 0; (i
< nset
); i
++)
411 close_trx(fp_out
[i
]);
415 int gmx_trjcat(int argc
, char *argv
[])
419 "[THISMODULE] concatenates several input trajectory files in sorted order. ",
420 "In case of double time frames the one in the later file is used. ",
421 "By specifying [TT]-settime[tt] you will be asked for the start time ",
422 "of each file. The input files are taken from the command line, ",
423 "such that a command like [TT]gmx trjcat -f *.trr -o fixed.trr[tt] should do ",
424 "the trick. Using [TT]-cat[tt], you can simply paste several files ",
425 "together without removal of frames with identical time stamps.[PAR]",
426 "One important option is inferred when the output file is amongst the",
427 "input files. In that case that particular file will be appended to",
428 "which implies you do not need to store double the amount of data.",
429 "Obviously the file to append to has to be the one with lowest starting",
430 "time since one can only append at the end of a file.[PAR]",
431 "If the [TT]-demux[tt] option is given, the N trajectories that are",
432 "read, are written in another order as specified in the [TT].xvg[tt] file.",
433 "The [TT].xvg[tt] file should contain something like:[PAR]",
434 "[TT]0 0 1 2 3 4 5[BR]",
435 "2 1 0 2 3 5 4[tt][BR]",
436 "Where the first number is the time, and subsequent numbers point to",
437 "trajectory indices.",
438 "The frames corresponding to the numbers present at the first line",
439 "are collected into the output trajectory. If the number of frames in",
440 "the trajectory does not match that in the [TT].xvg[tt] file then the program",
441 "tries to be smart. Beware."
443 static gmx_bool bVels
= TRUE
;
444 static gmx_bool bCat
= FALSE
;
445 static gmx_bool bSort
= TRUE
;
446 static gmx_bool bKeepLast
= FALSE
;
447 static gmx_bool bKeepLastAppend
= FALSE
;
448 static gmx_bool bOverwrite
= FALSE
;
449 static gmx_bool bSetTime
= FALSE
;
450 static gmx_bool bDeMux
;
451 static real begin
= -1;
452 static real end
= -1;
458 { "-b", FALSE
, etTIME
,
459 { &begin
}, "First time to use (%t)" },
460 { "-e", FALSE
, etTIME
,
461 { &end
}, "Last time to use (%t)" },
462 { "-dt", FALSE
, etTIME
,
463 { &dt
}, "Only write frame when t MOD dt = first time (%t)" },
464 { "-vel", FALSE
, etBOOL
,
465 { &bVels
}, "Read and write velocities if possible" },
466 { "-settime", FALSE
, etBOOL
,
467 { &bSetTime
}, "Change starting time interactively" },
468 { "-sort", FALSE
, etBOOL
,
469 { &bSort
}, "Sort trajectory files (not frames)" },
470 { "-keeplast", FALSE
, etBOOL
,
471 { &bKeepLast
}, "Keep overlapping frames at end of trajectory" },
472 { "-overwrite", FALSE
, etBOOL
,
473 { &bOverwrite
}, "Overwrite overlapping frames during appending" },
474 { "-cat", FALSE
, etBOOL
,
475 { &bCat
}, "Do not discard double time frames" }
477 #define npargs asize(pa)
478 int ftpin
, i
, frame
, frame_out
, step
= 0, trjout
= 0;
479 t_trxstatus
*status
, *trxout
= NULL
;
482 t_trxframe fr
, frout
;
483 char **fnms
, **fnms_out
, *in_file
, *out_file
;
485 gmx_bool bNewFile
, bIndex
, bWrite
;
486 int earliersteps
, nfile_in
, nfile_out
, *cont_type
, last_ok_step
;
487 real
*readtime
, *timest
, *settime
;
488 real first_time
= 0, lasttime
= NOTSET
, last_ok_t
= -1, timestep
;
489 real last_frame_time
, searchtime
;
491 atom_id
*index
= NULL
, imax
;
493 real
**val
= NULL
, *t
= NULL
, dt_remd
;
494 int n
, nset
, ftpout
= -1, prevEndStep
= 0, filetype
;
500 { efTRX
, "-f", NULL
, ffRDMULT
},
501 { efTRO
, "-o", NULL
, ffWRMULT
},
502 { efNDX
, "-n", "index", ffOPTRD
},
503 { efXVG
, "-demux", "remd", ffOPTRD
}
506 #define NFILE asize(fnm)
508 if (!parse_common_args(&argc
, argv
, PCA_BE_NICE
| PCA_TIME_UNIT
, NFILE
, fnm
,
509 asize(pa
), pa
, asize(desc
), desc
, 0, NULL
, &oenv
))
514 bIndex
= ftp2bSet(efNDX
, NFILE
, fnm
);
515 bDeMux
= ftp2bSet(efXVG
, NFILE
, fnm
);
516 bSort
= bSort
&& !bDeMux
;
521 printf("Select group for output\n");
522 rd_index(ftp2fn(efNDX
, NFILE
, fnm
), 1, &isize
, &index
, &grpname
);
525 for (i
= 1; i
< isize
; i
++)
527 imax
= max(imax
, index
[i
]);
534 val
= read_xvg_time(opt2fn("-demux", NFILE
, fnm
), TRUE
,
535 opt2parg_bSet("-b", npargs
, pa
), begin
,
536 opt2parg_bSet("-e", npargs
, pa
), end
, 1, &nset
, &n
,
538 printf("Read %d sets of %d points, dt = %g\n\n", nset
, n
, dt_remd
);
541 fprintf(debug
, "Dump of replica_index.xvg\n");
542 for (i
= 0; (i
< n
); i
++)
544 fprintf(debug
, "%10g", t
[i
]);
545 for (j
= 0; (j
< nset
); j
++)
547 fprintf(debug
, " %3d", gmx_nint(val
[j
][i
]));
549 fprintf(debug
, "\n");
554 nfile_in
= opt2fns(&fnms
, "-f", NFILE
, fnm
);
557 gmx_fatal(FARGS
, "No input files!" );
560 if (bDeMux
&& (nfile_in
!= nset
))
562 gmx_fatal(FARGS
, "You have specified %d files and %d entries in the demux table", nfile_in
, nset
);
565 ftpin
= fn2ftp(fnms
[0]);
567 for (i
= 1; i
< nfile_in
; i
++)
569 if (ftpin
!= fn2ftp(fnms
[i
]))
571 gmx_fatal(FARGS
, "All input files must be of the same format");
575 nfile_out
= opt2fns(&fnms_out
, "-o", NFILE
, fnm
);
578 gmx_fatal(FARGS
, "No output files!");
580 if ((nfile_out
> 1) && !bDeMux
)
582 gmx_fatal(FARGS
, "Don't know what to do with more than 1 output file if not demultiplexing");
584 else if (bDeMux
&& (nfile_out
!= nset
) && (nfile_out
!= 1))
586 gmx_fatal(FARGS
, "Number of output files should be 1 or %d (#input files), not %d", nset
, nfile_out
);
590 if (nfile_out
!= nset
)
592 char *buf
= strdup(fnms_out
[0]);
593 snew(fnms_out
, nset
);
594 for (i
= 0; (i
< nset
); i
++)
596 snew(fnms_out
[i
], strlen(buf
)+32);
597 sprintf(fnms_out
[i
], "%d_%s", i
, buf
);
601 do_demux(nfile_in
, fnms
, fnms_out
, n
, val
, t
, dt_remd
, isize
, index
, dt
, oenv
);
605 snew(readtime
, nfile_in
+1);
606 snew(timest
, nfile_in
+1);
607 scan_trj_files(fnms
, nfile_in
, readtime
, timest
, imax
, oenv
);
609 snew(settime
, nfile_in
+1);
610 snew(cont_type
, nfile_in
+1);
611 edit_files(fnms
, nfile_in
, readtime
, timest
, settime
, cont_type
, bSetTime
, bSort
,
614 /* Check whether the output file is amongst the input files
615 * This has to be done after sorting etc.
617 out_file
= fnms_out
[0];
618 ftpout
= fn2ftp(out_file
);
620 for (i
= 0; ((i
< nfile_in
) && (n_append
== -1)); i
++)
622 if (strcmp(fnms
[i
], out_file
) == 0)
629 fprintf(stderr
, "Will append to %s rather than creating a new file\n",
632 else if (n_append
!= -1)
634 gmx_fatal(FARGS
, "Can only append to the first file which is %s (not %s)",
639 /* Not checking input format, could be dangerous :-) */
640 /* Not checking output format, equally dangerous :-) */
644 /* the default is not to change the time at all,
645 * but this is overridden by the edit_files routine
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 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
))
704 bKeepLastAppend
= TRUE
;
706 trxout
= open_trx(out_file
, "a");
710 if (gmx_fio_getftp(stfio
) != efXTC
)
712 gmx_fatal(FARGS
, "Overwrite only supported for XTC." );
714 last_frame_time
= trx_get_time_of_final_frame(status
);
716 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
717 * or when seek time = 0 */
718 if (nfile_in
> 1 && settime
[1] < last_frame_time
+timest
[0]*0.5)
720 /* Jump to one time-frame before the start of next
722 searchtime
= settime
[1]-timest
[0]*1.25;
726 searchtime
= last_frame_time
;
728 if (xtc_seek_time(stfio
, searchtime
, fr
.natoms
, TRUE
))
730 gmx_fatal(FARGS
, "Error seeking to append position.");
732 read_next_frame(oenv
, status
, &fr
);
733 if (fabs(searchtime
- fr
.time
) > timest
[0]*0.5)
735 gmx_fatal(FARGS
, "Error seeking: attempted to seek to %f but got %f.",
736 searchtime
, fr
.time
);
739 fpos
= gmx_fio_ftell(stfio
);
741 trxout
= open_trx(out_file
, "r+");
742 if (gmx_fio_seek(trx_get_fileio(trxout
), fpos
))
744 gmx_fatal(FARGS
, "Error seeking to append position.");
747 printf("\n Will append after %f \n", lasttime
);
750 /* Lets stitch up some files */
751 timestep
= timest
[0];
752 for (i
= n_append
+1; (i
< nfile_in
); i
++)
756 /* set the next time from the last frame in previous file */
759 /* When writing TNG the step determine which frame to write. Use an
760 * offset to be able to increase steps properly when changing files. */
763 prevEndStep
= frout
.step
;
768 if (cont_type
[i
] == TIME_CONTINUE
)
771 begin
+= 0.5*timestep
;
772 settime
[i
] = frout
.time
;
773 cont_type
[i
] = TIME_EXPLICIT
;
775 else if (cont_type
[i
] == TIME_LAST
)
778 begin
+= 0.5*timestep
;
780 /* Or, if the time in the next part should be changed by the
781 * same amount, start at half a timestep from the last time
782 * so we dont repeat frames.
784 /* I don't understand the comment above, but for all the cases
785 * I tried the code seems to work properly. B. Hess 2008-4-2.
788 /* Or, if time is set explicitly, we check for overlap/gap */
789 if (cont_type
[i
] == TIME_EXPLICIT
)
791 if ( ( i
< nfile_in
) &&
792 ( frout
.time
< settime
[i
]-1.5*timestep
) )
794 fprintf(stderr
, "WARNING: Frames around t=%f %s have a different "
795 "spacing than the rest,\n"
796 "might be a gap or overlap that couldn't be corrected "
797 "automatically.\n", output_env_conv_time(oenv
, frout
.time
),
798 output_env_get_time_unit(oenv
));
803 /* if we don't have a timestep in the current file, use the old one */
806 timestep
= timest
[i
];
808 read_first_frame(oenv
, &status
, fnms
[i
], &fr
, FLAGS
);
812 fprintf(stderr
, "\nWARNING: Couldn't find a time in the frame.\n");
815 if (cont_type
[i
] == TIME_EXPLICIT
)
817 t_corr
= settime
[i
]-fr
.time
;
819 /* t_corr is the amount we want to change the time.
820 * If the user has chosen not to change the time for
821 * this part of the trajectory t_corr remains at
822 * the value it had in the last part, changing this
823 * by the same amount.
824 * If no value was given for the first trajectory part
825 * we let the time start at zero, see the edit_files routine.
831 if (lasttime
!= NOTSET
)
833 printf("lasttime %g\n", lasttime
);
838 /* copy the input frame to the output frame */
840 /* set the new time by adding the correct calculated above */
841 frout
.time
+= t_corr
;
844 frout
.step
+= prevEndStep
;
846 /* quit if we have reached the end of what should be written */
847 if ((end
> 0) && (frout
.time
> end
+GMX_REAL_EPS
))
853 /* determine if we should write this frame (dt is handled elsewhere) */
854 if (bCat
) /* write all frames of all files */
858 else if (bKeepLast
|| (bKeepLastAppend
&& i
== 1))
859 /* write till last frame of this traj
860 and skip first frame(s) of next traj */
862 bWrite
= ( frout
.time
> lasttime
+0.5*timestep
);
864 else /* write till first frame of next traj */
866 bWrite
= ( frout
.time
< settime
[i
+1]-0.5*timestep
);
869 if (bWrite
&& (frout
.time
>= begin
) )
874 first_time
= frout
.time
;
876 lasttime
= frout
.time
;
877 if (dt
== 0 || bRmod(frout
.time
, first_time
, dt
))
880 last_ok_t
= frout
.time
;
883 fprintf(stderr
, "\nContinue writing frames from %s t=%g %s, "
885 fnms
[i
], output_env_conv_time(oenv
, frout
.time
), output_env_get_time_unit(oenv
),
892 write_trxframe_indexed(trxout
, &frout
, isize
, index
,
897 write_trxframe(trxout
, &frout
, NULL
);
899 if ( ((frame
% 10) == 0) || (frame
< 10) )
901 fprintf(stderr
, " -> frame %6d time %8.3f %s \r",
902 frame_out
, output_env_conv_time(oenv
, frout
.time
), output_env_get_time_unit(oenv
));
907 while (read_next_frame(oenv
, status
, &fr
));
911 earliersteps
+= step
;
917 fprintf(stderr
, "\nLast frame written was %d, time %f %s\n",
918 frame
, output_env_conv_time(oenv
, last_ok_t
), output_env_get_time_unit(oenv
));