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