Cleaned up memory usage in gmx traj and trjconv
[gromacs.git] / src / gromacs / gmxana / gmx_trjcat.cpp
blob950c87fee7e88a9d2fd23eb36660fef4aa3ac27e
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,2016,2017, 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>
44 #include <string>
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/fileio/gmxfio.h"
49 #include "gromacs/fileio/pdbio.h"
50 #include "gromacs/fileio/tngio.h"
51 #include "gromacs/fileio/trxio.h"
52 #include "gromacs/fileio/xtcio.h"
53 #include "gromacs/fileio/xvgr.h"
54 #include "gromacs/gmxana/gmx_ana.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdtypes/md_enums.h"
57 #include "gromacs/topology/index.h"
58 #include "gromacs/trajectory/trajectoryframe.h"
59 #include "gromacs/utility/arraysize.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
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)
72 static void scan_trj_files(char **fnms, int nfiles, real *readtime,
73 real *timestep, int imax,
74 const gmx_output_env_t *oenv)
76 /* Check start time of all files */
77 int i, natoms = 0;
78 t_trxstatus *status;
79 t_trxframe fr;
80 gmx_bool ok;
82 for (i = 0; i < nfiles; i++)
84 ok = read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
86 if (!ok)
88 gmx_fatal(FARGS, "\nCouldn't read frame from file." );
90 if (fr.bTime)
92 readtime[i] = fr.time;
94 else
96 readtime[i] = 0;
97 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
100 if (i == 0)
102 natoms = fr.natoms;
104 else
106 if (imax == -1)
108 if (natoms != fr.natoms)
110 gmx_fatal(FARGS, "\nDifferent numbers of atoms (%d/%d) in files",
111 natoms, fr.natoms);
114 else
116 if (fr.natoms <= imax)
118 gmx_fatal(FARGS, "\nNot enough atoms (%d) for index group (%d)",
119 fr.natoms, imax);
123 ok = read_next_frame(oenv, status, &fr);
124 if (ok && fr.bTime)
126 timestep[i] = fr.time - readtime[i];
128 else
130 timestep[i] = 0;
133 close_trx(status);
134 if (fr.bX)
136 sfree(fr.x);
138 if (fr.bV)
140 sfree(fr.v);
142 if (fr.bF)
144 sfree(fr.f);
147 fprintf(stderr, "\n");
151 static void sort_files(char **fnms, real *settime, int nfile)
153 int i, j, minidx;
154 real timeswap;
155 char *chptr;
157 for (i = 0; i < nfile; i++)
159 minidx = i;
160 for (j = i + 1; j < nfile; j++)
162 if (settime[j] < settime[minidx])
164 minidx = j;
167 if (minidx != i)
169 timeswap = settime[i];
170 settime[i] = settime[minidx];
171 settime[minidx] = timeswap;
172 chptr = fnms[i];
173 fnms[i] = fnms[minidx];
174 fnms[minidx] = chptr;
179 static void edit_files(char **fnms, int nfiles, real *readtime, real *timestep,
180 real *settime, int *cont_type, gmx_bool bSetTime,
181 gmx_bool bSort, const gmx_output_env_t *oenv)
183 int i;
184 gmx_bool ok;
185 char inputstring[STRLEN], *chptr;
187 auto timeUnit = output_env_get_time_unit(oenv);
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 timeUnit.c_str());
201 fprintf(
202 stderr,
203 " File Current start (%s) New start (%s)\n"
204 "---------------------------------------------------------\n",
205 timeUnit.c_str(), timeUnit.c_str());
207 for (i = 0; i < nfiles; i++)
209 fprintf(stderr, "%25s %10.3f %s ", fnms[i],
210 output_env_conv_time(oenv, readtime[i]), timeUnit.c_str());
211 ok = FALSE;
214 if (nullptr == 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]), timeUnit.c_str(),
287 output_env_conv_time(oenv, timestep[i]), timeUnit.c_str());
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 int 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, nullptr);
390 else
392 write_trxframe(fp_out[j], &trx[i], nullptr);
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 bCat = FALSE;
443 static gmx_bool bSort = TRUE;
444 static gmx_bool bKeepLast = FALSE;
445 static gmx_bool bKeepLastAppend = FALSE;
446 static gmx_bool bOverwrite = FALSE;
447 static gmx_bool bSetTime = FALSE;
448 static gmx_bool bDeMux;
449 static real begin = -1;
450 static real end = -1;
451 static real dt = 0;
453 t_pargs
454 pa[] =
456 { "-b", FALSE, etTIME,
457 { &begin }, "First time to use (%t)" },
458 { "-e", FALSE, etTIME,
459 { &end }, "Last time to use (%t)" },
460 { "-dt", FALSE, etTIME,
461 { &dt }, "Only write frame when t MOD dt = first time (%t)" },
462 { "-settime", FALSE, etBOOL,
463 { &bSetTime }, "Change starting time interactively" },
464 { "-sort", FALSE, etBOOL,
465 { &bSort }, "Sort trajectory files (not frames)" },
466 { "-keeplast", FALSE, etBOOL,
467 { &bKeepLast }, "Keep overlapping frames at end of trajectory" },
468 { "-overwrite", FALSE, etBOOL,
469 { &bOverwrite }, "Overwrite overlapping frames during appending" },
470 { "-cat", FALSE, etBOOL,
471 { &bCat }, "Do not discard double time frames" }
473 #define npargs asize(pa)
474 int ftpin, i, frame, frame_out;
475 t_trxstatus *status, *trxout = nullptr;
476 real t_corr;
477 t_trxframe fr, frout;
478 char **fnms, **fnms_out, *out_file;
479 int n_append;
480 gmx_bool bNewFile, bIndex, bWrite;
481 int nfile_in, nfile_out, *cont_type;
482 real *readtime, *timest, *settime;
483 real first_time = 0, lasttime, last_ok_t = -1, timestep;
484 gmx_bool lastTimeSet = FALSE;
485 real last_frame_time, searchtime;
486 int isize = 0, j;
487 int *index = nullptr, imax;
488 char *grpname;
489 real **val = nullptr, *t = nullptr, dt_remd;
490 int n, nset, ftpout = -1, prevEndStep = 0, filetype;
491 gmx_off_t fpos;
492 gmx_output_env_t *oenv;
493 t_filenm fnm[] =
495 { efTRX, "-f", nullptr, ffRDMULT },
496 { efTRO, "-o", nullptr, ffWRMULT },
497 { efNDX, "-n", "index", ffOPTRD },
498 { efXVG, "-demux", "remd", ffOPTRD }
501 #define NFILE asize(fnm)
503 if (!parse_common_args(&argc, argv, PCA_TIME_UNIT, NFILE, fnm,
504 asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
506 return 0;
508 auto timeUnit = output_env_get_time_unit(oenv);
510 bIndex = ftp2bSet(efNDX, NFILE, fnm);
511 bDeMux = ftp2bSet(efXVG, NFILE, fnm);
512 bSort = bSort && !bDeMux;
514 imax = -1;
515 if (bIndex)
517 printf("Select group for output\n");
518 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
519 /* scan index */
520 imax = index[0];
521 for (i = 1; i < isize; i++)
523 imax = std::max(imax, index[i]);
526 if (bDeMux)
528 nset = 0;
529 dt_remd = 0;
530 val = read_xvg_time(opt2fn("-demux", NFILE, fnm), TRUE,
531 opt2parg_bSet("-b", npargs, pa), begin,
532 opt2parg_bSet("-e", npargs, pa), end, 1, &nset, &n,
533 &dt_remd, &t);
534 printf("Read %d sets of %d points, dt = %g\n\n", nset, n, dt_remd);
535 if (debug)
537 fprintf(debug, "Dump of replica_index.xvg\n");
538 for (i = 0; (i < n); i++)
540 fprintf(debug, "%10g", t[i]);
541 for (j = 0; (j < nset); j++)
543 fprintf(debug, " %3d", static_cast<int>(std::round(val[j][i])));
545 fprintf(debug, "\n");
550 nfile_in = opt2fns(&fnms, "-f", NFILE, fnm);
551 if (!nfile_in)
553 gmx_fatal(FARGS, "No input files!" );
556 if (bDeMux && (nfile_in != nset))
558 gmx_fatal(FARGS, "You have specified %d files and %d entries in the demux table", nfile_in, nset);
561 ftpin = fn2ftp(fnms[0]);
563 for (i = 1; i < nfile_in; i++)
565 if (ftpin != fn2ftp(fnms[i]))
567 gmx_fatal(FARGS, "All input files must be of the same format");
571 nfile_out = opt2fns(&fnms_out, "-o", NFILE, fnm);
572 if (!nfile_out)
574 gmx_fatal(FARGS, "No output files!");
576 if ((nfile_out > 1) && !bDeMux)
578 gmx_fatal(FARGS, "Don't know what to do with more than 1 output file if not demultiplexing");
580 else if (bDeMux && (nfile_out != nset) && (nfile_out != 1))
582 gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %d", nset, nfile_out);
584 if (bDeMux)
586 if (nfile_out != nset)
588 char *buf = gmx_strdup(fnms_out[0]);
589 snew(fnms_out, nset);
590 for (i = 0; (i < nset); i++)
592 snew(fnms_out[i], std::strlen(buf)+32);
593 sprintf(fnms_out[i], "%d_%s", i, buf);
595 sfree(buf);
597 do_demux(nfile_in, fnms, fnms_out, n, val, t, dt_remd, isize, index, dt, oenv);
599 else
601 snew(readtime, nfile_in+1);
602 snew(timest, nfile_in+1);
603 scan_trj_files(fnms, nfile_in, readtime, timest, imax, oenv);
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);
610 /* Check whether the output file is amongst the input files
611 * This has to be done after sorting etc.
613 out_file = fnms_out[0];
614 ftpout = fn2ftp(out_file);
615 n_append = -1;
616 for (i = 0; ((i < nfile_in) && (n_append == -1)); i++)
618 if (std::strcmp(fnms[i], out_file) == 0)
620 n_append = i;
623 if (n_append == 0)
625 fprintf(stderr, "Will append to %s rather than creating a new file\n",
626 out_file);
628 else if (n_append != -1)
630 gmx_fatal(FARGS, "Can only append to the first file which is %s (not %s)",
631 fnms[0], out_file);
634 /* Not checking input format, could be dangerous :-) */
635 /* Not checking output format, equally dangerous :-) */
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
642 t_corr = 0;
644 if (n_append == -1)
646 if (ftpout == efTNG)
648 if (ftpout != ftpin)
650 gmx_fatal(FARGS, "When writing TNG the input file format must also be TNG");
652 if (bIndex)
654 trjtools_gmx_prepare_tng_writing(out_file, 'w', nullptr, &trxout,
655 fnms[0], isize, nullptr, index, grpname);
657 else
659 trjtools_gmx_prepare_tng_writing(out_file, 'w', nullptr, &trxout,
660 fnms[0], -1, nullptr, nullptr, nullptr);
663 else
665 trxout = open_trx(out_file, "w");
667 std::memset(&frout, 0, sizeof(frout));
669 else
671 t_fileio *stfio;
673 if (!read_first_frame(oenv, &status, out_file, &fr, FLAGS))
675 gmx_fatal(FARGS, "Reading first frame from %s", out_file);
678 stfio = trx_get_fileio(status);
679 if (!bKeepLast && !bOverwrite)
681 fprintf(stderr, "\n\nWARNING: Appending without -overwrite implies -keeplast "
682 "between the first two files. \n"
683 "If the trajectories have an overlap and have not been written binary \n"
684 "reproducible this will produce an incorrect trajectory!\n\n");
686 filetype = gmx_fio_getftp(stfio);
687 /* Fails if last frame is incomplete
688 * We can't do anything about it without overwriting
689 * */
690 if (filetype == efXTC || filetype == efTNG)
692 lasttime = trx_get_time_of_final_frame(status);
693 fr.time = lasttime;
695 else
697 while (read_next_frame(oenv, status, &fr))
701 lasttime = fr.time;
703 lastTimeSet = TRUE;
704 bKeepLastAppend = TRUE;
705 close_trx(status);
706 trxout = open_trx(out_file, "a");
708 else if (bOverwrite)
710 if (gmx_fio_getftp(stfio) != efXTC)
712 gmx_fatal(FARGS, "Overwrite only supported for XTC." );
714 last_frame_time = trx_get_time_of_final_frame(status);
716 /* xtc_seek_time broken for trajectories containing only 1 or 2 frames
717 * or when seek time = 0 */
718 if (nfile_in > 1 && settime[1] < last_frame_time+timest[0]*0.5)
720 /* Jump to one time-frame before the start of next
721 * trajectory file */
722 searchtime = settime[1]-timest[0]*1.25;
724 else
726 searchtime = last_frame_time;
728 if (xtc_seek_time(stfio, searchtime, fr.natoms, TRUE))
730 gmx_fatal(FARGS, "Error seeking to append position.");
732 read_next_frame(oenv, status, &fr);
733 if (std::abs(searchtime - fr.time) > timest[0]*0.5)
735 gmx_fatal(FARGS, "Error seeking: attempted to seek to %f but got %f.",
736 searchtime, fr.time);
738 lasttime = fr.time;
739 lastTimeSet = TRUE;
740 fpos = gmx_fio_ftell(stfio);
741 close_trx(status);
742 trxout = open_trx(out_file, "r+");
743 if (gmx_fio_seek(trx_get_fileio(trxout), fpos))
745 gmx_fatal(FARGS, "Error seeking to append position.");
748 if (lastTimeSet)
750 printf("\n Will append after %f \n", lasttime);
752 frout = fr;
754 /* Lets stitch up some files */
755 timestep = timest[0];
756 for (i = n_append+1; (i < nfile_in); i++)
758 /* Open next file */
760 /* set the next time from the last frame in previous file */
761 if (i > 0)
763 /* When writing TNG the step determine which frame to write. Use an
764 * offset to be able to increase steps properly when changing files. */
765 if (ftpout == efTNG)
767 prevEndStep = frout.step;
770 if (frame_out >= 0)
772 if (cont_type[i] == TIME_CONTINUE)
774 begin = frout.time;
775 begin += 0.5*timestep;
776 settime[i] = frout.time;
777 cont_type[i] = TIME_EXPLICIT;
779 else if (cont_type[i] == TIME_LAST)
781 begin = frout.time;
782 begin += 0.5*timestep;
784 /* Or, if the time in the next part should be changed by the
785 * same amount, start at half a timestep from the last time
786 * so we dont repeat frames.
788 /* I don't understand the comment above, but for all the cases
789 * I tried the code seems to work properly. B. Hess 2008-4-2.
792 /* Or, if time is set explicitly, we check for overlap/gap */
793 if (cont_type[i] == TIME_EXPLICIT)
795 if ( ( i < nfile_in ) &&
796 ( frout.time < settime[i]-1.5*timestep ) )
798 fprintf(stderr, "WARNING: Frames around t=%f %s have a different "
799 "spacing than the rest,\n"
800 "might be a gap or overlap that couldn't be corrected "
801 "automatically.\n", output_env_conv_time(oenv, frout.time),
802 timeUnit.c_str());
807 /* if we don't have a timestep in the current file, use the old one */
808 if (timest[i] != 0)
810 timestep = timest[i];
812 read_first_frame(oenv, &status, fnms[i], &fr, FLAGS);
813 if (!fr.bTime)
815 fr.time = 0;
816 fprintf(stderr, "\nWARNING: Couldn't find a time in the frame.\n");
819 if (cont_type[i] == TIME_EXPLICIT)
821 t_corr = settime[i]-fr.time;
823 /* t_corr is the amount we want to change the time.
824 * If the user has chosen not to change the time for
825 * this part of the trajectory t_corr remains at
826 * the value it had in the last part, changing this
827 * by the same amount.
828 * If no value was given for the first trajectory part
829 * we let the time start at zero, see the edit_files routine.
832 bNewFile = TRUE;
834 if (!lastTimeSet)
836 lasttime = 0;
837 lastTimeSet = true;
839 printf("\n");
840 printf("lasttime %g\n", lasttime);
844 /* copy the input frame to the output frame */
845 frout = fr;
846 /* set the new time by adding the correct calculated above */
847 frout.time += t_corr;
848 if (ftpout == efTNG)
850 frout.step += prevEndStep;
852 /* quit if we have reached the end of what should be written */
853 if ((end > 0) && (frout.time > end+GMX_REAL_EPS))
855 i = nfile_in;
856 break;
859 /* determine if we should write this frame (dt is handled elsewhere) */
860 if (bCat) /* write all frames of all files */
862 bWrite = TRUE;
864 else if (bKeepLast || (bKeepLastAppend && i == 1))
865 /* write till last frame of this traj
866 and skip first frame(s) of next traj */
868 bWrite = ( frout.time > lasttime+0.5*timestep );
870 else /* write till first frame of next traj */
872 bWrite = ( frout.time < settime[i+1]-0.5*timestep );
875 if (bWrite && (frout.time >= begin) )
877 frame++;
878 if (frame_out == -1)
880 first_time = frout.time;
882 lasttime = frout.time;
883 lastTimeSet = TRUE;
884 if (dt == 0 || bRmod(frout.time, first_time, dt))
886 frame_out++;
887 last_ok_t = frout.time;
888 if (bNewFile)
890 fprintf(stderr, "\nContinue writing frames from %s t=%g %s, "
891 "frame=%d \n",
892 fnms[i], output_env_conv_time(oenv, frout.time), timeUnit.c_str(),
893 frame);
894 bNewFile = FALSE;
897 if (bIndex)
899 write_trxframe_indexed(trxout, &frout, isize, index,
900 nullptr);
902 else
904 write_trxframe(trxout, &frout, nullptr);
906 if ( ((frame % 10) == 0) || (frame < 10) )
908 fprintf(stderr, " -> frame %6d time %8.3f %s \r",
909 frame_out, output_env_conv_time(oenv, frout.time), timeUnit.c_str());
910 fflush(stderr);
915 while (read_next_frame(oenv, status, &fr));
917 close_trx(status);
919 if (trxout)
921 close_trx(trxout);
923 fprintf(stderr, "\nLast frame written was %d, time %f %s\n",
924 frame, output_env_conv_time(oenv, last_ok_t), timeUnit.c_str());
927 return 0;