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