Replaced NOTSET from typedefs.h by local solutions.
[gromacs.git] / src / gromacs / gmxana / gmx_trjcat.cpp
blob30971680fdc8910e13d6e5df910698fab2112ea3
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdlib>
41 #include <cstring>
43 #include <algorithm>
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/fileio/gmxfio.h"
48 #include "gromacs/fileio/pdbio.h"
49 #include "gromacs/fileio/tngio.h"
50 #include "gromacs/fileio/tngio_for_tools.h"
51 #include "gromacs/fileio/trx.h"
52 #include "gromacs/fileio/trxio.h"
53 #include "gromacs/fileio/xtcio.h"
54 #include "gromacs/fileio/xvgr.h"
55 #include "gromacs/gmxana/gmx_ana.h"
56 #include "gromacs/legacyheaders/names.h"
57 #include "gromacs/legacyheaders/typedefs.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/topology/index.h"
60 #include "gromacs/utility/arraysize.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/futil.h"
63 #include "gromacs/utility/smalloc.h"
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)
73 static void scan_trj_files(char **fnms, int nfiles, real *readtime,
74 real *timestep, atom_id imax,
75 const gmx_output_env_t *oenv)
77 /* Check start time of all files */
78 int i, natoms = 0;
79 t_trxstatus *status;
80 t_trxframe fr;
81 gmx_bool ok;
83 for (i = 0; i < nfiles; i++)
85 ok = read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
87 if (!ok)
89 gmx_fatal(FARGS, "\nCouldn't read frame from file." );
91 if (fr.bTime)
93 readtime[i] = fr.time;
95 else
97 readtime[i] = 0;
98 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
101 if (i == 0)
103 natoms = fr.natoms;
105 else
107 if (imax == NO_ATID)
109 if (natoms != fr.natoms)
111 gmx_fatal(FARGS, "\nDifferent numbers of atoms (%d/%d) in files",
112 natoms, fr.natoms);
115 else
117 if (fr.natoms <= imax)
119 gmx_fatal(FARGS, "\nNot enough atoms (%d) for index group (%d)",
120 fr.natoms, imax);
124 ok = read_next_frame(oenv, status, &fr);
125 if (ok && fr.bTime)
127 timestep[i] = fr.time - readtime[i];
129 else
131 timestep[i] = 0;
134 close_trj(status);
135 if (fr.bX)
137 sfree(fr.x);
139 if (fr.bV)
141 sfree(fr.v);
143 if (fr.bF)
145 sfree(fr.f);
148 fprintf(stderr, "\n");
152 static void sort_files(char **fnms, real *settime, int nfile)
154 int i, j, minidx;
155 real timeswap;
156 char *chptr;
158 for (i = 0; i < nfile; i++)
160 minidx = i;
161 for (j = i + 1; j < nfile; j++)
163 if (settime[j] < settime[minidx])
165 minidx = j;
168 if (minidx != i)
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;
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 gmx_output_env_t *oenv)
184 int i;
185 gmx_bool ok;
186 char inputstring[STRLEN], *chptr;
188 if (bSetTime)
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));
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));
207 for (i = 0; i < nfiles; i++)
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;
214 if (NULL == fgets(inputstring, STRLEN - 1, stdin))
216 gmx_fatal(FARGS, "Error reading user input" );
219 inputstring[std::strlen(inputstring)-1] = 0;
221 if (inputstring[0] == 'c' || inputstring[0] == 'C')
223 cont_type[i] = TIME_CONTINUE;
224 bSort = FALSE;
225 ok = TRUE;
226 settime[i] = FLT_MAX;
228 else if (inputstring[0] == 'l' ||
229 inputstring[0] == 'L')
231 cont_type[i] = TIME_LAST;
232 bSort = FALSE;
233 ok = TRUE;
234 settime[i] = FLT_MAX;
236 else
238 settime[i] = strtod(inputstring, &chptr)*
239 output_env_get_time_invfactor(oenv);
240 if (chptr == inputstring)
242 fprintf(stderr, "'%s' not recognized as a floating point number, 'c' or 'l'. "
243 "Try again: ", inputstring);
245 else
247 cont_type[i] = TIME_EXPLICIT;
248 ok = TRUE;
252 while (!ok);
254 if (cont_type[0] != TIME_EXPLICIT)
256 cont_type[0] = TIME_EXPLICIT;
257 settime[0] = 0;
260 else
262 for (i = 0; i < nfiles; i++)
264 settime[i] = readtime[i];
267 if (!bSort)
269 fprintf(stderr, "Sorting disabled.\n");
271 else
273 sort_files(fnms, settime, nfiles);
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++)
281 switch (cont_type[i])
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])
291 fprintf(stderr, " WARNING: same Start time as previous");
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;
304 fprintf(stderr, "\n");
306 settime[nfiles] = FLT_MAX;
307 cont_type[nfiles] = TIME_EXPLICIT;
308 readtime[nfiles] = FLT_MAX;
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 gmx_output_env_t *oenv)
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;
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++)
328 nnn = read_first_frame(oenv, &(fp_in[i]), fnms[i], &(trx[i]),
329 TRX_NEED_X);
330 if (natoms == -1)
332 natoms = nnn;
333 first_time = trx[i].time;
335 else if (natoms != nnn)
337 gmx_fatal(FARGS, "Trajectory file %s has %d atoms while previous trajs had %d atoms", fnms[i], nnn, natoms);
339 if (t == -1)
341 t = trx[i].time;
343 else if (t != trx[i].time)
345 gmx_fatal(FARGS, "Trajectory file %s has time %f while previous trajs had time %f", fnms[i], trx[i].time, t);
349 snew(fp_out, nset);
350 for (i = 0; (i < nset); i++)
352 fp_out[i] = open_trx(fnms_out[i], "w");
354 k = 0;
355 if (std::round(time[k] - t) != 0)
357 gmx_fatal(FARGS, "First time in demuxing table does not match trajectories");
361 while ((k+1 < nval) && ((trx[0].time - time[k+1]) > dt_remd*0.1))
363 k++;
365 if (debug)
367 fprintf(debug, "trx[0].time = %g, time[k] = %g\n", trx[0].time, time[k]);
369 for (i = 0; (i < nset); i++)
371 bSet[i] = FALSE;
373 for (i = 0; (i < nset); i++)
375 j = std::round(value[i][k]);
376 range_check(j, 0, nset);
377 if (bSet[j])
379 gmx_fatal(FARGS, "Demuxing the same replica %d twice at time %f",
380 j, trx[0].time);
382 bSet[j] = TRUE;
384 if (dt == 0 || bRmod(trx[i].time, first_time, dt))
386 if (index)
388 write_trxframe_indexed(fp_out[j], &trx[i], isize, index, NULL);
390 else
392 write_trxframe(fp_out[j], &trx[i], NULL);
397 bCont = (k < nval);
398 for (i = 0; (i < nset); i++)
400 bCont = bCont && read_next_frame(oenv, fp_in[i], &trx[i]);
403 while (bCont);
405 for (i = 0; (i < nset); i++)
407 close_trx(fp_in[i]);
408 close_trx(fp_out[i]);
412 int gmx_trjcat(int argc, char *argv[])
414 const char *desc[] =
416 "[THISMODULE] concatenates several input trajectory files in sorted order. ",
417 "In case of double time frames the one in the later file is used. ",
418 "By specifying [TT]-settime[tt] you will be asked for the start time ",
419 "of each file. The input files are taken from the command line, ",
420 "such that a command like [TT]gmx trjcat -f *.trr -o fixed.trr[tt] should do ",
421 "the trick. Using [TT]-cat[tt], you can simply paste several files ",
422 "together without removal of frames with identical time stamps.[PAR]",
423 "One important option is inferred when the output file is amongst the",
424 "input files. In that case that particular file will be appended to",
425 "which implies you do not need to store double the amount of data.",
426 "Obviously the file to append to has to be the one with lowest starting",
427 "time since one can only append at the end of a file.[PAR]",
428 "If the [TT]-demux[tt] option is given, the N trajectories that are",
429 "read, are written in another order as specified in the [REF].xvg[ref] file.",
430 "The [REF].xvg[ref] file should contain something like::",
432 " 0 0 1 2 3 4 5",
433 " 2 1 0 2 3 5 4",
435 "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 [REF].xvg[ref] file then the program",
440 "tries to be smart. Beware."
442 static gmx_bool bVels = TRUE;
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;
454 t_pargs
455 pa[] =
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 { "-vel", FALSE, etBOOL,
464 { &bVels }, "Read and write velocities if possible" },
465 { "-settime", FALSE, etBOOL,
466 { &bSetTime }, "Change starting time interactively" },
467 { "-sort", FALSE, etBOOL,
468 { &bSort }, "Sort trajectory files (not frames)" },
469 { "-keeplast", FALSE, etBOOL,
470 { &bKeepLast }, "Keep overlapping frames at end of trajectory" },
471 { "-overwrite", FALSE, etBOOL,
472 { &bOverwrite }, "Overwrite overlapping frames during appending" },
473 { "-cat", FALSE, etBOOL,
474 { &bCat }, "Do not discard double time frames" }
476 #define npargs asize(pa)
477 int ftpin, i, frame, frame_out;
478 t_trxstatus *status, *trxout = NULL;
479 real t_corr;
480 t_trxframe fr, frout;
481 char **fnms, **fnms_out, *out_file;
482 int n_append;
483 gmx_bool bNewFile, bIndex, bWrite;
484 int nfile_in, nfile_out, *cont_type;
485 real *readtime, *timest, *settime;
486 real first_time = 0, lasttime, last_ok_t = -1, timestep;
487 gmx_bool lastTimeSet = FALSE;
488 real last_frame_time, searchtime;
489 int isize = 0, j;
490 atom_id *index = NULL, imax;
491 char *grpname;
492 real **val = NULL, *t = NULL, dt_remd;
493 int n, nset, ftpout = -1, prevEndStep = 0, filetype;
494 gmx_off_t fpos;
495 gmx_output_env_t *oenv;
496 t_filenm fnm[] =
498 { efTRX, "-f", NULL, ffRDMULT },
499 { efTRO, "-o", NULL, ffWRMULT },
500 { efNDX, "-n", "index", ffOPTRD },
501 { efXVG, "-demux", "remd", ffOPTRD }
504 #define NFILE asize(fnm)
506 if (!parse_common_args(&argc, argv, PCA_TIME_UNIT, NFILE, fnm,
507 asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
509 return 0;
512 bIndex = ftp2bSet(efNDX, NFILE, fnm);
513 bDeMux = ftp2bSet(efXVG, NFILE, fnm);
514 bSort = bSort && !bDeMux;
516 imax = NO_ATID;
517 if (bIndex)
519 printf("Select group for output\n");
520 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
521 /* scan index */
522 imax = index[0];
523 for (i = 1; i < isize; i++)
525 imax = std::max(imax, index[i]);
528 if (bDeMux)
530 nset = 0;
531 dt_remd = 0;
532 val = read_xvg_time(opt2fn("-demux", NFILE, fnm), TRUE,
533 opt2parg_bSet("-b", npargs, pa), begin,
534 opt2parg_bSet("-e", npargs, pa), end, 1, &nset, &n,
535 &dt_remd, &t);
536 printf("Read %d sets of %d points, dt = %g\n\n", nset, n, dt_remd);
537 if (debug)
539 fprintf(debug, "Dump of replica_index.xvg\n");
540 for (i = 0; (i < n); i++)
542 fprintf(debug, "%10g", t[i]);
543 for (j = 0; (j < nset); j++)
545 fprintf(debug, " %3d", static_cast<int>(std::round(val[j][i])));
547 fprintf(debug, "\n");
552 nfile_in = opt2fns(&fnms, "-f", NFILE, fnm);
553 if (!nfile_in)
555 gmx_fatal(FARGS, "No input files!" );
558 if (bDeMux && (nfile_in != nset))
560 gmx_fatal(FARGS, "You have specified %d files and %d entries in the demux table", nfile_in, nset);
563 ftpin = fn2ftp(fnms[0]);
565 for (i = 1; i < nfile_in; i++)
567 if (ftpin != fn2ftp(fnms[i]))
569 gmx_fatal(FARGS, "All input files must be of the same format");
573 nfile_out = opt2fns(&fnms_out, "-o", NFILE, fnm);
574 if (!nfile_out)
576 gmx_fatal(FARGS, "No output files!");
578 if ((nfile_out > 1) && !bDeMux)
580 gmx_fatal(FARGS, "Don't know what to do with more than 1 output file if not demultiplexing");
582 else if (bDeMux && (nfile_out != nset) && (nfile_out != 1))
584 gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %d", nset, nfile_out);
586 if (bDeMux)
588 if (nfile_out != nset)
590 char *buf = gmx_strdup(fnms_out[0]);
591 snew(fnms_out, nset);
592 for (i = 0; (i < nset); i++)
594 snew(fnms_out[i], std::strlen(buf)+32);
595 sprintf(fnms_out[i], "%d_%s", i, buf);
597 sfree(buf);
599 do_demux(nfile_in, fnms, fnms_out, n, val, t, dt_remd, isize, index, dt, oenv);
601 else
603 snew(readtime, nfile_in+1);
604 snew(timest, nfile_in+1);
605 scan_trj_files(fnms, nfile_in, readtime, timest, imax, oenv);
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);
612 /* Check whether the output file is amongst the input files
613 * This has to be done after sorting etc.
615 out_file = fnms_out[0];
616 ftpout = fn2ftp(out_file);
617 n_append = -1;
618 for (i = 0; ((i < nfile_in) && (n_append == -1)); i++)
620 if (std::strcmp(fnms[i], out_file) == 0)
622 n_append = i;
625 if (n_append == 0)
627 fprintf(stderr, "Will append to %s rather than creating a new file\n",
628 out_file);
630 else if (n_append != -1)
632 gmx_fatal(FARGS, "Can only append to the first file which is %s (not %s)",
633 fnms[0], out_file);
636 /* Not checking input format, could be dangerous :-) */
637 /* Not checking output format, equally dangerous :-) */
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
644 t_corr = 0;
646 if (n_append == -1)
648 if (ftpout == efTNG)
650 if (ftpout != ftpin)
652 gmx_fatal(FARGS, "When writing TNG the input file format must also be TNG");
654 if (bIndex)
656 trjtools_gmx_prepare_tng_writing(out_file, 'w', NULL, &trxout,
657 fnms[0], isize, NULL, index, grpname);
659 else
661 trjtools_gmx_prepare_tng_writing(out_file, 'w', NULL, &trxout,
662 fnms[0], -1, NULL, NULL, NULL);
665 else
667 trxout = open_trx(out_file, "w");
669 std::memset(&frout, 0, sizeof(frout));
671 else
673 t_fileio *stfio;
675 if (!read_first_frame(oenv, &status, out_file, &fr, FLAGS))
677 gmx_fatal(FARGS, "Reading first frame from %s", out_file);
680 stfio = trx_get_fileio(status);
681 if (!bKeepLast && !bOverwrite)
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");
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 || filetype == efTNG)
694 lasttime = trx_get_time_of_final_frame(status);
695 fr.time = lasttime;
697 else
699 while (read_next_frame(oenv, status, &fr))
703 lasttime = fr.time;
705 lastTimeSet = TRUE;
706 bKeepLastAppend = TRUE;
707 close_trj(status);
708 trxout = open_trx(out_file, "a");
710 else if (bOverwrite)
712 if (gmx_fio_getftp(stfio) != efXTC)
714 gmx_fatal(FARGS, "Overwrite only supported for XTC." );
716 last_frame_time = trx_get_time_of_final_frame(status);
718 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
719 * or when seek time = 0 */
720 if (nfile_in > 1 && settime[1] < last_frame_time+timest[0]*0.5)
722 /* Jump to one time-frame before the start of next
723 * trajectory file */
724 searchtime = settime[1]-timest[0]*1.25;
726 else
728 searchtime = last_frame_time;
730 if (xtc_seek_time(stfio, searchtime, fr.natoms, TRUE))
732 gmx_fatal(FARGS, "Error seeking to append position.");
734 read_next_frame(oenv, status, &fr);
735 if (std::abs(searchtime - fr.time) > timest[0]*0.5)
737 gmx_fatal(FARGS, "Error seeking: attempted to seek to %f but got %f.",
738 searchtime, fr.time);
740 lasttime = fr.time;
741 lastTimeSet = TRUE;
742 fpos = gmx_fio_ftell(stfio);
743 close_trj(status);
744 trxout = open_trx(out_file, "r+");
745 if (gmx_fio_seek(trx_get_fileio(trxout), fpos))
747 gmx_fatal(FARGS, "Error seeking to append position.");
750 if (lastTimeSet)
752 printf("\n Will append after %f \n", lasttime);
754 frout = fr;
756 /* Lets stitch up some files */
757 timestep = timest[0];
758 for (i = n_append+1; (i < nfile_in); i++)
760 /* Open next file */
762 /* set the next time from the last frame in previous file */
763 if (i > 0)
765 /* When writing TNG the step determine which frame to write. Use an
766 * offset to be able to increase steps properly when changing files. */
767 if (ftpout == efTNG)
769 prevEndStep = frout.step;
772 if (frame_out >= 0)
774 if (cont_type[i] == TIME_CONTINUE)
776 begin = frout.time;
777 begin += 0.5*timestep;
778 settime[i] = frout.time;
779 cont_type[i] = TIME_EXPLICIT;
781 else if (cont_type[i] == TIME_LAST)
783 begin = frout.time;
784 begin += 0.5*timestep;
786 /* Or, if the time in the next part should be changed by the
787 * same amount, start at half a timestep from the last time
788 * so we dont repeat frames.
790 /* I don't understand the comment above, but for all the cases
791 * I tried the code seems to work properly. B. Hess 2008-4-2.
794 /* Or, if time is set explicitly, we check for overlap/gap */
795 if (cont_type[i] == TIME_EXPLICIT)
797 if ( ( i < nfile_in ) &&
798 ( frout.time < settime[i]-1.5*timestep ) )
800 fprintf(stderr, "WARNING: Frames around t=%f %s have a different "
801 "spacing than the rest,\n"
802 "might be a gap or overlap that couldn't be corrected "
803 "automatically.\n", output_env_conv_time(oenv, frout.time),
804 output_env_get_time_unit(oenv));
809 /* if we don't have a timestep in the current file, use the old one */
810 if (timest[i] != 0)
812 timestep = timest[i];
814 read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
815 if (!fr.bTime)
817 fr.time = 0;
818 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
821 if (cont_type[i] == TIME_EXPLICIT)
823 t_corr = settime[i]-fr.time;
825 /* t_corr is the amount we want to change the time.
826 * If the user has chosen not to change the time for
827 * this part of the trajectory t_corr remains at
828 * the value it had in the last part, changing this
829 * by the same amount.
830 * If no value was given for the first trajectory part
831 * we let the time start at zero, see the edit_files routine.
834 bNewFile = TRUE;
836 if (!lastTimeSet)
838 lasttime = 0;
839 lastTimeSet = true;
841 printf("\n");
842 printf("lasttime %g\n", lasttime);
846 /* copy the input frame to the output frame */
847 frout = fr;
848 /* set the new time by adding the correct calculated above */
849 frout.time += t_corr;
850 if (ftpout == efTNG)
852 frout.step += prevEndStep;
854 /* quit if we have reached the end of what should be written */
855 if ((end > 0) && (frout.time > end+GMX_REAL_EPS))
857 i = nfile_in;
858 break;
861 /* determine if we should write this frame (dt is handled elsewhere) */
862 if (bCat) /* write all frames of all files */
864 bWrite = TRUE;
866 else if (bKeepLast || (bKeepLastAppend && i == 1))
867 /* write till last frame of this traj
868 and skip first frame(s) of next traj */
870 bWrite = ( frout.time > lasttime+0.5*timestep );
872 else /* write till first frame of next traj */
874 bWrite = ( frout.time < settime[i+1]-0.5*timestep );
877 if (bWrite && (frout.time >= begin) )
879 frame++;
880 if (frame_out == -1)
882 first_time = frout.time;
884 lasttime = frout.time;
885 lastTimeSet = TRUE;
886 if (dt == 0 || bRmod(frout.time, first_time, dt))
888 frame_out++;
889 last_ok_t = frout.time;
890 if (bNewFile)
892 fprintf(stderr, "\nContinue writing frames from %s t=%g %s, "
893 "frame=%d \n",
894 fnms[i], output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv),
895 frame);
896 bNewFile = FALSE;
899 if (bIndex)
901 write_trxframe_indexed(trxout, &frout, isize, index,
902 NULL);
904 else
906 write_trxframe(trxout, &frout, NULL);
908 if ( ((frame % 10) == 0) || (frame < 10) )
910 fprintf(stderr, " -> frame %6d time %8.3f %s \r",
911 frame_out, output_env_conv_time(oenv, frout.time), output_env_get_time_unit(oenv));
916 while (read_next_frame(oenv, status, &fr));
918 close_trj(status);
920 if (trxout)
922 close_trx(trxout);
924 fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
925 frame, output_env_conv_time(oenv, last_ok_t), output_env_get_time_unit(oenv));
928 return 0;