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, 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.
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/commandline/pargs.h"
52 #include "gromacs/fileio/futil.h"
53 #include "gromacs/fileio/pdbio.h"
54 #include "gromacs/fileio/confio.h"
58 #include "gromacs/fileio/xtcio.h"
61 #include "gromacs/fileio/g87io.h"
64 #include "gromacs/fileio/xdrf.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;
482 t_trxframe fr, frout;
483 char **fnms, **fnms_out, *in_file, *out_file;
485 t_trxstatus *trxout = NULL;
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;
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 nfile_out = opt2fns(&fnms_out, "-o", NFILE, fnm);
569 gmx_fatal(FARGS, "No output files!");
571 if ((nfile_out > 1) && !bDeMux)
573 gmx_fatal(FARGS, "Don't know what to do with more than 1 output file if not demultiplexing");
575 else if (bDeMux && (nfile_out != nset) && (nfile_out != 1))
577 gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %d", nset, nfile_out);
581 if (nfile_out != nset)
583 char *buf = strdup(fnms_out[0]);
584 snew(fnms_out, nset);
585 for (i = 0; (i < nset); i++)
587 snew(fnms_out[i], strlen(buf)+32);
588 sprintf(fnms_out[i], "%d_%s", i, buf);
592 do_demux(nfile_in, fnms, fnms_out, n, val, t, dt_remd, isize, index, dt, oenv);
596 snew(readtime, nfile_in+1);
597 snew(timest, nfile_in+1);
598 scan_trj_files(fnms, nfile_in, readtime, timest, imax, oenv);
600 snew(settime, nfile_in+1);
601 snew(cont_type, nfile_in+1);
602 edit_files(fnms, nfile_in, readtime, timest, settime, cont_type, bSetTime, bSort,
605 /* Check whether the output file is amongst the input files
606 * This has to be done after sorting etc.
608 out_file = fnms_out[0];
610 for (i = 0; ((i < nfile_in) && (n_append == -1)); i++)
612 if (strcmp(fnms[i], out_file) == 0)
619 fprintf(stderr, "Will append to %s rather than creating a new file\n",
622 else if (n_append != -1)
624 gmx_fatal(FARGS, "Can only append to the first file which is %s (not %s)",
629 /* Not checking input format, could be dangerous :-) */
630 /* Not checking output format, equally dangerous :-) */
634 /* the default is not to change the time at all,
635 * but this is overridden by the edit_files routine
641 trxout = open_trx(out_file, "w");
642 memset(&frout, 0, sizeof(frout));
648 if (!read_first_frame(oenv, &status, out_file, &fr, FLAGS))
650 gmx_fatal(FARGS, "Reading first frame from %s", out_file);
653 stfio = trx_get_fileio(status);
654 if (!bKeepLast && !bOverwrite)
656 fprintf(stderr, "\n\nWARNING: Appending without -overwrite implies -keeplast "
657 "between the first two files. \n"
658 "If the trajectories have an overlap and have not been written binary \n"
659 "reproducible this will produce an incorrect trajectory!\n\n");
661 /* Fails if last frame is incomplete
662 * We can't do anything about it without overwriting
664 if (gmx_fio_getftp(stfio) == efXTC)
667 xdr_xtc_get_last_frame_time(gmx_fio_getfp(stfio),
668 gmx_fio_getxdr(stfio),
673 gmx_fatal(FARGS, "Error reading last frame. Maybe seek not supported." );
678 while (read_next_frame(oenv, status, &fr))
684 bKeepLastAppend = TRUE;
686 trxout = open_trx(out_file, "a");
690 if (gmx_fio_getftp(stfio) != efXTC)
692 gmx_fatal(FARGS, "Overwrite only supported for XTC." );
695 xdr_xtc_get_last_frame_time(gmx_fio_getfp(stfio),
696 gmx_fio_getxdr(stfio),
700 gmx_fatal(FARGS, "Error reading last frame. Maybe seek not supported." );
702 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
703 * or when seek time = 0 */
704 if (nfile_in > 1 && settime[1] < last_frame_time+timest[0]*0.5)
706 /* Jump to one time-frame before the start of next
708 searchtime = settime[1]-timest[0]*1.25;
712 searchtime = last_frame_time;
714 if (xtc_seek_time(stfio, searchtime, fr.natoms, TRUE))
716 gmx_fatal(FARGS, "Error seeking to append position.");
718 read_next_frame(oenv, status, &fr);
719 if (fabs(searchtime - fr.time) > timest[0]*0.5)
721 gmx_fatal(FARGS, "Error seeking: attempted to seek to %f but got %f.",
722 searchtime, fr.time);
725 fpos = gmx_fio_ftell(stfio);
727 trxout = open_trx(out_file, "r+");
728 if (gmx_fio_seek(trx_get_fileio(trxout), fpos))
730 gmx_fatal(FARGS, "Error seeking to append position.");
733 printf("\n Will append after %f \n", lasttime);
736 /* Lets stitch up some files */
737 timestep = timest[0];
738 for (i = n_append+1; (i < nfile_in); i++)
742 /* set the next time from the last frame in previous file */
747 if (cont_type[i] == TIME_CONTINUE)
750 begin += 0.5*timestep;
751 settime[i] = frout.time;
752 cont_type[i] = TIME_EXPLICIT;
754 else if (cont_type[i] == TIME_LAST)
757 begin += 0.5*timestep;
759 /* Or, if the time in the next part should be changed by the
760 * same amount, start at half a timestep from the last time
761 * so we dont repeat frames.
763 /* I don't understand the comment above, but for all the cases
764 * I tried the code seems to work properly. B. Hess 2008-4-2.
767 /* Or, if time is set explicitly, we check for overlap/gap */
768 if (cont_type[i] == TIME_EXPLICIT)
770 if ( ( i < nfile_in ) &&
771 ( frout.time < settime[i]-1.5*timestep ) )
773 fprintf(stderr, "WARNING: Frames around t=%f %s have a different "
774 "spacing than the rest,\n"
775 "might be a gap or overlap that couldn't be corrected "
776 "automatically.\n", output_env_conv_time(oenv, frout.time),
777 output_env_get_time_unit(oenv));
782 /* if we don't have a timestep in the current file, use the old one */
785 timestep = timest[i];
787 read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
791 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
794 if (cont_type[i] == TIME_EXPLICIT)
796 t_corr = settime[i]-fr.time;
798 /* t_corr is the amount we want to change the time.
799 * If the user has chosen not to change the time for
800 * this part of the trajectory t_corr remains at
801 * the value it had in the last part, changing this
802 * by the same amount.
803 * If no value was given for the first trajectory part
804 * we let the time start at zero, see the edit_files routine.
810 if (lasttime != NOTSET)
812 printf("lasttime %g\n", lasttime);
817 /* copy the input frame to the output frame */
819 /* set the new time by adding the correct calculated above */
820 frout.time += t_corr;
821 /* quit if we have reached the end of what should be written */
822 if ((end > 0) && (frout.time > end+GMX_REAL_EPS))
828 /* determine if we should write this frame (dt is handled elsewhere) */
829 if (bCat) /* write all frames of all files */
833 else if (bKeepLast || (bKeepLastAppend && i == 1))
834 /* write till last frame of this traj
835 and skip first frame(s) of next traj */
837 bWrite = ( frout.time > lasttime+0.5*timestep );
839 else /* write till first frame of next traj */
841 bWrite = ( frout.time < settime[i+1]-0.5*timestep );
844 if (bWrite && (frout.time >= begin) )
849 first_time = frout.time;
851 lasttime = frout.time;
852 if (dt == 0 || bRmod(frout.time, first_time, dt))
855 last_ok_t = frout.time;
858 fprintf(stderr, "\nContinue writing frames from %s t=%g %s, "
860 fnms[i], output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv),
867 write_trxframe_indexed(trxout, &frout, isize, index,
872 write_trxframe(trxout, &frout, NULL);
874 if ( ((frame % 10) == 0) || (frame < 10) )
876 fprintf(stderr, " -> frame %6d time %8.3f %s \r",
877 frame_out, output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv));
882 while (read_next_frame(oenv, status, &fr));
886 earliersteps += step;
892 fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
893 frame, output_env_conv_time(oenv, last_ok_t), output_env_get_time_unit(oenv));