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.
45 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/fileio/gmxfio.h"
48 #include "gromacs/fileio/tpxio.h"
49 #include "gromacs/fileio/trxio.h"
50 #include "gromacs/fileio/trnio.h"
51 #include "gromacs/fileio/tngio.h"
52 #include "gromacs/fileio/tngio_for_tools.h"
53 #include "gromacs/commandline/pargs.h"
54 #include "gromacs/fileio/futil.h"
55 #include "gromacs/fileio/pdbio.h"
56 #include "gromacs/fileio/confio.h"
60 #include "gromacs/fileio/xtcio.h"
66 #include "gmx_fatal.h"
68 #define TIME_EXPLICIT 0
69 #define TIME_CONTINUE 1
74 #define FLAGS (TRX_READ_X | TRX_READ_V | TRX_READ_F)
76 static void scan_trj_files(char **fnms, int nfiles, real *readtime,
77 real *timestep, atom_id imax,
78 const output_env_t oenv)
80 /* Check start time of all files */
87 for (i = 0; i < nfiles; i++)
89 ok = read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
93 gmx_fatal(FARGS, "\nCouldn't read frame from file." );
97 readtime[i] = fr.time;
102 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
113 if (natoms != fr.natoms)
115 gmx_fatal(FARGS, "\nDifferent numbers of atoms (%d/%d) in files",
121 if (fr.natoms <= imax)
123 gmx_fatal(FARGS, "\nNot enough atoms (%d) for index group (%d)",
128 ok = read_next_frame(oenv, status, &fr);
131 timestep[i] = fr.time - readtime[i];
152 fprintf(stderr, "\n");
156 static void sort_files(char **fnms, real *settime, int nfile)
162 for (i = 0; i < nfile; i++)
165 for (j = i + 1; j < nfile; j++)
167 if (settime[j] < settime[minidx])
174 timeswap = settime[i];
175 settime[i] = settime[minidx];
176 settime[minidx] = timeswap;
178 fnms[i] = fnms[minidx];
179 fnms[minidx] = chptr;
184 static void edit_files(char **fnms, int nfiles, real *readtime, real *timestep,
185 real *settime, int *cont_type, gmx_bool bSetTime,
186 gmx_bool bSort, const output_env_t oenv)
190 char inputstring[STRLEN], *chptr;
194 fprintf(stderr, "\n\nEnter the new start time (%s) for each file.\n"
195 "There are two special options, both disable sorting:\n\n"
196 "c (continue) - The start time is taken from the end\n"
197 "of the previous file. Use it when your continuation run\n"
198 "restarts with t=0.\n\n"
199 "l (last) - The time in this file will be changed the\n"
200 "same amount as in the previous. Use it when the time in the\n"
201 "new run continues from the end of the previous one,\n"
202 "since this takes possible overlap into account.\n\n",
203 output_env_get_time_unit(oenv));
207 " File Current start (%s) New start (%s)\n"
208 "---------------------------------------------------------\n",
209 output_env_get_time_unit(oenv), output_env_get_time_unit(oenv));
211 for (i = 0; i < nfiles; i++)
213 fprintf(stderr, "%25s %10.3f %s ", fnms[i],
214 output_env_conv_time(oenv, readtime[i]), output_env_get_time_unit(oenv));
218 if (NULL == fgets(inputstring, STRLEN - 1, stdin))
220 gmx_fatal(FARGS, "Error reading user input" );
223 inputstring[strlen(inputstring)-1] = 0;
225 if (inputstring[0] == 'c' || inputstring[0] == 'C')
227 cont_type[i] = TIME_CONTINUE;
230 settime[i] = FLT_MAX;
232 else if (inputstring[0] == 'l' ||
233 inputstring[0] == 'L')
235 cont_type[i] = TIME_LAST;
238 settime[i] = FLT_MAX;
242 settime[i] = strtod(inputstring, &chptr)*
243 output_env_get_time_invfactor(oenv);
244 if (chptr == inputstring)
246 fprintf(stderr, "'%s' not recognized as a floating point number, 'c' or 'l'. "
247 "Try again: ", inputstring);
251 cont_type[i] = TIME_EXPLICIT;
258 if (cont_type[0] != TIME_EXPLICIT)
260 cont_type[0] = TIME_EXPLICIT;
266 for (i = 0; i < nfiles; i++)
268 settime[i] = readtime[i];
273 fprintf(stderr, "Sorting disabled.\n");
277 sort_files(fnms, settime, nfiles);
279 /* Write out the new order and start times */
280 fprintf(stderr, "\nSummary of files and start times used:\n\n"
281 " File Start time Time step\n"
282 "---------------------------------------------------------\n");
283 for (i = 0; i < nfiles; i++)
285 switch (cont_type[i])
288 fprintf(stderr, "%25s %10.3f %s %10.3f %s",
290 output_env_conv_time(oenv, settime[i]), output_env_get_time_unit(oenv),
291 output_env_conv_time(oenv, timestep[i]), output_env_get_time_unit(oenv));
293 cont_type[i-1] == TIME_EXPLICIT && settime[i] == settime[i-1])
295 fprintf(stderr, " WARNING: same Start time as previous");
297 fprintf(stderr, "\n");
300 fprintf(stderr, "%25s Continue from last file\n", fnms[i]);
303 fprintf(stderr, "%25s Change by same amount as last file\n",
308 fprintf(stderr, "\n");
310 settime[nfiles] = FLT_MAX;
311 cont_type[nfiles] = TIME_EXPLICIT;
312 readtime[nfiles] = FLT_MAX;
315 static void do_demux(int nset, char *fnms[], char *fnms_out[], int nval,
316 real **value, real *time, real dt_remd, int isize,
317 atom_id index[], real dt, const output_env_t oenv)
319 int i, j, k, natoms, nnn;
320 t_trxstatus **fp_in, **fp_out;
321 gmx_bool bCont, *bSet;
322 real t, first_time = 0;
330 for (i = 0; (i < nset); i++)
332 nnn = read_first_frame(oenv, &(fp_in[i]), fnms[i], &(trx[i]),
337 first_time = trx[i].time;
339 else if (natoms != nnn)
341 gmx_fatal(FARGS, "Trajectory file %s has %d atoms while previous trajs had %d atoms", fnms[i], nnn, natoms);
347 else if (t != trx[i].time)
349 gmx_fatal(FARGS, "Trajectory file %s has time %f while previous trajs had time %f", fnms[i], trx[i].time, t);
354 for (i = 0; (i < nset); i++)
356 fp_out[i] = open_trx(fnms_out[i], "w");
359 if (gmx_nint(time[k] - t) != 0)
361 gmx_fatal(FARGS, "First time in demuxing table does not match trajectories");
365 while ((k+1 < nval) && ((trx[0].time - time[k+1]) > dt_remd*0.1))
371 fprintf(debug, "trx[0].time = %g, time[k] = %g\n", trx[0].time, time[k]);
373 for (i = 0; (i < nset); i++)
377 for (i = 0; (i < nset); i++)
379 j = gmx_nint(value[i][k]);
380 range_check(j, 0, nset);
383 gmx_fatal(FARGS, "Demuxing the same replica %d twice at time %f",
388 if (dt == 0 || bRmod(trx[i].time, first_time, dt))
392 write_trxframe_indexed(fp_out[j], &trx[i], isize, index, NULL);
396 write_trxframe(fp_out[j], &trx[i], NULL);
402 for (i = 0; (i < nset); i++)
404 bCont = bCont && read_next_frame(oenv, fp_in[i], &trx[i]);
409 for (i = 0; (i < nset); i++)
412 close_trx(fp_out[i]);
416 int gmx_trjcat(int argc, char *argv[])
420 "[THISMODULE] concatenates several input trajectory files in sorted order. ",
421 "In case of double time frames the one in the later file is used. ",
422 "By specifying [TT]-settime[tt] you will be asked for the start time ",
423 "of each file. The input files are taken from the command line, ",
424 "such that a command like [TT]gmx trjcat -f *.trr -o fixed.trr[tt] should do ",
425 "the trick. Using [TT]-cat[tt], you can simply paste several files ",
426 "together without removal of frames with identical time stamps.[PAR]",
427 "One important option is inferred when the output file is amongst the",
428 "input files. In that case that particular file will be appended to",
429 "which implies you do not need to store double the amount of data.",
430 "Obviously the file to append to has to be the one with lowest starting",
431 "time since one can only append at the end of a file.[PAR]",
432 "If the [TT]-demux[tt] option is given, the N trajectories that are",
433 "read, are written in another order as specified in the [TT].xvg[tt] file.",
434 "The [TT].xvg[tt] file should contain something like:[PAR]",
435 "[TT]0 0 1 2 3 4 5[BR]",
436 "2 1 0 2 3 5 4[tt][BR]",
437 "Where the first number is the time, and subsequent numbers point to",
438 "trajectory indices.",
439 "The frames corresponding to the numbers present at the first line",
440 "are collected into the output trajectory. If the number of frames in",
441 "the trajectory does not match that in the [TT].xvg[tt] file then the program",
442 "tries to be smart. Beware."
444 static gmx_bool bVels = TRUE;
445 static gmx_bool bCat = FALSE;
446 static gmx_bool bSort = TRUE;
447 static gmx_bool bKeepLast = FALSE;
448 static gmx_bool bKeepLastAppend = FALSE;
449 static gmx_bool bOverwrite = FALSE;
450 static gmx_bool bSetTime = FALSE;
451 static gmx_bool bDeMux;
452 static real begin = -1;
453 static real end = -1;
459 { "-b", FALSE, etTIME,
460 { &begin }, "First time to use (%t)" },
461 { "-e", FALSE, etTIME,
462 { &end }, "Last time to use (%t)" },
463 { "-dt", FALSE, etTIME,
464 { &dt }, "Only write frame when t MOD dt = first time (%t)" },
465 { "-vel", FALSE, etBOOL,
466 { &bVels }, "Read and write velocities if possible" },
467 { "-settime", FALSE, etBOOL,
468 { &bSetTime }, "Change starting time interactively" },
469 { "-sort", FALSE, etBOOL,
470 { &bSort }, "Sort trajectory files (not frames)" },
471 { "-keeplast", FALSE, etBOOL,
472 { &bKeepLast }, "Keep overlapping frames at end of trajectory" },
473 { "-overwrite", FALSE, etBOOL,
474 { &bOverwrite }, "Overwrite overlapping frames during appending" },
475 { "-cat", FALSE, etBOOL,
476 { &bCat }, "Do not discard double time frames" }
478 #define npargs asize(pa)
479 int ftpin, i, frame, frame_out, step = 0, trjout = 0;
480 t_trxstatus *status, *trxout = NULL;
483 t_trxframe fr, frout;
484 char **fnms, **fnms_out, *in_file, *out_file;
486 gmx_bool bNewFile, bIndex, bWrite;
487 int earliersteps, nfile_in, nfile_out, *cont_type, last_ok_step;
488 real *readtime, *timest, *settime;
489 real first_time = 0, lasttime = NOTSET, last_ok_t = -1, timestep;
490 real last_frame_time, searchtime;
492 atom_id *index = NULL, imax;
494 real **val = NULL, *t = NULL, dt_remd;
495 int n, nset, ftpout = -1, prevEndStep = 0, filetype;
501 { efTRX, "-f", NULL, ffRDMULT },
502 { efTRO, "-o", NULL, ffWRMULT },
503 { efNDX, "-n", "index", ffOPTRD },
504 { efXVG, "-demux", "remd", ffOPTRD }
507 #define NFILE asize(fnm)
509 if (!parse_common_args(&argc, argv, PCA_BE_NICE | PCA_TIME_UNIT, NFILE, fnm,
510 asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
515 bIndex = ftp2bSet(efNDX, NFILE, fnm);
516 bDeMux = ftp2bSet(efXVG, NFILE, fnm);
517 bSort = bSort && !bDeMux;
522 printf("Select group for output\n");
523 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
526 for (i = 1; i < isize; i++)
528 imax = max(imax, index[i]);
535 val = read_xvg_time(opt2fn("-demux", NFILE, fnm), TRUE,
536 opt2parg_bSet("-b", npargs, pa), begin,
537 opt2parg_bSet("-e", npargs, pa), end, 1, &nset, &n,
539 printf("Read %d sets of %d points, dt = %g\n\n", nset, n, dt_remd);
542 fprintf(debug, "Dump of replica_index.xvg\n");
543 for (i = 0; (i < n); i++)
545 fprintf(debug, "%10g", t[i]);
546 for (j = 0; (j < nset); j++)
548 fprintf(debug, " %3d", gmx_nint(val[j][i]));
550 fprintf(debug, "\n");
555 nfile_in = opt2fns(&fnms, "-f", NFILE, fnm);
558 gmx_fatal(FARGS, "No input files!" );
561 if (bDeMux && (nfile_in != nset))
563 gmx_fatal(FARGS, "You have specified %d files and %d entries in the demux table", nfile_in, nset);
566 ftpin = fn2ftp(fnms[0]);
568 for (i = 1; i < nfile_in; i++)
570 if (ftpin != fn2ftp(fnms[i]))
572 gmx_fatal(FARGS, "All input files must be of the same format");
576 nfile_out = opt2fns(&fnms_out, "-o", NFILE, fnm);
579 gmx_fatal(FARGS, "No output files!");
581 if ((nfile_out > 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 && (nfile_out != nset) && (nfile_out != 1))
587 gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %d", nset, nfile_out);
591 if (nfile_out != nset)
593 char *buf = strdup(fnms_out[0]);
594 snew(fnms_out, nset);
595 for (i = 0; (i < nset); i++)
597 snew(fnms_out[i], strlen(buf)+32);
598 sprintf(fnms_out[i], "%d_%s", i, buf);
602 do_demux(nfile_in, fnms, fnms_out, n, val, t, dt_remd, isize, index, dt, oenv);
606 snew(readtime, nfile_in+1);
607 snew(timest, nfile_in+1);
608 scan_trj_files(fnms, nfile_in, readtime, timest, imax, oenv);
610 snew(settime, nfile_in+1);
611 snew(cont_type, nfile_in+1);
612 edit_files(fnms, nfile_in, 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 out_file = fnms_out[0];
619 ftpout = fn2ftp(out_file);
621 for (i = 0; ((i < nfile_in) && (n_append == -1)); i++)
623 if (strcmp(fnms[i], 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)",
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 trjtools_gmx_prepare_tng_writing(out_file, 'w', NULL, &trxout,
657 fnms[0], isize, NULL, index, grpname);
661 trjtools_gmx_prepare_tng_writing(out_file, 'w', NULL, &trxout,
662 fnms[0], -1, NULL, NULL, NULL);
667 trxout = open_trx(out_file, "w");
669 memset(&frout, 0, sizeof(frout));
675 if (!read_first_frame(oenv, &status, out_file, &fr, FLAGS))
677 gmx_fatal(FARGS, "Reading first frame from %s", out_file);
680 stfio = trx_get_fileio(status);
681 if (!bKeepLast && !bOverwrite)
683 fprintf(stderr, "\n\nWARNING: Appending without -overwrite implies -keeplast "
684 "between the first two files. \n"
685 "If the trajectories have an overlap and have not been written binary \n"
686 "reproducible this will produce an incorrect trajectory!\n\n");
688 filetype = gmx_fio_getftp(stfio);
689 /* Fails if last frame is incomplete
690 * We can't do anything about it without overwriting
692 if (filetype == efXTC || filetype == efTNG)
694 lasttime = trx_get_time_of_final_frame(status);
699 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 (fabs(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);
740 fpos = gmx_fio_ftell(stfio);
742 trxout = open_trx(out_file, "r+");
743 if (gmx_fio_seek(trx_get_fileio(trxout), fpos))
745 gmx_fatal(FARGS, "Error seeking to append position.");
748 printf("\n Will append after %f \n", lasttime);
751 /* Lets stitch up some files */
752 timestep = timest[0];
753 for (i = n_append+1; (i < nfile_in); i++)
757 /* set the next time from the last frame in previous file */
760 /* When writing TNG the step determine which frame to write. Use an
761 * offset to be able to increase steps properly when changing files. */
764 prevEndStep = frout.step;
769 if (cont_type[i] == TIME_CONTINUE)
772 begin += 0.5*timestep;
773 settime[i] = frout.time;
774 cont_type[i] = TIME_EXPLICIT;
776 else if (cont_type[i] == TIME_LAST)
779 begin += 0.5*timestep;
781 /* Or, if the time in the next part should be changed by the
782 * same amount, start at half a timestep from the last time
783 * so we dont repeat frames.
785 /* I don't understand the comment above, but for all the cases
786 * I tried the code seems to work properly. B. Hess 2008-4-2.
789 /* Or, if time is set explicitly, we check for overlap/gap */
790 if (cont_type[i] == TIME_EXPLICIT)
792 if ( ( i < nfile_in ) &&
793 ( frout.time < settime[i]-1.5*timestep ) )
795 fprintf(stderr, "WARNING: Frames around t=%f %s have a different "
796 "spacing than the rest,\n"
797 "might be a gap or overlap that couldn't be corrected "
798 "automatically.\n", output_env_conv_time(oenv, frout.time),
799 output_env_get_time_unit(oenv));
804 /* if we don't have a timestep in the current file, use the old one */
807 timestep = timest[i];
809 read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
813 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
816 if (cont_type[i] == TIME_EXPLICIT)
818 t_corr = settime[i]-fr.time;
820 /* t_corr is the amount we want to change the time.
821 * If the user has chosen not to change the time for
822 * this part of the trajectory t_corr remains at
823 * the value it had in the last part, changing this
824 * by the same amount.
825 * If no value was given for the first trajectory part
826 * we let the time start at zero, see the edit_files routine.
832 if (lasttime != NOTSET)
834 printf("lasttime %g\n", lasttime);
839 /* copy the input frame to the output frame */
841 /* set the new time by adding the correct calculated above */
842 frout.time += t_corr;
845 frout.step += prevEndStep;
847 /* quit if we have reached the end of what should be written */
848 if ((end > 0) && (frout.time > end+GMX_REAL_EPS))
854 /* determine if we should write this frame (dt is handled elsewhere) */
855 if (bCat) /* write all frames of all files */
859 else if (bKeepLast || (bKeepLastAppend && i == 1))
860 /* write till last frame of this traj
861 and skip first frame(s) of next traj */
863 bWrite = ( frout.time > lasttime+0.5*timestep );
865 else /* write till first frame of next traj */
867 bWrite = ( frout.time < settime[i+1]-0.5*timestep );
870 if (bWrite && (frout.time >= begin) )
875 first_time = frout.time;
877 lasttime = frout.time;
878 if (dt == 0 || bRmod(frout.time, first_time, dt))
881 last_ok_t = frout.time;
884 fprintf(stderr, "\nContinue writing frames from %s t=%g %s, "
886 fnms[i], output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv),
893 write_trxframe_indexed(trxout, &frout, isize, index,
898 write_trxframe(trxout, &frout, NULL);
900 if ( ((frame % 10) == 0) || (frame < 10) )
902 fprintf(stderr, " -> frame %6d time %8.3f %s \r",
903 frame_out, output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv));
908 while (read_next_frame(oenv, status, &fr));
912 earliersteps += step;
918 fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
919 frame, output_env_conv_time(oenv, last_ok_t), output_env_get_time_unit(oenv));