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