2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/fileio/gmxfio.h"
49 #include "gromacs/fileio/pdbio.h"
50 #include "gromacs/fileio/tngio.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/arrayref.h"
60 #include "gromacs/utility/arraysize.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/futil.h"
63 #include "gromacs/utility/smalloc.h"
64 #include "gromacs/utility/stringutil.h"
66 #define TIME_EXPLICIT 0
67 #define TIME_CONTINUE 1
72 #define FLAGS (TRX_READ_X | TRX_READ_V | TRX_READ_F)
74 static void scan_trj_files(gmx::ArrayRef<const std::string> files,
76 real *timestep, int imax,
77 const gmx_output_env_t *oenv)
79 /* Check start time of all files */
85 for (gmx::index i = 0; i < files.ssize(); i++)
87 ok = read_first_frame(oenv, &status, files[i].c_str(), &fr, FLAGS);
91 gmx_fatal(FARGS, "\nCouldn't read frame from file." );
95 readtime[i] = fr.time;
100 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
111 if (natoms != fr.natoms)
113 gmx_fatal(FARGS, "\nDifferent numbers of atoms (%d/%d) in files",
119 if (fr.natoms <= imax)
121 gmx_fatal(FARGS, "\nNot enough atoms (%d) for index group (%d)",
126 ok = read_next_frame(oenv, status, &fr);
129 timestep[i] = fr.time - readtime[i];
150 fprintf(stderr, "\n");
154 static void sort_files(gmx::ArrayRef<std::string> files, real *settime)
156 for (gmx::index i = 0; i < files.ssize(); i++)
158 gmx::index minidx = i;
159 for (gmx::index j = i + 1; j < files.ssize(); j++)
161 if (settime[j] < settime[minidx])
168 real timeswap = settime[i];
169 settime[i] = settime[minidx];
170 settime[minidx] = timeswap;
171 std::swap(files[i], files[minidx]);
176 static void edit_files(gmx::ArrayRef<std::string> files,
177 real *readtime, real *timestep,
178 real *settime, int *cont_type, gmx_bool bSetTime,
179 gmx_bool bSort, const gmx_output_env_t *oenv)
182 char inputstring[STRLEN], *chptr;
184 auto timeUnit = output_env_get_time_unit(oenv);
187 fprintf(stderr, "\n\nEnter the new start time (%s) for each file.\n"
188 "There are two special options, both disable sorting:\n\n"
189 "c (continue) - The start time is taken from the end\n"
190 "of the previous file. Use it when your continuation run\n"
191 "restarts with t=0.\n\n"
192 "l (last) - The time in this file will be changed the\n"
193 "same amount as in the previous. Use it when the time in the\n"
194 "new run continues from the end of the previous one,\n"
195 "since this takes possible overlap into account.\n\n",
200 " File Current start (%s) New start (%s)\n"
201 "---------------------------------------------------------\n",
202 timeUnit.c_str(), timeUnit.c_str());
204 for (gmx::index i = 0; i < files.ssize(); i++)
206 fprintf(stderr, "%25s %10.3f %s ", files[i].c_str(),
207 output_env_conv_time(oenv, readtime[i]), timeUnit.c_str());
211 if (nullptr == fgets(inputstring, STRLEN - 1, stdin))
213 gmx_fatal(FARGS, "Error reading user input" );
216 inputstring[std::strlen(inputstring)-1] = 0;
218 if (inputstring[0] == 'c' || inputstring[0] == 'C')
220 cont_type[i] = TIME_CONTINUE;
223 settime[i] = FLT_MAX;
225 else if (inputstring[0] == 'l' ||
226 inputstring[0] == 'L')
228 cont_type[i] = TIME_LAST;
231 settime[i] = FLT_MAX;
235 settime[i] = strtod(inputstring, &chptr)*
236 output_env_get_time_invfactor(oenv);
237 if (chptr == inputstring)
239 fprintf(stderr, "'%s' not recognized as a floating point number, 'c' or 'l'. "
240 "Try again: ", inputstring);
244 cont_type[i] = TIME_EXPLICIT;
251 if (cont_type[0] != TIME_EXPLICIT)
253 cont_type[0] = TIME_EXPLICIT;
259 for (gmx::index i = 0; i < files.ssize(); i++)
261 settime[i] = readtime[i];
266 fprintf(stderr, "Sorting disabled.\n");
270 sort_files(files, settime);
272 /* Write out the new order and start times */
273 fprintf(stderr, "\nSummary of files and start times used:\n\n"
274 " File Start time Time step\n"
275 "---------------------------------------------------------\n");
276 for (gmx::index i = 0; i < files.ssize(); i++)
278 switch (cont_type[i])
281 fprintf(stderr, "%25s %10.3f %s %10.3f %s",
283 output_env_conv_time(oenv, settime[i]), timeUnit.c_str(),
284 output_env_conv_time(oenv, timestep[i]), timeUnit.c_str());
286 cont_type[i-1] == TIME_EXPLICIT && settime[i] == settime[i-1])
288 fprintf(stderr, " WARNING: same Start time as previous");
290 fprintf(stderr, "\n");
293 fprintf(stderr, "%25s Continue from last file\n", files[i].c_str());
296 fprintf(stderr, "%25s Change by same amount as last file\n",
301 fprintf(stderr, "\n");
303 settime[files.size()] = FLT_MAX;
304 cont_type[files.size()] = TIME_EXPLICIT;
305 readtime[files.size()] = FLT_MAX;
308 static void do_demux(gmx::ArrayRef<const std::string> inFiles,
309 gmx::ArrayRef<const std::string> outFiles, int nval,
310 real **value, real *time, real dt_remd, int isize,
311 int index[], real dt, const gmx_output_env_t *oenv)
314 t_trxstatus **fp_in, **fp_out;
315 gmx_bool bCont, *bSet;
316 real t, first_time = 0;
319 snew(fp_in, inFiles.size());
320 snew(trx, inFiles.size());
321 snew(bSet, inFiles.size());
324 for (gmx::index i = 0; i < inFiles.ssize(); i++)
326 read_first_frame(oenv, &(fp_in[i]), inFiles[i].c_str(), &(trx[i]),
330 natoms = trx[i].natoms;
331 first_time = trx[i].time;
333 else if (natoms != trx[i].natoms)
335 gmx_fatal(FARGS, "Trajectory file %s has %d atoms while previous trajs had %d atoms", inFiles[i].c_str(), trx[i].natoms, natoms);
341 else if (t != trx[i].time)
343 gmx_fatal(FARGS, "Trajectory file %s has time %f while previous trajs had time %f", inFiles[i].c_str(), trx[i].time, t);
347 snew(fp_out, inFiles.size());
348 for (gmx::index i = 0; i < inFiles.ssize(); i++)
350 fp_out[i] = open_trx(outFiles[i].c_str(), "w");
353 if (std::round(time[k] - t) != 0)
355 gmx_fatal(FARGS, "First time in demuxing table does not match trajectories");
359 while ((k+1 < nval) && ((trx[0].time - time[k+1]) > dt_remd*0.1))
365 fprintf(debug, "trx[0].time = %g, time[k] = %g\n", trx[0].time, time[k]);
367 for (gmx::index i = 0; i < inFiles.ssize(); i++)
371 for (gmx::index i = 0; i < inFiles.ssize(); i++)
373 int j = gmx::roundToInt(value[i][k]);
374 range_check(j, 0, inFiles.size());
377 gmx_fatal(FARGS, "Demuxing the same replica %d twice at time %f",
382 if (dt == 0 || bRmod(trx[i].time, first_time, dt))
386 write_trxframe_indexed(fp_out[j], &trx[i], isize, index, nullptr);
390 write_trxframe(fp_out[j], &trx[i], nullptr);
396 for (gmx::index i = 0; i < inFiles.ssize(); i++)
398 bCont = bCont && read_next_frame(oenv, fp_in[i], &trx[i]);
403 for (gmx::index i = 0; i < inFiles.ssize(); i++)
406 close_trx(fp_out[i]);
410 int gmx_trjcat(int argc, char *argv[])
414 "[THISMODULE] concatenates several input trajectory files in sorted order. ",
415 "In case of double time frames the one in the later file is used. ",
416 "By specifying [TT]-settime[tt] you will be asked for the start time ",
417 "of each file. The input files are taken from the command line, ",
418 "such that a command like [TT]gmx trjcat -f *.trr -o fixed.trr[tt] should do ",
419 "the trick. Using [TT]-cat[tt], you can simply paste several files ",
420 "together without removal of frames with identical time stamps.[PAR]",
421 "One important option is inferred when the output file is amongst the",
422 "input files. In that case that particular file will be appended to",
423 "which implies you do not need to store double the amount of data.",
424 "Obviously the file to append to has to be the one with lowest starting",
425 "time since one can only append at the end of a file.[PAR]",
426 "If the [TT]-demux[tt] option is given, the N trajectories that are",
427 "read, are written in another order as specified in the [REF].xvg[ref] file.",
428 "The [REF].xvg[ref] file should contain something like::",
433 "The first number is the time, and subsequent numbers point to",
434 "trajectory indices.",
435 "The frames corresponding to the numbers present at the first line",
436 "are collected into the output trajectory. If the number of frames in",
437 "the trajectory does not match that in the [REF].xvg[ref] file then the program",
438 "tries to be smart. Beware."
440 static gmx_bool bCat = FALSE;
441 static gmx_bool bSort = TRUE;
442 static gmx_bool bKeepLast = FALSE;
443 static gmx_bool bKeepLastAppend = FALSE;
444 static gmx_bool bOverwrite = FALSE;
445 static gmx_bool bSetTime = FALSE;
446 static gmx_bool bDeMux;
447 static real begin = -1;
448 static real end = -1;
454 { "-b", FALSE, etTIME,
455 { &begin }, "First time to use (%t)" },
456 { "-e", FALSE, etTIME,
457 { &end }, "Last time to use (%t)" },
458 { "-dt", FALSE, etTIME,
459 { &dt }, "Only write frame when t MOD dt = first time (%t)" },
460 { "-settime", FALSE, etBOOL,
461 { &bSetTime }, "Change starting time interactively" },
462 { "-sort", FALSE, etBOOL,
463 { &bSort }, "Sort trajectory files (not frames)" },
464 { "-keeplast", FALSE, etBOOL,
465 { &bKeepLast }, "Keep overlapping frames at end of trajectory" },
466 { "-overwrite", FALSE, etBOOL,
467 { &bOverwrite }, "Overwrite overlapping frames during appending" },
468 { "-cat", FALSE, etBOOL,
469 { &bCat }, "Do not discard double time frames" }
471 #define npargs asize(pa)
472 int ftpin, i, frame, frame_out;
473 t_trxstatus *status, *trxout = nullptr;
475 t_trxframe fr, frout;
477 gmx_bool bNewFile, bIndex, bWrite;
479 real *readtime, *timest, *settime;
480 real first_time = 0, lasttime = 0, last_ok_t = -1, timestep;
481 gmx_bool lastTimeSet = FALSE;
482 real last_frame_time, searchtime;
484 int *index = nullptr, imax;
486 real **val = nullptr, *t = nullptr, dt_remd;
487 int n, nset, ftpout = -1, prevEndStep = 0, filetype;
489 gmx_output_env_t *oenv;
492 { efTRX, "-f", nullptr, ffRDMULT },
493 { efTRO, "-o", nullptr, ffWRMULT },
494 { efNDX, "-n", "index", ffOPTRD },
495 { efXVG, "-demux", "remd", ffOPTRD }
498 #define NFILE asize(fnm)
500 if (!parse_common_args(&argc, argv, PCA_TIME_UNIT, NFILE, fnm,
501 asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
505 fprintf(stdout, "Note that major changes are planned in future for "
506 "trjcat, to improve usability and utility.");
508 auto timeUnit = output_env_get_time_unit(oenv);
510 bIndex = ftp2bSet(efNDX, NFILE, fnm);
511 bDeMux = ftp2bSet(efXVG, NFILE, fnm);
512 bSort = bSort && !bDeMux;
517 printf("Select group for output\n");
518 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
521 for (i = 1; i < isize; i++)
523 imax = std::max(imax, index[i]);
530 val = read_xvg_time(opt2fn("-demux", NFILE, fnm), TRUE,
531 opt2parg_bSet("-b", npargs, pa), begin,
532 opt2parg_bSet("-e", npargs, pa), end, 1, &nset, &n,
534 printf("Read %d sets of %d points, dt = %g\n\n", nset, n, dt_remd);
537 fprintf(debug, "Dump of replica_index.xvg\n");
538 for (i = 0; (i < n); i++)
540 fprintf(debug, "%10g", t[i]);
541 for (j = 0; (j < nset); j++)
543 fprintf(debug, " %3d", static_cast<int>(std::round(val[j][i])));
545 fprintf(debug, "\n");
550 gmx::ArrayRef<const std::string> inFiles = opt2fns("-f", NFILE, fnm);
553 gmx_fatal(FARGS, "No input files!" );
556 if (bDeMux && ssize(inFiles) != nset)
558 gmx_fatal(FARGS, "You have specified %td files and %d entries in the demux table", inFiles.ssize(), nset);
561 ftpin = fn2ftp(inFiles[0].c_str());
563 if (ftpin != efTRR && ftpin != efXTC && ftpin != efTNG)
565 gmx_fatal(FARGS, "gmx trjcat can only handle binary trajectory formats (trr, xtc, tng)");
568 for (const std::string &inFile : inFiles)
570 if (ftpin != fn2ftp(inFile.c_str()))
572 gmx_fatal(FARGS, "All input files must be of the same (trr, xtc or tng) format");
576 gmx::ArrayRef<const std::string> outFiles = opt2fns("-o", NFILE, fnm);
577 if (outFiles.empty())
579 gmx_fatal(FARGS, "No output files!");
581 if ((outFiles.size() > 1) && !bDeMux)
583 gmx_fatal(FARGS, "Don't know what to do with more than 1 output file if not demultiplexing");
585 else if (bDeMux && ssize(outFiles) != nset && outFiles.size() != 1)
587 gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %td", nset, outFiles.ssize());
591 auto outFilesDemux = gmx::copyOf(outFiles);
592 if (gmx::ssize(outFilesDemux) != nset)
594 std::string name = outFilesDemux[0];
595 outFilesDemux.resize(nset);
596 for (i = 0; (i < nset); i++)
598 outFilesDemux[0] = gmx::formatString("%d_%s", i, name.c_str());
601 do_demux(inFiles, outFilesDemux, n, val, t, dt_remd, isize, index, dt, oenv);
605 snew(readtime, inFiles.size() + 1);
606 snew(timest, inFiles.size() + 1);
607 scan_trj_files(inFiles, readtime, timest, imax, oenv);
609 snew(settime, inFiles.size() + 1);
610 snew(cont_type, inFiles.size() + 1);
611 auto inFilesEdited = gmx::copyOf(inFiles);
612 edit_files(inFilesEdited, readtime, timest, settime, cont_type, bSetTime, bSort,
615 /* Check whether the output file is amongst the input files
616 * This has to be done after sorting etc.
618 const char *out_file = outFiles[0].c_str();
619 ftpout = fn2ftp(out_file);
621 for (size_t i = 0; i < inFilesEdited.size() && n_append == -1; i++)
623 if (std::strcmp(inFilesEdited[i].c_str(), out_file) == 0)
630 fprintf(stderr, "Will append to %s rather than creating a new file\n",
633 else if (n_append != -1)
635 gmx_fatal(FARGS, "Can only append to the first file which is %s (not %s)",
636 inFilesEdited[0].c_str(), out_file);
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 gmx_fatal(FARGS, "When writing TNG the input file format must also be TNG");
659 trxout = trjtools_gmx_prepare_tng_writing(out_file, 'w', nullptr,
660 inFilesEdited[0].c_str(), isize, nullptr, gmx::arrayRefFromArray(index, isize), grpname);
664 trxout = trjtools_gmx_prepare_tng_writing(out_file, 'w', nullptr,
665 inFilesEdited[0].c_str(), -1, nullptr, gmx::EmptyArrayRef(), nullptr);
670 trxout = open_trx(out_file, "w");
672 std::memset(&frout, 0, sizeof(frout));
678 if (!read_first_frame(oenv, &status, out_file, &fr, FLAGS))
680 gmx_fatal(FARGS, "Reading first frame from %s", out_file);
683 stfio = trx_get_fileio(status);
684 if (!bKeepLast && !bOverwrite)
686 fprintf(stderr, "\n\nWARNING: Appending without -overwrite implies -keeplast "
687 "between the first two files. \n"
688 "If the trajectories have an overlap and have not been written binary \n"
689 "reproducible this will produce an incorrect trajectory!\n\n");
691 filetype = gmx_fio_getftp(stfio);
692 /* Fails if last frame is incomplete
693 * We can't do anything about it without overwriting
695 if (filetype == efXTC || filetype == efTNG)
697 lasttime = trx_get_time_of_final_frame(status);
702 while (read_next_frame(oenv, status, &fr))
709 bKeepLastAppend = TRUE;
711 trxout = open_trx(out_file, "a");
715 if (gmx_fio_getftp(stfio) != efXTC)
717 gmx_fatal(FARGS, "Overwrite only supported for XTC." );
719 last_frame_time = trx_get_time_of_final_frame(status);
721 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
722 * or when seek time = 0 */
723 if (inFilesEdited.size() > 1 && settime[1] < last_frame_time+timest[0]*0.5)
725 /* Jump to one time-frame before the start of next
727 searchtime = settime[1]-timest[0]*1.25;
731 searchtime = last_frame_time;
733 if (xtc_seek_time(stfio, searchtime, fr.natoms, TRUE))
735 gmx_fatal(FARGS, "Error seeking to append position.");
737 read_next_frame(oenv, status, &fr);
738 if (std::abs(searchtime - fr.time) > timest[0]*0.5)
740 gmx_fatal(FARGS, "Error seeking: attempted to seek to %f but got %f.",
741 searchtime, fr.time);
745 fpos = gmx_fio_ftell(stfio);
747 trxout = open_trx(out_file, "r+");
748 if (gmx_fio_seek(trx_get_fileio(trxout), fpos))
750 gmx_fatal(FARGS, "Error seeking to append position.");
755 printf("\n Will append after %f \n", lasttime);
759 /* Lets stitch up some files */
760 timestep = timest[0];
761 for (size_t i = n_append + 1; i < inFilesEdited.size(); i++)
765 /* set the next time from the last frame in previous file */
768 /* When writing TNG the step determine which frame to write. Use an
769 * offset to be able to increase steps properly when changing files. */
772 prevEndStep = frout.step;
777 if (cont_type[i] == TIME_CONTINUE)
780 begin += 0.5*timestep;
781 settime[i] = frout.time;
782 cont_type[i] = TIME_EXPLICIT;
784 else if (cont_type[i] == TIME_LAST)
787 begin += 0.5*timestep;
789 /* Or, if the time in the next part should be changed by the
790 * same amount, start at half a timestep from the last time
791 * so we dont repeat frames.
793 /* I don't understand the comment above, but for all the cases
794 * I tried the code seems to work properly. B. Hess 2008-4-2.
797 /* Or, if time is set explicitly, we check for overlap/gap */
798 if (cont_type[i] == TIME_EXPLICIT)
800 if (i < inFilesEdited.size() &&
801 frout.time < settime[i] - 1.5*timestep)
803 fprintf(stderr, "WARNING: Frames around t=%f %s have a different "
804 "spacing than the rest,\n"
805 "might be a gap or overlap that couldn't be corrected "
806 "automatically.\n", output_env_conv_time(oenv, frout.time),
812 /* if we don't have a timestep in the current file, use the old one */
815 timestep = timest[i];
817 read_first_frame(oenv, &status, inFilesEdited[i].c_str(), &fr, FLAGS);
821 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
824 if (cont_type[i] == TIME_EXPLICIT)
826 t_corr = settime[i]-fr.time;
828 /* t_corr is the amount we want to change the time.
829 * If the user has chosen not to change the time for
830 * this part of the trajectory t_corr remains at
831 * the value it had in the last part, changing this
832 * by the same amount.
833 * If no value was given for the first trajectory part
834 * we let the time start at zero, see the edit_files routine.
845 printf("lasttime %g\n", lasttime);
849 /* copy the input frame to the output frame */
851 /* set the new time by adding the correct calculated above */
852 frout.time += t_corr;
855 frout.step += prevEndStep;
857 /* quit if we have reached the end of what should be written */
858 if ((end > 0) && (frout.time > end+GMX_REAL_EPS))
860 i = inFilesEdited.size();
864 /* determine if we should write this frame (dt is handled elsewhere) */
865 if (bCat) /* write all frames of all files */
869 else if (bKeepLast || (bKeepLastAppend && i == 1))
870 /* write till last frame of this traj
871 and skip first frame(s) of next traj */
873 bWrite = ( frout.time > lasttime+0.5*timestep );
875 else /* write till first frame of next traj */
877 bWrite = ( frout.time < settime[i+1]-0.5*timestep );
880 if (bWrite && (frout.time >= begin) )
885 first_time = frout.time;
887 lasttime = frout.time;
889 if (dt == 0 || bRmod(frout.time, first_time, dt))
892 last_ok_t = frout.time;
895 fprintf(stderr, "\nContinue writing frames from %s t=%g %s, "
897 inFilesEdited[i].c_str(),
898 output_env_conv_time(oenv, frout.time), timeUnit.c_str(),
905 write_trxframe_indexed(trxout, &frout, isize, index,
910 write_trxframe(trxout, &frout, nullptr);
912 if ( ((frame % 10) == 0) || (frame < 10) )
914 fprintf(stderr, " -> frame %6d time %8.3f %s \r",
915 frame_out, output_env_conv_time(oenv, frout.time), timeUnit.c_str());
921 while (read_next_frame(oenv, status, &fr));
929 fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
930 frame, output_env_conv_time(oenv, last_ok_t), timeUnit.c_str());