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