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