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