a651272a226dc4d5ff68fdc86fc150074797cb69
[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,2015, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
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 #include "gmxpre.h"
38
39 #include <math.h>
40 #include <stdlib.h>
41 #include <string.h>
42
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"
63
64 #define TIME_EXPLICIT 0
65 #define TIME_CONTINUE 1
66 #define TIME_LAST     2
67 #ifndef FLT_MAX
68 #define FLT_MAX 1e36
69 #endif
70 #define FLAGS (TRX_READ_X | TRX_READ_V | TRX_READ_F)
71
72 static void scan_trj_files(char **fnms, int nfiles, real *readtime,
73                            real *timestep, atom_id imax,
74                            const output_env_t oenv)
75 {
76     /* Check start time of all files */
77     int          i, natoms = 0;
78     t_trxstatus *status;
79     real         t;
80     t_trxframe   fr;
81     gmx_bool     ok;
82
83     for (i = 0; i < nfiles; i++)
84     {
85         ok = read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
86
87         if (!ok)
88         {
89             gmx_fatal(FARGS, "\nCouldn't read frame from file." );
90         }
91         if (fr.bTime)
92         {
93             readtime[i] = fr.time;
94         }
95         else
96         {
97             readtime[i] = 0;
98             fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
99         }
100
101         if (i == 0)
102         {
103             natoms = fr.natoms;
104         }
105         else
106         {
107             if (imax == NO_ATID)
108             {
109                 if (natoms != fr.natoms)
110                 {
111                     gmx_fatal(FARGS, "\nDifferent numbers of atoms (%d/%d) in files",
112                               natoms, fr.natoms);
113                 }
114             }
115             else
116             {
117                 if (fr.natoms <= imax)
118                 {
119                     gmx_fatal(FARGS, "\nNot enough atoms (%d) for index group (%d)",
120                               fr.natoms, imax);
121                 }
122             }
123         }
124         ok = read_next_frame(oenv, status, &fr);
125         if (ok && fr.bTime)
126         {
127             timestep[i] = fr.time - readtime[i];
128         }
129         else
130         {
131             timestep[i] = 0;
132         }
133
134         close_trj(status);
135         if (fr.bX)
136         {
137             sfree(fr.x);
138         }
139         if (fr.bV)
140         {
141             sfree(fr.v);
142         }
143         if (fr.bF)
144         {
145             sfree(fr.f);
146         }
147     }
148     fprintf(stderr, "\n");
149
150 }
151
152 static void sort_files(char **fnms, real *settime, int nfile)
153 {
154     int   i, j, minidx;
155     real  timeswap;
156     char *chptr;
157
158     for (i = 0; i < nfile; i++)
159     {
160         minidx = i;
161         for (j = i + 1; j < nfile; j++)
162         {
163             if (settime[j] < settime[minidx])
164             {
165                 minidx = j;
166             }
167         }
168         if (minidx != i)
169         {
170             timeswap        = settime[i];
171             settime[i]      = settime[minidx];
172             settime[minidx] = timeswap;
173             chptr           = fnms[i];
174             fnms[i]         = fnms[minidx];
175             fnms[minidx]    = chptr;
176         }
177     }
178 }
179
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)
183 {
184     int      i;
185     gmx_bool ok;
186     char     inputstring[STRLEN], *chptr;
187
188     if (bSetTime)
189     {
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));
200
201         fprintf(
202                 stderr,
203                 "          File             Current start (%s)  New start (%s)\n"
204                 "---------------------------------------------------------\n",
205                 output_env_get_time_unit(oenv), output_env_get_time_unit(oenv));
206
207         for (i = 0; i < nfiles; i++)
208         {
209             fprintf(stderr, "%25s   %10.3f %s          ", fnms[i],
210                     output_env_conv_time(oenv, readtime[i]), output_env_get_time_unit(oenv));
211             ok = FALSE;
212             do
213             {
214                 if (NULL == fgets(inputstring, STRLEN - 1, stdin))
215                 {
216                     gmx_fatal(FARGS, "Error reading user input" );
217                 }
218
219                 inputstring[strlen(inputstring)-1] = 0;
220
221                 if (inputstring[0] == 'c' || inputstring[0] == 'C')
222                 {
223                     cont_type[i] = TIME_CONTINUE;
224                     bSort        = FALSE;
225                     ok           = TRUE;
226                     settime[i]   = FLT_MAX;
227                 }
228                 else if (inputstring[0] == 'l' ||
229                          inputstring[0] == 'L')
230                 {
231                     cont_type[i] = TIME_LAST;
232                     bSort        = FALSE;
233                     ok           = TRUE;
234                     settime[i]   = FLT_MAX;
235                 }
236                 else
237                 {
238                     settime[i] = strtod(inputstring, &chptr)*
239                         output_env_get_time_invfactor(oenv);
240                     if (chptr == inputstring)
241                     {
242                         fprintf(stderr, "'%s' not recognized as a floating point number, 'c' or 'l'. "
243                                 "Try again: ", inputstring);
244                     }
245                     else
246                     {
247                         cont_type[i] = TIME_EXPLICIT;
248                         ok           = TRUE;
249                     }
250                 }
251             }
252             while (!ok);
253         }
254         if (cont_type[0] != TIME_EXPLICIT)
255         {
256             cont_type[0] = TIME_EXPLICIT;
257             settime[0]   = 0;
258         }
259     }
260     else
261     {
262         for (i = 0; i < nfiles; i++)
263         {
264             settime[i] = readtime[i];
265         }
266     }
267     if (!bSort)
268     {
269         fprintf(stderr, "Sorting disabled.\n");
270     }
271     else
272     {
273         sort_files(fnms, settime, nfiles);
274     }
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++)
280     {
281         switch (cont_type[i])
282         {
283             case TIME_EXPLICIT:
284                 fprintf(stderr, "%25s   %10.3f %s   %10.3f %s",
285                         fnms[i],
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));
288                 if (i > 0 &&
289                     cont_type[i-1] == TIME_EXPLICIT && settime[i] == settime[i-1])
290                 {
291                     fprintf(stderr, " WARNING: same Start time as previous");
292                 }
293                 fprintf(stderr, "\n");
294                 break;
295             case TIME_CONTINUE:
296                 fprintf(stderr, "%25s        Continue from last file\n", fnms[i]);
297                 break;
298             case TIME_LAST:
299                 fprintf(stderr, "%25s        Change by same amount as last file\n",
300                         fnms[i]);
301                 break;
302         }
303     }
304     fprintf(stderr, "\n");
305
306     settime[nfiles]   = FLT_MAX;
307     cont_type[nfiles] = TIME_EXPLICIT;
308     readtime[nfiles]  = FLT_MAX;
309 }
310
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)
314 {
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;
319     t_trxframe   *trx;
320
321     snew(fp_in, nset);
322     snew(trx, nset);
323     snew(bSet, nset);
324     natoms = -1;
325     t      = -1;
326     for (i = 0; (i < nset); i++)
327     {
328         nnn = read_first_frame(oenv, &(fp_in[i]), fnms[i], &(trx[i]),
329                                TRX_NEED_X);
330         if (natoms == -1)
331         {
332             natoms     = nnn;
333             first_time = trx[i].time;
334         }
335         else if (natoms != nnn)
336         {
337             gmx_fatal(FARGS, "Trajectory file %s has %d atoms while previous trajs had %d atoms", fnms[i], nnn, natoms);
338         }
339         if (t == -1)
340         {
341             t = trx[i].time;
342         }
343         else if (t != trx[i].time)
344         {
345             gmx_fatal(FARGS, "Trajectory file %s has time %f while previous trajs had time %f", fnms[i], trx[i].time, t);
346         }
347     }
348
349     snew(fp_out, nset);
350     for (i = 0; (i < nset); i++)
351     {
352         fp_out[i] = open_trx(fnms_out[i], "w");
353     }
354     k = 0;
355     if (gmx_nint(time[k] - t) != 0)
356     {
357         gmx_fatal(FARGS, "First time in demuxing table does not match trajectories");
358     }
359     do
360     {
361         while ((k+1 < nval) && ((trx[0].time - time[k+1]) > dt_remd*0.1))
362         {
363             k++;
364         }
365         if (debug)
366         {
367             fprintf(debug, "trx[0].time = %g, time[k] = %g\n", trx[0].time, time[k]);
368         }
369         for (i = 0; (i < nset); i++)
370         {
371             bSet[i] = FALSE;
372         }
373         for (i = 0; (i < nset); i++)
374         {
375             j = gmx_nint(value[i][k]);
376             range_check(j, 0, nset);
377             if (bSet[j])
378             {
379                 gmx_fatal(FARGS, "Demuxing the same replica %d twice at time %f",
380                           j, trx[0].time);
381             }
382             bSet[j] = TRUE;
383
384             if (dt == 0 || bRmod(trx[i].time, first_time, dt))
385             {
386                 if (index)
387                 {
388                     write_trxframe_indexed(fp_out[j], &trx[i], isize, index, NULL);
389                 }
390                 else
391                 {
392                     write_trxframe(fp_out[j], &trx[i], NULL);
393                 }
394             }
395         }
396
397         bCont = (k < nval);
398         for (i = 0; (i < nset); i++)
399         {
400             bCont = bCont && read_next_frame(oenv, fp_in[i], &trx[i]);
401         }
402     }
403     while (bCont);
404
405     for (i = 0; (i < nset); i++)
406     {
407         close_trx(fp_in[i]);
408         close_trx(fp_out[i]);
409     }
410 }
411
412 int gmx_trjcat(int argc, char *argv[])
413 {
414     const char     *desc[] =
415     {
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 [REF].xvg[ref] file.",
430         "The [REF].xvg[ref] file should contain something like::",
431         "",
432         "    0  0  1  2  3  4  5",
433         "    2  1  0  2  3  5  4",
434         "",
435         "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 [REF].xvg[ref] 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_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 = gmx_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 (ftpout != ftpin)
653                 {
654                     gmx_fatal(FARGS, "When writing TNG the input file format must also be TNG");
655                 }
656                 if (bIndex)
657                 {
658                     trjtools_gmx_prepare_tng_writing(out_file, 'w', NULL, &trxout,
659                                                      fnms[0], isize, NULL, index, grpname);
660                 }
661                 else
662                 {
663                     trjtools_gmx_prepare_tng_writing(out_file, 'w', NULL, &trxout,
664                                                      fnms[0], -1, NULL, NULL, NULL);
665                 }
666             }
667             else
668             {
669                 trxout = open_trx(out_file, "w");
670             }
671             memset(&frout, 0, sizeof(frout));
672         }
673         else
674         {
675             t_fileio *stfio;
676
677             if (!read_first_frame(oenv, &status, out_file, &fr, FLAGS))
678             {
679                 gmx_fatal(FARGS, "Reading first frame from %s", out_file);
680             }
681
682             stfio = trx_get_fileio(status);
683             if (!bKeepLast && !bOverwrite)
684             {
685                 fprintf(stderr, "\n\nWARNING: Appending without -overwrite implies -keeplast "
686                         "between the first two files. \n"
687                         "If the trajectories have an overlap and have not been written binary \n"
688                         "reproducible this will produce an incorrect trajectory!\n\n");
689
690                 filetype = gmx_fio_getftp(stfio);
691                 /* Fails if last frame is incomplete
692                  * We can't do anything about it without overwriting
693                  * */
694                 if (filetype == efXTC || filetype == efTNG)
695                 {
696                     lasttime = trx_get_time_of_final_frame(status);
697                     fr.time  = lasttime;
698                 }
699                 else
700                 {
701                     while (read_next_frame(oenv, status, &fr))
702                     {
703                         ;
704                     }
705                     lasttime = fr.time;
706                 }
707                 bKeepLastAppend = TRUE;
708                 close_trj(status);
709                 trxout = open_trx(out_file, "a");
710             }
711             else if (bOverwrite)
712             {
713                 if (gmx_fio_getftp(stfio) != efXTC)
714                 {
715                     gmx_fatal(FARGS, "Overwrite only supported for XTC." );
716                 }
717                 last_frame_time = trx_get_time_of_final_frame(status);
718
719                 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
720                  *     or when seek time = 0 */
721                 if (nfile_in > 1 && settime[1] < last_frame_time+timest[0]*0.5)
722                 {
723                     /* Jump to one time-frame before the start of next
724                      *  trajectory file */
725                     searchtime = settime[1]-timest[0]*1.25;
726                 }
727                 else
728                 {
729                     searchtime = last_frame_time;
730                 }
731                 if (xtc_seek_time(stfio, searchtime, fr.natoms, TRUE))
732                 {
733                     gmx_fatal(FARGS, "Error seeking to append position.");
734                 }
735                 read_next_frame(oenv, status, &fr);
736                 if (fabs(searchtime - fr.time) > timest[0]*0.5)
737                 {
738                     gmx_fatal(FARGS, "Error seeking: attempted to seek to %f but got %f.",
739                               searchtime, fr.time);
740                 }
741                 lasttime = fr.time;
742                 fpos     = gmx_fio_ftell(stfio);
743                 close_trj(status);
744                 trxout = open_trx(out_file, "r+");
745                 if (gmx_fio_seek(trx_get_fileio(trxout), fpos))
746                 {
747                     gmx_fatal(FARGS, "Error seeking to append position.");
748                 }
749             }
750             printf("\n Will append after %f \n", lasttime);
751             frout = fr;
752         }
753         /* Lets stitch up some files */
754         timestep = timest[0];
755         for (i = n_append+1; (i < nfile_in); i++)
756         {
757             /* Open next file */
758
759             /* set the next time from the last frame in previous file */
760             if (i > 0)
761             {
762                 /* When writing TNG the step determine which frame to write. Use an
763                  * offset to be able to increase steps properly when changing files. */
764                 if (ftpout == efTNG)
765                 {
766                     prevEndStep = frout.step;
767                 }
768
769                 if (frame_out >= 0)
770                 {
771                     if (cont_type[i] == TIME_CONTINUE)
772                     {
773                         begin        = frout.time;
774                         begin       += 0.5*timestep;
775                         settime[i]   = frout.time;
776                         cont_type[i] = TIME_EXPLICIT;
777                     }
778                     else if (cont_type[i] == TIME_LAST)
779                     {
780                         begin  = frout.time;
781                         begin += 0.5*timestep;
782                     }
783                     /* Or, if the time in the next part should be changed by the
784                      * same amount, start at half a timestep from the last time
785                      * so we dont repeat frames.
786                      */
787                     /* I don't understand the comment above, but for all the cases
788                      * I tried the code seems to work properly. B. Hess 2008-4-2.
789                      */
790                 }
791                 /* Or, if time is set explicitly, we check for overlap/gap */
792                 if (cont_type[i] == TIME_EXPLICIT)
793                 {
794                     if ( ( i < nfile_in ) &&
795                          ( frout.time < settime[i]-1.5*timestep ) )
796                     {
797                         fprintf(stderr, "WARNING: Frames around t=%f %s have a different "
798                                 "spacing than the rest,\n"
799                                 "might be a gap or overlap that couldn't be corrected "
800                                 "automatically.\n", output_env_conv_time(oenv, frout.time),
801                                 output_env_get_time_unit(oenv));
802                     }
803                 }
804             }
805
806             /* if we don't have a timestep in the current file, use the old one */
807             if (timest[i] != 0)
808             {
809                 timestep = timest[i];
810             }
811             read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
812             if (!fr.bTime)
813             {
814                 fr.time = 0;
815                 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
816             }
817
818             if (cont_type[i] == TIME_EXPLICIT)
819             {
820                 t_corr = settime[i]-fr.time;
821             }
822             /* t_corr is the amount we want to change the time.
823              * If the user has chosen not to change the time for
824              * this part of the trajectory t_corr remains at
825              * the value it had in the last part, changing this
826              * by the same amount.
827              * If no value was given for the first trajectory part
828              * we let the time start at zero, see the edit_files routine.
829              */
830
831             bNewFile = TRUE;
832
833             printf("\n");
834             if (lasttime != NOTSET)
835             {
836                 printf("lasttime %g\n", lasttime);
837             }
838
839             do
840             {
841                 /* copy the input frame to the output frame */
842                 frout = fr;
843                 /* set the new time by adding the correct calculated above */
844                 frout.time += t_corr;
845                 if (ftpout == efTNG)
846                 {
847                     frout.step += prevEndStep;
848                 }
849                 /* quit if we have reached the end of what should be written */
850                 if ((end > 0) && (frout.time > end+GMX_REAL_EPS))
851                 {
852                     i = nfile_in;
853                     break;
854                 }
855
856                 /* determine if we should write this frame (dt is handled elsewhere) */
857                 if (bCat) /* write all frames of all files */
858                 {
859                     bWrite = TRUE;
860                 }
861                 else if (bKeepLast || (bKeepLastAppend && i == 1))
862                 /* write till last frame of this traj
863                    and skip first frame(s) of next traj */
864                 {
865                     bWrite = ( frout.time > lasttime+0.5*timestep );
866                 }
867                 else /* write till first frame of next traj */
868                 {
869                     bWrite = ( frout.time < settime[i+1]-0.5*timestep );
870                 }
871
872                 if (bWrite && (frout.time >= begin) )
873                 {
874                     frame++;
875                     if (frame_out == -1)
876                     {
877                         first_time = frout.time;
878                     }
879                     lasttime = frout.time;
880                     if (dt == 0 || bRmod(frout.time, first_time, dt))
881                     {
882                         frame_out++;
883                         last_ok_t = frout.time;
884                         if (bNewFile)
885                         {
886                             fprintf(stderr, "\nContinue writing frames from %s t=%g %s, "
887                                     "frame=%d      \n",
888                                     fnms[i], output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv),
889                                     frame);
890                             bNewFile = FALSE;
891                         }
892
893                         if (bIndex)
894                         {
895                             write_trxframe_indexed(trxout, &frout, isize, index,
896                                                    NULL);
897                         }
898                         else
899                         {
900                             write_trxframe(trxout, &frout, NULL);
901                         }
902                         if ( ((frame % 10) == 0) || (frame < 10) )
903                         {
904                             fprintf(stderr, " ->  frame %6d time %8.3f %s     \r",
905                                     frame_out, output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv));
906                         }
907                     }
908                 }
909             }
910             while (read_next_frame(oenv, status, &fr));
911
912             close_trj(status);
913
914             earliersteps += step;
915         }
916         if (trxout)
917         {
918             close_trx(trxout);
919         }
920         fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
921                 frame, output_env_conv_time(oenv, last_ok_t), output_env_get_time_unit(oenv));
922     }
923
924     return 0;
925 }