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"
61 #include "gromacs/fileio/xtcio.h"
64 #include "gromacs/fileio/xvgr.h"
67 #include "gromacs/utility/fatalerror.h"
69 #define TIME_EXPLICIT 0
70 #define TIME_CONTINUE 1
75 #define FLAGS (TRX_READ_X | TRX_READ_V | TRX_READ_F)
77 static void scan_trj_files(char **fnms, int nfiles, real *readtime,
78 real *timestep, atom_id imax,
79 const output_env_t oenv)
81 /* Check start time of all files */
88 for (i = 0; i < nfiles; i++)
90 ok = read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
94 gmx_fatal(FARGS, "\nCouldn't read frame from file." );
98 readtime[i] = fr.time;
103 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
114 if (natoms != fr.natoms)
116 gmx_fatal(FARGS, "\nDifferent numbers of atoms (%d/%d) in files",
122 if (fr.natoms <= imax)
124 gmx_fatal(FARGS, "\nNot enough atoms (%d) for index group (%d)",
129 ok = read_next_frame(oenv, status, &fr);
132 timestep[i] = fr.time - readtime[i];
153 fprintf(stderr, "\n");
157 static void sort_files(char **fnms, real *settime, int nfile)
163 for (i = 0; i < nfile; i++)
166 for (j = i + 1; j < nfile; j++)
168 if (settime[j] < settime[minidx])
175 timeswap = settime[i];
176 settime[i] = settime[minidx];
177 settime[minidx] = timeswap;
179 fnms[i] = fnms[minidx];
180 fnms[minidx] = chptr;
185 static void edit_files(char **fnms, int nfiles, real *readtime, real *timestep,
186 real *settime, int *cont_type, gmx_bool bSetTime,
187 gmx_bool bSort, const output_env_t oenv)
191 char inputstring[STRLEN], *chptr;
195 fprintf(stderr, "\n\nEnter the new start time (%s) for each file.\n"
196 "There are two special options, both disable sorting:\n\n"
197 "c (continue) - The start time is taken from the end\n"
198 "of the previous file. Use it when your continuation run\n"
199 "restarts with t=0.\n\n"
200 "l (last) - The time in this file will be changed the\n"
201 "same amount as in the previous. Use it when the time in the\n"
202 "new run continues from the end of the previous one,\n"
203 "since this takes possible overlap into account.\n\n",
204 output_env_get_time_unit(oenv));
208 " File Current start (%s) New start (%s)\n"
209 "---------------------------------------------------------\n",
210 output_env_get_time_unit(oenv), output_env_get_time_unit(oenv));
212 for (i = 0; i < nfiles; i++)
214 fprintf(stderr, "%25s %10.3f %s ", fnms[i],
215 output_env_conv_time(oenv, readtime[i]), output_env_get_time_unit(oenv));
219 if (NULL == fgets(inputstring, STRLEN - 1, stdin))
221 gmx_fatal(FARGS, "Error reading user input" );
224 inputstring[strlen(inputstring)-1] = 0;
226 if (inputstring[0] == 'c' || inputstring[0] == 'C')
228 cont_type[i] = TIME_CONTINUE;
231 settime[i] = FLT_MAX;
233 else if (inputstring[0] == 'l' ||
234 inputstring[0] == 'L')
236 cont_type[i] = TIME_LAST;
239 settime[i] = FLT_MAX;
243 settime[i] = strtod(inputstring, &chptr)*
244 output_env_get_time_invfactor(oenv);
245 if (chptr == inputstring)
247 fprintf(stderr, "'%s' not recognized as a floating point number, 'c' or 'l'. "
248 "Try again: ", inputstring);
252 cont_type[i] = TIME_EXPLICIT;
259 if (cont_type[0] != TIME_EXPLICIT)
261 cont_type[0] = TIME_EXPLICIT;
267 for (i = 0; i < nfiles; i++)
269 settime[i] = readtime[i];
274 fprintf(stderr, "Sorting disabled.\n");
278 sort_files(fnms, settime, nfiles);
280 /* Write out the new order and start times */
281 fprintf(stderr, "\nSummary of files and start times used:\n\n"
282 " File Start time Time step\n"
283 "---------------------------------------------------------\n");
284 for (i = 0; i < nfiles; i++)
286 switch (cont_type[i])
289 fprintf(stderr, "%25s %10.3f %s %10.3f %s",
291 output_env_conv_time(oenv, settime[i]), output_env_get_time_unit(oenv),
292 output_env_conv_time(oenv, timestep[i]), output_env_get_time_unit(oenv));
294 cont_type[i-1] == TIME_EXPLICIT && settime[i] == settime[i-1])
296 fprintf(stderr, " WARNING: same Start time as previous");
298 fprintf(stderr, "\n");
301 fprintf(stderr, "%25s Continue from last file\n", fnms[i]);
304 fprintf(stderr, "%25s Change by same amount as last file\n",
309 fprintf(stderr, "\n");
311 settime[nfiles] = FLT_MAX;
312 cont_type[nfiles] = TIME_EXPLICIT;
313 readtime[nfiles] = FLT_MAX;
316 static void do_demux(int nset, char *fnms[], char *fnms_out[], int nval,
317 real **value, real *time, real dt_remd, int isize,
318 atom_id index[], real dt, const output_env_t oenv)
320 int i, j, k, natoms, nnn;
321 t_trxstatus **fp_in, **fp_out;
322 gmx_bool bCont, *bSet;
323 real t, first_time = 0;
331 for (i = 0; (i < nset); i++)
333 nnn = read_first_frame(oenv, &(fp_in[i]), fnms[i], &(trx[i]),
338 first_time = trx[i].time;
340 else if (natoms != nnn)
342 gmx_fatal(FARGS, "Trajectory file %s has %d atoms while previous trajs had %d atoms", fnms[i], nnn, natoms);
348 else if (t != trx[i].time)
350 gmx_fatal(FARGS, "Trajectory file %s has time %f while previous trajs had time %f", fnms[i], trx[i].time, t);
355 for (i = 0; (i < nset); i++)
357 fp_out[i] = open_trx(fnms_out[i], "w");
360 if (gmx_nint(time[k] - t) != 0)
362 gmx_fatal(FARGS, "First time in demuxing table does not match trajectories");
366 while ((k+1 < nval) && ((trx[0].time - time[k+1]) > dt_remd*0.1))
372 fprintf(debug, "trx[0].time = %g, time[k] = %g\n", trx[0].time, time[k]);
374 for (i = 0; (i < nset); i++)
378 for (i = 0; (i < nset); i++)
380 j = gmx_nint(value[i][k]);
381 range_check(j, 0, nset);
384 gmx_fatal(FARGS, "Demuxing the same replica %d twice at time %f",
389 if (dt == 0 || bRmod(trx[i].time, first_time, dt))
393 write_trxframe_indexed(fp_out[j], &trx[i], isize, index, NULL);
397 write_trxframe(fp_out[j], &trx[i], NULL);
403 for (i = 0; (i < nset); i++)
405 bCont = bCont && read_next_frame(oenv, fp_in[i], &trx[i]);
410 for (i = 0; (i < nset); i++)
413 close_trx(fp_out[i]);
417 int gmx_trjcat(int argc, char *argv[])
421 "[THISMODULE] concatenates several input trajectory files in sorted order. ",
422 "In case of double time frames the one in the later file is used. ",
423 "By specifying [TT]-settime[tt] you will be asked for the start time ",
424 "of each file. The input files are taken from the command line, ",
425 "such that a command like [TT]gmx trjcat -f *.trr -o fixed.trr[tt] should do ",
426 "the trick. Using [TT]-cat[tt], you can simply paste several files ",
427 "together without removal of frames with identical time stamps.[PAR]",
428 "One important option is inferred when the output file is amongst the",
429 "input files. In that case that particular file will be appended to",
430 "which implies you do not need to store double the amount of data.",
431 "Obviously the file to append to has to be the one with lowest starting",
432 "time since one can only append at the end of a file.[PAR]",
433 "If the [TT]-demux[tt] option is given, the N trajectories that are",
434 "read, are written in another order as specified in the [TT].xvg[tt] file.",
435 "The [TT].xvg[tt] file should contain something like:[PAR]",
436 "[TT]0 0 1 2 3 4 5[BR]",
437 "2 1 0 2 3 5 4[tt][BR]",
438 "Where the first number is the time, and subsequent numbers point to",
439 "trajectory indices.",
440 "The frames corresponding to the numbers present at the first line",
441 "are collected into the output trajectory. If the number of frames in",
442 "the trajectory does not match that in the [TT].xvg[tt] file then the program",
443 "tries to be smart. Beware."
445 static gmx_bool bVels = TRUE;
446 static gmx_bool bCat = FALSE;
447 static gmx_bool bSort = TRUE;
448 static gmx_bool bKeepLast = FALSE;
449 static gmx_bool bKeepLastAppend = FALSE;
450 static gmx_bool bOverwrite = FALSE;
451 static gmx_bool bSetTime = FALSE;
452 static gmx_bool bDeMux;
453 static real begin = -1;
454 static real end = -1;
460 { "-b", FALSE, etTIME,
461 { &begin }, "First time to use (%t)" },
462 { "-e", FALSE, etTIME,
463 { &end }, "Last time to use (%t)" },
464 { "-dt", FALSE, etTIME,
465 { &dt }, "Only write frame when t MOD dt = first time (%t)" },
466 { "-vel", FALSE, etBOOL,
467 { &bVels }, "Read and write velocities if possible" },
468 { "-settime", FALSE, etBOOL,
469 { &bSetTime }, "Change starting time interactively" },
470 { "-sort", FALSE, etBOOL,
471 { &bSort }, "Sort trajectory files (not frames)" },
472 { "-keeplast", FALSE, etBOOL,
473 { &bKeepLast }, "Keep overlapping frames at end of trajectory" },
474 { "-overwrite", FALSE, etBOOL,
475 { &bOverwrite }, "Overwrite overlapping frames during appending" },
476 { "-cat", FALSE, etBOOL,
477 { &bCat }, "Do not discard double time frames" }
479 #define npargs asize(pa)
480 int ftpin, i, frame, frame_out, step = 0, trjout = 0;
481 t_trxstatus *status, *trxout = NULL;
484 t_trxframe fr, frout;
485 char **fnms, **fnms_out, *in_file, *out_file;
487 gmx_bool bNewFile, bIndex, bWrite;
488 int earliersteps, nfile_in, nfile_out, *cont_type, last_ok_step;
489 real *readtime, *timest, *settime;
490 real first_time = 0, lasttime = NOTSET, last_ok_t = -1, timestep;
491 real last_frame_time, searchtime;
493 atom_id *index = NULL, imax;
495 real **val = NULL, *t = NULL, dt_remd;
496 int n, nset, ftpout = -1, prevEndStep = 0, filetype;
502 { efTRX, "-f", NULL, ffRDMULT },
503 { efTRO, "-o", NULL, ffWRMULT },
504 { efNDX, "-n", "index", ffOPTRD },
505 { efXVG, "-demux", "remd", ffOPTRD }
508 #define NFILE asize(fnm)
510 if (!parse_common_args(&argc, argv, PCA_BE_NICE | PCA_TIME_UNIT, NFILE, fnm,
511 asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
516 bIndex = ftp2bSet(efNDX, NFILE, fnm);
517 bDeMux = ftp2bSet(efXVG, NFILE, fnm);
518 bSort = bSort && !bDeMux;
523 printf("Select group for output\n");
524 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
527 for (i = 1; i < isize; i++)
529 imax = max(imax, index[i]);
536 val = read_xvg_time(opt2fn("-demux", NFILE, fnm), TRUE,
537 opt2parg_bSet("-b", npargs, pa), begin,
538 opt2parg_bSet("-e", npargs, pa), end, 1, &nset, &n,
540 printf("Read %d sets of %d points, dt = %g\n\n", nset, n, dt_remd);
543 fprintf(debug, "Dump of replica_index.xvg\n");
544 for (i = 0; (i < n); i++)
546 fprintf(debug, "%10g", t[i]);
547 for (j = 0; (j < nset); j++)
549 fprintf(debug, " %3d", gmx_nint(val[j][i]));
551 fprintf(debug, "\n");
556 nfile_in = opt2fns(&fnms, "-f", NFILE, fnm);
559 gmx_fatal(FARGS, "No input files!" );
562 if (bDeMux && (nfile_in != nset))
564 gmx_fatal(FARGS, "You have specified %d files and %d entries in the demux table", nfile_in, nset);
567 ftpin = fn2ftp(fnms[0]);
569 for (i = 1; i < nfile_in; i++)
571 if (ftpin != fn2ftp(fnms[i]))
573 gmx_fatal(FARGS, "All input files must be of the same format");
577 nfile_out = opt2fns(&fnms_out, "-o", NFILE, fnm);
580 gmx_fatal(FARGS, "No output files!");
582 if ((nfile_out > 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 && (nfile_out != nset) && (nfile_out != 1))
588 gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %d", nset, nfile_out);
592 if (nfile_out != nset)
594 char *buf = strdup(fnms_out[0]);
595 snew(fnms_out, nset);
596 for (i = 0; (i < nset); i++)
598 snew(fnms_out[i], strlen(buf)+32);
599 sprintf(fnms_out[i], "%d_%s", i, buf);
603 do_demux(nfile_in, fnms, fnms_out, n, val, t, dt_remd, isize, index, dt, oenv);
607 snew(readtime, nfile_in+1);
608 snew(timest, nfile_in+1);
609 scan_trj_files(fnms, nfile_in, readtime, timest, imax, oenv);
611 snew(settime, nfile_in+1);
612 snew(cont_type, nfile_in+1);
613 edit_files(fnms, nfile_in, 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 out_file = fnms_out[0];
620 ftpout = fn2ftp(out_file);
622 for (i = 0; ((i < nfile_in) && (n_append == -1)); i++)
624 if (strcmp(fnms[i], 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)",
641 /* Not checking input format, could be dangerous :-) */
642 /* Not checking output format, equally dangerous :-) */
646 /* the default is not to change the time at all,
647 * but this is overridden by the edit_files routine
657 trjtools_gmx_prepare_tng_writing(out_file, 'w', NULL, &trxout,
658 fnms[0], isize, NULL, index, grpname);
662 trjtools_gmx_prepare_tng_writing(out_file, 'w', NULL, &trxout,
663 fnms[0], -1, NULL, NULL, NULL);
668 trxout = open_trx(out_file, "w");
670 memset(&frout, 0, sizeof(frout));
676 if (!read_first_frame(oenv, &status, out_file, &fr, FLAGS))
678 gmx_fatal(FARGS, "Reading first frame from %s", out_file);
681 stfio = trx_get_fileio(status);
682 if (!bKeepLast && !bOverwrite)
684 fprintf(stderr, "\n\nWARNING: Appending without -overwrite implies -keeplast "
685 "between the first two files. \n"
686 "If the trajectories have an overlap and have not been written binary \n"
687 "reproducible this will produce an incorrect trajectory!\n\n");
689 filetype = gmx_fio_getftp(stfio);
690 /* Fails if last frame is incomplete
691 * We can't do anything about it without overwriting
693 if (filetype == efXTC || filetype == efTNG)
695 lasttime = trx_get_time_of_final_frame(status);
700 while (read_next_frame(oenv, status, &fr))
706 bKeepLastAppend = TRUE;
708 trxout = open_trx(out_file, "a");
712 if (gmx_fio_getftp(stfio) != efXTC)
714 gmx_fatal(FARGS, "Overwrite only supported for XTC." );
716 last_frame_time = trx_get_time_of_final_frame(status);
718 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
719 * or when seek time = 0 */
720 if (nfile_in > 1 && settime[1] < last_frame_time+timest[0]*0.5)
722 /* Jump to one time-frame before the start of next
724 searchtime = settime[1]-timest[0]*1.25;
728 searchtime = last_frame_time;
730 if (xtc_seek_time(stfio, searchtime, fr.natoms, TRUE))
732 gmx_fatal(FARGS, "Error seeking to append position.");
734 read_next_frame(oenv, status, &fr);
735 if (fabs(searchtime - fr.time) > timest[0]*0.5)
737 gmx_fatal(FARGS, "Error seeking: attempted to seek to %f but got %f.",
738 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.");
749 printf("\n Will append after %f \n", lasttime);
752 /* Lets stitch up some files */
753 timestep = timest[0];
754 for (i = n_append+1; (i < nfile_in); i++)
758 /* set the next time from the last frame in previous file */
761 /* When writing TNG the step determine which frame to write. Use an
762 * offset to be able to increase steps properly when changing files. */
765 prevEndStep = frout.step;
770 if (cont_type[i] == TIME_CONTINUE)
773 begin += 0.5*timestep;
774 settime[i] = frout.time;
775 cont_type[i] = TIME_EXPLICIT;
777 else if (cont_type[i] == TIME_LAST)
780 begin += 0.5*timestep;
782 /* Or, if the time in the next part should be changed by the
783 * same amount, start at half a timestep from the last time
784 * so we dont repeat frames.
786 /* I don't understand the comment above, but for all the cases
787 * I tried the code seems to work properly. B. Hess 2008-4-2.
790 /* Or, if time is set explicitly, we check for overlap/gap */
791 if (cont_type[i] == TIME_EXPLICIT)
793 if ( ( i < nfile_in ) &&
794 ( frout.time < settime[i]-1.5*timestep ) )
796 fprintf(stderr, "WARNING: Frames around t=%f %s have a different "
797 "spacing than the rest,\n"
798 "might be a gap or overlap that couldn't be corrected "
799 "automatically.\n", output_env_conv_time(oenv, frout.time),
800 output_env_get_time_unit(oenv));
805 /* if we don't have a timestep in the current file, use the old one */
808 timestep = timest[i];
810 read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
814 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
817 if (cont_type[i] == TIME_EXPLICIT)
819 t_corr = settime[i]-fr.time;
821 /* t_corr is the amount we want to change the time.
822 * If the user has chosen not to change the time for
823 * this part of the trajectory t_corr remains at
824 * the value it had in the last part, changing this
825 * by the same amount.
826 * If no value was given for the first trajectory part
827 * we let the time start at zero, see the edit_files routine.
833 if (lasttime != NOTSET)
835 printf("lasttime %g\n", lasttime);
840 /* copy the input frame to the output frame */
842 /* set the new time by adding the correct calculated above */
843 frout.time += t_corr;
846 frout.step += prevEndStep;
848 /* quit if we have reached the end of what should be written */
849 if ((end > 0) && (frout.time > end+GMX_REAL_EPS))
855 /* determine if we should write this frame (dt is handled elsewhere) */
856 if (bCat) /* write all frames of all files */
860 else if (bKeepLast || (bKeepLastAppend && i == 1))
861 /* write till last frame of this traj
862 and skip first frame(s) of next traj */
864 bWrite = ( frout.time > lasttime+0.5*timestep );
866 else /* write till first frame of next traj */
868 bWrite = ( frout.time < settime[i+1]-0.5*timestep );
871 if (bWrite && (frout.time >= begin) )
876 first_time = frout.time;
878 lasttime = frout.time;
879 if (dt == 0 || bRmod(frout.time, first_time, dt))
882 last_ok_t = frout.time;
885 fprintf(stderr, "\nContinue writing frames from %s t=%g %s, "
887 fnms[i], output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv),
894 write_trxframe_indexed(trxout, &frout, isize, index,
899 write_trxframe(trxout, &frout, NULL);
901 if ( ((frame % 10) == 0) || (frame < 10) )
903 fprintf(stderr, " -> frame %6d time %8.3f %s \r",
904 frame_out, output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv));
909 while (read_next_frame(oenv, status, &fr));
913 earliersteps += step;
919 fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
920 frame, output_env_conv_time(oenv, last_ok_t), output_env_get_time_unit(oenv));