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.
48 #include "gromacs/commandline/pargs.h"
49 #include "gromacs/fileio/confio.h"
50 #include "gromacs/fileio/gmxfio.h"
51 #include "gromacs/fileio/pdbio.h"
52 #include "gromacs/fileio/tngio.h"
53 #include "gromacs/fileio/trxio.h"
54 #include "gromacs/fileio/xtcio.h"
55 #include "gromacs/fileio/xvgr.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/topology/index.h"
59 #include "gromacs/trajectory/trajectoryframe.h"
60 #include "gromacs/utility/arrayref.h"
61 #include "gromacs/utility/arraysize.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/futil.h"
64 #include "gromacs/utility/smalloc.h"
65 #include "gromacs/utility/stringutil.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(gmx::ArrayRef<const std::string> files,
77 real *timestep, int imax,
78 const gmx_output_env_t *oenv)
80 /* Check start time of all files */
86 for (gmx::index i = 0; i < files.ssize(); i++)
88 ok = read_first_frame(oenv, &status, files[i].c_str(), &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(gmx::ArrayRef<std::string> files, real *settime)
157 for (gmx::index i = 0; i < files.ssize(); i++)
159 gmx::index minidx = i;
160 for (gmx::index j = i + 1; j < files.ssize(); j++)
162 if (settime[j] < settime[minidx])
169 real timeswap = settime[i];
170 settime[i] = settime[minidx];
171 settime[minidx] = timeswap;
172 std::swap(files[i], files[minidx]);
177 static void edit_files(gmx::ArrayRef<std::string> files,
178 real *readtime, real *timestep,
179 real *settime, int *cont_type, gmx_bool bSetTime,
180 gmx_bool bSort, const gmx_output_env_t *oenv)
183 char inputstring[STRLEN], *chptr;
185 auto timeUnit = output_env_get_time_unit(oenv);
188 fprintf(stderr, "\n\nEnter the new start time (%s) for each file.\n"
189 "There are two special options, both disable sorting:\n\n"
190 "c (continue) - The start time is taken from the end\n"
191 "of the previous file. Use it when your continuation run\n"
192 "restarts with t=0.\n\n"
193 "l (last) - The time in this file will be changed the\n"
194 "same amount as in the previous. Use it when the time in the\n"
195 "new run continues from the end of the previous one,\n"
196 "since this takes possible overlap into account.\n\n",
201 " File Current start (%s) New start (%s)\n"
202 "---------------------------------------------------------\n",
203 timeUnit.c_str(), timeUnit.c_str());
205 for (gmx::index i = 0; i < files.ssize(); i++)
207 fprintf(stderr, "%25s %10.3f %s ", files[i].c_str(),
208 output_env_conv_time(oenv, readtime[i]), timeUnit.c_str());
212 if (nullptr == fgets(inputstring, STRLEN - 1, stdin))
214 gmx_fatal(FARGS, "Error reading user input" );
217 inputstring[std::strlen(inputstring)-1] = 0;
219 if (inputstring[0] == 'c' || inputstring[0] == 'C')
221 cont_type[i] = TIME_CONTINUE;
224 settime[i] = FLT_MAX;
226 else if (inputstring[0] == 'l' ||
227 inputstring[0] == 'L')
229 cont_type[i] = TIME_LAST;
232 settime[i] = FLT_MAX;
236 settime[i] = strtod(inputstring, &chptr)*
237 output_env_get_time_invfactor(oenv);
238 if (chptr == inputstring)
240 fprintf(stderr, "'%s' not recognized as a floating point number, 'c' or 'l'. "
241 "Try again: ", inputstring);
245 cont_type[i] = TIME_EXPLICIT;
252 if (cont_type[0] != TIME_EXPLICIT)
254 cont_type[0] = TIME_EXPLICIT;
260 for (gmx::index i = 0; i < files.ssize(); i++)
262 settime[i] = readtime[i];
267 fprintf(stderr, "Sorting disabled.\n");
271 sort_files(files, settime);
273 /* Write out the new order and start times */
274 fprintf(stderr, "\nSummary of files and start times used:\n\n"
275 " File Start time Time step\n"
276 "---------------------------------------------------------\n");
277 for (gmx::index i = 0; i < files.ssize(); i++)
279 switch (cont_type[i])
282 fprintf(stderr, "%25s %10.3f %s %10.3f %s",
284 output_env_conv_time(oenv, settime[i]), timeUnit.c_str(),
285 output_env_conv_time(oenv, timestep[i]), timeUnit.c_str());
287 cont_type[i-1] == TIME_EXPLICIT && settime[i] == settime[i-1])
289 fprintf(stderr, " WARNING: same Start time as previous");
291 fprintf(stderr, "\n");
294 fprintf(stderr, "%25s Continue from last file\n", files[i].c_str());
297 fprintf(stderr, "%25s Change by same amount as last file\n",
302 fprintf(stderr, "\n");
304 settime[files.size()] = FLT_MAX;
305 cont_type[files.size()] = TIME_EXPLICIT;
306 readtime[files.size()] = FLT_MAX;
309 static void do_demux(gmx::ArrayRef<const std::string> inFiles,
310 gmx::ArrayRef<const std::string> outFiles, int nval,
311 real **value, real *time, real dt_remd, int isize,
312 int index[], real dt, const gmx_output_env_t *oenv)
315 t_trxstatus **fp_in, **fp_out;
316 gmx_bool bCont, *bSet;
317 real t, first_time = 0;
320 snew(fp_in, inFiles.size());
321 snew(trx, inFiles.size());
322 snew(bSet, inFiles.size());
325 for (gmx::index i = 0; i < inFiles.ssize(); i++)
327 read_first_frame(oenv, &(fp_in[i]), inFiles[i].c_str(), &(trx[i]),
331 natoms = trx[i].natoms;
332 first_time = trx[i].time;
334 else if (natoms != trx[i].natoms)
336 gmx_fatal(FARGS, "Trajectory file %s has %d atoms while previous trajs had %d atoms", inFiles[i].c_str(), trx[i].natoms, natoms);
342 else if (t != trx[i].time)
344 gmx_fatal(FARGS, "Trajectory file %s has time %f while previous trajs had time %f", inFiles[i].c_str(), trx[i].time, t);
348 snew(fp_out, inFiles.size());
349 for (gmx::index i = 0; i < inFiles.ssize(); i++)
351 fp_out[i] = open_trx(outFiles[i].c_str(), "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 (gmx::index i = 0; i < inFiles.ssize(); i++)
372 for (gmx::index i = 0; i < inFiles.ssize(); i++)
374 int j = gmx::roundToInt(value[i][k]);
375 range_check(j, 0, inFiles.size());
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, nullptr);
391 write_trxframe(fp_out[j], &trx[i], nullptr);
397 for (gmx::index i = 0; i < inFiles.ssize(); i++)
399 bCont = bCont && read_next_frame(oenv, fp_in[i], &trx[i]);
404 for (gmx::index i = 0; i < inFiles.ssize(); 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 bCat = FALSE;
442 static gmx_bool bSort = TRUE;
443 static gmx_bool bKeepLast = FALSE;
444 static gmx_bool bKeepLastAppend = FALSE;
445 static gmx_bool bOverwrite = FALSE;
446 static gmx_bool bSetTime = FALSE;
447 static gmx_bool bDeMux;
448 static real begin = -1;
449 static real end = -1;
455 { "-b", FALSE, etTIME,
456 { &begin }, "First time to use (%t)" },
457 { "-e", FALSE, etTIME,
458 { &end }, "Last time to use (%t)" },
459 { "-dt", FALSE, etTIME,
460 { &dt }, "Only write frame when t MOD dt = first time (%t)" },
461 { "-settime", FALSE, etBOOL,
462 { &bSetTime }, "Change starting time interactively" },
463 { "-sort", FALSE, etBOOL,
464 { &bSort }, "Sort trajectory files (not frames)" },
465 { "-keeplast", FALSE, etBOOL,
466 { &bKeepLast }, "Keep overlapping frames at end of trajectory" },
467 { "-overwrite", FALSE, etBOOL,
468 { &bOverwrite }, "Overwrite overlapping frames during appending" },
469 { "-cat", FALSE, etBOOL,
470 { &bCat }, "Do not discard double time frames" }
472 #define npargs asize(pa)
473 int ftpin, i, frame, frame_out;
474 t_trxstatus *status, *trxout = nullptr;
476 t_trxframe fr, frout;
478 gmx_bool bNewFile, bIndex, bWrite;
480 real *readtime, *timest, *settime;
481 real first_time = 0, lasttime = 0, last_ok_t = -1, timestep;
482 gmx_bool lastTimeSet = FALSE;
483 real last_frame_time, searchtime;
485 int *index = nullptr, imax;
487 real **val = nullptr, *t = nullptr, dt_remd;
488 int n, nset, ftpout = -1, prevEndStep = 0, filetype;
490 gmx_output_env_t *oenv;
493 { efTRX, "-f", nullptr, ffRDMULT },
494 { efTRO, "-o", nullptr, ffWRMULT },
495 { efNDX, "-n", "index", ffOPTRD },
496 { efXVG, "-demux", "remd", ffOPTRD }
499 #define NFILE asize(fnm)
501 if (!parse_common_args(&argc, argv, PCA_TIME_UNIT, NFILE, fnm,
502 asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
506 fprintf(stdout, "Note that major changes are planned in future for "
507 "trjcat, to improve usability and utility.");
509 auto timeUnit = output_env_get_time_unit(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 gmx::ArrayRef<const std::string> inFiles = opt2fns("-f", NFILE, fnm);
554 gmx_fatal(FARGS, "No input files!" );
557 if (bDeMux && ssize(inFiles) != nset)
559 gmx_fatal(FARGS, "You have specified %td files and %d entries in the demux table", inFiles.ssize(), nset);
562 ftpin = fn2ftp(inFiles[0].c_str());
564 if (ftpin != efTRR && ftpin != efXTC && ftpin != efTNG)
566 gmx_fatal(FARGS, "gmx trjcat can only handle binary trajectory formats (trr, xtc, tng)");
569 for (const std::string &inFile : inFiles)
571 if (ftpin != fn2ftp(inFile.c_str()))
573 gmx_fatal(FARGS, "All input files must be of the same (trr, xtc or tng) format");
577 gmx::ArrayRef<const std::string> outFiles = opt2fns("-o", NFILE, fnm);
578 if (outFiles.empty())
580 gmx_fatal(FARGS, "No output files!");
582 if ((outFiles.size() > 1) && !bDeMux)
584 gmx_fatal(FARGS, "Don't know what to do with more than 1 output file if not demultiplexing");
586 else if (bDeMux && ssize(outFiles) != nset && outFiles.size() != 1)
588 gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %td", nset, outFiles.ssize());
592 auto outFilesDemux = gmx::copyOf(outFiles);
593 if (gmx::ssize(outFilesDemux) != nset)
595 std::string name = outFilesDemux[0];
596 outFilesDemux.resize(nset);
597 for (i = 0; (i < nset); i++)
599 outFilesDemux[0] = gmx::formatString("%d_%s", i, name.c_str());
602 do_demux(inFiles, outFilesDemux, n, val, t, dt_remd, isize, index, dt, oenv);
606 snew(readtime, inFiles.size() + 1);
607 snew(timest, inFiles.size() + 1);
608 scan_trj_files(inFiles, readtime, timest, imax, oenv);
610 snew(settime, inFiles.size() + 1);
611 snew(cont_type, inFiles.size() + 1);
612 auto inFilesEdited = gmx::copyOf(inFiles);
613 edit_files(inFilesEdited, readtime, timest, settime, cont_type, bSetTime, bSort,
616 /* Check whether the output file is amongst the input files
617 * This has to be done after sorting etc.
619 const char *out_file = outFiles[0].c_str();
620 ftpout = fn2ftp(out_file);
622 for (size_t i = 0; i < inFilesEdited.size() && n_append == -1; i++)
624 if (std::strcmp(inFilesEdited[i].c_str(), out_file) == 0)
631 fprintf(stderr, "Will append to %s rather than creating a new file\n",
634 else if (n_append != -1)
636 gmx_fatal(FARGS, "Can only append to the first file which is %s (not %s)",
637 inFilesEdited[0].c_str(), out_file);
640 /* Not checking input format, could be dangerous :-) */
641 /* Not checking output format, equally dangerous :-) */
645 /* the default is not to change the time at all,
646 * but this is overridden by the edit_files routine
656 gmx_fatal(FARGS, "When writing TNG the input file format must also be TNG");
660 trxout = trjtools_gmx_prepare_tng_writing(out_file, 'w', nullptr,
661 inFilesEdited[0].c_str(), isize, nullptr, gmx::arrayRefFromArray(index, isize), grpname);
665 trxout = trjtools_gmx_prepare_tng_writing(out_file, 'w', nullptr,
666 inFilesEdited[0].c_str(), -1, nullptr, {}, nullptr);
671 trxout = open_trx(out_file, "w");
673 std::memset(&frout, 0, sizeof(frout));
679 if (!read_first_frame(oenv, &status, out_file, &fr, FLAGS))
681 gmx_fatal(FARGS, "Reading first frame from %s", out_file);
684 stfio = trx_get_fileio(status);
685 if (!bKeepLast && !bOverwrite)
687 fprintf(stderr, "\n\nWARNING: Appending without -overwrite implies -keeplast "
688 "between the first two files. \n"
689 "If the trajectories have an overlap and have not been written binary \n"
690 "reproducible this will produce an incorrect trajectory!\n\n");
692 filetype = gmx_fio_getftp(stfio);
693 /* Fails if last frame is incomplete
694 * We can't do anything about it without overwriting
696 if (filetype == efXTC || filetype == efTNG)
698 lasttime = trx_get_time_of_final_frame(status);
703 while (read_next_frame(oenv, status, &fr))
710 bKeepLastAppend = TRUE;
712 trxout = open_trx(out_file, "a");
716 if (gmx_fio_getftp(stfio) != efXTC)
718 gmx_fatal(FARGS, "Overwrite only supported for XTC." );
720 last_frame_time = trx_get_time_of_final_frame(status);
722 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
723 * or when seek time = 0 */
724 if (inFilesEdited.size() > 1 && settime[1] < last_frame_time+timest[0]*0.5)
726 /* Jump to one time-frame before the start of next
728 searchtime = settime[1]-timest[0]*1.25;
732 searchtime = last_frame_time;
734 if (xtc_seek_time(stfio, searchtime, fr.natoms, TRUE))
736 gmx_fatal(FARGS, "Error seeking to append position.");
738 read_next_frame(oenv, status, &fr);
739 if (std::abs(searchtime - fr.time) > timest[0]*0.5)
741 gmx_fatal(FARGS, "Error seeking: attempted to seek to %f but got %f.",
742 searchtime, fr.time);
746 fpos = gmx_fio_ftell(stfio);
748 trxout = open_trx(out_file, "r+");
749 if (gmx_fio_seek(trx_get_fileio(trxout), fpos))
751 gmx_fatal(FARGS, "Error seeking to append position.");
756 printf("\n Will append after %f \n", lasttime);
760 /* Lets stitch up some files */
761 timestep = timest[0];
762 for (size_t i = n_append + 1; i < inFilesEdited.size(); i++)
766 /* set the next time from the last frame in previous file */
769 /* When writing TNG the step determine which frame to write. Use an
770 * offset to be able to increase steps properly when changing files. */
773 prevEndStep = frout.step;
778 if (cont_type[i] == TIME_CONTINUE)
781 begin += 0.5*timestep;
782 settime[i] = frout.time;
783 cont_type[i] = TIME_EXPLICIT;
785 else if (cont_type[i] == TIME_LAST)
788 begin += 0.5*timestep;
790 /* Or, if the time in the next part should be changed by the
791 * same amount, start at half a timestep from the last time
792 * so we dont repeat frames.
794 /* I don't understand the comment above, but for all the cases
795 * I tried the code seems to work properly. B. Hess 2008-4-2.
798 /* Or, if time is set explicitly, we check for overlap/gap */
799 if (cont_type[i] == TIME_EXPLICIT)
801 if (i < inFilesEdited.size() &&
802 frout.time < settime[i] - 1.5*timestep)
804 fprintf(stderr, "WARNING: Frames around t=%f %s have a different "
805 "spacing than the rest,\n"
806 "might be a gap or overlap that couldn't be corrected "
807 "automatically.\n", output_env_conv_time(oenv, frout.time),
813 /* if we don't have a timestep in the current file, use the old one */
816 timestep = timest[i];
818 read_first_frame(oenv, &status, inFilesEdited[i].c_str(), &fr, FLAGS);
822 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
825 if (cont_type[i] == TIME_EXPLICIT)
827 t_corr = settime[i]-fr.time;
829 /* t_corr is the amount we want to change the time.
830 * If the user has chosen not to change the time for
831 * this part of the trajectory t_corr remains at
832 * the value it had in the last part, changing this
833 * by the same amount.
834 * If no value was given for the first trajectory part
835 * we let the time start at zero, see the edit_files routine.
846 printf("lasttime %g\n", lasttime);
850 /* copy the input frame to the output frame */
852 /* set the new time by adding the correct calculated above */
853 frout.time += t_corr;
856 frout.step += prevEndStep;
858 /* quit if we have reached the end of what should be written */
859 if ((end > 0) && (frout.time > end+GMX_REAL_EPS))
861 i = inFilesEdited.size();
865 /* determine if we should write this frame (dt is handled elsewhere) */
866 if (bCat) /* write all frames of all files */
870 else if (bKeepLast || (bKeepLastAppend && i == 1))
871 /* write till last frame of this traj
872 and skip first frame(s) of next traj */
874 bWrite = ( frout.time > lasttime+0.5*timestep );
876 else /* write till first frame of next traj */
878 bWrite = ( frout.time < settime[i+1]-0.5*timestep );
881 if (bWrite && (frout.time >= begin) )
886 first_time = frout.time;
888 lasttime = frout.time;
890 if (dt == 0 || bRmod(frout.time, first_time, dt))
893 last_ok_t = frout.time;
896 fprintf(stderr, "\nContinue writing frames from %s t=%g %s, "
898 inFilesEdited[i].c_str(),
899 output_env_conv_time(oenv, frout.time), timeUnit.c_str(),
906 write_trxframe_indexed(trxout, &frout, isize, index,
911 write_trxframe(trxout, &frout, nullptr);
913 if ( ((frame % 10) == 0) || (frame < 10) )
915 fprintf(stderr, " -> frame %6d time %8.3f %s \r",
916 frame_out, output_env_conv_time(oenv, frout.time), timeUnit.c_str());
922 while (read_next_frame(oenv, status, &fr));
930 fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
931 frame, output_env_conv_time(oenv, last_ok_t), timeUnit.c_str());