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.
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/tpxio.h"
50 #include "gromacs/fileio/trnio.h"
51 #include "gromacs/fileio/trxio.h"
52 #include "gromacs/fileio/xtcio.h"
53 #include "gromacs/fileio/xvgr.h"
54 #include "gromacs/gmxana/gmx_ana.h"
55 #include "gromacs/legacyheaders/macros.h"
56 #include "gromacs/legacyheaders/names.h"
57 #include "gromacs/legacyheaders/typedefs.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/topology/index.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
64 #define TIME_EXPLICIT 0
65 #define TIME_CONTINUE 1
70 #define FLAGS (TRX_READ_X | TRX_READ_V | TRX_READ_F)
72 static void scan_trj_files(char **fnms, int nfiles, real *readtime,
73 real *timestep, atom_id imax,
74 const output_env_t oenv)
76 /* Check start time of all files */
83 for (i = 0; i < nfiles; i++)
85 ok = read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
89 gmx_fatal(FARGS, "\nCouldn't read frame from file." );
93 readtime[i] = fr.time;
98 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
109 if (natoms != fr.natoms)
111 gmx_fatal(FARGS, "\nDifferent numbers of atoms (%d/%d) in files",
117 if (fr.natoms <= imax)
119 gmx_fatal(FARGS, "\nNot enough atoms (%d) for index group (%d)",
124 ok = read_next_frame(oenv, status, &fr);
127 timestep[i] = fr.time - readtime[i];
148 fprintf(stderr, "\n");
152 static void sort_files(char **fnms, real *settime, int nfile)
158 for (i = 0; i < nfile; i++)
161 for (j = i + 1; j < nfile; j++)
163 if (settime[j] < settime[minidx])
170 timeswap = settime[i];
171 settime[i] = settime[minidx];
172 settime[minidx] = timeswap;
174 fnms[i] = fnms[minidx];
175 fnms[minidx] = chptr;
180 static void edit_files(char **fnms, int nfiles, real *readtime, real *timestep,
181 real *settime, int *cont_type, gmx_bool bSetTime,
182 gmx_bool bSort, const output_env_t oenv)
186 char inputstring[STRLEN], *chptr;
190 fprintf(stderr, "\n\nEnter the new start time (%s) for each file.\n"
191 "There are two special options, both disable sorting:\n\n"
192 "c (continue) - The start time is taken from the end\n"
193 "of the previous file. Use it when your continuation run\n"
194 "restarts with t=0.\n\n"
195 "l (last) - The time in this file will be changed the\n"
196 "same amount as in the previous. Use it when the time in the\n"
197 "new run continues from the end of the previous one,\n"
198 "since this takes possible overlap into account.\n\n",
199 output_env_get_time_unit(oenv));
203 " File Current start (%s) New start (%s)\n"
204 "---------------------------------------------------------\n",
205 output_env_get_time_unit(oenv), output_env_get_time_unit(oenv));
207 for (i = 0; i < nfiles; i++)
209 fprintf(stderr, "%25s %10.3f %s ", fnms[i],
210 output_env_conv_time(oenv, readtime[i]), output_env_get_time_unit(oenv));
214 if (NULL == fgets(inputstring, STRLEN - 1, stdin))
216 gmx_fatal(FARGS, "Error reading user input" );
219 inputstring[strlen(inputstring)-1] = 0;
221 if (inputstring[0] == 'c' || inputstring[0] == 'C')
223 cont_type[i] = TIME_CONTINUE;
226 settime[i] = FLT_MAX;
228 else if (inputstring[0] == 'l' ||
229 inputstring[0] == 'L')
231 cont_type[i] = TIME_LAST;
234 settime[i] = FLT_MAX;
238 settime[i] = strtod(inputstring, &chptr)*
239 output_env_get_time_invfactor(oenv);
240 if (chptr == inputstring)
242 fprintf(stderr, "'%s' not recognized as a floating point number, 'c' or 'l'. "
243 "Try again: ", inputstring);
247 cont_type[i] = TIME_EXPLICIT;
254 if (cont_type[0] != TIME_EXPLICIT)
256 cont_type[0] = TIME_EXPLICIT;
262 for (i = 0; i < nfiles; i++)
264 settime[i] = readtime[i];
269 fprintf(stderr, "Sorting disabled.\n");
273 sort_files(fnms, settime, nfiles);
275 /* Write out the new order and start times */
276 fprintf(stderr, "\nSummary of files and start times used:\n\n"
277 " File Start time Time step\n"
278 "---------------------------------------------------------\n");
279 for (i = 0; i < nfiles; i++)
281 switch (cont_type[i])
284 fprintf(stderr, "%25s %10.3f %s %10.3f %s",
286 output_env_conv_time(oenv, settime[i]), output_env_get_time_unit(oenv),
287 output_env_conv_time(oenv, timestep[i]), output_env_get_time_unit(oenv));
289 cont_type[i-1] == TIME_EXPLICIT && settime[i] == settime[i-1])
291 fprintf(stderr, " WARNING: same Start time as previous");
293 fprintf(stderr, "\n");
296 fprintf(stderr, "%25s Continue from last file\n", fnms[i]);
299 fprintf(stderr, "%25s Change by same amount as last file\n",
304 fprintf(stderr, "\n");
306 settime[nfiles] = FLT_MAX;
307 cont_type[nfiles] = TIME_EXPLICIT;
308 readtime[nfiles] = FLT_MAX;
311 static void do_demux(int nset, char *fnms[], char *fnms_out[], int nval,
312 real **value, real *time, real dt_remd, int isize,
313 atom_id index[], real dt, const output_env_t oenv)
315 int i, j, k, natoms, nnn;
316 t_trxstatus **fp_in, **fp_out;
317 gmx_bool bCont, *bSet;
318 real t, first_time = 0;
326 for (i = 0; (i < nset); i++)
328 nnn = read_first_frame(oenv, &(fp_in[i]), fnms[i], &(trx[i]),
333 first_time = trx[i].time;
335 else if (natoms != nnn)
337 gmx_fatal(FARGS, "Trajectory file %s has %d atoms while previous trajs had %d atoms", fnms[i], nnn, natoms);
343 else if (t != trx[i].time)
345 gmx_fatal(FARGS, "Trajectory file %s has time %f while previous trajs had time %f", fnms[i], trx[i].time, t);
350 for (i = 0; (i < nset); i++)
352 fp_out[i] = open_trx(fnms_out[i], "w");
355 if (gmx_nint(time[k] - t) != 0)
357 gmx_fatal(FARGS, "First time in demuxing table does not match trajectories");
361 while ((k+1 < nval) && ((trx[0].time - time[k+1]) > dt_remd*0.1))
367 fprintf(debug, "trx[0].time = %g, time[k] = %g\n", trx[0].time, time[k]);
369 for (i = 0; (i < nset); i++)
373 for (i = 0; (i < nset); i++)
375 j = gmx_nint(value[i][k]);
376 range_check(j, 0, nset);
379 gmx_fatal(FARGS, "Demuxing the same replica %d twice at time %f",
384 if (dt == 0 || bRmod(trx[i].time, first_time, dt))
388 write_trxframe_indexed(fp_out[j], &trx[i], isize, index, NULL);
392 write_trxframe(fp_out[j], &trx[i], NULL);
398 for (i = 0; (i < nset); i++)
400 bCont = bCont && read_next_frame(oenv, fp_in[i], &trx[i]);
405 for (i = 0; (i < nset); i++)
408 close_trx(fp_out[i]);
412 int gmx_trjcat(int argc, char *argv[])
416 "[THISMODULE] concatenates several input trajectory files in sorted order. ",
417 "In case of double time frames the one in the later file is used. ",
418 "By specifying [TT]-settime[tt] you will be asked for the start time ",
419 "of each file. The input files are taken from the command line, ",
420 "such that a command like [TT]gmx trjcat -f *.trr -o fixed.trr[tt] should do ",
421 "the trick. Using [TT]-cat[tt], you can simply paste several files ",
422 "together without removal of frames with identical time stamps.[PAR]",
423 "One important option is inferred when the output file is amongst the",
424 "input files. In that case that particular file will be appended to",
425 "which implies you do not need to store double the amount of data.",
426 "Obviously the file to append to has to be the one with lowest starting",
427 "time since one can only append at the end of a file.[PAR]",
428 "If the [TT]-demux[tt] option is given, the N trajectories that are",
429 "read, are written in another order as specified in the [TT].xvg[tt] file.",
430 "The [TT].xvg[tt] file should contain something like:[PAR]",
431 "[TT]0 0 1 2 3 4 5[BR]",
432 "2 1 0 2 3 5 4[tt][BR]",
433 "Where 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 [TT].xvg[tt] 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 trjtools_gmx_prepare_tng_writing(out_file, 'w', NULL, &trxout,
653 fnms[0], isize, NULL, index, grpname);
657 trjtools_gmx_prepare_tng_writing(out_file, 'w', NULL, &trxout,
658 fnms[0], -1, NULL, NULL, NULL);
663 trxout = open_trx(out_file, "w");
665 memset(&frout, 0, sizeof(frout));
671 if (!read_first_frame(oenv, &status, out_file, &fr, FLAGS))
673 gmx_fatal(FARGS, "Reading first frame from %s", out_file);
676 stfio = trx_get_fileio(status);
677 if (!bKeepLast && !bOverwrite)
679 fprintf(stderr, "\n\nWARNING: Appending without -overwrite implies -keeplast "
680 "between the first two files. \n"
681 "If the trajectories have an overlap and have not been written binary \n"
682 "reproducible this will produce an incorrect trajectory!\n\n");
684 filetype = gmx_fio_getftp(stfio);
685 /* Fails if last frame is incomplete
686 * We can't do anything about it without overwriting
688 if (filetype == efXTC || filetype == efTNG)
690 lasttime = trx_get_time_of_final_frame(status);
695 while (read_next_frame(oenv, status, &fr))
701 bKeepLastAppend = TRUE;
703 trxout = open_trx(out_file, "a");
707 if (gmx_fio_getftp(stfio) != efXTC)
709 gmx_fatal(FARGS, "Overwrite only supported for XTC." );
711 last_frame_time = trx_get_time_of_final_frame(status);
713 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
714 * or when seek time = 0 */
715 if (nfile_in > 1 && settime[1] < last_frame_time+timest[0]*0.5)
717 /* Jump to one time-frame before the start of next
719 searchtime = settime[1]-timest[0]*1.25;
723 searchtime = last_frame_time;
725 if (xtc_seek_time(stfio, searchtime, fr.natoms, TRUE))
727 gmx_fatal(FARGS, "Error seeking to append position.");
729 read_next_frame(oenv, status, &fr);
730 if (fabs(searchtime - fr.time) > timest[0]*0.5)
732 gmx_fatal(FARGS, "Error seeking: attempted to seek to %f but got %f.",
733 searchtime, fr.time);
736 fpos = gmx_fio_ftell(stfio);
738 trxout = open_trx(out_file, "r+");
739 if (gmx_fio_seek(trx_get_fileio(trxout), fpos))
741 gmx_fatal(FARGS, "Error seeking to append position.");
744 printf("\n Will append after %f \n", lasttime);
747 /* Lets stitch up some files */
748 timestep = timest[0];
749 for (i = n_append+1; (i < nfile_in); i++)
753 /* set the next time from the last frame in previous file */
756 /* When writing TNG the step determine which frame to write. Use an
757 * offset to be able to increase steps properly when changing files. */
760 prevEndStep = frout.step;
765 if (cont_type[i] == TIME_CONTINUE)
768 begin += 0.5*timestep;
769 settime[i] = frout.time;
770 cont_type[i] = TIME_EXPLICIT;
772 else if (cont_type[i] == TIME_LAST)
775 begin += 0.5*timestep;
777 /* Or, if the time in the next part should be changed by the
778 * same amount, start at half a timestep from the last time
779 * so we dont repeat frames.
781 /* I don't understand the comment above, but for all the cases
782 * I tried the code seems to work properly. B. Hess 2008-4-2.
785 /* Or, if time is set explicitly, we check for overlap/gap */
786 if (cont_type[i] == TIME_EXPLICIT)
788 if ( ( i < nfile_in ) &&
789 ( frout.time < settime[i]-1.5*timestep ) )
791 fprintf(stderr, "WARNING: Frames around t=%f %s have a different "
792 "spacing than the rest,\n"
793 "might be a gap or overlap that couldn't be corrected "
794 "automatically.\n", output_env_conv_time(oenv, frout.time),
795 output_env_get_time_unit(oenv));
800 /* if we don't have a timestep in the current file, use the old one */
803 timestep = timest[i];
805 read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
809 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
812 if (cont_type[i] == TIME_EXPLICIT)
814 t_corr = settime[i]-fr.time;
816 /* t_corr is the amount we want to change the time.
817 * If the user has chosen not to change the time for
818 * this part of the trajectory t_corr remains at
819 * the value it had in the last part, changing this
820 * by the same amount.
821 * If no value was given for the first trajectory part
822 * we let the time start at zero, see the edit_files routine.
828 if (lasttime != NOTSET)
830 printf("lasttime %g\n", lasttime);
835 /* copy the input frame to the output frame */
837 /* set the new time by adding the correct calculated above */
838 frout.time += t_corr;
841 frout.step += prevEndStep;
843 /* quit if we have reached the end of what should be written */
844 if ((end > 0) && (frout.time > end+GMX_REAL_EPS))
850 /* determine if we should write this frame (dt is handled elsewhere) */
851 if (bCat) /* write all frames of all files */
855 else if (bKeepLast || (bKeepLastAppend && i == 1))
856 /* write till last frame of this traj
857 and skip first frame(s) of next traj */
859 bWrite = ( frout.time > lasttime+0.5*timestep );
861 else /* write till first frame of next traj */
863 bWrite = ( frout.time < settime[i+1]-0.5*timestep );
866 if (bWrite && (frout.time >= begin) )
871 first_time = frout.time;
873 lasttime = frout.time;
874 if (dt == 0 || bRmod(frout.time, first_time, dt))
877 last_ok_t = frout.time;
880 fprintf(stderr, "\nContinue writing frames from %s t=%g %s, "
882 fnms[i], output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv),
889 write_trxframe_indexed(trxout, &frout, isize, index,
894 write_trxframe(trxout, &frout, NULL);
896 if ( ((frame % 10) == 0) || (frame < 10) )
898 fprintf(stderr, " -> frame %6d time %8.3f %s \r",
899 frame_out, output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv));
904 while (read_next_frame(oenv, status, &fr));
908 earliersteps += step;
914 fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
915 frame, output_env_conv_time(oenv, last_ok_t), output_env_get_time_unit(oenv));