Merge branch 'release-2019' into master
[alexxy/gromacs.git] / src / gromacs / tools / trjconv.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "trjconv.h"
40
41 #include <cmath>
42 #include <cstdlib>
43 #include <cstring>
44
45 #include <algorithm>
46 #include <memory>
47
48 #include "gromacs/commandline/pargs.h"
49 #include "gromacs/commandline/viewit.h"
50 #include "gromacs/fileio/confio.h"
51 #include "gromacs/fileio/g96io.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/groio.h"
54 #include "gromacs/fileio/pdbio.h"
55 #include "gromacs/fileio/tngio.h"
56 #include "gromacs/fileio/tpxio.h"
57 #include "gromacs/fileio/trrio.h"
58 #include "gromacs/fileio/trxio.h"
59 #include "gromacs/fileio/xtcio.h"
60 #include "gromacs/fileio/xvgr.h"
61 #include "gromacs/math/do_fit.h"
62 #include "gromacs/math/functions.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/pbcutil/pbcmethods.h"
67 #include "gromacs/pbcutil/rmpbc.h"
68 #include "gromacs/topology/index.h"
69 #include "gromacs/topology/topology.h"
70 #include "gromacs/trajectory/trajectoryframe.h"
71 #include "gromacs/utility/arrayref.h"
72 #include "gromacs/utility/arraysize.h"
73 #include "gromacs/utility/fatalerror.h"
74 #include "gromacs/utility/futil.h"
75 #include "gromacs/utility/smalloc.h"
76
77 static void mk_filenm(char *base, const char *ext, int ndigit, int file_nr,
78                       char out_file[])
79 {
80     char nbuf[128];
81     int  nd = 0, fnr;
82
83     std::strcpy(out_file, base);
84     fnr = file_nr;
85     do
86     {
87         fnr /= 10;
88         nd++;
89     }
90     while (fnr > 0);
91
92     if (nd < ndigit)
93     {
94         std::strncat(out_file, "00000000000", ndigit-nd);
95     }
96     sprintf(nbuf, "%d.", file_nr);
97     std::strcat(out_file, nbuf);
98     std::strcat(out_file, ext);
99 }
100
101 static void check_trr(const char *fn)
102 {
103     if (fn2ftp(fn) != efTRR)
104     {
105         gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
106     }
107 }
108
109 static void do_trunc(const char *fn, real t0)
110 {
111     t_fileio        *in;
112     FILE            *fp;
113     gmx_bool         bStop, bOK;
114     gmx_trr_header_t sh;
115     gmx_off_t        fpos;
116     char             yesno[256];
117     int              j;
118     real             t = 0;
119
120     if (t0 == -1)
121     {
122         gmx_fatal(FARGS, "You forgot to set the truncation time");
123     }
124
125     /* Check whether this is a .trr file */
126     check_trr(fn);
127
128     in   = gmx_trr_open(fn, "r");
129     fp   = gmx_fio_getfp(in);
130     if (fp == nullptr)
131     {
132         fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
133         gmx_trr_close(in);
134     }
135     else
136     {
137         j     = 0;
138         fpos  = gmx_fio_ftell(in);
139         bStop = FALSE;
140         while (!bStop && gmx_trr_read_frame_header(in, &sh, &bOK))
141         {
142             gmx_trr_read_frame_data(in, &sh, nullptr, nullptr, nullptr, nullptr);
143             fpos = gmx_ftell(fp);
144             t    = sh.t;
145             if (t >= t0)
146             {
147                 gmx_fseek(fp, fpos, SEEK_SET);
148                 bStop = TRUE;
149             }
150         }
151         if (bStop)
152         {
153             fprintf(stderr, "Do you REALLY want to truncate this trajectory (%s) at:\n"
154                     "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
155                     fn, j, t, static_cast<long int>(fpos));
156             if (1 != scanf("%s", yesno))
157             {
158                 gmx_fatal(FARGS, "Error reading user input");
159             }
160             if (std::strcmp(yesno, "YES") == 0)
161             {
162                 fprintf(stderr, "Once again, I'm gonna DO this...\n");
163                 gmx_trr_close(in);
164                 if (0 != gmx_truncate(fn, fpos))
165                 {
166                     gmx_fatal(FARGS, "Error truncating file %s", fn);
167                 }
168             }
169             else
170             {
171                 fprintf(stderr, "Ok, I'll forget about it\n");
172             }
173         }
174         else
175         {
176             fprintf(stderr, "Already at end of file (t=%g)...\n", t);
177             gmx_trr_close(in);
178         }
179     }
180 }
181
182 /*! \brief Read a full molecular topology if useful and available.
183  *
184  * If the input trajectory file is not in TNG format, and the output
185  * file is in TNG format, then we want to try to read a full topology
186  * (if available), so that we can write molecule information to the
187  * output file. The full topology provides better molecule information
188  * than is available from the normal t_topology data used by GROMACS
189  * tools.
190  *
191  * Also, the t_topology is only read under (different) particular
192  * conditions. If both apply, then a .tpr file might be read
193  * twice. Trying to fix this redundancy while trjconv is still an
194  * all-purpose tool does not seem worthwhile.
195  *
196  * Because of the way gmx_prepare_tng_writing is implemented, the case
197  * where the input TNG file has no molecule information will never
198  * lead to an output TNG file having molecule information. Since
199  * molecule information will generally be present if the input TNG
200  * file was written by a GROMACS tool, this seems like reasonable
201  * behaviour. */
202 static std::unique_ptr<gmx_mtop_t>
203 read_mtop_for_tng(const char *tps_file,
204                   const char *input_file,
205                   const char *output_file)
206 {
207     std::unique_ptr<gmx_mtop_t> mtop;
208
209     if (fn2bTPX(tps_file) &&
210         efTNG != fn2ftp(input_file) &&
211         efTNG == fn2ftp(output_file))
212     {
213         int temp_natoms = -1;
214         mtop = std::make_unique<gmx_mtop_t>();
215         read_tpx(tps_file, nullptr, nullptr, &temp_natoms,
216                  nullptr, nullptr, mtop.get());
217     }
218
219     return mtop;
220 }
221
222 int gmx_trjconv(int argc, char *argv[])
223 {
224     const char *desc[] = {
225         "[THISMODULE] can convert trajectory files in many ways:",
226         "",
227         "* from one format to another",
228         "* select a subset of atoms",
229         "* change the periodicity representation",
230         "* keep multimeric molecules together",
231         "* center atoms in the box",
232         "* fit atoms to reference structure",
233         "* reduce the number of frames",
234         "* change the timestamps of the frames ([TT]-t0[tt] and [TT]-timestep[tt])",
235         "* select frames within a certain range of a quantity given",
236         "  in an [REF].xvg[ref] file.",
237         "",
238         "The option to write subtrajectories (-sub) based on the information obtained from",
239         "cluster analysis has been removed from [THISMODULE] and is now part of",
240         "[gmx extract-cluster]",
241         "",
242         "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
243         "[PAR]",
244
245         "The following formats are supported for input and output:",
246         "[REF].xtc[ref], [REF].trr[ref], [REF].gro[ref], [TT].g96[tt]",
247         "and [REF].pdb[ref].",
248         "The file formats are detected from the file extension.",
249         "The precision of the [REF].xtc[ref] output is taken from the",
250         "input file for [REF].xtc[ref], [REF].gro[ref] and [REF].pdb[ref],",
251         "and from the [TT]-ndec[tt] option for other input formats. The precision",
252         "is always taken from [TT]-ndec[tt], when this option is set.",
253         "All other formats have fixed precision. [REF].trr[ref]",
254         "output can be single or double precision, depending on the precision",
255         "of the [THISMODULE] binary.",
256         "Note that velocities are only supported in",
257         "[REF].trr[ref], [REF].gro[ref] and [TT].g96[tt] files.[PAR]",
258
259         "Option [TT]-sep[tt] can be used to write every frame to a separate",
260         "[TT].gro, .g96[tt] or [REF].pdb[ref] file. By default, all frames all written to one file.",
261         "[REF].pdb[ref] files with all frames concatenated can be viewed with",
262         "[TT]rasmol -nmrpdb[tt].[PAR]",
263
264         "It is possible to select part of your trajectory and write it out",
265         "to a new trajectory file in order to save disk space, e.g. for leaving",
266         "out the water from a trajectory of a protein in water.",
267         "[BB]ALWAYS[bb] put the original trajectory on tape!",
268         "We recommend to use the portable [REF].xtc[ref] format for your analysis",
269         "to save disk space and to have portable files.[PAR]",
270
271         "There are two options for fitting the trajectory to a reference",
272         "either for essential dynamics analysis, etc.",
273         "The first option is just plain fitting to a reference structure",
274         "in the structure file. The second option is a progressive fit",
275         "in which the first timeframe is fitted to the reference structure ",
276         "in the structure file to obtain and each subsequent timeframe is ",
277         "fitted to the previously fitted structure. This way a continuous",
278         "trajectory is generated, which might not be the case when using the",
279         "regular fit method, e.g. when your protein undergoes large",
280         "conformational transitions.[PAR]",
281
282         "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
283         "treatment:",
284         "",
285         " * [TT]mol[tt] puts the center of mass of molecules in the box,",
286         "   and requires a run input file to be supplied with [TT]-s[tt].",
287         " * [TT]res[tt] puts the center of mass of residues in the box.",
288         " * [TT]atom[tt] puts all the atoms in the box.",
289         " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
290         "   them back. This has the effect that all molecules",
291         "   will remain whole (provided they were whole in the initial",
292         "   conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
293         "   molecules may diffuse out of the box. The starting configuration",
294         "   for this procedure is taken from the structure file, if one is",
295         "   supplied, otherwise it is the first frame.",
296         " * [TT]cluster[tt] clusters all the atoms in the selected index",
297         "   such that they are all closest to the center of mass of the cluster,",
298         "   which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
299         "   results if you in fact have a cluster. Luckily that can be checked",
300         "   afterwards using a trajectory viewer. Note also that if your molecules",
301         "   are broken this will not work either.",
302         " * [TT]whole[tt] only makes broken molecules whole.",
303         "",
304
305         "Option [TT]-ur[tt] sets the unit cell representation for options",
306         "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
307         "All three options give different results for triclinic boxes and",
308         "identical results for rectangular boxes.",
309         "[TT]rect[tt] is the ordinary brick shape.",
310         "[TT]tric[tt] is the triclinic unit cell.",
311         "[TT]compact[tt] puts all atoms at the closest distance from the center",
312         "of the box. This can be useful for visualizing e.g. truncated octahedra",
313         "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
314         "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
315         "is set differently.[PAR]",
316
317         "Option [TT]-center[tt] centers the system in the box. The user can",
318         "select the group which is used to determine the geometrical center.",
319         "Option [TT]-boxcenter[tt] sets the location of the center of the box",
320         "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
321         "[TT]tric[tt]: half of the sum of the box vectors,",
322         "[TT]rect[tt]: half of the box diagonal,",
323         "[TT]zero[tt]: zero.",
324         "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
325         "want all molecules in the box after the centering.[PAR]",
326
327         "Option [TT]-box[tt] sets the size of the new box. This option only works",
328         "for leading dimensions and is thus generally only useful for rectangular boxes.",
329         "If you want to modify only some of the dimensions, e.g. when reading from",
330         "a trajectory, you can use -1 for those dimensions that should stay the same",
331
332         "It is not always possible to use combinations of [TT]-pbc[tt],",
333         "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
334         "you want in one call to [THISMODULE]. Consider using multiple",
335         "calls, and check out the GROMACS website for suggestions.[PAR]",
336
337         "With [TT]-dt[tt], it is possible to reduce the number of ",
338         "frames in the output. This option relies on the accuracy of the times",
339         "in your input trajectory, so if these are inaccurate use the",
340         "[TT]-timestep[tt] option to modify the time (this can be done",
341         "simultaneously). For making smooth movies, the program [gmx-filter]",
342         "can reduce the number of frames while using low-pass frequency",
343         "filtering, this reduces aliasing of high frequency motions.[PAR]",
344
345         "Using [TT]-trunc[tt] [THISMODULE] can truncate [REF].trr[ref] in place, i.e.",
346         "without copying the file. This is useful when a run has crashed",
347         "during disk I/O (i.e. full disk), or when two contiguous",
348         "trajectories must be concatenated without having double frames.[PAR]",
349
350         "Option [TT]-dump[tt] can be used to extract a frame at or near",
351         "one specific time from your trajectory, but only works reliably",
352         "if the time interval between frames is uniform.[PAR]",
353
354         "Option [TT]-drop[tt] reads an [REF].xvg[ref] file with times and values.",
355         "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
356         "frames with a value below and above the value of the respective options",
357         "will not be written."
358     };
359
360     int         pbc_enum;
361     enum
362     {
363         epSel,
364         epNone,
365         epComMol,
366         epComRes,
367         epComAtom,
368         epNojump,
369         epCluster,
370         epWhole,
371         epNR
372     };
373     const char *pbc_opt[epNR + 1] =
374     {
375         nullptr, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
376         nullptr
377     };
378
379     int         unitcell_enum;
380     const char *unitcell_opt[euNR+1] =
381     { nullptr, "rect", "tric", "compact", nullptr };
382
383     enum
384     {
385         ecSel, ecTric, ecRect, ecZero, ecNR
386     };
387     const char *center_opt[ecNR+1] =
388     { nullptr, "tric", "rect", "zero", nullptr };
389     int         ecenter;
390
391     int         fit_enum;
392     enum
393     {
394         efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
395     };
396     const char *fit[efNR + 1] =
397     {
398         nullptr, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
399         "progressive", nullptr
400     };
401
402     static gmx_bool  bSeparate     = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
403     static gmx_bool  bCenter       = FALSE;
404     static int       skip_nr       = 1, ndec = 3, nzero = 0;
405     static real      tzero         = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
406     static rvec      newbox        = {0, 0, 0}, shift = {0, 0, 0}, trans = {0, 0, 0};
407     static char     *exec_command  = nullptr;
408     static real      dropunder     = 0, dropover = 0;
409     static gmx_bool  bRound        = FALSE;
410
411     t_pargs
412         pa[] =
413     {
414         { "-skip", FALSE, etINT,
415           { &skip_nr }, "Only write every nr-th frame" },
416         { "-dt", FALSE, etTIME,
417           { &delta_t },
418           "Only write frame when t MOD dt = first time (%t)" },
419         { "-round", FALSE, etBOOL,
420           { &bRound }, "Round measurements to nearest picosecond"},
421         { "-dump", FALSE, etTIME,
422           { &tdump }, "Dump frame nearest specified time (%t)" },
423         { "-t0", FALSE, etTIME,
424           { &tzero },
425           "Starting time (%t) (default: don't change)" },
426         { "-timestep", FALSE, etTIME,
427           { &timestep },
428           "Change time step between input frames (%t)" },
429         { "-pbc", FALSE, etENUM,
430           { pbc_opt },
431           "PBC treatment (see help text for full description)" },
432         { "-ur", FALSE, etENUM,
433           { unitcell_opt }, "Unit-cell representation" },
434         { "-center", FALSE, etBOOL,
435           { &bCenter }, "Center atoms in box" },
436         { "-boxcenter", FALSE, etENUM,
437           { center_opt }, "Center for -pbc and -center" },
438         { "-box", FALSE, etRVEC,
439           { newbox },
440           "Size for new cubic box (default: read from input)" },
441         { "-trans", FALSE, etRVEC,
442           { trans },
443           "All coordinates will be translated by trans. This "
444           "can advantageously be combined with -pbc mol -ur "
445           "compact." },
446         { "-shift", FALSE, etRVEC,
447           { shift },
448           "All coordinates will be shifted by framenr*shift" },
449         { "-fit", FALSE, etENUM,
450           { fit },
451           "Fit molecule to ref structure in the structure file" },
452         { "-ndec", FALSE, etINT,
453           { &ndec },
454           "Number of decimal places to write to .xtc output" },
455         { "-vel", FALSE, etBOOL,
456           { &bVels }, "Read and write velocities if possible" },
457         { "-force", FALSE, etBOOL,
458           { &bForce }, "Read and write forces if possible" },
459         { "-trunc", FALSE, etTIME,
460           { &ttrunc },
461           "Truncate input trajectory file after this time (%t)" },
462         { "-exec", FALSE, etSTR,
463           { &exec_command },
464           "Execute command for every output frame with the "
465           "frame number as argument" },
466         { "-split", FALSE, etTIME,
467           { &split_t },
468           "Start writing new file when t MOD split = first "
469           "time (%t)" },
470         { "-sep", FALSE, etBOOL,
471           { &bSeparate },
472           "Write each frame to a separate .gro, .g96 or .pdb "
473           "file" },
474         { "-nzero", FALSE, etINT,
475           { &nzero },
476           "If the -sep flag is set, use these many digits "
477           "for the file numbers and prepend zeros as needed" },
478         { "-dropunder", FALSE, etREAL,
479           { &dropunder }, "Drop all frames below this value" },
480         { "-dropover", FALSE, etREAL,
481           { &dropover }, "Drop all frames above this value" },
482         { "-conect", FALSE, etBOOL,
483           { &bCONECT },
484           "Add conect records when writing [REF].pdb[ref] files. Useful "
485           "for visualization of non-standard molecules, e.g. "
486           "coarse grained ones" }
487     };
488 #define NPA asize(pa)
489
490     FILE             *out    = nullptr;
491     t_trxstatus      *trxout = nullptr;
492     t_trxstatus      *trxin;
493     int               file_nr;
494     t_trxframe        fr, frout;
495     int               flags;
496     rvec             *xmem  = nullptr, *vmem = nullptr, *fmem = nullptr;
497     rvec             *xp    = nullptr, x_shift, hbox;
498     real             *w_rls = nullptr;
499     int               m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
500 #define SKIP 10
501     t_topology       *top   = nullptr;
502     gmx_conect        gc    = nullptr;
503     int               ePBC  = -1;
504     t_atoms          *atoms = nullptr, useatoms;
505     matrix            top_box;
506     int              *index = nullptr, *cindex = nullptr;
507     char             *grpnm = nullptr;
508     int              *frindex, nrfri;
509     char             *frname;
510     int               ifit;
511     int              *ind_fit;
512     char             *gn_fit;
513     int               ndrop           = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
514     double          **dropval;
515     real              tshift = 0, dt = -1, prec;
516     gmx_bool          bFit, bPFit, bReset;
517     int               nfitdim;
518     gmx_rmpbc_t       gpbc = nullptr;
519     gmx_bool          bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
520     gmx_bool          bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
521     gmx_bool          bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetXtcPrec, bNeedPrec;
522     gmx_bool          bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
523     gmx_bool          bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
524     gmx_bool          bWriteFrame, bSplitHere;
525     const char       *top_file, *in_file, *out_file = nullptr;
526     char              out_file2[256], *charpt;
527     char             *outf_base = nullptr;
528     const char       *outf_ext  = nullptr;
529     char              top_title[256], timestr[32], stepstr[32], filemode[5];
530     gmx_output_env_t *oenv;
531
532     t_filenm          fnm[] = {
533         { efTRX, "-f",   nullptr,      ffREAD  },
534         { efTRO, "-o",   nullptr,      ffWRITE },
535         { efTPS, nullptr,   nullptr,      ffOPTRD },
536         { efNDX, nullptr,   nullptr,      ffOPTRD },
537         { efNDX, "-fr",  "frames",  ffOPTRD },
538         { efNDX, "-sub", "cluster", ffOPTRD },
539         { efXVG, "-drop", "drop",    ffOPTRD }
540     };
541 #define NFILE asize(fnm)
542
543     if (!parse_common_args(&argc, argv,
544                            PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW |
545                            PCA_TIME_UNIT,
546                            NFILE, fnm, NPA, pa, asize(desc), desc,
547                            0, nullptr, &oenv))
548     {
549         return 0;
550     }
551     fprintf(stdout, "Note that major changes are planned in future for "
552             "trjconv, to improve usability and utility.\n");
553
554     top_file = ftp2fn(efTPS, NFILE, fnm);
555
556     /* Check command line */
557     in_file = opt2fn("-f", NFILE, fnm);
558
559     if (ttrunc != -1)
560     {
561         do_trunc(in_file, ttrunc);
562     }
563     else
564     {
565         /* mark active cmdline options */
566         bSetBox     = opt2parg_bSet("-box", NPA, pa);
567         bSetTime    = opt2parg_bSet("-t0", NPA, pa);
568         bSetXtcPrec = opt2parg_bSet("-ndec", NPA, pa);
569         bSetUR      = opt2parg_bSet("-ur", NPA, pa);
570         bExec       = opt2parg_bSet("-exec", NPA, pa);
571         bTimeStep   = opt2parg_bSet("-timestep", NPA, pa);
572         bTDump      = opt2parg_bSet("-dump", NPA, pa);
573         bDropUnder  = opt2parg_bSet("-dropunder", NPA, pa);
574         bDropOver   = opt2parg_bSet("-dropover", NPA, pa);
575         bTrans      = opt2parg_bSet("-trans", NPA, pa);
576         bSplit      = (split_t != 0);
577
578         /* parse enum options */
579         fit_enum      = nenum(fit);
580         bFit          = (fit_enum == efFit || fit_enum == efFitXY);
581         bReset        = (fit_enum == efReset || fit_enum == efResetXY);
582         bPFit         = fit_enum == efPFit;
583         pbc_enum      = nenum(pbc_opt);
584         bPBCWhole     = pbc_enum == epWhole;
585         bPBCcomRes    = pbc_enum == epComRes;
586         bPBCcomMol    = pbc_enum == epComMol;
587         bPBCcomAtom   = pbc_enum == epComAtom;
588         bNoJump       = pbc_enum == epNojump;
589         bCluster      = pbc_enum == epCluster;
590         bPBC          = pbc_enum != epNone;
591         unitcell_enum = nenum(unitcell_opt);
592         ecenter       = nenum(center_opt) - ecTric;
593
594         /* set and check option dependencies */
595         if (bPFit)
596         {
597             bFit = TRUE;        /* for pfit, fit *must* be set */
598         }
599         if (bFit)
600         {
601             bReset = TRUE;       /* for fit, reset *must* be set */
602         }
603         nfitdim = 0;
604         if (bFit || bReset)
605         {
606             nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
607         }
608         bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
609
610         if (bSetUR)
611         {
612             if (!(bPBCcomRes || bPBCcomMol ||  bPBCcomAtom))
613             {
614                 fprintf(stderr,
615                         "WARNING: Option for unitcell representation (-ur %s)\n"
616                         "         only has effect in combination with -pbc %s, %s or %s.\n"
617                         "         Ingoring unitcell representation.\n\n",
618                         unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
619             }
620         }
621         if (bFit && bPBC)
622         {
623             gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
624                       "Please do the PBC condition treatment first and then run trjconv in a second step\n"
625                       "for the rotational fit.\n"
626                       "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
627                       "results!");
628         }
629
630         /* ndec for XTC writing is in nr of decimal places, prec is a multiplication factor: */
631         prec = 1;
632         for (i = 0; i < ndec; i++)
633         {
634             prec *= 10;
635         }
636
637         bIndex = ftp2bSet(efNDX, NFILE, fnm);
638
639
640         /* Determine output type */
641         out_file = opt2fn("-o", NFILE, fnm);
642         int ftp  = fn2ftp(out_file);
643         fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
644         bNeedPrec = (ftp == efXTC);
645         int ftpin = fn2ftp(in_file);
646         if (bVels)
647         {
648             /* check if velocities are possible in input and output files */
649             bVels = (ftp == efTRR || ftp == efGRO ||
650                      ftp == efG96 || ftp == efTNG)
651                 && (ftpin == efTRR || ftpin == efGRO ||
652                     ftpin == efG96 || ftpin == efTNG || ftpin == efCPT);
653         }
654         if (bSeparate || bSplit)
655         {
656             outf_ext = std::strrchr(out_file, '.');
657             if (outf_ext == nullptr)
658             {
659                 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
660             }
661             outf_base = gmx_strdup(out_file);
662             outf_base[outf_ext - out_file] = '\0';
663         }
664
665         bool bSubTraj = opt2bSet("-sub", NFILE, fnm);
666         if (bSubTraj)
667         {
668             gmx_fatal(FARGS, "The -sub option has been removed from gmx trjconv and is now part\n"
669                       "of gmx extract-cluster and does nothing here\n");
670         }
671
672         /* skipping */
673         if (skip_nr <= 0)
674         {
675             gmx_fatal(FARGS, "Argument for -skip (%d) needs to be greater or equal to 1.", skip_nr);
676         }
677
678         std::unique_ptr<gmx_mtop_t> mtop = read_mtop_for_tng(top_file, in_file, out_file);
679
680         /* Determine whether to read a topology */
681         bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
682                 bRmPBC || bReset || bPBCcomMol || bCluster ||
683                 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
684
685         /* Determine if when can read index groups */
686         bIndex = (bIndex || bTPS);
687
688         if (bTPS)
689         {
690             snew(top, 1);
691             read_tps_conf(top_file, top, &ePBC, &xp, nullptr, top_box,
692                           bReset || bPBCcomRes);
693             std::strncpy(top_title, *top->name, 255);
694             top_title[255] = '\0';
695             atoms          = &top->atoms;
696
697             if (0 == top->mols.nr && (bCluster || bPBCcomMol))
698             {
699                 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
700             }
701
702             /* top_title is only used for gro and pdb,
703              * the header in such a file is top_title, followed by
704              * t= ... and/or step= ...
705              * to prevent double t= or step=, remove it from top_title.
706              * From GROMACS-2018 we only write t/step when the frame actually
707              * has a valid time/step, so we need to check for both separately.
708              */
709             if ((charpt = std::strstr(top_title, " t= ")))
710             {
711                 charpt[0] = '\0';
712             }
713             if ((charpt = std::strstr(top_title, " step= ")))
714             {
715                 charpt[0] = '\0';
716             }
717
718             if (bCONECT)
719             {
720                 gc = gmx_conect_generate(top);
721             }
722             if (bRmPBC)
723             {
724                 gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
725             }
726         }
727
728         /* get frame number index */
729         frindex = nullptr;
730         if (opt2bSet("-fr", NFILE, fnm))
731         {
732             printf("Select groups of frame number indices:\n");
733             rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, &frindex, &frname);
734             if (debug)
735             {
736                 for (i = 0; i < nrfri; i++)
737                 {
738                     fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
739                 }
740             }
741         }
742
743         /* get index groups etc. */
744         if (bReset)
745         {
746             printf("Select group for %s fit\n",
747                    bFit ? "least squares" : "translational");
748             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
749                       1, &ifit, &ind_fit, &gn_fit);
750
751             if (bFit)
752             {
753                 if (ifit < 2)
754                 {
755                     gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
756                 }
757                 else if (ifit == 3)
758                 {
759                     fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
760                 }
761             }
762         }
763         else if (bCluster)
764         {
765             printf("Select group for clustering\n");
766             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
767                       1, &ifit, &ind_fit, &gn_fit);
768         }
769
770         if (bIndex)
771         {
772             if (bCenter)
773             {
774                 printf("Select group for centering\n");
775                 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
776                           1, &ncent, &cindex, &grpnm);
777             }
778             printf("Select group for output\n");
779             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
780                       1, &nout, &index, &grpnm);
781         }
782         else
783         {
784             /* no index file, so read natoms from TRX */
785             if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
786             {
787                 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
788             }
789             natoms = fr.natoms;
790             close_trx(trxin);
791             sfree(fr.x);
792             snew(index, natoms);
793             for (i = 0; i < natoms; i++)
794             {
795                 index[i] = i;
796             }
797             nout = natoms;
798             if (bCenter)
799             {
800                 ncent  = nout;
801                 cindex = index;
802             }
803         }
804
805         if (bReset)
806         {
807             snew(w_rls, atoms->nr);
808             for (i = 0; (i < ifit); i++)
809             {
810                 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
811             }
812
813             /* Restore reference structure and set to origin,
814                store original location (to put structure back) */
815             if (bRmPBC)
816             {
817                 gmx_rmpbc(gpbc, top->atoms.nr, top_box, xp);
818             }
819             copy_rvec(xp[index[0]], x_shift);
820             reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, nullptr, xp, w_rls);
821             rvec_dec(x_shift, xp[index[0]]);
822         }
823         else
824         {
825             clear_rvec(x_shift);
826         }
827
828         if (bDropUnder || bDropOver)
829         {
830             /* Read the .xvg file with the drop values */
831             fprintf(stderr, "\nReading drop file ...");
832             ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
833             fprintf(stderr, " %d time points\n", ndrop);
834             if (ndrop == 0 || ncol < 2)
835             {
836                 gmx_fatal(FARGS, "Found no data points in %s",
837                           opt2fn("-drop", NFILE, fnm));
838             }
839             drop0 = 0;
840             drop1 = 0;
841         }
842
843         /* Make atoms struct for output in GRO or PDB files */
844         if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
845         {
846             /* get memory for stuff to go in .pdb file, and initialize
847              * the pdbinfo structure part if the input has it.
848              */
849             init_t_atoms(&useatoms, atoms->nr, atoms->havePdbInfo);
850             sfree(useatoms.resinfo);
851             useatoms.resinfo = atoms->resinfo;
852             for (i = 0; (i < nout); i++)
853             {
854                 useatoms.atomname[i] = atoms->atomname[index[i]];
855                 useatoms.atom[i]     = atoms->atom[index[i]];
856                 if (atoms->havePdbInfo)
857                 {
858                     useatoms.pdbinfo[i]  = atoms->pdbinfo[index[i]];
859                 }
860                 useatoms.nres        = std::max(useatoms.nres, useatoms.atom[i].resind+1);
861             }
862             useatoms.nr = nout;
863         }
864         /* select what to read */
865         if (ftp == efTRR)
866         {
867             flags = TRX_READ_X;
868         }
869         else
870         {
871             flags = TRX_NEED_X;
872         }
873         if (bVels)
874         {
875             flags = flags | TRX_READ_V;
876         }
877         if (bForce)
878         {
879             flags = flags | TRX_READ_F;
880         }
881
882         /* open trx file for reading */
883         bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
884         if (fr.bPrec)
885         {
886             fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
887         }
888         if (bNeedPrec)
889         {
890             if (bSetXtcPrec || !fr.bPrec)
891             {
892                 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
893             }
894             else
895             {
896                 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
897             }
898         }
899
900         if (bHaveFirstFrame)
901         {
902             if (bTDump)
903             {
904                 // Determine timestep (assuming constant spacing for now) if we
905                 // need to dump frames based on time. This is required so we do not
906                 // skip the first frame if that was the one that should have been dumped
907                 double firstFrameTime = fr.time;
908                 if (read_next_frame(oenv, trxin, &fr))
909                 {
910                     dt     = fr.time - firstFrameTime;
911                     bDTset = TRUE;
912                     if (dt <= 0)
913                     {
914                         fprintf(stderr, "Warning: Frame times are not incrementing - will dump first frame.\n");
915                     }
916                 }
917                 // Now close and reopen so we are at first frame again
918                 close_trx(trxin);
919                 done_frame(&fr);
920                 // Reopen at first frame (We already know it exists if we got here)
921                 read_first_frame(oenv, &trxin, in_file, &fr, flags);
922             }
923
924             set_trxframe_ePBC(&fr, ePBC);
925             natoms = fr.natoms;
926
927             if (bSetTime)
928             {
929                 tshift = tzero-fr.time;
930             }
931             else
932             {
933                 tzero = fr.time;
934             }
935
936             bCopy = FALSE;
937             if (bIndex)
938             {
939                 /* check if index is meaningful */
940                 for (i = 0; i < nout; i++)
941                 {
942                     if (index[i] >= natoms)
943                     {
944                         gmx_fatal(FARGS,
945                                   "Index[%d] %d is larger than the number of atoms in the\n"
946                                   "trajectory file (%d). There is a mismatch in the contents\n"
947                                   "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
948                     }
949                     bCopy = bCopy || (i != index[i]);
950                 }
951             }
952
953             /* open output for writing */
954             std::strcpy(filemode, "w");
955             switch (ftp)
956             {
957                 case efTNG:
958                     trxout = trjtools_gmx_prepare_tng_writing(out_file,
959                                                               filemode[0],
960                                                               trxin,
961                                                               nullptr,
962                                                               nout,
963                                                               mtop.get(),
964                                                               gmx::arrayRefFromArray(index, nout),
965                                                               grpnm);
966                     break;
967                 case efXTC:
968                 case efTRR:
969                     out = nullptr;
970                     if (!bSplit)
971                     {
972                         trxout = open_trx(out_file, filemode);
973                     }
974                     break;
975                 case efGRO:
976                 case efG96:
977                 case efPDB:
978                     if (!bSeparate && !bSplit)
979                     {
980                         out = gmx_ffopen(out_file, filemode);
981                     }
982                     break;
983                 default:
984                     gmx_incons("Illegal output file format");
985             }
986
987             if (bCopy)
988             {
989                 snew(xmem, nout);
990                 if (bVels)
991                 {
992                     snew(vmem, nout);
993                 }
994                 if (bForce)
995                 {
996                     snew(fmem, nout);
997                 }
998             }
999
1000             /* Start the big loop over frames */
1001             file_nr  =  0;
1002             frame    =  0;
1003             outframe =  0;
1004             model_nr =  0;
1005
1006             /* Main loop over frames */
1007             do
1008             {
1009                 if (!fr.bStep)
1010                 {
1011                     /* set the step */
1012                     fr.step = newstep;
1013                     newstep++;
1014                 }
1015
1016                 if (bSetBox)
1017                 {
1018                     /* generate new box */
1019                     if (!fr.bBox)
1020                     {
1021                         clear_mat(fr.box);
1022                     }
1023                     for (m = 0; m < DIM; m++)
1024                     {
1025                         if (newbox[m] >= 0)
1026                         {
1027                             fr.box[m][m] = newbox[m];
1028                         }
1029                         else
1030                         {
1031                             if (!fr.bBox)
1032                             {
1033                                 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1034                             }
1035                         }
1036                     }
1037                 }
1038
1039                 if (bTrans)
1040                 {
1041                     for (i = 0; i < natoms; i++)
1042                     {
1043                         rvec_inc(fr.x[i], trans);
1044                     }
1045                 }
1046
1047                 if (bTDump)
1048                 {
1049                     // If we could not read two frames or times are not incrementing
1050                     // we have almost no idea what to do,
1051                     // but dump the first frame so output is not broken.
1052                     if (dt <= 0 || !bDTset)
1053                     {
1054                         bDumpFrame = true;
1055                     }
1056                     else
1057                     {
1058                         // Dump the frame if we are less than half a frame time
1059                         // below it. This will also ensure we at least dump a
1060                         // somewhat reasonable frame if the spacing is unequal
1061                         // and we have overrun the frame time. Once we dump one
1062                         // frame based on time we quit, so it does not matter
1063                         // that this might be true for all subsequent frames too.
1064                         bDumpFrame = (fr.time > tdump-0.5*dt);
1065                     }
1066                 }
1067                 else
1068                 {
1069                     bDumpFrame = FALSE;
1070                 }
1071
1072                 /* determine if an atom jumped across the box and reset it if so */
1073                 if (bNoJump && (bTPS || frame != 0))
1074                 {
1075                     for (d = 0; d < DIM; d++)
1076                     {
1077                         hbox[d] = 0.5*fr.box[d][d];
1078                     }
1079                     for (i = 0; i < natoms; i++)
1080                     {
1081                         if (bReset)
1082                         {
1083                             rvec_dec(fr.x[i], x_shift);
1084                         }
1085                         for (m = DIM-1; m >= 0; m--)
1086                         {
1087                             if (hbox[m] > 0)
1088                             {
1089                                 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1090                                 {
1091                                     for (d = 0; d <= m; d++)
1092                                     {
1093                                         fr.x[i][d] += fr.box[m][d];
1094                                     }
1095                                 }
1096                                 while (fr.x[i][m]-xp[i][m] > hbox[m])
1097                                 {
1098                                     for (d = 0; d <= m; d++)
1099                                     {
1100                                         fr.x[i][d] -= fr.box[m][d];
1101                                     }
1102                                 }
1103                             }
1104                         }
1105                     }
1106                 }
1107                 else if (bCluster)
1108                 {
1109                     calc_pbc_cluster(ecenter, ifit, top, ePBC, fr.x, ind_fit, fr.box);
1110                 }
1111
1112                 if (bPFit)
1113                 {
1114                     /* Now modify the coords according to the flags,
1115                        for normal fit, this is only done for output frames */
1116                     if (bRmPBC)
1117                     {
1118                         gmx_rmpbc_trxfr(gpbc, &fr);
1119                     }
1120
1121                     reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1122                     do_fit(natoms, w_rls, xp, fr.x);
1123                 }
1124
1125                 /* store this set of coordinates for future use */
1126                 if (bPFit || bNoJump)
1127                 {
1128                     if (xp == nullptr)
1129                     {
1130                         snew(xp, natoms);
1131                     }
1132                     for (i = 0; (i < natoms); i++)
1133                     {
1134                         copy_rvec(fr.x[i], xp[i]);
1135                         rvec_inc(fr.x[i], x_shift);
1136                     }
1137                 }
1138
1139                 if (frindex)
1140                 {
1141                     /* see if we have a frame from the frame index group */
1142                     for (i = 0; i < nrfri && !bDumpFrame; i++)
1143                     {
1144                         bDumpFrame = frame == frindex[i];
1145                     }
1146                 }
1147                 if (debug && bDumpFrame)
1148                 {
1149                     fprintf(debug, "dumping %d\n", frame);
1150                 }
1151
1152                 bWriteFrame =
1153                     ( ( !bTDump && (frindex == nullptr) && frame % skip_nr == 0 ) || bDumpFrame );
1154
1155                 if (bWriteFrame && (bDropUnder || bDropOver))
1156                 {
1157                     while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1158                     {
1159                         drop0 = drop1;
1160                         drop1++;
1161                     }
1162                     if (std::abs(dropval[0][drop0] - fr.time)
1163                         < std::abs(dropval[0][drop1] - fr.time))
1164                     {
1165                         dropuse = drop0;
1166                     }
1167                     else
1168                     {
1169                         dropuse = drop1;
1170                     }
1171                     if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1172                         (bDropOver && dropval[1][dropuse] > dropover))
1173                     {
1174                         bWriteFrame = FALSE;
1175                     }
1176                 }
1177
1178                 if (bWriteFrame)
1179                 {
1180                     /* We should avoid modifying the input frame,
1181                      * but since here we don't have the output frame yet,
1182                      * we introduce a temporary output frame time variable.
1183                      */
1184                     real frout_time;
1185
1186                     frout_time = fr.time;
1187
1188                     /* calc new time */
1189                     if (bTimeStep)
1190                     {
1191                         frout_time = tzero + frame*timestep;
1192                     }
1193                     else if (bSetTime)
1194                     {
1195                         frout_time += tshift;
1196                     }
1197
1198                     if (bTDump)
1199                     {
1200                         fprintf(stderr, "\nDumping frame at t= %g %s\n",
1201                                 output_env_conv_time(oenv, frout_time), output_env_get_time_unit(oenv).c_str());
1202                     }
1203
1204                     /* check for writing at each delta_t */
1205                     bDoIt = (delta_t == 0);
1206                     if (!bDoIt)
1207                     {
1208                         if (!bRound)
1209                         {
1210                             bDoIt = bRmod(frout_time, tzero, delta_t);
1211                         }
1212                         else
1213                         {
1214                             /* round() is not C89 compatible, so we do this:  */
1215                             bDoIt = bRmod(std::floor(frout_time+0.5), std::floor(tzero+0.5),
1216                                           std::floor(delta_t+0.5));
1217                         }
1218                     }
1219
1220                     if (bDoIt || bTDump)
1221                     {
1222                         /* print sometimes */
1223                         if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1224                         {
1225                             fprintf(stderr, " ->  frame %6d time %8.3f      \r",
1226                                     outframe, output_env_conv_time(oenv, frout_time));
1227                             fflush(stderr);
1228                         }
1229
1230                         if (!bPFit)
1231                         {
1232                             /* Now modify the coords according to the flags,
1233                                for PFit we did this already! */
1234
1235                             if (bRmPBC)
1236                             {
1237                                 gmx_rmpbc_trxfr(gpbc, &fr);
1238                             }
1239
1240                             if (bReset)
1241                             {
1242                                 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1243                                 if (bFit)
1244                                 {
1245                                     do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1246                                 }
1247                                 if (!bCenter)
1248                                 {
1249                                     for (i = 0; i < natoms; i++)
1250                                     {
1251                                         rvec_inc(fr.x[i], x_shift);
1252                                     }
1253                                 }
1254                             }
1255
1256                             if (bCenter)
1257                             {
1258                                 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1259                             }
1260                         }
1261
1262                         auto positionsArrayRef = gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec *>(fr.x), natoms);
1263                         if (bPBCcomAtom)
1264                         {
1265                             switch (unitcell_enum)
1266                             {
1267                                 case euRect:
1268                                     put_atoms_in_box(ePBC, fr.box, positionsArrayRef);
1269                                     break;
1270                                 case euTric:
1271                                     put_atoms_in_triclinic_unitcell(ecenter, fr.box, positionsArrayRef);
1272                                     break;
1273                                 case euCompact:
1274                                     put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1275                                                                   positionsArrayRef);
1276                                     break;
1277                             }
1278                         }
1279                         if (bPBCcomRes)
1280                         {
1281                             put_residue_com_in_box(unitcell_enum, ecenter,
1282                                                    natoms, atoms->atom, ePBC, fr.box, fr.x);
1283                         }
1284                         if (bPBCcomMol)
1285                         {
1286                             put_molecule_com_in_box(unitcell_enum, ecenter,
1287                                                     &top->mols,
1288                                                     natoms, atoms->atom, ePBC, fr.box, fr.x);
1289                         }
1290                         /* Copy the input trxframe struct to the output trxframe struct */
1291                         frout        = fr;
1292                         frout.time   = frout_time;
1293                         frout.bV     = (frout.bV && bVels);
1294                         frout.bF     = (frout.bF && bForce);
1295                         frout.natoms = nout;
1296                         if (bNeedPrec && (bSetXtcPrec || !fr.bPrec))
1297                         {
1298                             frout.bPrec = TRUE;
1299                             frout.prec  = prec;
1300                         }
1301                         if (bCopy)
1302                         {
1303                             frout.x = xmem;
1304                             if (frout.bV)
1305                             {
1306                                 frout.v = vmem;
1307                             }
1308                             if (frout.bF)
1309                             {
1310                                 frout.f = fmem;
1311                             }
1312                             for (i = 0; i < nout; i++)
1313                             {
1314                                 copy_rvec(fr.x[index[i]], frout.x[i]);
1315                                 if (frout.bV)
1316                                 {
1317                                     copy_rvec(fr.v[index[i]], frout.v[i]);
1318                                 }
1319                                 if (frout.bF)
1320                                 {
1321                                     copy_rvec(fr.f[index[i]], frout.f[i]);
1322                                 }
1323                             }
1324                         }
1325
1326                         if (opt2parg_bSet("-shift", NPA, pa))
1327                         {
1328                             for (i = 0; i < nout; i++)
1329                             {
1330                                 for (d = 0; d < DIM; d++)
1331                                 {
1332                                     frout.x[i][d] += outframe*shift[d];
1333                                 }
1334                             }
1335                         }
1336
1337                         if (!bRound)
1338                         {
1339                             bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1340                         }
1341                         else
1342                         {
1343                             /* round() is not C89 compatible, so we do this: */
1344                             bSplitHere = bSplit && bRmod(std::floor(frout.time+0.5),
1345                                                          std::floor(tzero+0.5),
1346                                                          std::floor(split_t+0.5));
1347                         }
1348                         if (bSeparate || bSplitHere)
1349                         {
1350                             mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1351                         }
1352
1353                         std::string title;
1354                         switch (ftp)
1355                         {
1356                             case efTNG:
1357                                 write_tng_frame(trxout, &frout);
1358                                 // TODO when trjconv behaves better: work how to read and write lambda
1359                                 break;
1360                             case efTRR:
1361                             case efXTC:
1362                                 if (bSplitHere)
1363                                 {
1364                                     if (trxout)
1365                                     {
1366                                         close_trx(trxout);
1367                                     }
1368                                     trxout = open_trx(out_file2, filemode);
1369                                 }
1370                                 write_trxframe(trxout, &frout, gc);
1371                                 break;
1372                             case efGRO:
1373                             case efG96:
1374                             case efPDB:
1375                                 // Only add a generator statement if title is empty,
1376                                 // to avoid multiple generated-by statements from various programs
1377                                 if (std::strlen(top_title) == 0)
1378                                 {
1379                                     sprintf(top_title, "Generated by trjconv");
1380                                 }
1381                                 if (frout.bTime)
1382                                 {
1383                                     sprintf(timestr, " t= %9.5f", frout.time);
1384                                 }
1385                                 else
1386                                 {
1387                                     std::strcpy(timestr, "");
1388                                 }
1389                                 if (frout.bStep)
1390                                 {
1391                                     sprintf(stepstr, " step= %" PRId64, frout.step);
1392                                 }
1393                                 else
1394                                 {
1395                                     std::strcpy(stepstr, "");
1396                                 }
1397                                 title = gmx::formatString("%s%s%s", top_title, timestr, stepstr);
1398                                 if (bSeparate || bSplitHere)
1399                                 {
1400                                     out = gmx_ffopen(out_file2, "w");
1401                                 }
1402                                 switch (ftp)
1403                                 {
1404                                     case efGRO:
1405                                         write_hconf_p(out, title.c_str(), &useatoms,
1406                                                       frout.x, frout.bV ? frout.v : nullptr, frout.box);
1407                                         break;
1408                                     case efPDB:
1409                                         fprintf(out, "REMARK    GENERATED BY TRJCONV\n");
1410                                         /* if reading from pdb, we want to keep the original
1411                                            model numbering else we write the output frame
1412                                            number plus one, because model 0 is not allowed in pdb */
1413                                         if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1414                                         {
1415                                             model_nr = fr.step;
1416                                         }
1417                                         else
1418                                         {
1419                                             model_nr++;
1420                                         }
1421                                         write_pdbfile(out, title.c_str(), &useatoms, frout.x,
1422                                                       frout.ePBC, frout.box, ' ', model_nr, gc);
1423                                         break;
1424                                     case efG96:
1425                                         const char *outputTitle = "";
1426                                         if (bSeparate || bTDump)
1427                                         {
1428                                             outputTitle = title.c_str();
1429                                             if (bTPS)
1430                                             {
1431                                                 frout.bAtoms = TRUE;
1432                                             }
1433                                             frout.atoms  = &useatoms;
1434                                             frout.bStep  = FALSE;
1435                                             frout.bTime  = FALSE;
1436                                         }
1437                                         else
1438                                         {
1439                                             if (outframe == 0)
1440                                             {
1441                                                 outputTitle = title.c_str();
1442                                             }
1443                                             frout.bAtoms = FALSE;
1444                                             frout.bStep  = TRUE;
1445                                             frout.bTime  = TRUE;
1446                                         }
1447                                         write_g96_conf(out, outputTitle, &frout, -1, nullptr);
1448                                 }
1449                                 if (bSeparate || bSplitHere)
1450                                 {
1451                                     gmx_ffclose(out);
1452                                     out = nullptr;
1453                                 }
1454                                 break;
1455                             default:
1456                                 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1457                         }
1458                         if (bSeparate || bSplitHere)
1459                         {
1460                             file_nr++;
1461                         }
1462
1463                         /* execute command */
1464                         if (bExec)
1465                         {
1466                             char c[255];
1467                             sprintf(c, "%s  %d", exec_command, file_nr-1);
1468                             /*fprintf(stderr,"Executing '%s'\n",c);*/
1469                             if (0 != system(c))
1470                             {
1471                                 gmx_fatal(FARGS, "Error executing command: %s", c);
1472                             }
1473                         }
1474                         outframe++;
1475                     }
1476                 }
1477                 frame++;
1478                 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1479             }
1480             while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1481         }
1482
1483         if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1484         {
1485             fprintf(stderr, "\nWARNING no output, "
1486                     "last frame read at t=%g\n", fr.time);
1487         }
1488         fprintf(stderr, "\n");
1489
1490         close_trx(trxin);
1491         sfree(outf_base);
1492
1493         if (bRmPBC)
1494         {
1495             gmx_rmpbc_done(gpbc);
1496         }
1497
1498         if (trxout)
1499         {
1500             close_trx(trxout);
1501         }
1502         else if (out != nullptr)
1503         {
1504             gmx_ffclose(out);
1505         }
1506     }
1507
1508     if (bTPS)
1509     {
1510         done_top(top);
1511         sfree(top);
1512     }
1513     sfree(xp);
1514     sfree(xmem);
1515     sfree(vmem);
1516     sfree(fmem);
1517     sfree(grpnm);
1518     sfree(index);
1519     sfree(cindex);
1520     done_frame(&fr);
1521
1522     do_view(oenv, out_file, nullptr);
1523
1524     output_env_done(oenv);
1525     return 0;
1526 }