3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
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[])
417 "[TT]trjcat[tt] concatenates several input trajectory files in sorted order. ",
418 "In case of double time frames the one in the later file is used. ",
419 "By specifying [TT]-settime[tt] you will be asked for the start time ",
420 "of each file. The input files are taken from the command line, ",
421 "such that a command like [TT]trjcat -f *.trr -o fixed.trr[tt] should do ",
422 "the trick. Using [TT]-cat[tt], you can simply paste several files ",
423 "together without removal of frames with identical time stamps.[PAR]",
424 "One important option is inferred when the output file is amongst the",
425 "input files. In that case that particular file will be appended to",
426 "which implies you do not need to store double the amount of data.",
427 "Obviously the file to append to has to be the one with lowest starting",
428 "time since one can only append at the end of a file.[PAR]",
429 "If the [TT]-demux[tt] option is given, the N trajectories that are",
430 "read, are written in another order as specified in the [TT].xvg[tt] file.",
431 "The [TT].xvg[tt] file should contain something like:[PAR]",
432 "[TT]0 0 1 2 3 4 5[BR]",
433 "2 1 0 2 3 5 4[tt][BR]",
434 "Where the first number is the time, and subsequent numbers point to",
435 "trajectory indices.",
436 "The frames corresponding to the numbers present at the first line",
437 "are collected into the output trajectory. If the number of frames in",
438 "the trajectory does not match that in the [TT].xvg[tt] file then the program",
439 "tries to be smart. Beware."
441 static gmx_bool bVels = TRUE;
443 static gmx_bool bCat = FALSE;
444 static gmx_bool bSort = TRUE;
445 static gmx_bool bKeepLast = FALSE;
446 static gmx_bool bKeepLastAppend = FALSE;
447 static gmx_bool bOverwrite = FALSE;
448 static gmx_bool bSetTime = FALSE;
449 static gmx_bool bDeMux;
450 static real begin = -1;
451 static real end = -1;
457 { "-b", FALSE, etTIME,
458 { &begin }, "First time to use (%t)" },
459 { "-e", FALSE, etTIME,
460 { &end }, "Last time to use (%t)" },
461 { "-dt", FALSE, etTIME,
462 { &dt }, "Only write frame when t MOD dt = first time (%t)" },
463 { "-prec", FALSE, etINT,
464 { &prec }, "Precision for [TT].xtc[tt] and [TT].gro[tt] writing in number of decimal places" },
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;
483 t_trxframe fr, frout;
484 char **fnms, **fnms_out, *in_file, *out_file;
486 t_trxstatus *trxout = NULL;
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;
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 parse_common_args(&argc, argv, PCA_BE_NICE | PCA_TIME_UNIT, NFILE, fnm,
511 asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
513 bIndex = ftp2bSet(efNDX, NFILE, fnm);
514 bDeMux = ftp2bSet(efXVG, NFILE, fnm);
515 bSort = bSort && !bDeMux;
520 printf("Select group for output\n");
521 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
524 for (i = 1; i < isize; i++)
526 imax = max(imax, index[i]);
533 val = read_xvg_time(opt2fn("-demux", NFILE, fnm), TRUE,
534 opt2parg_bSet("-b", npargs, pa), begin,
535 opt2parg_bSet("-e", npargs, pa), end, 1, &nset, &n,
537 printf("Read %d sets of %d points, dt = %g\n\n", nset, n, dt_remd);
540 fprintf(debug, "Dump of replica_index.xvg\n");
541 for (i = 0; (i < n); i++)
543 fprintf(debug, "%10g", t[i]);
544 for (j = 0; (j < nset); j++)
546 fprintf(debug, " %3d", gmx_nint(val[j][i]));
548 fprintf(debug, "\n");
552 /* prec is in nr of decimal places, xtcprec is a multiplication factor: */
554 for (i = 0; i < prec; i++)
559 nfile_in = opt2fns(&fnms, "-f", NFILE, fnm);
562 gmx_fatal(FARGS, "No input files!" );
565 if (bDeMux && (nfile_in != nset))
567 gmx_fatal(FARGS, "You have specified %d files and %d entries in the demux table", nfile_in, nset);
570 nfile_out = opt2fns(&fnms_out, "-o", NFILE, fnm);
573 gmx_fatal(FARGS, "No output files!");
575 if ((nfile_out > 1) && !bDeMux)
577 gmx_fatal(FARGS, "Don't know what to do with more than 1 output file if not demultiplexing");
579 else if (bDeMux && (nfile_out != nset) && (nfile_out != 1))
581 gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %d", nset, nfile_out);
585 if (nfile_out != nset)
587 char *buf = strdup(fnms_out[0]);
588 snew(fnms_out, nset);
589 for (i = 0; (i < nset); i++)
591 snew(fnms_out[i], strlen(buf)+32);
592 sprintf(fnms_out[i], "%d_%s", i, buf);
596 do_demux(nfile_in, fnms, fnms_out, n, val, t, dt_remd, isize, index, dt, oenv);
600 snew(readtime, nfile_in+1);
601 snew(timest, nfile_in+1);
602 scan_trj_files(fnms, nfile_in, readtime, timest, imax, oenv);
604 snew(settime, nfile_in+1);
605 snew(cont_type, nfile_in+1);
606 edit_files(fnms, nfile_in, readtime, timest, settime, cont_type, bSetTime, bSort,
609 /* Check whether the output file is amongst the input files
610 * This has to be done after sorting etc.
612 out_file = fnms_out[0];
614 for (i = 0; ((i < nfile_in) && (n_append == -1)); i++)
616 if (strcmp(fnms[i], out_file) == 0)
623 fprintf(stderr, "Will append to %s rather than creating a new file\n",
626 else if (n_append != -1)
628 gmx_fatal(FARGS, "Can only append to the first file which is %s (not %s)",
633 /* Not checking input format, could be dangerous :-) */
634 /* Not checking output format, equally dangerous :-) */
638 /* the default is not to change the time at all,
639 * but this is overridden by the edit_files routine
645 trxout = open_trx(out_file, "w");
646 memset(&frout, 0, sizeof(frout));
652 if (!read_first_frame(oenv, &status, out_file, &fr, FLAGS))
654 gmx_fatal(FARGS, "Reading first frame from %s", out_file);
657 stfio = trx_get_fileio(status);
658 if (!bKeepLast && !bOverwrite)
660 fprintf(stderr, "\n\nWARNING: Appending without -overwrite implies -keeplast "
661 "between the first two files. \n"
662 "If the trajectories have an overlap and have not been written binary \n"
663 "reproducible this will produce an incorrect trajectory!\n\n");
665 /* Fails if last frame is incomplete
666 * We can't do anything about it without overwriting
668 if (gmx_fio_getftp(stfio) == efXTC)
671 xdr_xtc_get_last_frame_time(gmx_fio_getfp(stfio),
672 gmx_fio_getxdr(stfio),
677 gmx_fatal(FARGS, "Error reading last frame. Maybe seek not supported." );
682 while (read_next_frame(oenv, status, &fr))
688 bKeepLastAppend = TRUE;
690 trxout = open_trx(out_file, "a");
694 if (gmx_fio_getftp(stfio) != efXTC)
696 gmx_fatal(FARGS, "Overwrite only supported for XTC." );
699 xdr_xtc_get_last_frame_time(gmx_fio_getfp(stfio),
700 gmx_fio_getxdr(stfio),
704 gmx_fatal(FARGS, "Error reading last frame. Maybe seek not supported." );
706 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
707 * or when seek time = 0 */
708 if (nfile_in > 1 && settime[1] < last_frame_time+timest[0]*0.5)
710 /* Jump to one time-frame before the start of next
712 searchtime = settime[1]-timest[0]*1.25;
716 searchtime = last_frame_time;
718 if (xtc_seek_time(stfio, searchtime, fr.natoms, TRUE))
720 gmx_fatal(FARGS, "Error seeking to append position.");
722 read_next_frame(oenv, status, &fr);
723 if (fabs(searchtime - fr.time) > timest[0]*0.5)
725 gmx_fatal(FARGS, "Error seeking: attempted to seek to %f but got %f.",
726 searchtime, fr.time);
729 fpos = gmx_fio_ftell(stfio);
731 trxout = open_trx(out_file, "r+");
732 if (gmx_fio_seek(trx_get_fileio(trxout), fpos))
734 gmx_fatal(FARGS, "Error seeking to append position.");
737 printf("\n Will append after %f \n", lasttime);
740 /* Lets stitch up some files */
741 timestep = timest[0];
742 for (i = n_append+1; (i < nfile_in); i++)
746 /* set the next time from the last frame in previous file */
751 if (cont_type[i] == TIME_CONTINUE)
754 begin += 0.5*timestep;
755 settime[i] = frout.time;
756 cont_type[i] = TIME_EXPLICIT;
758 else if (cont_type[i] == TIME_LAST)
761 begin += 0.5*timestep;
763 /* Or, if the time in the next part should be changed by the
764 * same amount, start at half a timestep from the last time
765 * so we dont repeat frames.
767 /* I don't understand the comment above, but for all the cases
768 * I tried the code seems to work properly. B. Hess 2008-4-2.
771 /* Or, if time is set explicitly, we check for overlap/gap */
772 if (cont_type[i] == TIME_EXPLICIT)
774 if ( ( i < nfile_in ) &&
775 ( frout.time < settime[i]-1.5*timestep ) )
777 fprintf(stderr, "WARNING: Frames around t=%f %s have a different "
778 "spacing than the rest,\n"
779 "might be a gap or overlap that couldn't be corrected "
780 "automatically.\n", output_env_conv_time(oenv, frout.time),
781 output_env_get_time_unit(oenv));
786 /* if we don't have a timestep in the current file, use the old one */
789 timestep = timest[i];
791 read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
795 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
798 if (cont_type[i] == TIME_EXPLICIT)
800 t_corr = settime[i]-fr.time;
802 /* t_corr is the amount we want to change the time.
803 * If the user has chosen not to change the time for
804 * this part of the trajectory t_corr remains at
805 * the value it had in the last part, changing this
806 * by the same amount.
807 * If no value was given for the first trajectory part
808 * we let the time start at zero, see the edit_files routine.
814 if (lasttime != NOTSET)
816 printf("lasttime %g\n", lasttime);
821 /* copy the input frame to the output frame */
823 /* set the new time by adding the correct calculated above */
824 frout.time += t_corr;
825 /* quit if we have reached the end of what should be written */
826 if ((end > 0) && (frout.time > end+GMX_REAL_EPS))
832 /* determine if we should write this frame (dt is handled elsewhere) */
833 if (bCat) /* write all frames of all files */
837 else if (bKeepLast || (bKeepLastAppend && i == 1))
838 /* write till last frame of this traj
839 and skip first frame(s) of next traj */
841 bWrite = ( frout.time > lasttime+0.5*timestep );
843 else /* write till first frame of next traj */
845 bWrite = ( frout.time < settime[i+1]-0.5*timestep );
848 if (bWrite && (frout.time >= begin) )
853 first_time = frout.time;
855 lasttime = frout.time;
856 if (dt == 0 || bRmod(frout.time, first_time, dt))
859 last_ok_t = frout.time;
862 fprintf(stderr, "\nContinue writing frames from %s t=%g %s, "
864 fnms[i], output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv),
871 write_trxframe_indexed(trxout, &frout, isize, index,
876 write_trxframe(trxout, &frout, NULL);
878 if ( ((frame % 10) == 0) || (frame < 10) )
880 fprintf(stderr, " -> frame %6d time %8.3f %s \r",
881 frame_out, output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv));
886 while (read_next_frame(oenv, status, &fr));
890 earliersteps += step;
896 fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
897 frame, output_env_conv_time(oenv, last_ok_t), output_env_get_time_unit(oenv));