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, 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.
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/fileio/gmxfio.h"
46 #include "gromacs/fileio/pdbio.h"
47 #include "gromacs/fileio/tngio.h"
48 #include "gromacs/fileio/tngio_for_tools.h"
49 #include "gromacs/fileio/trxio.h"
50 #include "gromacs/fileio/xtcio.h"
51 #include "gromacs/fileio/xvgr.h"
52 #include "gromacs/gmxana/gmx_ana.h"
53 #include "gromacs/legacyheaders/macros.h"
54 #include "gromacs/legacyheaders/names.h"
55 #include "gromacs/legacyheaders/typedefs.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/topology/index.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/smalloc.h"
62 #define TIME_EXPLICIT 0
63 #define TIME_CONTINUE 1
68 #define FLAGS (TRX_READ_X | TRX_READ_V | TRX_READ_F)
70 static void scan_trj_files(char **fnms, int nfiles, real *readtime,
71 real *timestep, atom_id imax,
72 const output_env_t oenv)
74 /* Check start time of all files */
81 for (i = 0; i < nfiles; i++)
83 ok = read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
87 gmx_fatal(FARGS, "\nCouldn't read frame from file." );
91 readtime[i] = fr.time;
96 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
107 if (natoms != fr.natoms)
109 gmx_fatal(FARGS, "\nDifferent numbers of atoms (%d/%d) in files",
115 if (fr.natoms <= imax)
117 gmx_fatal(FARGS, "\nNot enough atoms (%d) for index group (%d)",
122 ok = read_next_frame(oenv, status, &fr);
125 timestep[i] = fr.time - readtime[i];
146 fprintf(stderr, "\n");
150 static void sort_files(char **fnms, real *settime, int nfile)
156 for (i = 0; i < nfile; i++)
159 for (j = i + 1; j < nfile; j++)
161 if (settime[j] < settime[minidx])
168 timeswap = settime[i];
169 settime[i] = settime[minidx];
170 settime[minidx] = timeswap;
172 fnms[i] = fnms[minidx];
173 fnms[minidx] = chptr;
178 static void edit_files(char **fnms, int nfiles, real *readtime, real *timestep,
179 real *settime, int *cont_type, gmx_bool bSetTime,
180 gmx_bool bSort, const output_env_t oenv)
184 char inputstring[STRLEN], *chptr;
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",
197 output_env_get_time_unit(oenv));
201 " File Current start (%s) New start (%s)\n"
202 "---------------------------------------------------------\n",
203 output_env_get_time_unit(oenv), output_env_get_time_unit(oenv));
205 for (i = 0; i < nfiles; i++)
207 fprintf(stderr, "%25s %10.3f %s ", fnms[i],
208 output_env_conv_time(oenv, readtime[i]), output_env_get_time_unit(oenv));
212 if (NULL == fgets(inputstring, STRLEN - 1, stdin))
214 gmx_fatal(FARGS, "Error reading user input" );
217 inputstring[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 (i = 0; i < nfiles; i++)
262 settime[i] = readtime[i];
267 fprintf(stderr, "Sorting disabled.\n");
271 sort_files(fnms, settime, nfiles);
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 (i = 0; i < nfiles; i++)
279 switch (cont_type[i])
282 fprintf(stderr, "%25s %10.3f %s %10.3f %s",
284 output_env_conv_time(oenv, settime[i]), output_env_get_time_unit(oenv),
285 output_env_conv_time(oenv, timestep[i]), output_env_get_time_unit(oenv));
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", fnms[i]);
297 fprintf(stderr, "%25s Change by same amount as last file\n",
302 fprintf(stderr, "\n");
304 settime[nfiles] = FLT_MAX;
305 cont_type[nfiles] = TIME_EXPLICIT;
306 readtime[nfiles] = FLT_MAX;
309 static void do_demux(int nset, char *fnms[], char *fnms_out[], int nval,
310 real **value, real *time, real dt_remd, int isize,
311 atom_id index[], real dt, const output_env_t oenv)
313 int i, j, k, natoms, nnn;
314 t_trxstatus **fp_in, **fp_out;
315 gmx_bool bCont, *bSet;
316 real t, first_time = 0;
324 for (i = 0; (i < nset); i++)
326 nnn = read_first_frame(oenv, &(fp_in[i]), fnms[i], &(trx[i]),
331 first_time = trx[i].time;
333 else if (natoms != nnn)
335 gmx_fatal(FARGS, "Trajectory file %s has %d atoms while previous trajs had %d atoms", fnms[i], nnn, natoms);
341 else if (t != trx[i].time)
343 gmx_fatal(FARGS, "Trajectory file %s has time %f while previous trajs had time %f", fnms[i], trx[i].time, t);
348 for (i = 0; (i < nset); i++)
350 fp_out[i] = open_trx(fnms_out[i], "w");
353 if (gmx_nint(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 (i = 0; (i < nset); i++)
371 for (i = 0; (i < nset); i++)
373 j = gmx_nint(value[i][k]);
374 range_check(j, 0, nset);
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, NULL);
390 write_trxframe(fp_out[j], &trx[i], NULL);
396 for (i = 0; (i < nset); i++)
398 bCont = bCont && read_next_frame(oenv, fp_in[i], &trx[i]);
403 for (i = 0; (i < nset); 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 bVels = TRUE;
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 { "-vel", FALSE, etBOOL,
462 { &bVels }, "Read and write velocities if possible" },
463 { "-settime", FALSE, etBOOL,
464 { &bSetTime }, "Change starting time interactively" },
465 { "-sort", FALSE, etBOOL,
466 { &bSort }, "Sort trajectory files (not frames)" },
467 { "-keeplast", FALSE, etBOOL,
468 { &bKeepLast }, "Keep overlapping frames at end of trajectory" },
469 { "-overwrite", FALSE, etBOOL,
470 { &bOverwrite }, "Overwrite overlapping frames during appending" },
471 { "-cat", FALSE, etBOOL,
472 { &bCat }, "Do not discard double time frames" }
474 #define npargs asize(pa)
475 int ftpin, i, frame, frame_out, step = 0, trjout = 0;
476 t_trxstatus *status, *trxout = NULL;
479 t_trxframe fr, frout;
480 char **fnms, **fnms_out, *in_file, *out_file;
482 gmx_bool bNewFile, bIndex, bWrite;
483 int earliersteps, nfile_in, nfile_out, *cont_type, last_ok_step;
484 real *readtime, *timest, *settime;
485 real first_time = 0, lasttime = NOTSET, last_ok_t = -1, timestep;
486 real last_frame_time, searchtime;
488 atom_id *index = NULL, imax;
490 real **val = NULL, *t = NULL, dt_remd;
491 int n, nset, ftpout = -1, prevEndStep = 0, filetype;
497 { efTRX, "-f", NULL, ffRDMULT },
498 { efTRO, "-o", NULL, ffWRMULT },
499 { efNDX, "-n", "index", ffOPTRD },
500 { efXVG, "-demux", "remd", ffOPTRD }
503 #define NFILE asize(fnm)
505 if (!parse_common_args(&argc, argv, PCA_TIME_UNIT, NFILE, fnm,
506 asize(pa), pa, asize(desc), desc, 0, NULL, &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 = 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", gmx_nint(val[j][i]));
546 fprintf(debug, "\n");
551 nfile_in = opt2fns(&fnms, "-f", NFILE, fnm);
554 gmx_fatal(FARGS, "No input files!" );
557 if (bDeMux && (nfile_in != nset))
559 gmx_fatal(FARGS, "You have specified %d files and %d entries in the demux table", nfile_in, nset);
562 ftpin = fn2ftp(fnms[0]);
564 for (i = 1; i < nfile_in; i++)
566 if (ftpin != fn2ftp(fnms[i]))
568 gmx_fatal(FARGS, "All input files must be of the same format");
572 nfile_out = opt2fns(&fnms_out, "-o", NFILE, fnm);
575 gmx_fatal(FARGS, "No output files!");
577 if ((nfile_out > 1) && !bDeMux)
579 gmx_fatal(FARGS, "Don't know what to do with more than 1 output file if not demultiplexing");
581 else if (bDeMux && (nfile_out != nset) && (nfile_out != 1))
583 gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %d", nset, nfile_out);
587 if (nfile_out != nset)
589 char *buf = gmx_strdup(fnms_out[0]);
590 snew(fnms_out, nset);
591 for (i = 0; (i < nset); i++)
593 snew(fnms_out[i], strlen(buf)+32);
594 sprintf(fnms_out[i], "%d_%s", i, buf);
598 do_demux(nfile_in, fnms, fnms_out, n, val, t, dt_remd, isize, index, dt, oenv);
602 snew(readtime, nfile_in+1);
603 snew(timest, nfile_in+1);
604 scan_trj_files(fnms, nfile_in, readtime, timest, imax, oenv);
606 snew(settime, nfile_in+1);
607 snew(cont_type, nfile_in+1);
608 edit_files(fnms, nfile_in, readtime, timest, settime, cont_type, bSetTime, bSort,
611 /* Check whether the output file is amongst the input files
612 * This has to be done after sorting etc.
614 out_file = fnms_out[0];
615 ftpout = fn2ftp(out_file);
617 for (i = 0; ((i < nfile_in) && (n_append == -1)); i++)
619 if (strcmp(fnms[i], out_file) == 0)
626 fprintf(stderr, "Will append to %s rather than creating a new file\n",
629 else if (n_append != -1)
631 gmx_fatal(FARGS, "Can only append to the first file which is %s (not %s)",
636 /* Not checking input format, could be dangerous :-) */
637 /* Not checking output format, equally dangerous :-) */
641 /* the default is not to change the time at all,
642 * but this is overridden by the edit_files routine
652 gmx_fatal(FARGS, "When writing TNG the input file format must also be TNG");
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));