Python detection consolidation.
[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,2018,2019, 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 #include <string>
45
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/fileio/gmxfio.h"
49 #include "gromacs/fileio/pdbio.h"
50 #include "gromacs/fileio/tngio.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/math/vec.h"
56 #include "gromacs/mdtypes/md_enums.h"
57 #include "gromacs/topology/index.h"
58 #include "gromacs/trajectory/trajectoryframe.h"
59 #include "gromacs/utility/arrayref.h"
60 #include "gromacs/utility/arraysize.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/futil.h"
63 #include "gromacs/utility/smalloc.h"
64 #include "gromacs/utility/stringutil.h"
65
66 #define TIME_EXPLICIT 0
67 #define TIME_CONTINUE 1
68 #define TIME_LAST     2
69 #ifndef FLT_MAX
70 #define FLT_MAX 1e36
71 #endif
72 #define FLAGS (TRX_READ_X | TRX_READ_V | TRX_READ_F)
73
74 static void scan_trj_files(gmx::ArrayRef<const std::string> files,
75                            real *readtime,
76                            real *timestep, int imax,
77                            const gmx_output_env_t *oenv)
78 {
79     /* Check start time of all files */
80     int          natoms = 0;
81     t_trxstatus *status;
82     t_trxframe   fr;
83     bool         ok;
84
85     for (gmx::index i = 0; i < files.ssize(); i++)
86     {
87         ok = read_first_frame(oenv, &status, files[i].c_str(), &fr, FLAGS);
88
89         if (!ok)
90         {
91             gmx_fatal(FARGS, "\nCouldn't read frame from file." );
92         }
93         if (fr.bTime)
94         {
95             readtime[i] = fr.time;
96         }
97         else
98         {
99             readtime[i] = 0;
100             fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
101         }
102
103         if (i == 0)
104         {
105             natoms = fr.natoms;
106         }
107         else
108         {
109             if (imax == -1)
110             {
111                 if (natoms != fr.natoms)
112                 {
113                     gmx_fatal(FARGS, "\nDifferent numbers of atoms (%d/%d) in files",
114                               natoms, fr.natoms);
115                 }
116             }
117             else
118             {
119                 if (fr.natoms <= imax)
120                 {
121                     gmx_fatal(FARGS, "\nNot enough atoms (%d) for index group (%d)",
122                               fr.natoms, imax);
123                 }
124             }
125         }
126         ok = read_next_frame(oenv, status, &fr);
127         if (ok && fr.bTime)
128         {
129             timestep[i] = fr.time - readtime[i];
130         }
131         else
132         {
133             timestep[i] = 0;
134         }
135
136         close_trx(status);
137         if (fr.bX)
138         {
139             sfree(fr.x);
140         }
141         if (fr.bV)
142         {
143             sfree(fr.v);
144         }
145         if (fr.bF)
146         {
147             sfree(fr.f);
148         }
149     }
150     fprintf(stderr, "\n");
151
152 }
153
154 static void sort_files(gmx::ArrayRef<std::string> files, real *settime)
155 {
156     for (gmx::index i = 0; i < files.ssize(); i++)
157     {
158         gmx::index minidx = i;
159         for (gmx::index j = i + 1; j < files.ssize(); j++)
160         {
161             if (settime[j] < settime[minidx])
162             {
163                 minidx = j;
164             }
165         }
166         if (minidx != i)
167         {
168             real timeswap   = settime[i];
169             settime[i]      = settime[minidx];
170             settime[minidx] = timeswap;
171             std::swap(files[i], files[minidx]);
172         }
173     }
174 }
175
176 static void edit_files(gmx::ArrayRef<std::string> files,
177                        real *readtime, real *timestep,
178                        real *settime, int *cont_type, gmx_bool bSetTime,
179                        gmx_bool bSort, const gmx_output_env_t *oenv)
180 {
181     gmx_bool ok;
182     char     inputstring[STRLEN], *chptr;
183
184     auto     timeUnit = output_env_get_time_unit(oenv);
185     if (bSetTime)
186     {
187         fprintf(stderr, "\n\nEnter the new start time (%s) for each file.\n"
188                 "There are two special options, both disable sorting:\n\n"
189                 "c (continue) - The start time is taken from the end\n"
190                 "of the previous file. Use it when your continuation run\n"
191                 "restarts with t=0.\n\n"
192                 "l (last) - The time in this file will be changed the\n"
193                 "same amount as in the previous. Use it when the time in the\n"
194                 "new run continues from the end of the previous one,\n"
195                 "since this takes possible overlap into account.\n\n",
196                 timeUnit.c_str());
197
198         fprintf(
199                 stderr,
200                 "          File             Current start (%s)  New start (%s)\n"
201                 "---------------------------------------------------------\n",
202                 timeUnit.c_str(), timeUnit.c_str());
203
204         for (gmx::index i = 0; i < files.ssize(); i++)
205         {
206             fprintf(stderr, "%25s   %10.3f %s          ", files[i].c_str(),
207                     output_env_conv_time(oenv, readtime[i]), timeUnit.c_str());
208             ok = FALSE;
209             do
210             {
211                 if (nullptr == fgets(inputstring, STRLEN - 1, stdin))
212                 {
213                     gmx_fatal(FARGS, "Error reading user input" );
214                 }
215
216                 inputstring[std::strlen(inputstring)-1] = 0;
217
218                 if (inputstring[0] == 'c' || inputstring[0] == 'C')
219                 {
220                     cont_type[i] = TIME_CONTINUE;
221                     bSort        = FALSE;
222                     ok           = TRUE;
223                     settime[i]   = FLT_MAX;
224                 }
225                 else if (inputstring[0] == 'l' ||
226                          inputstring[0] == 'L')
227                 {
228                     cont_type[i] = TIME_LAST;
229                     bSort        = FALSE;
230                     ok           = TRUE;
231                     settime[i]   = FLT_MAX;
232                 }
233                 else
234                 {
235                     settime[i] = strtod(inputstring, &chptr)*
236                         output_env_get_time_invfactor(oenv);
237                     if (chptr == inputstring)
238                     {
239                         fprintf(stderr, "'%s' not recognized as a floating point number, 'c' or 'l'. "
240                                 "Try again: ", inputstring);
241                     }
242                     else
243                     {
244                         cont_type[i] = TIME_EXPLICIT;
245                         ok           = TRUE;
246                     }
247                 }
248             }
249             while (!ok);
250         }
251         if (cont_type[0] != TIME_EXPLICIT)
252         {
253             cont_type[0] = TIME_EXPLICIT;
254             settime[0]   = 0;
255         }
256     }
257     else
258     {
259         for (gmx::index i = 0; i < files.ssize(); i++)
260         {
261             settime[i] = readtime[i];
262         }
263     }
264     if (!bSort)
265     {
266         fprintf(stderr, "Sorting disabled.\n");
267     }
268     else
269     {
270         sort_files(files, settime);
271     }
272     /* Write out the new order and start times */
273     fprintf(stderr, "\nSummary of files and start times used:\n\n"
274             "          File                Start time       Time step\n"
275             "---------------------------------------------------------\n");
276     for (gmx::index i = 0; i < files.ssize(); i++)
277     {
278         switch (cont_type[i])
279         {
280             case TIME_EXPLICIT:
281                 fprintf(stderr, "%25s   %10.3f %s   %10.3f %s",
282                         files[i].c_str(),
283                         output_env_conv_time(oenv, settime[i]), timeUnit.c_str(),
284                         output_env_conv_time(oenv, timestep[i]), timeUnit.c_str());
285                 if (i > 0 &&
286                     cont_type[i-1] == TIME_EXPLICIT && settime[i] == settime[i-1])
287                 {
288                     fprintf(stderr, " WARNING: same Start time as previous");
289                 }
290                 fprintf(stderr, "\n");
291                 break;
292             case TIME_CONTINUE:
293                 fprintf(stderr, "%25s        Continue from last file\n", files[i].c_str());
294                 break;
295             case TIME_LAST:
296                 fprintf(stderr, "%25s        Change by same amount as last file\n",
297                         files[i].c_str());
298                 break;
299         }
300     }
301     fprintf(stderr, "\n");
302
303     settime[files.size()]   = FLT_MAX;
304     cont_type[files.size()] = TIME_EXPLICIT;
305     readtime[files.size()]  = FLT_MAX;
306 }
307
308 static void do_demux(gmx::ArrayRef<const std::string> inFiles,
309                      gmx::ArrayRef<const std::string> outFiles, 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           k, natoms;
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, inFiles.size());
320     snew(trx, inFiles.size());
321     snew(bSet, inFiles.size());
322     natoms = -1;
323     t      = -1;
324     for (gmx::index i = 0; i < inFiles.ssize(); i++)
325     {
326         read_first_frame(oenv, &(fp_in[i]), inFiles[i].c_str(), &(trx[i]),
327                          TRX_NEED_X);
328         if (natoms == -1)
329         {
330             natoms     = trx[i].natoms;
331             first_time = trx[i].time;
332         }
333         else if (natoms != trx[i].natoms)
334         {
335             gmx_fatal(FARGS, "Trajectory file %s has %d atoms while previous trajs had %d atoms", inFiles[i].c_str(), trx[i].natoms, 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", inFiles[i].c_str(), trx[i].time, t);
344         }
345     }
346
347     snew(fp_out, inFiles.size());
348     for (gmx::index i = 0; i < inFiles.ssize(); i++)
349     {
350         fp_out[i] = open_trx(outFiles[i].c_str(), "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 (gmx::index i = 0; i < inFiles.ssize(); i++)
368         {
369             bSet[i] = FALSE;
370         }
371         for (gmx::index i = 0; i < inFiles.ssize(); i++)
372         {
373             int j = gmx::roundToInt(value[i][k]);
374             range_check(j, 0, inFiles.size());
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 (gmx::index i = 0; i < inFiles.ssize(); i++)
397         {
398             bCont = bCont && read_next_frame(oenv, fp_in[i], &trx[i]);
399         }
400     }
401     while (bCont);
402
403     for (gmx::index i = 0; i < inFiles.ssize(); 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     int               n_append;
477     gmx_bool          bNewFile, bIndex, bWrite;
478     int              *cont_type;
479     real             *readtime, *timest, *settime;
480     real              first_time  = 0, lasttime = 0, last_ok_t = -1, timestep;
481     gmx_bool          lastTimeSet = FALSE;
482     real              last_frame_time, searchtime;
483     int               isize = 0, j;
484     int              *index = nullptr, imax;
485     char             *grpname;
486     real            **val = nullptr, *t = nullptr, dt_remd;
487     int               n, nset, ftpout = -1, prevEndStep = 0, filetype;
488     gmx_off_t         fpos;
489     gmx_output_env_t *oenv;
490     t_filenm          fnm[] =
491     {
492         { efTRX, "-f", nullptr, ffRDMULT },
493         { efTRO, "-o", nullptr, ffWRMULT },
494         { efNDX, "-n", "index", ffOPTRD },
495         { efXVG, "-demux", "remd", ffOPTRD }
496     };
497
498 #define NFILE asize(fnm)
499
500     if (!parse_common_args(&argc, argv, PCA_TIME_UNIT, NFILE, fnm,
501                            asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
502     {
503         return 0;
504     }
505     fprintf(stdout, "Note that major changes are planned in future for "
506             "trjcat, to improve usability and utility.");
507
508     auto timeUnit = output_env_get_time_unit(oenv);
509
510     bIndex = ftp2bSet(efNDX, NFILE, fnm);
511     bDeMux = ftp2bSet(efXVG, NFILE, fnm);
512     bSort  = bSort && !bDeMux;
513
514     imax = -1;
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", static_cast<int>(std::round(val[j][i])));
544                 }
545                 fprintf(debug, "\n");
546             }
547         }
548     }
549
550     gmx::ArrayRef<const std::string> inFiles = opt2fns("-f", NFILE, fnm);
551     if (inFiles.empty())
552     {
553         gmx_fatal(FARGS, "No input files!" );
554     }
555
556     if (bDeMux && ssize(inFiles) != nset)
557     {
558         gmx_fatal(FARGS, "You have specified %td files and %d entries in the demux table", inFiles.ssize(), nset);
559     }
560
561     ftpin = fn2ftp(inFiles[0].c_str());
562
563     if (ftpin != efTRR && ftpin != efXTC && ftpin != efTNG)
564     {
565         gmx_fatal(FARGS, "gmx trjcat can only handle binary trajectory formats (trr, xtc, tng)");
566     }
567
568     for (const std::string &inFile : inFiles)
569     {
570         if (ftpin != fn2ftp(inFile.c_str()))
571         {
572             gmx_fatal(FARGS, "All input files must be of the same (trr, xtc or tng) format");
573         }
574     }
575
576     gmx::ArrayRef<const std::string> outFiles = opt2fns("-o", NFILE, fnm);
577     if (outFiles.empty())
578     {
579         gmx_fatal(FARGS, "No output files!");
580     }
581     if ((outFiles.size() > 1) && !bDeMux)
582     {
583         gmx_fatal(FARGS, "Don't know what to do with more than 1 output file if  not demultiplexing");
584     }
585     else if (bDeMux && ssize(outFiles) != nset && outFiles.size() != 1)
586     {
587         gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %td", nset, outFiles.ssize());
588     }
589     if (bDeMux)
590     {
591         auto outFilesDemux = gmx::copyOf(outFiles);
592         if (gmx::ssize(outFilesDemux) != nset)
593         {
594             std::string name = outFilesDemux[0];
595             outFilesDemux.resize(nset);
596             for (i = 0; (i < nset); i++)
597             {
598                 outFilesDemux[0] = gmx::formatString("%d_%s", i, name.c_str());
599             }
600         }
601         do_demux(inFiles, outFilesDemux, n, val, t, dt_remd, isize, index, dt, oenv);
602     }
603     else
604     {
605         snew(readtime, inFiles.size() + 1);
606         snew(timest, inFiles.size() + 1);
607         scan_trj_files(inFiles, readtime, timest, imax, oenv);
608
609         snew(settime, inFiles.size() + 1);
610         snew(cont_type, inFiles.size() + 1);
611         auto inFilesEdited = gmx::copyOf(inFiles);
612         edit_files(inFilesEdited, readtime, timest, settime, cont_type, bSetTime, bSort,
613                    oenv);
614
615         /* Check whether the output file is amongst the input files
616          * This has to be done after sorting etc.
617          */
618         const char *out_file = outFiles[0].c_str();
619         ftpout   = fn2ftp(out_file);
620         n_append = -1;
621         for (size_t i = 0; i < inFilesEdited.size() && n_append == -1; i++)
622         {
623             if (std::strcmp(inFilesEdited[i].c_str(), out_file) == 0)
624             {
625                 n_append = i;
626             }
627         }
628         if (n_append == 0)
629         {
630             fprintf(stderr, "Will append to %s rather than creating a new file\n",
631                     out_file);
632         }
633         else if (n_append != -1)
634         {
635             gmx_fatal(FARGS, "Can only append to the first file which is %s (not %s)",
636                       inFilesEdited[0].c_str(), out_file);
637         }
638
639         /* Not checking input format, could be dangerous :-) */
640         /* Not checking output format, equally dangerous :-) */
641
642         frame     = -1;
643         frame_out = -1;
644         /* the default is not to change the time at all,
645          * but this is overridden by the edit_files routine
646          */
647         t_corr = 0;
648
649         if (n_append == -1)
650         {
651             if (ftpout == efTNG)
652             {
653                 if (ftpout != ftpin)
654                 {
655                     gmx_fatal(FARGS, "When writing TNG the input file format must also be TNG");
656                 }
657                 if (bIndex)
658                 {
659                     trxout = trjtools_gmx_prepare_tng_writing(out_file, 'w', nullptr,
660                                                               inFilesEdited[0].c_str(), isize, nullptr, gmx::arrayRefFromArray(index, isize), grpname);
661                 }
662                 else
663                 {
664                     trxout = trjtools_gmx_prepare_tng_writing(out_file, 'w', nullptr,
665                                                               inFilesEdited[0].c_str(), -1, nullptr, {}, nullptr);
666                 }
667             }
668             else
669             {
670                 trxout = open_trx(out_file, "w");
671             }
672             std::memset(&frout, 0, sizeof(frout));
673         }
674         else
675         {
676             t_fileio *stfio;
677
678             if (!read_first_frame(oenv, &status, out_file, &fr, FLAGS))
679             {
680                 gmx_fatal(FARGS, "Reading first frame from %s", out_file);
681             }
682
683             stfio = trx_get_fileio(status);
684             if (!bKeepLast && !bOverwrite)
685             {
686                 fprintf(stderr, "\n\nWARNING: Appending without -overwrite implies -keeplast "
687                         "between the first two files. \n"
688                         "If the trajectories have an overlap and have not been written binary \n"
689                         "reproducible this will produce an incorrect trajectory!\n\n");
690
691                 filetype = gmx_fio_getftp(stfio);
692                 /* Fails if last frame is incomplete
693                  * We can't do anything about it without overwriting
694                  * */
695                 if (filetype == efXTC || filetype == efTNG)
696                 {
697                     lasttime = trx_get_time_of_final_frame(status);
698                     fr.time  = lasttime;
699                 }
700                 else
701                 {
702                     while (read_next_frame(oenv, status, &fr))
703                     {
704                         ;
705                     }
706                     lasttime = fr.time;
707                 }
708                 lastTimeSet     = TRUE;
709                 bKeepLastAppend = TRUE;
710                 close_trx(status);
711                 trxout = open_trx(out_file, "a");
712             }
713             else if (bOverwrite)
714             {
715                 if (gmx_fio_getftp(stfio) != efXTC)
716                 {
717                     gmx_fatal(FARGS, "Overwrite only supported for XTC." );
718                 }
719                 last_frame_time = trx_get_time_of_final_frame(status);
720
721                 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
722                  *     or when seek time = 0 */
723                 if (inFilesEdited.size() > 1 && settime[1] < last_frame_time+timest[0]*0.5)
724                 {
725                     /* Jump to one time-frame before the start of next
726                      *  trajectory file */
727                     searchtime = settime[1]-timest[0]*1.25;
728                 }
729                 else
730                 {
731                     searchtime = last_frame_time;
732                 }
733                 if (xtc_seek_time(stfio, searchtime, fr.natoms, TRUE))
734                 {
735                     gmx_fatal(FARGS, "Error seeking to append position.");
736                 }
737                 read_next_frame(oenv, status, &fr);
738                 if (std::abs(searchtime - fr.time) > timest[0]*0.5)
739                 {
740                     gmx_fatal(FARGS, "Error seeking: attempted to seek to %f but got %f.",
741                               searchtime, fr.time);
742                 }
743                 lasttime    = fr.time;
744                 lastTimeSet = TRUE;
745                 fpos        = gmx_fio_ftell(stfio);
746                 close_trx(status);
747                 trxout = open_trx(out_file, "r+");
748                 if (gmx_fio_seek(trx_get_fileio(trxout), fpos))
749                 {
750                     gmx_fatal(FARGS, "Error seeking to append position.");
751                 }
752             }
753             if (lastTimeSet)
754             {
755                 printf("\n Will append after %f \n", lasttime);
756             }
757             frout = fr;
758         }
759         /* Lets stitch up some files */
760         timestep = timest[0];
761         for (size_t i = n_append + 1; i < inFilesEdited.size(); i++)
762         {
763             /* Open next file */
764
765             /* set the next time from the last frame in previous file */
766             if (i > 0)
767             {
768                 /* When writing TNG the step determine which frame to write. Use an
769                  * offset to be able to increase steps properly when changing files. */
770                 if (ftpout == efTNG)
771                 {
772                     prevEndStep = frout.step;
773                 }
774
775                 if (frame_out >= 0)
776                 {
777                     if (cont_type[i] == TIME_CONTINUE)
778                     {
779                         begin        = frout.time;
780                         begin       += 0.5*timestep;
781                         settime[i]   = frout.time;
782                         cont_type[i] = TIME_EXPLICIT;
783                     }
784                     else if (cont_type[i] == TIME_LAST)
785                     {
786                         begin  = frout.time;
787                         begin += 0.5*timestep;
788                     }
789                     /* Or, if the time in the next part should be changed by the
790                      * same amount, start at half a timestep from the last time
791                      * so we dont repeat frames.
792                      */
793                     /* I don't understand the comment above, but for all the cases
794                      * I tried the code seems to work properly. B. Hess 2008-4-2.
795                      */
796                 }
797                 /* Or, if time is set explicitly, we check for overlap/gap */
798                 if (cont_type[i] == TIME_EXPLICIT)
799                 {
800                     if (i < inFilesEdited.size() &&
801                         frout.time < settime[i] - 1.5*timestep)
802                     {
803                         fprintf(stderr, "WARNING: Frames around t=%f %s have a different "
804                                 "spacing than the rest,\n"
805                                 "might be a gap or overlap that couldn't be corrected "
806                                 "automatically.\n", output_env_conv_time(oenv, frout.time),
807                                 timeUnit.c_str());
808                     }
809                 }
810             }
811
812             /* if we don't have a timestep in the current file, use the old one */
813             if (timest[i] != 0)
814             {
815                 timestep = timest[i];
816             }
817             read_first_frame(oenv, &status, inFilesEdited[i].c_str(), &fr, FLAGS);
818             if (!fr.bTime)
819             {
820                 fr.time = 0;
821                 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
822             }
823
824             if (cont_type[i] == TIME_EXPLICIT)
825             {
826                 t_corr = settime[i]-fr.time;
827             }
828             /* t_corr is the amount we want to change the time.
829              * If the user has chosen not to change the time for
830              * this part of the trajectory t_corr remains at
831              * the value it had in the last part, changing this
832              * by the same amount.
833              * If no value was given for the first trajectory part
834              * we let the time start at zero, see the edit_files routine.
835              */
836
837             bNewFile = TRUE;
838
839             if (!lastTimeSet)
840             {
841                 lasttime    = 0;
842                 lastTimeSet = true;
843             }
844             printf("\n");
845             printf("lasttime %g\n", lasttime);
846
847             do
848             {
849                 /* copy the input frame to the output frame */
850                 frout = fr;
851                 /* set the new time by adding the correct calculated above */
852                 frout.time += t_corr;
853                 if (ftpout == efTNG)
854                 {
855                     frout.step += prevEndStep;
856                 }
857                 /* quit if we have reached the end of what should be written */
858                 if ((end > 0) && (frout.time > end+GMX_REAL_EPS))
859                 {
860                     i = inFilesEdited.size();
861                     break;
862                 }
863
864                 /* determine if we should write this frame (dt is handled elsewhere) */
865                 if (bCat) /* write all frames of all files */
866                 {
867                     bWrite = TRUE;
868                 }
869                 else if (bKeepLast || (bKeepLastAppend && i == 1))
870                 /* write till last frame of this traj
871                    and skip first frame(s) of next traj */
872                 {
873                     bWrite = ( frout.time > lasttime+0.5*timestep );
874                 }
875                 else /* write till first frame of next traj */
876                 {
877                     bWrite = ( frout.time < settime[i+1]-0.5*timestep );
878                 }
879
880                 if (bWrite && (frout.time >= begin) )
881                 {
882                     frame++;
883                     if (frame_out == -1)
884                     {
885                         first_time = frout.time;
886                     }
887                     lasttime    = frout.time;
888                     lastTimeSet = TRUE;
889                     if (dt == 0 || bRmod(frout.time, first_time, dt))
890                     {
891                         frame_out++;
892                         last_ok_t = frout.time;
893                         if (bNewFile)
894                         {
895                             fprintf(stderr, "\nContinue writing frames from %s t=%g %s, "
896                                     "frame=%d      \n",
897                                     inFilesEdited[i].c_str(),
898                                     output_env_conv_time(oenv, frout.time), timeUnit.c_str(),
899                                     frame);
900                             bNewFile = FALSE;
901                         }
902
903                         if (bIndex)
904                         {
905                             write_trxframe_indexed(trxout, &frout, isize, index,
906                                                    nullptr);
907                         }
908                         else
909                         {
910                             write_trxframe(trxout, &frout, nullptr);
911                         }
912                         if ( ((frame % 10) == 0) || (frame < 10) )
913                         {
914                             fprintf(stderr, " ->  frame %6d time %8.3f %s     \r",
915                                     frame_out, output_env_conv_time(oenv, frout.time), timeUnit.c_str());
916                             fflush(stderr);
917                         }
918                     }
919                 }
920             }
921             while (read_next_frame(oenv, status, &fr));
922
923             close_trx(status);
924         }
925         if (trxout)
926         {
927             close_trx(trxout);
928         }
929         fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
930                 frame, output_env_conv_time(oenv, last_ok_t), timeUnit.c_str());
931     }
932
933     return 0;
934 }