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