Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_trjcat.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <string.h>
42 #include <math.h>
43 #include "macros.h"
44 #include "sysstuff.h"
45 #include "smalloc.h"
46 #include "typedefs.h"
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"
55 #include "names.h"
56 #include "index.h"
57 #include "vec.h"
58 #include "gromacs/fileio/xtcio.h"
59 #include "do_fit.h"
60 #include "rmpbc.h"
61 #include "gromacs/fileio/g87io.h"
62 #include "pbc.h"
63 #include "xvgr.h"
64 #include "gromacs/fileio/xdrf.h"
65 #include "gmx_ana.h"
66
67 #define TIME_EXPLICIT 0
68 #define TIME_CONTINUE 1
69 #define TIME_LAST     2
70 #ifndef FLT_MAX
71 #define FLT_MAX 1e36
72 #endif
73 #define FLAGS (TRX_READ_X | TRX_READ_V | TRX_READ_F)
74
75 static void scan_trj_files(char **fnms, int nfiles, real *readtime,
76                            real *timestep, atom_id imax,
77                            const output_env_t oenv)
78 {
79     /* Check start time of all files */
80     int          i, natoms = 0;
81     t_trxstatus *status;
82     real         t;
83     t_trxframe   fr;
84     gmx_bool     ok;
85
86     for (i = 0; i < nfiles; i++)
87     {
88         ok = read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
89
90         if (!ok)
91         {
92             gmx_fatal(FARGS, "\nCouldn't read frame from file." );
93         }
94         if (fr.bTime)
95         {
96             readtime[i] = fr.time;
97         }
98         else
99         {
100             readtime[i] = 0;
101             fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
102         }
103
104         if (i == 0)
105         {
106             natoms = fr.natoms;
107         }
108         else
109         {
110             if (imax == NO_ATID)
111             {
112                 if (natoms != fr.natoms)
113                 {
114                     gmx_fatal(FARGS, "\nDifferent numbers of atoms (%d/%d) in files",
115                               natoms, fr.natoms);
116                 }
117             }
118             else
119             {
120                 if (fr.natoms <= imax)
121                 {
122                     gmx_fatal(FARGS, "\nNot enough atoms (%d) for index group (%d)",
123                               fr.natoms, imax);
124                 }
125             }
126         }
127         ok = read_next_frame(oenv, status, &fr);
128         if (ok && fr.bTime)
129         {
130             timestep[i] = fr.time - readtime[i];
131         }
132         else
133         {
134             timestep[i] = 0;
135         }
136
137         close_trj(status);
138         if (fr.bX)
139         {
140             sfree(fr.x);
141         }
142         if (fr.bV)
143         {
144             sfree(fr.v);
145         }
146         if (fr.bF)
147         {
148             sfree(fr.f);
149         }
150     }
151     fprintf(stderr, "\n");
152
153 }
154
155 static void sort_files(char **fnms, real *settime, int nfile)
156 {
157     int   i, j, minidx;
158     real  timeswap;
159     char *chptr;
160
161     for (i = 0; i < nfile; i++)
162     {
163         minidx = i;
164         for (j = i + 1; j < nfile; j++)
165         {
166             if (settime[j] < settime[minidx])
167             {
168                 minidx = j;
169             }
170         }
171         if (minidx != i)
172         {
173             timeswap        = settime[i];
174             settime[i]      = settime[minidx];
175             settime[minidx] = timeswap;
176             chptr           = fnms[i];
177             fnms[i]         = fnms[minidx];
178             fnms[minidx]    = chptr;
179         }
180     }
181 }
182
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)
186 {
187     int      i;
188     gmx_bool ok;
189     char     inputstring[STRLEN], *chptr;
190
191     if (bSetTime)
192     {
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));
203
204         fprintf(
205                 stderr,
206                 "          File             Current start (%s)  New start (%s)\n"
207                 "---------------------------------------------------------\n",
208                 output_env_get_time_unit(oenv), output_env_get_time_unit(oenv));
209
210         for (i = 0; i < nfiles; i++)
211         {
212             fprintf(stderr, "%25s   %10.3f %s          ", fnms[i],
213                     output_env_conv_time(oenv, readtime[i]), output_env_get_time_unit(oenv));
214             ok = FALSE;
215             do
216             {
217                 if (NULL == fgets(inputstring, STRLEN - 1, stdin))
218                 {
219                     gmx_fatal(FARGS, "Error reading user input" );
220                 }
221
222                 inputstring[strlen(inputstring)-1] = 0;
223
224                 if (inputstring[0] == 'c' || inputstring[0] == 'C')
225                 {
226                     cont_type[i] = TIME_CONTINUE;
227                     bSort        = FALSE;
228                     ok           = TRUE;
229                     settime[i]   = FLT_MAX;
230                 }
231                 else if (inputstring[0] == 'l' ||
232                          inputstring[0] == 'L')
233                 {
234                     cont_type[i] = TIME_LAST;
235                     bSort        = FALSE;
236                     ok           = TRUE;
237                     settime[i]   = FLT_MAX;
238                 }
239                 else
240                 {
241                     settime[i] = strtod(inputstring, &chptr)*
242                         output_env_get_time_invfactor(oenv);
243                     if (chptr == inputstring)
244                     {
245                         fprintf(stderr, "'%s' not recognized as a floating point number, 'c' or 'l'. "
246                                 "Try again: ", inputstring);
247                     }
248                     else
249                     {
250                         cont_type[i] = TIME_EXPLICIT;
251                         ok           = TRUE;
252                     }
253                 }
254             }
255             while (!ok);
256         }
257         if (cont_type[0] != TIME_EXPLICIT)
258         {
259             cont_type[0] = TIME_EXPLICIT;
260             settime[0]   = 0;
261         }
262     }
263     else
264     {
265         for (i = 0; i < nfiles; i++)
266         {
267             settime[i] = readtime[i];
268         }
269     }
270     if (!bSort)
271     {
272         fprintf(stderr, "Sorting disabled.\n");
273     }
274     else
275     {
276         sort_files(fnms, settime, nfiles);
277     }
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++)
283     {
284         switch (cont_type[i])
285         {
286             case TIME_EXPLICIT:
287                 fprintf(stderr, "%25s   %10.3f %s   %10.3f %s",
288                         fnms[i],
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));
291                 if (i > 0 &&
292                     cont_type[i-1] == TIME_EXPLICIT && settime[i] == settime[i-1])
293                 {
294                     fprintf(stderr, " WARNING: same Start time as previous");
295                 }
296                 fprintf(stderr, "\n");
297                 break;
298             case TIME_CONTINUE:
299                 fprintf(stderr, "%25s        Continue from last file\n", fnms[i]);
300                 break;
301             case TIME_LAST:
302                 fprintf(stderr, "%25s        Change by same amount as last file\n",
303                         fnms[i]);
304                 break;
305         }
306     }
307     fprintf(stderr, "\n");
308
309     settime[nfiles]   = FLT_MAX;
310     cont_type[nfiles] = TIME_EXPLICIT;
311     readtime[nfiles]  = FLT_MAX;
312 }
313
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)
317 {
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;
322     t_trxframe   *trx;
323
324     snew(fp_in, nset);
325     snew(trx, nset);
326     snew(bSet, nset);
327     natoms = -1;
328     t      = -1;
329     for (i = 0; (i < nset); i++)
330     {
331         nnn = read_first_frame(oenv, &(fp_in[i]), fnms[i], &(trx[i]),
332                                TRX_NEED_X);
333         if (natoms == -1)
334         {
335             natoms     = nnn;
336             first_time = trx[i].time;
337         }
338         else if (natoms != nnn)
339         {
340             gmx_fatal(FARGS, "Trajectory file %s has %d atoms while previous trajs had %d atoms", fnms[i], nnn, natoms);
341         }
342         if (t == -1)
343         {
344             t = trx[i].time;
345         }
346         else if (t != trx[i].time)
347         {
348             gmx_fatal(FARGS, "Trajectory file %s has time %f while previous trajs had time %f", fnms[i], trx[i].time, t);
349         }
350     }
351
352     snew(fp_out, nset);
353     for (i = 0; (i < nset); i++)
354     {
355         fp_out[i] = open_trx(fnms_out[i], "w");
356     }
357     k = 0;
358     if (gmx_nint(time[k] - t) != 0)
359     {
360         gmx_fatal(FARGS, "First time in demuxing table does not match trajectories");
361     }
362     do
363     {
364         while ((k+1 < nval) && ((trx[0].time - time[k+1]) > dt_remd*0.1))
365         {
366             k++;
367         }
368         if (debug)
369         {
370             fprintf(debug, "trx[0].time = %g, time[k] = %g\n", trx[0].time, time[k]);
371         }
372         for (i = 0; (i < nset); i++)
373         {
374             bSet[i] = FALSE;
375         }
376         for (i = 0; (i < nset); i++)
377         {
378             j = gmx_nint(value[i][k]);
379             range_check(j, 0, nset);
380             if (bSet[j])
381             {
382                 gmx_fatal(FARGS, "Demuxing the same replica %d twice at time %f",
383                           j, trx[0].time);
384             }
385             bSet[j] = TRUE;
386
387             if (dt == 0 || bRmod(trx[i].time, first_time, dt))
388             {
389                 if (index)
390                 {
391                     write_trxframe_indexed(fp_out[j], &trx[i], isize, index, NULL);
392                 }
393                 else
394                 {
395                     write_trxframe(fp_out[j], &trx[i], NULL);
396                 }
397             }
398         }
399
400         bCont = (k < nval);
401         for (i = 0; (i < nset); i++)
402         {
403             bCont = bCont && read_next_frame(oenv, fp_in[i], &trx[i]);
404         }
405     }
406     while (bCont);
407
408     for (i = 0; (i < nset); i++)
409     {
410         close_trx(fp_in[i]);
411         close_trx(fp_out[i]);
412     }
413 }
414
415 int gmx_trjcat(int argc, char *argv[])
416 {
417     const char     *desc[] =
418     {
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."
442     };
443     static gmx_bool bVels           = TRUE;
444     static int      prec            = 3;
445     static gmx_bool bCat            = FALSE;
446     static gmx_bool bSort           = TRUE;
447     static gmx_bool bKeepLast       = FALSE;
448     static gmx_bool bKeepLastAppend = FALSE;
449     static gmx_bool bOverwrite      = FALSE;
450     static gmx_bool bSetTime        = FALSE;
451     static gmx_bool bDeMux;
452     static real     begin = -1;
453     static real     end   = -1;
454     static real     dt    = 0;
455
456     t_pargs
457         pa[] =
458     {
459         { "-b", FALSE, etTIME,
460           { &begin }, "First time to use (%t)" },
461         { "-e", FALSE, etTIME,
462           { &end }, "Last time to use (%t)" },
463         { "-dt", FALSE, etTIME,
464           { &dt }, "Only write frame when t MOD dt = first time (%t)" },
465         { "-prec", FALSE, etINT,
466           { &prec }, "Precision for [TT].xtc[tt] and [TT].gro[tt] writing in number of decimal places" },
467         { "-vel", FALSE, etBOOL,
468           { &bVels }, "Read and write velocities if possible" },
469         { "-settime", FALSE, etBOOL,
470           { &bSetTime }, "Change starting time interactively" },
471         { "-sort", FALSE, etBOOL,
472           { &bSort }, "Sort trajectory files (not frames)" },
473         { "-keeplast", FALSE, etBOOL,
474           { &bKeepLast }, "Keep overlapping frames at end of trajectory" },
475         { "-overwrite", FALSE, etBOOL,
476           { &bOverwrite }, "Overwrite overlapping frames during appending" },
477         { "-cat", FALSE, etBOOL,
478           { &bCat }, "Do not discard double time frames" }
479     };
480 #define npargs asize(pa)
481     int          ftpin, i, frame, frame_out, step = 0, trjout = 0;
482     t_trxstatus *status;
483     rvec        *x, *v;
484     real         xtcpr, t_corr;
485     t_trxframe   fr, frout;
486     char       **fnms, **fnms_out, *in_file, *out_file;
487     int          n_append;
488     t_trxstatus *trxout = NULL;
489     gmx_bool     bNewFile, bIndex, bWrite;
490     int          earliersteps, nfile_in, nfile_out, *cont_type, last_ok_step;
491     real        *readtime, *timest, *settime;
492     real         first_time = 0, lasttime = NOTSET, last_ok_t = -1, timestep;
493     real         last_frame_time, searchtime;
494     int          isize, j;
495     atom_id     *index = NULL, imax;
496     char        *grpname;
497     real       **val = NULL, *t = NULL, dt_remd;
498     int          n, nset;
499     gmx_bool     bOK;
500     gmx_off_t    fpos;
501     output_env_t oenv;
502     t_filenm     fnm[] =
503     {
504         { efTRX, "-f", NULL, ffRDMULT },
505         { efTRO, "-o", NULL, ffWRMULT },
506         { efNDX, "-n", "index", ffOPTRD },
507         { efXVG, "-demux", "remd", ffOPTRD }
508     };
509
510 #define NFILE asize(fnm)
511
512     if (!parse_common_args(&argc, argv, PCA_BE_NICE | PCA_TIME_UNIT, NFILE, fnm,
513                            asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
514     {
515         return 0;
516     }
517
518     bIndex = ftp2bSet(efNDX, NFILE, fnm);
519     bDeMux = ftp2bSet(efXVG, NFILE, fnm);
520     bSort  = bSort && !bDeMux;
521
522     imax = NO_ATID;
523     if (bIndex)
524     {
525         printf("Select group for output\n");
526         rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
527         /* scan index */
528         imax = index[0];
529         for (i = 1; i < isize; i++)
530         {
531             imax = max(imax, index[i]);
532         }
533     }
534     if (bDeMux)
535     {
536         nset    = 0;
537         dt_remd = 0;
538         val     = read_xvg_time(opt2fn("-demux", NFILE, fnm), TRUE,
539                                 opt2parg_bSet("-b", npargs, pa), begin,
540                                 opt2parg_bSet("-e", npargs, pa), end, 1, &nset, &n,
541                                 &dt_remd, &t);
542         printf("Read %d sets of %d points, dt = %g\n\n", nset, n, dt_remd);
543         if (debug)
544         {
545             fprintf(debug, "Dump of replica_index.xvg\n");
546             for (i = 0; (i < n); i++)
547             {
548                 fprintf(debug, "%10g", t[i]);
549                 for (j = 0; (j < nset); j++)
550                 {
551                     fprintf(debug, "  %3d", gmx_nint(val[j][i]));
552                 }
553                 fprintf(debug, "\n");
554             }
555         }
556     }
557     /* prec is in nr of decimal places, xtcprec is a multiplication factor: */
558     xtcpr = 1;
559     for (i = 0; i < prec; i++)
560     {
561         xtcpr *= 10;
562     }
563
564     nfile_in = opt2fns(&fnms, "-f", NFILE, fnm);
565     if (!nfile_in)
566     {
567         gmx_fatal(FARGS, "No input files!" );
568     }
569
570     if (bDeMux && (nfile_in != nset))
571     {
572         gmx_fatal(FARGS, "You have specified %d files and %d entries in the demux table", nfile_in, nset);
573     }
574
575     nfile_out = opt2fns(&fnms_out, "-o", NFILE, fnm);
576     if (!nfile_out)
577     {
578         gmx_fatal(FARGS, "No output files!");
579     }
580     if ((nfile_out > 1) && !bDeMux)
581     {
582         gmx_fatal(FARGS, "Don't know what to do with more than 1 output file if  not demultiplexing");
583     }
584     else if (bDeMux && (nfile_out != nset) && (nfile_out != 1))
585     {
586         gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %d", nset, nfile_out);
587     }
588     if (bDeMux)
589     {
590         if (nfile_out != nset)
591         {
592             char *buf = strdup(fnms_out[0]);
593             snew(fnms_out, nset);
594             for (i = 0; (i < nset); i++)
595             {
596                 snew(fnms_out[i], strlen(buf)+32);
597                 sprintf(fnms_out[i], "%d_%s", i, buf);
598             }
599             sfree(buf);
600         }
601         do_demux(nfile_in, fnms, fnms_out, n, val, t, dt_remd, isize, index, dt, oenv);
602     }
603     else
604     {
605         snew(readtime, nfile_in+1);
606         snew(timest, nfile_in+1);
607         scan_trj_files(fnms, nfile_in, readtime, timest, imax, oenv);
608
609         snew(settime, nfile_in+1);
610         snew(cont_type, nfile_in+1);
611         edit_files(fnms, nfile_in, readtime, timest, settime, cont_type, bSetTime, bSort,
612                    oenv);
613
614         /* Check whether the output file is amongst the input files
615          * This has to be done after sorting etc.
616          */
617         out_file = fnms_out[0];
618         n_append = -1;
619         for (i = 0; ((i < nfile_in) && (n_append == -1)); i++)
620         {
621             if (strcmp(fnms[i], out_file) == 0)
622             {
623                 n_append = i;
624             }
625         }
626         if (n_append == 0)
627         {
628             fprintf(stderr, "Will append to %s rather than creating a new file\n",
629                     out_file);
630         }
631         else if (n_append != -1)
632         {
633             gmx_fatal(FARGS, "Can only append to the first file which is %s (not %s)",
634                       fnms[0], out_file);
635         }
636         earliersteps = 0;
637
638         /* Not checking input format, could be dangerous :-) */
639         /* Not checking output format, equally dangerous :-) */
640
641         frame     = -1;
642         frame_out = -1;
643         /* the default is not to change the time at all,
644          * but this is overridden by the edit_files routine
645          */
646         t_corr = 0;
647
648         if (n_append == -1)
649         {
650             trxout = open_trx(out_file, "w");
651             memset(&frout, 0, sizeof(frout));
652         }
653         else
654         {
655             t_fileio *stfio;
656
657             if (!read_first_frame(oenv, &status, out_file, &fr, FLAGS))
658             {
659                 gmx_fatal(FARGS, "Reading first frame from %s", out_file);
660             }
661
662             stfio = trx_get_fileio(status);
663             if (!bKeepLast && !bOverwrite)
664             {
665                 fprintf(stderr, "\n\nWARNING: Appending without -overwrite implies -keeplast "
666                         "between the first two files. \n"
667                         "If the trajectories have an overlap and have not been written binary \n"
668                         "reproducible this will produce an incorrect trajectory!\n\n");
669
670                 /* Fails if last frame is incomplete
671                  * We can't do anything about it without overwriting
672                  * */
673                 if (gmx_fio_getftp(stfio) == efXTC)
674                 {
675                     lasttime =
676                         xdr_xtc_get_last_frame_time(gmx_fio_getfp(stfio),
677                                                     gmx_fio_getxdr(stfio),
678                                                     fr.natoms, &bOK);
679                     fr.time = lasttime;
680                     if (!bOK)
681                     {
682                         gmx_fatal(FARGS, "Error reading last frame. Maybe seek not supported." );
683                     }
684                 }
685                 else
686                 {
687                     while (read_next_frame(oenv, status, &fr))
688                     {
689                         ;
690                     }
691                     lasttime = fr.time;
692                 }
693                 bKeepLastAppend = TRUE;
694                 close_trj(status);
695                 trxout = open_trx(out_file, "a");
696             }
697             else if (bOverwrite)
698             {
699                 if (gmx_fio_getftp(stfio) != efXTC)
700                 {
701                     gmx_fatal(FARGS, "Overwrite only supported for XTC." );
702                 }
703                 last_frame_time =
704                     xdr_xtc_get_last_frame_time(gmx_fio_getfp(stfio),
705                                                 gmx_fio_getxdr(stfio),
706                                                 fr.natoms, &bOK);
707                 if (!bOK)
708                 {
709                     gmx_fatal(FARGS, "Error reading last frame. Maybe seek not supported." );
710                 }
711                 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
712                  *     or when seek time = 0 */
713                 if (nfile_in > 1 && settime[1] < last_frame_time+timest[0]*0.5)
714                 {
715                     /* Jump to one time-frame before the start of next
716                      *  trajectory file */
717                     searchtime = settime[1]-timest[0]*1.25;
718                 }
719                 else
720                 {
721                     searchtime = last_frame_time;
722                 }
723                 if (xtc_seek_time(stfio, searchtime, fr.natoms, TRUE))
724                 {
725                     gmx_fatal(FARGS, "Error seeking to append position.");
726                 }
727                 read_next_frame(oenv, status, &fr);
728                 if (fabs(searchtime - fr.time) > timest[0]*0.5)
729                 {
730                     gmx_fatal(FARGS, "Error seeking: attempted to seek to %f but got %f.",
731                               searchtime, fr.time);
732                 }
733                 lasttime = fr.time;
734                 fpos     = gmx_fio_ftell(stfio);
735                 close_trj(status);
736                 trxout = open_trx(out_file, "r+");
737                 if (gmx_fio_seek(trx_get_fileio(trxout), fpos))
738                 {
739                     gmx_fatal(FARGS, "Error seeking to append position.");
740                 }
741             }
742             printf("\n Will append after %f \n", lasttime);
743             frout = fr;
744         }
745         /* Lets stitch up some files */
746         timestep = timest[0];
747         for (i = n_append+1; (i < nfile_in); i++)
748         {
749             /* Open next file */
750
751             /* set the next time from the last frame in previous file */
752             if (i > 0)
753             {
754                 if (frame_out >= 0)
755                 {
756                     if (cont_type[i] == TIME_CONTINUE)
757                     {
758                         begin        = frout.time;
759                         begin       += 0.5*timestep;
760                         settime[i]   = frout.time;
761                         cont_type[i] = TIME_EXPLICIT;
762                     }
763                     else if (cont_type[i] == TIME_LAST)
764                     {
765                         begin  = frout.time;
766                         begin += 0.5*timestep;
767                     }
768                     /* Or, if the time in the next part should be changed by the
769                      * same amount, start at half a timestep from the last time
770                      * so we dont repeat frames.
771                      */
772                     /* I don't understand the comment above, but for all the cases
773                      * I tried the code seems to work properly. B. Hess 2008-4-2.
774                      */
775                 }
776                 /* Or, if time is set explicitly, we check for overlap/gap */
777                 if (cont_type[i] == TIME_EXPLICIT)
778                 {
779                     if ( ( i < nfile_in ) &&
780                          ( frout.time < settime[i]-1.5*timestep ) )
781                     {
782                         fprintf(stderr, "WARNING: Frames around t=%f %s have a different "
783                                 "spacing than the rest,\n"
784                                 "might be a gap or overlap that couldn't be corrected "
785                                 "automatically.\n", output_env_conv_time(oenv, frout.time),
786                                 output_env_get_time_unit(oenv));
787                     }
788                 }
789             }
790
791             /* if we don't have a timestep in the current file, use the old one */
792             if (timest[i] != 0)
793             {
794                 timestep = timest[i];
795             }
796             read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
797             if (!fr.bTime)
798             {
799                 fr.time = 0;
800                 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
801             }
802
803             if (cont_type[i] == TIME_EXPLICIT)
804             {
805                 t_corr = settime[i]-fr.time;
806             }
807             /* t_corr is the amount we want to change the time.
808              * If the user has chosen not to change the time for
809              * this part of the trajectory t_corr remains at
810              * the value it had in the last part, changing this
811              * by the same amount.
812              * If no value was given for the first trajectory part
813              * we let the time start at zero, see the edit_files routine.
814              */
815
816             bNewFile = TRUE;
817
818             printf("\n");
819             if (lasttime != NOTSET)
820             {
821                 printf("lasttime %g\n", lasttime);
822             }
823
824             do
825             {
826                 /* copy the input frame to the output frame */
827                 frout = fr;
828                 /* set the new time by adding the correct calculated above */
829                 frout.time += t_corr;
830                 /* quit if we have reached the end of what should be written */
831                 if ((end > 0) && (frout.time > end+GMX_REAL_EPS))
832                 {
833                     i = nfile_in;
834                     break;
835                 }
836
837                 /* determine if we should write this frame (dt is handled elsewhere) */
838                 if (bCat) /* write all frames of all files */
839                 {
840                     bWrite = TRUE;
841                 }
842                 else if (bKeepLast || (bKeepLastAppend && i == 1))
843                 /* write till last frame of this traj
844                    and skip first frame(s) of next traj */
845                 {
846                     bWrite = ( frout.time > lasttime+0.5*timestep );
847                 }
848                 else /* write till first frame of next traj */
849                 {
850                     bWrite = ( frout.time < settime[i+1]-0.5*timestep );
851                 }
852
853                 if (bWrite && (frout.time >= begin) )
854                 {
855                     frame++;
856                     if (frame_out == -1)
857                     {
858                         first_time = frout.time;
859                     }
860                     lasttime = frout.time;
861                     if (dt == 0 || bRmod(frout.time, first_time, dt))
862                     {
863                         frame_out++;
864                         last_ok_t = frout.time;
865                         if (bNewFile)
866                         {
867                             fprintf(stderr, "\nContinue writing frames from %s t=%g %s, "
868                                     "frame=%d      \n",
869                                     fnms[i], output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv),
870                                     frame);
871                             bNewFile = FALSE;
872                         }
873
874                         if (bIndex)
875                         {
876                             write_trxframe_indexed(trxout, &frout, isize, index,
877                                                    NULL);
878                         }
879                         else
880                         {
881                             write_trxframe(trxout, &frout, NULL);
882                         }
883                         if ( ((frame % 10) == 0) || (frame < 10) )
884                         {
885                             fprintf(stderr, " ->  frame %6d time %8.3f %s     \r",
886                                     frame_out, output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv));
887                         }
888                     }
889                 }
890             }
891             while (read_next_frame(oenv, status, &fr));
892
893             close_trj(status);
894
895             earliersteps += step;
896         }
897         if (trxout)
898         {
899             close_trx(trxout);
900         }
901         fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
902                 frame, output_env_conv_time(oenv, last_ok_t), output_env_get_time_unit(oenv));
903     }
904
905     return 0;
906 }