Moved second half of gmxana tools to C++
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_trjcat.cpp
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 <cmath>
40 #include <cstdlib>
41 #include <cstring>
42
43 #include <algorithm>
44
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/fileio/gmxfio.h"
48 #include "gromacs/fileio/pdbio.h"
49 #include "gromacs/fileio/tngio.h"
50 #include "gromacs/fileio/tngio_for_tools.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     t_trxframe   fr;
80     gmx_bool     ok;
81
82     for (i = 0; i < nfiles; i++)
83     {
84         ok = read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
85
86         if (!ok)
87         {
88             gmx_fatal(FARGS, "\nCouldn't read frame from file." );
89         }
90         if (fr.bTime)
91         {
92             readtime[i] = fr.time;
93         }
94         else
95         {
96             readtime[i] = 0;
97             fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
98         }
99
100         if (i == 0)
101         {
102             natoms = fr.natoms;
103         }
104         else
105         {
106             if (imax == NO_ATID)
107             {
108                 if (natoms != fr.natoms)
109                 {
110                     gmx_fatal(FARGS, "\nDifferent numbers of atoms (%d/%d) in files",
111                               natoms, fr.natoms);
112                 }
113             }
114             else
115             {
116                 if (fr.natoms <= imax)
117                 {
118                     gmx_fatal(FARGS, "\nNot enough atoms (%d) for index group (%d)",
119                               fr.natoms, imax);
120                 }
121             }
122         }
123         ok = read_next_frame(oenv, status, &fr);
124         if (ok && fr.bTime)
125         {
126             timestep[i] = fr.time - readtime[i];
127         }
128         else
129         {
130             timestep[i] = 0;
131         }
132
133         close_trj(status);
134         if (fr.bX)
135         {
136             sfree(fr.x);
137         }
138         if (fr.bV)
139         {
140             sfree(fr.v);
141         }
142         if (fr.bF)
143         {
144             sfree(fr.f);
145         }
146     }
147     fprintf(stderr, "\n");
148
149 }
150
151 static void sort_files(char **fnms, real *settime, int nfile)
152 {
153     int   i, j, minidx;
154     real  timeswap;
155     char *chptr;
156
157     for (i = 0; i < nfile; i++)
158     {
159         minidx = i;
160         for (j = i + 1; j < nfile; j++)
161         {
162             if (settime[j] < settime[minidx])
163             {
164                 minidx = j;
165             }
166         }
167         if (minidx != i)
168         {
169             timeswap        = settime[i];
170             settime[i]      = settime[minidx];
171             settime[minidx] = timeswap;
172             chptr           = fnms[i];
173             fnms[i]         = fnms[minidx];
174             fnms[minidx]    = chptr;
175         }
176     }
177 }
178
179 static void edit_files(char **fnms, int nfiles, real *readtime, real *timestep,
180                        real *settime, int *cont_type, gmx_bool bSetTime,
181                        gmx_bool bSort, const output_env_t oenv)
182 {
183     int      i;
184     gmx_bool ok;
185     char     inputstring[STRLEN], *chptr;
186
187     if (bSetTime)
188     {
189         fprintf(stderr, "\n\nEnter the new start time (%s) for each file.\n"
190                 "There are two special options, both disable sorting:\n\n"
191                 "c (continue) - The start time is taken from the end\n"
192                 "of the previous file. Use it when your continuation run\n"
193                 "restarts with t=0.\n\n"
194                 "l (last) - The time in this file will be changed the\n"
195                 "same amount as in the previous. Use it when the time in the\n"
196                 "new run continues from the end of the previous one,\n"
197                 "since this takes possible overlap into account.\n\n",
198                 output_env_get_time_unit(oenv));
199
200         fprintf(
201                 stderr,
202                 "          File             Current start (%s)  New start (%s)\n"
203                 "---------------------------------------------------------\n",
204                 output_env_get_time_unit(oenv), output_env_get_time_unit(oenv));
205
206         for (i = 0; i < nfiles; i++)
207         {
208             fprintf(stderr, "%25s   %10.3f %s          ", fnms[i],
209                     output_env_conv_time(oenv, readtime[i]), output_env_get_time_unit(oenv));
210             ok = FALSE;
211             do
212             {
213                 if (NULL == fgets(inputstring, STRLEN - 1, stdin))
214                 {
215                     gmx_fatal(FARGS, "Error reading user input" );
216                 }
217
218                 inputstring[std::strlen(inputstring)-1] = 0;
219
220                 if (inputstring[0] == 'c' || inputstring[0] == 'C')
221                 {
222                     cont_type[i] = TIME_CONTINUE;
223                     bSort        = FALSE;
224                     ok           = TRUE;
225                     settime[i]   = FLT_MAX;
226                 }
227                 else if (inputstring[0] == 'l' ||
228                          inputstring[0] == 'L')
229                 {
230                     cont_type[i] = TIME_LAST;
231                     bSort        = FALSE;
232                     ok           = TRUE;
233                     settime[i]   = FLT_MAX;
234                 }
235                 else
236                 {
237                     settime[i] = strtod(inputstring, &chptr)*
238                         output_env_get_time_invfactor(oenv);
239                     if (chptr == inputstring)
240                     {
241                         fprintf(stderr, "'%s' not recognized as a floating point number, 'c' or 'l'. "
242                                 "Try again: ", inputstring);
243                     }
244                     else
245                     {
246                         cont_type[i] = TIME_EXPLICIT;
247                         ok           = TRUE;
248                     }
249                 }
250             }
251             while (!ok);
252         }
253         if (cont_type[0] != TIME_EXPLICIT)
254         {
255             cont_type[0] = TIME_EXPLICIT;
256             settime[0]   = 0;
257         }
258     }
259     else
260     {
261         for (i = 0; i < nfiles; i++)
262         {
263             settime[i] = readtime[i];
264         }
265     }
266     if (!bSort)
267     {
268         fprintf(stderr, "Sorting disabled.\n");
269     }
270     else
271     {
272         sort_files(fnms, settime, nfiles);
273     }
274     /* Write out the new order and start times */
275     fprintf(stderr, "\nSummary of files and start times used:\n\n"
276             "          File                Start time       Time step\n"
277             "---------------------------------------------------------\n");
278     for (i = 0; i < nfiles; i++)
279     {
280         switch (cont_type[i])
281         {
282             case TIME_EXPLICIT:
283                 fprintf(stderr, "%25s   %10.3f %s   %10.3f %s",
284                         fnms[i],
285                         output_env_conv_time(oenv, settime[i]), output_env_get_time_unit(oenv),
286                         output_env_conv_time(oenv, timestep[i]), output_env_get_time_unit(oenv));
287                 if (i > 0 &&
288                     cont_type[i-1] == TIME_EXPLICIT && settime[i] == settime[i-1])
289                 {
290                     fprintf(stderr, " WARNING: same Start time as previous");
291                 }
292                 fprintf(stderr, "\n");
293                 break;
294             case TIME_CONTINUE:
295                 fprintf(stderr, "%25s        Continue from last file\n", fnms[i]);
296                 break;
297             case TIME_LAST:
298                 fprintf(stderr, "%25s        Change by same amount as last file\n",
299                         fnms[i]);
300                 break;
301         }
302     }
303     fprintf(stderr, "\n");
304
305     settime[nfiles]   = FLT_MAX;
306     cont_type[nfiles] = TIME_EXPLICIT;
307     readtime[nfiles]  = FLT_MAX;
308 }
309
310 static void do_demux(int nset, char *fnms[], char *fnms_out[], int nval,
311                      real **value, real *time, real dt_remd, int isize,
312                      atom_id index[], real dt, const output_env_t oenv)
313 {
314     int           i, j, k, natoms, nnn;
315     t_trxstatus **fp_in, **fp_out;
316     gmx_bool      bCont, *bSet;
317     real          t, first_time = 0;
318     t_trxframe   *trx;
319
320     snew(fp_in, nset);
321     snew(trx, nset);
322     snew(bSet, nset);
323     natoms = -1;
324     t      = -1;
325     for (i = 0; (i < nset); i++)
326     {
327         nnn = read_first_frame(oenv, &(fp_in[i]), fnms[i], &(trx[i]),
328                                TRX_NEED_X);
329         if (natoms == -1)
330         {
331             natoms     = nnn;
332             first_time = trx[i].time;
333         }
334         else if (natoms != nnn)
335         {
336             gmx_fatal(FARGS, "Trajectory file %s has %d atoms while previous trajs had %d atoms", fnms[i], nnn, natoms);
337         }
338         if (t == -1)
339         {
340             t = trx[i].time;
341         }
342         else if (t != trx[i].time)
343         {
344             gmx_fatal(FARGS, "Trajectory file %s has time %f while previous trajs had time %f", fnms[i], trx[i].time, t);
345         }
346     }
347
348     snew(fp_out, nset);
349     for (i = 0; (i < nset); i++)
350     {
351         fp_out[i] = open_trx(fnms_out[i], "w");
352     }
353     k = 0;
354     if (gmx_nint(time[k] - t) != 0)
355     {
356         gmx_fatal(FARGS, "First time in demuxing table does not match trajectories");
357     }
358     do
359     {
360         while ((k+1 < nval) && ((trx[0].time - time[k+1]) > dt_remd*0.1))
361         {
362             k++;
363         }
364         if (debug)
365         {
366             fprintf(debug, "trx[0].time = %g, time[k] = %g\n", trx[0].time, time[k]);
367         }
368         for (i = 0; (i < nset); i++)
369         {
370             bSet[i] = FALSE;
371         }
372         for (i = 0; (i < nset); i++)
373         {
374             j = gmx_nint(value[i][k]);
375             range_check(j, 0, nset);
376             if (bSet[j])
377             {
378                 gmx_fatal(FARGS, "Demuxing the same replica %d twice at time %f",
379                           j, trx[0].time);
380             }
381             bSet[j] = TRUE;
382
383             if (dt == 0 || bRmod(trx[i].time, first_time, dt))
384             {
385                 if (index)
386                 {
387                     write_trxframe_indexed(fp_out[j], &trx[i], isize, index, NULL);
388                 }
389                 else
390                 {
391                     write_trxframe(fp_out[j], &trx[i], NULL);
392                 }
393             }
394         }
395
396         bCont = (k < nval);
397         for (i = 0; (i < nset); i++)
398         {
399             bCont = bCont && read_next_frame(oenv, fp_in[i], &trx[i]);
400         }
401     }
402     while (bCont);
403
404     for (i = 0; (i < nset); i++)
405     {
406         close_trx(fp_in[i]);
407         close_trx(fp_out[i]);
408     }
409 }
410
411 int gmx_trjcat(int argc, char *argv[])
412 {
413     const char     *desc[] =
414     {
415         "[THISMODULE] concatenates several input trajectory files in sorted order. ",
416         "In case of double time frames the one in the later file is used. ",
417         "By specifying [TT]-settime[tt] you will be asked for the start time ",
418         "of each file. The input files are taken from the command line, ",
419         "such that a command like [TT]gmx trjcat -f *.trr -o fixed.trr[tt] should do ",
420         "the trick. Using [TT]-cat[tt], you can simply paste several files ",
421         "together without removal of frames with identical time stamps.[PAR]",
422         "One important option is inferred when the output file is amongst the",
423         "input files. In that case that particular file will be appended to",
424         "which implies you do not need to store double the amount of data.",
425         "Obviously the file to append to has to be the one with lowest starting",
426         "time since one can only append at the end of a file.[PAR]",
427         "If the [TT]-demux[tt] option is given, the N trajectories that are",
428         "read, are written in another order as specified in the [REF].xvg[ref] file.",
429         "The [REF].xvg[ref] file should contain something like::",
430         "",
431         "    0  0  1  2  3  4  5",
432         "    2  1  0  2  3  5  4",
433         "",
434         "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 [REF].xvg[ref] 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;
477     t_trxstatus *status, *trxout = NULL;
478     real         t_corr;
479     t_trxframe   fr, frout;
480     char       **fnms, **fnms_out, *out_file;
481     int          n_append;
482     gmx_bool     bNewFile, bIndex, bWrite;
483     int          nfile_in, nfile_out, *cont_type;
484     real        *readtime, *timest, *settime;
485     real         first_time = 0, lasttime = NOTSET, last_ok_t = -1, timestep;
486     real         last_frame_time, searchtime;
487     int          isize = 0, j;
488     atom_id     *index = NULL, imax;
489     char        *grpname;
490     real       **val = NULL, *t = NULL, dt_remd;
491     int          n, nset, ftpout = -1, prevEndStep = 0, filetype;
492     gmx_off_t    fpos;
493     output_env_t oenv;
494     t_filenm     fnm[] =
495     {
496         { efTRX, "-f", NULL, ffRDMULT },
497         { efTRO, "-o", NULL, ffWRMULT },
498         { efNDX, "-n", "index", ffOPTRD },
499         { efXVG, "-demux", "remd", ffOPTRD }
500     };
501
502 #define NFILE asize(fnm)
503
504     if (!parse_common_args(&argc, argv, PCA_TIME_UNIT, NFILE, fnm,
505                            asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
506     {
507         return 0;
508     }
509
510     bIndex = ftp2bSet(efNDX, NFILE, fnm);
511     bDeMux = ftp2bSet(efXVG, NFILE, fnm);
512     bSort  = bSort && !bDeMux;
513
514     imax = NO_ATID;
515     if (bIndex)
516     {
517         printf("Select group for output\n");
518         rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
519         /* scan index */
520         imax = index[0];
521         for (i = 1; i < isize; i++)
522         {
523             imax = std::max(imax, index[i]);
524         }
525     }
526     if (bDeMux)
527     {
528         nset    = 0;
529         dt_remd = 0;
530         val     = read_xvg_time(opt2fn("-demux", NFILE, fnm), TRUE,
531                                 opt2parg_bSet("-b", npargs, pa), begin,
532                                 opt2parg_bSet("-e", npargs, pa), end, 1, &nset, &n,
533                                 &dt_remd, &t);
534         printf("Read %d sets of %d points, dt = %g\n\n", nset, n, dt_remd);
535         if (debug)
536         {
537             fprintf(debug, "Dump of replica_index.xvg\n");
538             for (i = 0; (i < n); i++)
539             {
540                 fprintf(debug, "%10g", t[i]);
541                 for (j = 0; (j < nset); j++)
542                 {
543                     fprintf(debug, "  %3d", gmx_nint(val[j][i]));
544                 }
545                 fprintf(debug, "\n");
546             }
547         }
548     }
549
550     nfile_in = opt2fns(&fnms, "-f", NFILE, fnm);
551     if (!nfile_in)
552     {
553         gmx_fatal(FARGS, "No input files!" );
554     }
555
556     if (bDeMux && (nfile_in != nset))
557     {
558         gmx_fatal(FARGS, "You have specified %d files and %d entries in the demux table", nfile_in, nset);
559     }
560
561     ftpin = fn2ftp(fnms[0]);
562
563     for (i = 1; i < nfile_in; i++)
564     {
565         if (ftpin != fn2ftp(fnms[i]))
566         {
567             gmx_fatal(FARGS, "All input files must be of the same format");
568         }
569     }
570
571     nfile_out = opt2fns(&fnms_out, "-o", NFILE, fnm);
572     if (!nfile_out)
573     {
574         gmx_fatal(FARGS, "No output files!");
575     }
576     if ((nfile_out > 1) && !bDeMux)
577     {
578         gmx_fatal(FARGS, "Don't know what to do with more than 1 output file if  not demultiplexing");
579     }
580     else if (bDeMux && (nfile_out != nset) && (nfile_out != 1))
581     {
582         gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %d", nset, nfile_out);
583     }
584     if (bDeMux)
585     {
586         if (nfile_out != nset)
587         {
588             char *buf = gmx_strdup(fnms_out[0]);
589             snew(fnms_out, nset);
590             for (i = 0; (i < nset); i++)
591             {
592                 snew(fnms_out[i], std::strlen(buf)+32);
593                 sprintf(fnms_out[i], "%d_%s", i, buf);
594             }
595             sfree(buf);
596         }
597         do_demux(nfile_in, fnms, fnms_out, n, val, t, dt_remd, isize, index, dt, oenv);
598     }
599     else
600     {
601         snew(readtime, nfile_in+1);
602         snew(timest, nfile_in+1);
603         scan_trj_files(fnms, nfile_in, readtime, timest, imax, oenv);
604
605         snew(settime, nfile_in+1);
606         snew(cont_type, nfile_in+1);
607         edit_files(fnms, nfile_in, readtime, timest, settime, cont_type, bSetTime, bSort,
608                    oenv);
609
610         /* Check whether the output file is amongst the input files
611          * This has to be done after sorting etc.
612          */
613         out_file = fnms_out[0];
614         ftpout   = fn2ftp(out_file);
615         n_append = -1;
616         for (i = 0; ((i < nfile_in) && (n_append == -1)); i++)
617         {
618             if (std::strcmp(fnms[i], out_file) == 0)
619             {
620                 n_append = i;
621             }
622         }
623         if (n_append == 0)
624         {
625             fprintf(stderr, "Will append to %s rather than creating a new file\n",
626                     out_file);
627         }
628         else if (n_append != -1)
629         {
630             gmx_fatal(FARGS, "Can only append to the first file which is %s (not %s)",
631                       fnms[0], out_file);
632         }
633
634         /* Not checking input format, could be dangerous :-) */
635         /* Not checking output format, equally dangerous :-) */
636
637         frame     = -1;
638         frame_out = -1;
639         /* the default is not to change the time at all,
640          * but this is overridden by the edit_files routine
641          */
642         t_corr = 0;
643
644         if (n_append == -1)
645         {
646             if (ftpout == efTNG)
647             {
648                 if (ftpout != ftpin)
649                 {
650                     gmx_fatal(FARGS, "When writing TNG the input file format must also be TNG");
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             std::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 (std::abs(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         if (trxout)
911         {
912             close_trx(trxout);
913         }
914         fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
915                 frame, output_env_conv_time(oenv, last_ok_t), output_env_get_time_unit(oenv));
916     }
917
918     return 0;
919 }