Merge 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         "* cut the trajectory in small subtrajectories according",
236         "  to information in an index file. This allows subsequent analysis of",
237         "  the subtrajectories that could, for example, be the result of a",
238         "  cluster analysis. Use option [TT]-sub[tt].",
239         "  This assumes that the entries in the index file are frame numbers and",
240         "  dumps each group in the index file to a separate trajectory file.",
241         "* select frames within a certain range of a quantity given",
242         "  in an [REF].xvg[ref] file.",
243         "",
244         "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
245         "[PAR]",
246
247         "The following formats are supported for input and output:",
248         "[REF].xtc[ref], [REF].trr[ref], [REF].gro[ref], [TT].g96[tt]",
249         "and [REF].pdb[ref].",
250         "The file formats are detected from the file extension.",
251         "The precision of the [REF].xtc[ref] output is taken from the",
252         "input file for [REF].xtc[ref], [REF].gro[ref] and [REF].pdb[ref],",
253         "and from the [TT]-ndec[tt] option for other input formats. The precision",
254         "is always taken from [TT]-ndec[tt], when this option is set.",
255         "All other formats have fixed precision. [REF].trr[ref]",
256         "output can be single or double precision, depending on the precision",
257         "of the [THISMODULE] binary.",
258         "Note that velocities are only supported in",
259         "[REF].trr[ref], [REF].gro[ref] and [TT].g96[tt] files.[PAR]",
260
261         "Option [TT]-sep[tt] can be used to write every frame to a separate",
262         "[TT].gro, .g96[tt] or [REF].pdb[ref] file. By default, all frames all written to one file.",
263         "[REF].pdb[ref] files with all frames concatenated can be viewed with",
264         "[TT]rasmol -nmrpdb[tt].[PAR]",
265
266         "It is possible to select part of your trajectory and write it out",
267         "to a new trajectory file in order to save disk space, e.g. for leaving",
268         "out the water from a trajectory of a protein in water.",
269         "[BB]ALWAYS[bb] put the original trajectory on tape!",
270         "We recommend to use the portable [REF].xtc[ref] format for your analysis",
271         "to save disk space and to have portable files.[PAR]",
272
273         "There are two options for fitting the trajectory to a reference",
274         "either for essential dynamics analysis, etc.",
275         "The first option is just plain fitting to a reference structure",
276         "in the structure file. The second option is a progressive fit",
277         "in which the first timeframe is fitted to the reference structure ",
278         "in the structure file to obtain and each subsequent timeframe is ",
279         "fitted to the previously fitted structure. This way a continuous",
280         "trajectory is generated, which might not be the case when using the",
281         "regular fit method, e.g. when your protein undergoes large",
282         "conformational transitions.[PAR]",
283
284         "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
285         "treatment:",
286         "",
287         " * [TT]mol[tt] puts the center of mass of molecules in the box,",
288         "   and requires a run input file to be supplied with [TT]-s[tt].",
289         " * [TT]res[tt] puts the center of mass of residues in the box.",
290         " * [TT]atom[tt] puts all the atoms in the box.",
291         " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
292         "   them back. This has the effect that all molecules",
293         "   will remain whole (provided they were whole in the initial",
294         "   conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
295         "   molecules may diffuse out of the box. The starting configuration",
296         "   for this procedure is taken from the structure file, if one is",
297         "   supplied, otherwise it is the first frame.",
298         " * [TT]cluster[tt] clusters all the atoms in the selected index",
299         "   such that they are all closest to the center of mass of the cluster,",
300         "   which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
301         "   results if you in fact have a cluster. Luckily that can be checked",
302         "   afterwards using a trajectory viewer. Note also that if your molecules",
303         "   are broken this will not work either.",
304         " * [TT]whole[tt] only makes broken molecules whole.",
305         "",
306
307         "Option [TT]-ur[tt] sets the unit cell representation for options",
308         "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
309         "All three options give different results for triclinic boxes and",
310         "identical results for rectangular boxes.",
311         "[TT]rect[tt] is the ordinary brick shape.",
312         "[TT]tric[tt] is the triclinic unit cell.",
313         "[TT]compact[tt] puts all atoms at the closest distance from the center",
314         "of the box. This can be useful for visualizing e.g. truncated octahedra",
315         "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
316         "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
317         "is set differently.[PAR]",
318
319         "Option [TT]-center[tt] centers the system in the box. The user can",
320         "select the group which is used to determine the geometrical center.",
321         "Option [TT]-boxcenter[tt] sets the location of the center of the box",
322         "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
323         "[TT]tric[tt]: half of the sum of the box vectors,",
324         "[TT]rect[tt]: half of the box diagonal,",
325         "[TT]zero[tt]: zero.",
326         "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
327         "want all molecules in the box after the centering.[PAR]",
328
329         "Option [TT]-box[tt] sets the size of the new box. This option only works",
330         "for leading dimensions and is thus generally only useful for rectangular boxes.",
331         "If you want to modify only some of the dimensions, e.g. when reading from",
332         "a trajectory, you can use -1 for those dimensions that should stay the same",
333
334         "It is not always possible to use combinations of [TT]-pbc[tt],",
335         "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
336         "you want in one call to [THISMODULE]. Consider using multiple",
337         "calls, and check out the GROMACS website for suggestions.[PAR]",
338
339         "With [TT]-dt[tt], it is possible to reduce the number of ",
340         "frames in the output. This option relies on the accuracy of the times",
341         "in your input trajectory, so if these are inaccurate use the",
342         "[TT]-timestep[tt] option to modify the time (this can be done",
343         "simultaneously). For making smooth movies, the program [gmx-filter]",
344         "can reduce the number of frames while using low-pass frequency",
345         "filtering, this reduces aliasing of high frequency motions.[PAR]",
346
347         "Using [TT]-trunc[tt] [THISMODULE] can truncate [REF].trr[ref] in place, i.e.",
348         "without copying the file. This is useful when a run has crashed",
349         "during disk I/O (i.e. full disk), or when two contiguous",
350         "trajectories must be concatenated without having double frames.[PAR]",
351
352         "Option [TT]-dump[tt] can be used to extract a frame at or near",
353         "one specific time from your trajectory, but only works reliably",
354         "if the time interval between frames is uniform.[PAR]",
355
356         "Option [TT]-drop[tt] reads an [REF].xvg[ref] file with times and values.",
357         "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
358         "frames with a value below and above the value of the respective options",
359         "will not be written."
360     };
361
362     int         pbc_enum;
363     enum
364     {
365         epSel,
366         epNone,
367         epComMol,
368         epComRes,
369         epComAtom,
370         epNojump,
371         epCluster,
372         epWhole,
373         epNR
374     };
375     const char *pbc_opt[epNR + 1] =
376     {
377         nullptr, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
378         nullptr
379     };
380
381     int         unitcell_enum;
382     const char *unitcell_opt[euNR+1] =
383     { nullptr, "rect", "tric", "compact", nullptr };
384
385     enum
386     {
387         ecSel, ecTric, ecRect, ecZero, ecNR
388     };
389     const char *center_opt[ecNR+1] =
390     { nullptr, "tric", "rect", "zero", nullptr };
391     int         ecenter;
392
393     int         fit_enum;
394     enum
395     {
396         efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
397     };
398     const char *fit[efNR + 1] =
399     {
400         nullptr, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
401         "progressive", nullptr
402     };
403
404     static gmx_bool  bSeparate     = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
405     static gmx_bool  bCenter       = FALSE;
406     static int       skip_nr       = 1, ndec = 3, nzero = 0;
407     static real      tzero         = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
408     static rvec      newbox        = {0, 0, 0}, shift = {0, 0, 0}, trans = {0, 0, 0};
409     static char     *exec_command  = nullptr;
410     static real      dropunder     = 0, dropover = 0;
411     static gmx_bool  bRound        = FALSE;
412
413     t_pargs
414         pa[] =
415     {
416         { "-skip", FALSE, etINT,
417           { &skip_nr }, "Only write every nr-th frame" },
418         { "-dt", FALSE, etTIME,
419           { &delta_t },
420           "Only write frame when t MOD dt = first time (%t)" },
421         { "-round", FALSE, etBOOL,
422           { &bRound }, "Round measurements to nearest picosecond"},
423         { "-dump", FALSE, etTIME,
424           { &tdump }, "Dump frame nearest specified time (%t)" },
425         { "-t0", FALSE, etTIME,
426           { &tzero },
427           "Starting time (%t) (default: don't change)" },
428         { "-timestep", FALSE, etTIME,
429           { &timestep },
430           "Change time step between input frames (%t)" },
431         { "-pbc", FALSE, etENUM,
432           { pbc_opt },
433           "PBC treatment (see help text for full description)" },
434         { "-ur", FALSE, etENUM,
435           { unitcell_opt }, "Unit-cell representation" },
436         { "-center", FALSE, etBOOL,
437           { &bCenter }, "Center atoms in box" },
438         { "-boxcenter", FALSE, etENUM,
439           { center_opt }, "Center for -pbc and -center" },
440         { "-box", FALSE, etRVEC,
441           { newbox },
442           "Size for new cubic box (default: read from input)" },
443         { "-trans", FALSE, etRVEC,
444           { trans },
445           "All coordinates will be translated by trans. This "
446           "can advantageously be combined with -pbc mol -ur "
447           "compact." },
448         { "-shift", FALSE, etRVEC,
449           { shift },
450           "All coordinates will be shifted by framenr*shift" },
451         { "-fit", FALSE, etENUM,
452           { fit },
453           "Fit molecule to ref structure in the structure file" },
454         { "-ndec", FALSE, etINT,
455           { &ndec },
456           "Number of decimal places to write to .xtc output" },
457         { "-vel", FALSE, etBOOL,
458           { &bVels }, "Read and write velocities if possible" },
459         { "-force", FALSE, etBOOL,
460           { &bForce }, "Read and write forces if possible" },
461         { "-trunc", FALSE, etTIME,
462           { &ttrunc },
463           "Truncate input trajectory file after this time (%t)" },
464         { "-exec", FALSE, etSTR,
465           { &exec_command },
466           "Execute command for every output frame with the "
467           "frame number as argument" },
468         { "-split", FALSE, etTIME,
469           { &split_t },
470           "Start writing new file when t MOD split = first "
471           "time (%t)" },
472         { "-sep", FALSE, etBOOL,
473           { &bSeparate },
474           "Write each frame to a separate .gro, .g96 or .pdb "
475           "file" },
476         { "-nzero", FALSE, etINT,
477           { &nzero },
478           "If the -sep flag is set, use these many digits "
479           "for the file numbers and prepend zeros as needed" },
480         { "-dropunder", FALSE, etREAL,
481           { &dropunder }, "Drop all frames below this value" },
482         { "-dropover", FALSE, etREAL,
483           { &dropover }, "Drop all frames above this value" },
484         { "-conect", FALSE, etBOOL,
485           { &bCONECT },
486           "Add conect records when writing [REF].pdb[ref] files. Useful "
487           "for visualization of non-standard molecules, e.g. "
488           "coarse grained ones" }
489     };
490 #define NPA asize(pa)
491
492     FILE             *out    = nullptr;
493     t_trxstatus      *trxout = nullptr;
494     t_trxstatus      *trxin;
495     int               file_nr;
496     t_trxframe        fr, frout;
497     int               flags;
498     rvec             *xmem  = nullptr, *vmem = nullptr, *fmem = nullptr;
499     rvec             *xp    = nullptr, x_shift, hbox;
500     real             *w_rls = nullptr;
501     int               m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
502 #define SKIP 10
503     t_topology       *top   = nullptr;
504     gmx_conect        gc    = nullptr;
505     int               ePBC  = -1;
506     t_atoms          *atoms = nullptr, useatoms;
507     matrix            top_box;
508     int              *index = nullptr, *cindex = nullptr;
509     char             *grpnm = nullptr;
510     int              *frindex, nrfri;
511     char             *frname;
512     int               ifit, my_clust = -1;
513     int              *ind_fit;
514     char             *gn_fit;
515     t_cluster_ndx    *clust           = nullptr;
516     t_trxstatus     **clust_status    = nullptr;
517     int              *clust_status_id = nullptr;
518     int               ntrxopen        = 0;
519     int              *nfwritten       = nullptr;
520     int               ndrop           = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
521     double          **dropval;
522     real              tshift = 0, dt = -1, prec;
523     gmx_bool          bFit, bPFit, bReset;
524     int               nfitdim;
525     gmx_rmpbc_t       gpbc = nullptr;
526     gmx_bool          bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
527     gmx_bool          bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
528     gmx_bool          bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetXtcPrec, bNeedPrec;
529     gmx_bool          bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
530     gmx_bool          bSubTraj = FALSE, bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
531     gmx_bool          bWriteFrame, bSplitHere;
532     const char       *top_file, *in_file, *out_file = nullptr;
533     char              out_file2[256], *charpt;
534     char             *outf_base = nullptr;
535     const char       *outf_ext  = nullptr;
536     char              top_title[256], timestr[32], stepstr[32], filemode[5];
537     gmx_output_env_t *oenv;
538
539     t_filenm          fnm[] = {
540         { efTRX, "-f",   nullptr,      ffREAD  },
541         { efTRO, "-o",   nullptr,      ffWRITE },
542         { efTPS, nullptr,   nullptr,      ffOPTRD },
543         { efNDX, nullptr,   nullptr,      ffOPTRD },
544         { efNDX, "-fr",  "frames",  ffOPTRD },
545         { efNDX, "-sub", "cluster", ffOPTRD },
546         { efXVG, "-drop", "drop",    ffOPTRD }
547     };
548 #define NFILE asize(fnm)
549
550     if (!parse_common_args(&argc, argv,
551                            PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW |
552                            PCA_TIME_UNIT,
553                            NFILE, fnm, NPA, pa, asize(desc), desc,
554                            0, nullptr, &oenv))
555     {
556         return 0;
557     }
558     fprintf(stdout, "Note that major changes are planned in future for "
559             "trjconv, to improve usability and utility.\n");
560
561     top_file = ftp2fn(efTPS, NFILE, fnm);
562
563     /* Check command line */
564     in_file = opt2fn("-f", NFILE, fnm);
565
566     if (ttrunc != -1)
567     {
568         do_trunc(in_file, ttrunc);
569     }
570     else
571     {
572         /* mark active cmdline options */
573         bSetBox     = opt2parg_bSet("-box", NPA, pa);
574         bSetTime    = opt2parg_bSet("-t0", NPA, pa);
575         bSetXtcPrec = opt2parg_bSet("-ndec", NPA, pa);
576         bSetUR      = opt2parg_bSet("-ur", NPA, pa);
577         bExec       = opt2parg_bSet("-exec", NPA, pa);
578         bTimeStep   = opt2parg_bSet("-timestep", NPA, pa);
579         bTDump      = opt2parg_bSet("-dump", NPA, pa);
580         bDropUnder  = opt2parg_bSet("-dropunder", NPA, pa);
581         bDropOver   = opt2parg_bSet("-dropover", NPA, pa);
582         bTrans      = opt2parg_bSet("-trans", NPA, pa);
583         bSplit      = (split_t != 0);
584
585         /* parse enum options */
586         fit_enum      = nenum(fit);
587         bFit          = (fit_enum == efFit || fit_enum == efFitXY);
588         bReset        = (fit_enum == efReset || fit_enum == efResetXY);
589         bPFit         = fit_enum == efPFit;
590         pbc_enum      = nenum(pbc_opt);
591         bPBCWhole     = pbc_enum == epWhole;
592         bPBCcomRes    = pbc_enum == epComRes;
593         bPBCcomMol    = pbc_enum == epComMol;
594         bPBCcomAtom   = pbc_enum == epComAtom;
595         bNoJump       = pbc_enum == epNojump;
596         bCluster      = pbc_enum == epCluster;
597         bPBC          = pbc_enum != epNone;
598         unitcell_enum = nenum(unitcell_opt);
599         ecenter       = nenum(center_opt) - ecTric;
600
601         /* set and check option dependencies */
602         if (bPFit)
603         {
604             bFit = TRUE;        /* for pfit, fit *must* be set */
605         }
606         if (bFit)
607         {
608             bReset = TRUE;       /* for fit, reset *must* be set */
609         }
610         nfitdim = 0;
611         if (bFit || bReset)
612         {
613             nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
614         }
615         bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
616
617         if (bSetUR)
618         {
619             if (!(bPBCcomRes || bPBCcomMol ||  bPBCcomAtom))
620             {
621                 fprintf(stderr,
622                         "WARNING: Option for unitcell representation (-ur %s)\n"
623                         "         only has effect in combination with -pbc %s, %s or %s.\n"
624                         "         Ingoring unitcell representation.\n\n",
625                         unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
626             }
627         }
628         if (bFit && bPBC)
629         {
630             gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
631                       "Please do the PBC condition treatment first and then run trjconv in a second step\n"
632                       "for the rotational fit.\n"
633                       "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
634                       "results!");
635         }
636
637         /* ndec for XTC writing is in nr of decimal places, prec is a multiplication factor: */
638         prec = 1;
639         for (i = 0; i < ndec; i++)
640         {
641             prec *= 10;
642         }
643
644         bIndex = ftp2bSet(efNDX, NFILE, fnm);
645
646
647         /* Determine output type */
648         out_file = opt2fn("-o", NFILE, fnm);
649         int ftp  = fn2ftp(out_file);
650         fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
651         bNeedPrec = (ftp == efXTC);
652         int ftpin = fn2ftp(in_file);
653         if (bVels)
654         {
655             /* check if velocities are possible in input and output files */
656             bVels = (ftp == efTRR || ftp == efGRO ||
657                      ftp == efG96 || ftp == efTNG)
658                 && (ftpin == efTRR || ftpin == efGRO ||
659                     ftpin == efG96 || ftpin == efTNG || ftpin == efCPT);
660         }
661         if (bSeparate || bSplit)
662         {
663             outf_ext = std::strrchr(out_file, '.');
664             if (outf_ext == nullptr)
665             {
666                 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
667             }
668             outf_base = gmx_strdup(out_file);
669             outf_base[outf_ext - out_file] = '\0';
670         }
671
672         bSubTraj = opt2bSet("-sub", NFILE, fnm);
673         if (bSubTraj)
674         {
675             if ((ftp != efXTC) && (ftp != efTRR))
676             {
677                 /* It seems likely that other trajectory file types
678                  * could work here. */
679                 gmx_fatal(FARGS, "Can only use the sub option with output file types "
680                           "xtc and trr");
681             }
682             clust = cluster_index(nullptr, opt2fn("-sub", NFILE, fnm));
683
684             /* Check for number of files disabled, as FOPEN_MAX is not the correct
685              * number to check for. In my linux box it is only 16.
686              */
687             if (/* DISABLES CODE */ (false))
688             //if (clust->clust->nr > FOPEN_MAX-4)
689             {
690                 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
691                           " trajectories.\ntry splitting the index file in %d parts.\n"
692                           "FOPEN_MAX = %d",
693                           clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
694             }
695             gmx_warning("The -sub option could require as many open output files as there are\n"
696                         "index groups in the file (%d). If you get I/O errors opening new files,\n"
697                         "try reducing the number of index groups in the file, and perhaps\n"
698                         "using trjconv -sub several times on different chunks of your index file.\n",
699                         clust->clust->nr);
700
701             snew(clust_status, clust->clust->nr);
702             snew(clust_status_id, clust->clust->nr);
703             snew(nfwritten, clust->clust->nr);
704             for (i = 0; (i < clust->clust->nr); i++)
705             {
706                 clust_status[i]    = nullptr;
707                 clust_status_id[i] = -1;
708             }
709             bSeparate = bSplit = FALSE;
710         }
711         /* skipping */
712         if (skip_nr <= 0)
713         {
714             gmx_fatal(FARGS, "Argument for -skip (%d) needs to be greater or equal to 1.", skip_nr);
715         }
716
717         std::unique_ptr<gmx_mtop_t> mtop = read_mtop_for_tng(top_file, in_file, out_file);
718
719         /* Determine whether to read a topology */
720         bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
721                 bRmPBC || bReset || bPBCcomMol || bCluster ||
722                 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
723
724         /* Determine if when can read index groups */
725         bIndex = (bIndex || bTPS);
726
727         if (bTPS)
728         {
729             snew(top, 1);
730             read_tps_conf(top_file, top, &ePBC, &xp, nullptr, top_box,
731                           bReset || bPBCcomRes);
732             std::strncpy(top_title, *top->name, 255);
733             top_title[255] = '\0';
734             atoms          = &top->atoms;
735
736             if (0 == top->mols.nr && (bCluster || bPBCcomMol))
737             {
738                 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
739             }
740
741             /* top_title is only used for gro and pdb,
742              * the header in such a file is top_title, followed by
743              * t= ... and/or step= ...
744              * to prevent double t= or step=, remove it from top_title.
745              * From GROMACS-2018 we only write t/step when the frame actually
746              * has a valid time/step, so we need to check for both separately.
747              */
748             if ((charpt = std::strstr(top_title, " t= ")))
749             {
750                 charpt[0] = '\0';
751             }
752             if ((charpt = std::strstr(top_title, " step= ")))
753             {
754                 charpt[0] = '\0';
755             }
756
757             if (bCONECT)
758             {
759                 gc = gmx_conect_generate(top);
760             }
761             if (bRmPBC)
762             {
763                 gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
764             }
765         }
766
767         /* get frame number index */
768         frindex = nullptr;
769         if (opt2bSet("-fr", NFILE, fnm))
770         {
771             printf("Select groups of frame number indices:\n");
772             rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, &frindex, &frname);
773             if (debug)
774             {
775                 for (i = 0; i < nrfri; i++)
776                 {
777                     fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
778                 }
779             }
780         }
781
782         /* get index groups etc. */
783         if (bReset)
784         {
785             printf("Select group for %s fit\n",
786                    bFit ? "least squares" : "translational");
787             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
788                       1, &ifit, &ind_fit, &gn_fit);
789
790             if (bFit)
791             {
792                 if (ifit < 2)
793                 {
794                     gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
795                 }
796                 else if (ifit == 3)
797                 {
798                     fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
799                 }
800             }
801         }
802         else if (bCluster)
803         {
804             printf("Select group for clustering\n");
805             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
806                       1, &ifit, &ind_fit, &gn_fit);
807         }
808
809         if (bIndex)
810         {
811             if (bCenter)
812             {
813                 printf("Select group for centering\n");
814                 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
815                           1, &ncent, &cindex, &grpnm);
816             }
817             printf("Select group for output\n");
818             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
819                       1, &nout, &index, &grpnm);
820         }
821         else
822         {
823             /* no index file, so read natoms from TRX */
824             if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
825             {
826                 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
827             }
828             natoms = fr.natoms;
829             close_trx(trxin);
830             sfree(fr.x);
831             snew(index, natoms);
832             for (i = 0; i < natoms; i++)
833             {
834                 index[i] = i;
835             }
836             nout = natoms;
837             if (bCenter)
838             {
839                 ncent  = nout;
840                 cindex = index;
841             }
842         }
843
844         if (bReset)
845         {
846             snew(w_rls, atoms->nr);
847             for (i = 0; (i < ifit); i++)
848             {
849                 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
850             }
851
852             /* Restore reference structure and set to origin,
853                store original location (to put structure back) */
854             if (bRmPBC)
855             {
856                 gmx_rmpbc(gpbc, top->atoms.nr, top_box, xp);
857             }
858             copy_rvec(xp[index[0]], x_shift);
859             reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, nullptr, xp, w_rls);
860             rvec_dec(x_shift, xp[index[0]]);
861         }
862         else
863         {
864             clear_rvec(x_shift);
865         }
866
867         if (bDropUnder || bDropOver)
868         {
869             /* Read the .xvg file with the drop values */
870             fprintf(stderr, "\nReading drop file ...");
871             ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
872             fprintf(stderr, " %d time points\n", ndrop);
873             if (ndrop == 0 || ncol < 2)
874             {
875                 gmx_fatal(FARGS, "Found no data points in %s",
876                           opt2fn("-drop", NFILE, fnm));
877             }
878             drop0 = 0;
879             drop1 = 0;
880         }
881
882         /* Make atoms struct for output in GRO or PDB files */
883         if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
884         {
885             /* get memory for stuff to go in .pdb file, and initialize
886              * the pdbinfo structure part if the input has it.
887              */
888             init_t_atoms(&useatoms, atoms->nr, atoms->havePdbInfo);
889             sfree(useatoms.resinfo);
890             useatoms.resinfo = atoms->resinfo;
891             for (i = 0; (i < nout); i++)
892             {
893                 useatoms.atomname[i] = atoms->atomname[index[i]];
894                 useatoms.atom[i]     = atoms->atom[index[i]];
895                 if (atoms->havePdbInfo)
896                 {
897                     useatoms.pdbinfo[i]  = atoms->pdbinfo[index[i]];
898                 }
899                 useatoms.nres        = std::max(useatoms.nres, useatoms.atom[i].resind+1);
900             }
901             useatoms.nr = nout;
902         }
903         /* select what to read */
904         if (ftp == efTRR)
905         {
906             flags = TRX_READ_X;
907         }
908         else
909         {
910             flags = TRX_NEED_X;
911         }
912         if (bVels)
913         {
914             flags = flags | TRX_READ_V;
915         }
916         if (bForce)
917         {
918             flags = flags | TRX_READ_F;
919         }
920
921         /* open trx file for reading */
922         bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
923         if (fr.bPrec)
924         {
925             fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
926         }
927         if (bNeedPrec)
928         {
929             if (bSetXtcPrec || !fr.bPrec)
930             {
931                 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
932             }
933             else
934             {
935                 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
936             }
937         }
938
939         if (bHaveFirstFrame)
940         {
941             if (bTDump)
942             {
943                 // Determine timestep (assuming constant spacing for now) if we
944                 // need to dump frames based on time. This is required so we do not
945                 // skip the first frame if that was the one that should have been dumped
946                 double firstFrameTime = fr.time;
947                 if (read_next_frame(oenv, trxin, &fr))
948                 {
949                     dt     = fr.time - firstFrameTime;
950                     bDTset = TRUE;
951                     if (dt <= 0)
952                     {
953                         fprintf(stderr, "Warning: Frame times are not incrementing - will dump first frame.\n");
954                     }
955                 }
956                 // Now close and reopen so we are at first frame again
957                 close_trx(trxin);
958                 done_frame(&fr);
959                 // Reopen at first frame (We already know it exists if we got here)
960                 read_first_frame(oenv, &trxin, in_file, &fr, flags);
961             }
962
963             set_trxframe_ePBC(&fr, ePBC);
964             natoms = fr.natoms;
965
966             if (bSetTime)
967             {
968                 tshift = tzero-fr.time;
969             }
970             else
971             {
972                 tzero = fr.time;
973             }
974
975             bCopy = FALSE;
976             if (bIndex)
977             {
978                 /* check if index is meaningful */
979                 for (i = 0; i < nout; i++)
980                 {
981                     if (index[i] >= natoms)
982                     {
983                         gmx_fatal(FARGS,
984                                   "Index[%d] %d is larger than the number of atoms in the\n"
985                                   "trajectory file (%d). There is a mismatch in the contents\n"
986                                   "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
987                     }
988                     bCopy = bCopy || (i != index[i]);
989                 }
990             }
991
992             /* open output for writing */
993             std::strcpy(filemode, "w");
994             switch (ftp)
995             {
996                 case efTNG:
997                     trxout = trjtools_gmx_prepare_tng_writing(out_file,
998                                                               filemode[0],
999                                                               trxin,
1000                                                               nullptr,
1001                                                               nout,
1002                                                               mtop.get(),
1003                                                               gmx::arrayRefFromArray(index, nout),
1004                                                               grpnm);
1005                     break;
1006                 case efXTC:
1007                 case efTRR:
1008                     out = nullptr;
1009                     if (!bSplit && !bSubTraj)
1010                     {
1011                         trxout = open_trx(out_file, filemode);
1012                     }
1013                     break;
1014                 case efGRO:
1015                 case efG96:
1016                 case efPDB:
1017                     if (( !bSeparate && !bSplit ) && !bSubTraj)
1018                     {
1019                         out = gmx_ffopen(out_file, filemode);
1020                     }
1021                     break;
1022                 default:
1023                     gmx_incons("Illegal output file format");
1024             }
1025
1026             if (bCopy)
1027             {
1028                 snew(xmem, nout);
1029                 if (bVels)
1030                 {
1031                     snew(vmem, nout);
1032                 }
1033                 if (bForce)
1034                 {
1035                     snew(fmem, nout);
1036                 }
1037             }
1038
1039             /* Start the big loop over frames */
1040             file_nr  =  0;
1041             frame    =  0;
1042             outframe =  0;
1043             model_nr =  0;
1044
1045             /* Main loop over frames */
1046             do
1047             {
1048                 if (!fr.bStep)
1049                 {
1050                     /* set the step */
1051                     fr.step = newstep;
1052                     newstep++;
1053                 }
1054                 if (bSubTraj)
1055                 {
1056                     /*if (frame >= clust->clust->nra)
1057                        gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1058                     if (frame > clust->maxframe)
1059                     {
1060                         my_clust = -1;
1061                     }
1062                     else
1063                     {
1064                         my_clust = clust->inv_clust[frame];
1065                     }
1066                     if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1067                         (my_clust == -1))
1068                     {
1069                         my_clust = -1;
1070                     }
1071                 }
1072
1073                 if (bSetBox)
1074                 {
1075                     /* generate new box */
1076                     if (!fr.bBox)
1077                     {
1078                         clear_mat(fr.box);
1079                     }
1080                     for (m = 0; m < DIM; m++)
1081                     {
1082                         if (newbox[m] >= 0)
1083                         {
1084                             fr.box[m][m] = newbox[m];
1085                         }
1086                         else
1087                         {
1088                             if (!fr.bBox)
1089                             {
1090                                 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1091                             }
1092                         }
1093                     }
1094                 }
1095
1096                 if (bTrans)
1097                 {
1098                     for (i = 0; i < natoms; i++)
1099                     {
1100                         rvec_inc(fr.x[i], trans);
1101                     }
1102                 }
1103
1104                 if (bTDump)
1105                 {
1106                     // If we could not read two frames or times are not incrementing
1107                     // we have almost no idea what to do,
1108                     // but dump the first frame so output is not broken.
1109                     if (dt <= 0 || !bDTset)
1110                     {
1111                         bDumpFrame = true;
1112                     }
1113                     else
1114                     {
1115                         // Dump the frame if we are less than half a frame time
1116                         // below it. This will also ensure we at least dump a
1117                         // somewhat reasonable frame if the spacing is unequal
1118                         // and we have overrun the frame time. Once we dump one
1119                         // frame based on time we quit, so it does not matter
1120                         // that this might be true for all subsequent frames too.
1121                         bDumpFrame = (fr.time > tdump-0.5*dt);
1122                     }
1123                 }
1124                 else
1125                 {
1126                     bDumpFrame = FALSE;
1127                 }
1128
1129                 /* determine if an atom jumped across the box and reset it if so */
1130                 if (bNoJump && (bTPS || frame != 0))
1131                 {
1132                     for (d = 0; d < DIM; d++)
1133                     {
1134                         hbox[d] = 0.5*fr.box[d][d];
1135                     }
1136                     for (i = 0; i < natoms; i++)
1137                     {
1138                         if (bReset)
1139                         {
1140                             rvec_dec(fr.x[i], x_shift);
1141                         }
1142                         for (m = DIM-1; m >= 0; m--)
1143                         {
1144                             if (hbox[m] > 0)
1145                             {
1146                                 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1147                                 {
1148                                     for (d = 0; d <= m; d++)
1149                                     {
1150                                         fr.x[i][d] += fr.box[m][d];
1151                                     }
1152                                 }
1153                                 while (fr.x[i][m]-xp[i][m] > hbox[m])
1154                                 {
1155                                     for (d = 0; d <= m; d++)
1156                                     {
1157                                         fr.x[i][d] -= fr.box[m][d];
1158                                     }
1159                                 }
1160                             }
1161                         }
1162                     }
1163                 }
1164                 else if (bCluster)
1165                 {
1166                     calc_pbc_cluster(ecenter, ifit, top, ePBC, fr.x, ind_fit, fr.box);
1167                 }
1168
1169                 if (bPFit)
1170                 {
1171                     /* Now modify the coords according to the flags,
1172                        for normal fit, this is only done for output frames */
1173                     if (bRmPBC)
1174                     {
1175                         gmx_rmpbc_trxfr(gpbc, &fr);
1176                     }
1177
1178                     reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1179                     do_fit(natoms, w_rls, xp, fr.x);
1180                 }
1181
1182                 /* store this set of coordinates for future use */
1183                 if (bPFit || bNoJump)
1184                 {
1185                     if (xp == nullptr)
1186                     {
1187                         snew(xp, natoms);
1188                     }
1189                     for (i = 0; (i < natoms); i++)
1190                     {
1191                         copy_rvec(fr.x[i], xp[i]);
1192                         rvec_inc(fr.x[i], x_shift);
1193                     }
1194                 }
1195
1196                 if (frindex)
1197                 {
1198                     /* see if we have a frame from the frame index group */
1199                     for (i = 0; i < nrfri && !bDumpFrame; i++)
1200                     {
1201                         bDumpFrame = frame == frindex[i];
1202                     }
1203                 }
1204                 if (debug && bDumpFrame)
1205                 {
1206                     fprintf(debug, "dumping %d\n", frame);
1207                 }
1208
1209                 bWriteFrame =
1210                     ( ( !bTDump && (frindex == nullptr) && frame % skip_nr == 0 ) || bDumpFrame );
1211
1212                 if (bWriteFrame && (bDropUnder || bDropOver))
1213                 {
1214                     while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1215                     {
1216                         drop0 = drop1;
1217                         drop1++;
1218                     }
1219                     if (std::abs(dropval[0][drop0] - fr.time)
1220                         < std::abs(dropval[0][drop1] - fr.time))
1221                     {
1222                         dropuse = drop0;
1223                     }
1224                     else
1225                     {
1226                         dropuse = drop1;
1227                     }
1228                     if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1229                         (bDropOver && dropval[1][dropuse] > dropover))
1230                     {
1231                         bWriteFrame = FALSE;
1232                     }
1233                 }
1234
1235                 if (bWriteFrame)
1236                 {
1237                     /* We should avoid modifying the input frame,
1238                      * but since here we don't have the output frame yet,
1239                      * we introduce a temporary output frame time variable.
1240                      */
1241                     real frout_time;
1242
1243                     frout_time = fr.time;
1244
1245                     /* calc new time */
1246                     if (bTimeStep)
1247                     {
1248                         frout_time = tzero + frame*timestep;
1249                     }
1250                     else if (bSetTime)
1251                     {
1252                         frout_time += tshift;
1253                     }
1254
1255                     if (bTDump)
1256                     {
1257                         fprintf(stderr, "\nDumping frame at t= %g %s\n",
1258                                 output_env_conv_time(oenv, frout_time), output_env_get_time_unit(oenv).c_str());
1259                     }
1260
1261                     /* check for writing at each delta_t */
1262                     bDoIt = (delta_t == 0);
1263                     if (!bDoIt)
1264                     {
1265                         if (!bRound)
1266                         {
1267                             bDoIt = bRmod(frout_time, tzero, delta_t);
1268                         }
1269                         else
1270                         {
1271                             /* round() is not C89 compatible, so we do this:  */
1272                             bDoIt = bRmod(std::floor(frout_time+0.5), std::floor(tzero+0.5),
1273                                           std::floor(delta_t+0.5));
1274                         }
1275                     }
1276
1277                     if (bDoIt || bTDump)
1278                     {
1279                         /* print sometimes */
1280                         if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1281                         {
1282                             fprintf(stderr, " ->  frame %6d time %8.3f      \r",
1283                                     outframe, output_env_conv_time(oenv, frout_time));
1284                             fflush(stderr);
1285                         }
1286
1287                         if (!bPFit)
1288                         {
1289                             /* Now modify the coords according to the flags,
1290                                for PFit we did this already! */
1291
1292                             if (bRmPBC)
1293                             {
1294                                 gmx_rmpbc_trxfr(gpbc, &fr);
1295                             }
1296
1297                             if (bReset)
1298                             {
1299                                 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1300                                 if (bFit)
1301                                 {
1302                                     do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1303                                 }
1304                                 if (!bCenter)
1305                                 {
1306                                     for (i = 0; i < natoms; i++)
1307                                     {
1308                                         rvec_inc(fr.x[i], x_shift);
1309                                     }
1310                                 }
1311                             }
1312
1313                             if (bCenter)
1314                             {
1315                                 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1316                             }
1317                         }
1318
1319                         auto positionsArrayRef = gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec *>(fr.x), natoms);
1320                         if (bPBCcomAtom)
1321                         {
1322                             switch (unitcell_enum)
1323                             {
1324                                 case euRect:
1325                                     put_atoms_in_box(ePBC, fr.box, positionsArrayRef);
1326                                     break;
1327                                 case euTric:
1328                                     put_atoms_in_triclinic_unitcell(ecenter, fr.box, positionsArrayRef);
1329                                     break;
1330                                 case euCompact:
1331                                     put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1332                                                                   positionsArrayRef);
1333                                     break;
1334                             }
1335                         }
1336                         if (bPBCcomRes)
1337                         {
1338                             put_residue_com_in_box(unitcell_enum, ecenter,
1339                                                    natoms, atoms->atom, ePBC, fr.box, fr.x);
1340                         }
1341                         if (bPBCcomMol)
1342                         {
1343                             put_molecule_com_in_box(unitcell_enum, ecenter,
1344                                                     &top->mols,
1345                                                     natoms, atoms->atom, ePBC, fr.box, fr.x);
1346                         }
1347                         /* Copy the input trxframe struct to the output trxframe struct */
1348                         frout        = fr;
1349                         frout.time   = frout_time;
1350                         frout.bV     = (frout.bV && bVels);
1351                         frout.bF     = (frout.bF && bForce);
1352                         frout.natoms = nout;
1353                         if (bNeedPrec && (bSetXtcPrec || !fr.bPrec))
1354                         {
1355                             frout.bPrec = TRUE;
1356                             frout.prec  = prec;
1357                         }
1358                         if (bCopy)
1359                         {
1360                             frout.x = xmem;
1361                             if (frout.bV)
1362                             {
1363                                 frout.v = vmem;
1364                             }
1365                             if (frout.bF)
1366                             {
1367                                 frout.f = fmem;
1368                             }
1369                             for (i = 0; i < nout; i++)
1370                             {
1371                                 copy_rvec(fr.x[index[i]], frout.x[i]);
1372                                 if (frout.bV)
1373                                 {
1374                                     copy_rvec(fr.v[index[i]], frout.v[i]);
1375                                 }
1376                                 if (frout.bF)
1377                                 {
1378                                     copy_rvec(fr.f[index[i]], frout.f[i]);
1379                                 }
1380                             }
1381                         }
1382
1383                         if (opt2parg_bSet("-shift", NPA, pa))
1384                         {
1385                             for (i = 0; i < nout; i++)
1386                             {
1387                                 for (d = 0; d < DIM; d++)
1388                                 {
1389                                     frout.x[i][d] += outframe*shift[d];
1390                                 }
1391                             }
1392                         }
1393
1394                         if (!bRound)
1395                         {
1396                             bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1397                         }
1398                         else
1399                         {
1400                             /* round() is not C89 compatible, so we do this: */
1401                             bSplitHere = bSplit && bRmod(std::floor(frout.time+0.5),
1402                                                          std::floor(tzero+0.5),
1403                                                          std::floor(split_t+0.5));
1404                         }
1405                         if (bSeparate || bSplitHere)
1406                         {
1407                             mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1408                         }
1409
1410                         std::string title;
1411                         switch (ftp)
1412                         {
1413                             case efTNG:
1414                                 write_tng_frame(trxout, &frout);
1415                                 // TODO when trjconv behaves better: work how to read and write lambda
1416                                 break;
1417                             case efTRR:
1418                             case efXTC:
1419                                 if (bSplitHere)
1420                                 {
1421                                     if (trxout)
1422                                     {
1423                                         close_trx(trxout);
1424                                     }
1425                                     trxout = open_trx(out_file2, filemode);
1426                                 }
1427                                 if (bSubTraj)
1428                                 {
1429                                     if (my_clust != -1)
1430                                     {
1431                                         char buf[STRLEN];
1432                                         if (clust_status_id[my_clust] == -1)
1433                                         {
1434                                             sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1435                                             clust_status[my_clust]    = open_trx(buf, "w");
1436                                             clust_status_id[my_clust] = 1;
1437                                             ntrxopen++;
1438                                         }
1439                                         else if (clust_status_id[my_clust] == -2)
1440                                         {
1441                                             gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1442                                                       clust->grpname[my_clust], ntrxopen, frame,
1443                                                       my_clust);
1444                                         }
1445                                         write_trxframe(clust_status[my_clust], &frout, gc);
1446                                         nfwritten[my_clust]++;
1447                                         if (nfwritten[my_clust] ==
1448                                             (clust->clust->index[my_clust+1]-
1449                                              clust->clust->index[my_clust]))
1450                                         {
1451                                             close_trx(clust_status[my_clust]);
1452                                             clust_status[my_clust]    = nullptr;
1453                                             clust_status_id[my_clust] = -2;
1454                                             ntrxopen--;
1455                                             if (ntrxopen < 0)
1456                                             {
1457                                                 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1458                                             }
1459                                         }
1460                                     }
1461                                 }
1462                                 else
1463                                 {
1464                                     write_trxframe(trxout, &frout, gc);
1465                                 }
1466                                 break;
1467                             case efGRO:
1468                             case efG96:
1469                             case efPDB:
1470                                 // Only add a generator statement if title is empty,
1471                                 // to avoid multiple generated-by statements from various programs
1472                                 if (std::strlen(top_title) == 0)
1473                                 {
1474                                     sprintf(top_title, "Generated by trjconv");
1475                                 }
1476                                 if (frout.bTime)
1477                                 {
1478                                     sprintf(timestr, " t= %9.5f", frout.time);
1479                                 }
1480                                 else
1481                                 {
1482                                     std::strcpy(timestr, "");
1483                                 }
1484                                 if (frout.bStep)
1485                                 {
1486                                     sprintf(stepstr, " step= %" PRId64, frout.step);
1487                                 }
1488                                 else
1489                                 {
1490                                     std::strcpy(stepstr, "");
1491                                 }
1492                                 title = gmx::formatString("%s%s%s", top_title, timestr, stepstr);
1493                                 if (bSeparate || bSplitHere)
1494                                 {
1495                                     out = gmx_ffopen(out_file2, "w");
1496                                 }
1497                                 switch (ftp)
1498                                 {
1499                                     case efGRO:
1500                                         write_hconf_p(out, title.c_str(), &useatoms,
1501                                                       frout.x, frout.bV ? frout.v : nullptr, frout.box);
1502                                         break;
1503                                     case efPDB:
1504                                         fprintf(out, "REMARK    GENERATED BY TRJCONV\n");
1505                                         /* if reading from pdb, we want to keep the original
1506                                            model numbering else we write the output frame
1507                                            number plus one, because model 0 is not allowed in pdb */
1508                                         if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1509                                         {
1510                                             model_nr = fr.step;
1511                                         }
1512                                         else
1513                                         {
1514                                             model_nr++;
1515                                         }
1516                                         write_pdbfile(out, title.c_str(), &useatoms, frout.x,
1517                                                       frout.ePBC, frout.box, ' ', model_nr, gc);
1518                                         break;
1519                                     case efG96:
1520                                         const char *outputTitle = "";
1521                                         if (bSeparate || bTDump)
1522                                         {
1523                                             outputTitle = title.c_str();
1524                                             if (bTPS)
1525                                             {
1526                                                 frout.bAtoms = TRUE;
1527                                             }
1528                                             frout.atoms  = &useatoms;
1529                                             frout.bStep  = FALSE;
1530                                             frout.bTime  = FALSE;
1531                                         }
1532                                         else
1533                                         {
1534                                             if (outframe == 0)
1535                                             {
1536                                                 outputTitle = title.c_str();
1537                                             }
1538                                             frout.bAtoms = FALSE;
1539                                             frout.bStep  = TRUE;
1540                                             frout.bTime  = TRUE;
1541                                         }
1542                                         write_g96_conf(out, outputTitle, &frout, -1, nullptr);
1543                                 }
1544                                 if (bSeparate || bSplitHere)
1545                                 {
1546                                     gmx_ffclose(out);
1547                                     out = nullptr;
1548                                 }
1549                                 break;
1550                             default:
1551                                 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1552                         }
1553                         if (bSeparate || bSplitHere)
1554                         {
1555                             file_nr++;
1556                         }
1557
1558                         /* execute command */
1559                         if (bExec)
1560                         {
1561                             char c[255];
1562                             sprintf(c, "%s  %d", exec_command, file_nr-1);
1563                             /*fprintf(stderr,"Executing '%s'\n",c);*/
1564                             if (0 != system(c))
1565                             {
1566                                 gmx_fatal(FARGS, "Error executing command: %s", c);
1567                             }
1568                         }
1569                         outframe++;
1570                     }
1571                 }
1572                 frame++;
1573                 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1574             }
1575             while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1576         }
1577
1578         if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1579         {
1580             fprintf(stderr, "\nWARNING no output, "
1581                     "last frame read at t=%g\n", fr.time);
1582         }
1583         fprintf(stderr, "\n");
1584
1585         close_trx(trxin);
1586         sfree(outf_base);
1587
1588         if (bRmPBC)
1589         {
1590             gmx_rmpbc_done(gpbc);
1591         }
1592
1593         if (trxout)
1594         {
1595             close_trx(trxout);
1596         }
1597         else if (out != nullptr)
1598         {
1599             gmx_ffclose(out);
1600         }
1601         if (bSubTraj)
1602         {
1603             for (i = 0; (i < clust->clust->nr); i++)
1604             {
1605                 if (clust_status_id[i] >= 0)
1606                 {
1607                     close_trx(clust_status[i]);
1608                 }
1609             }
1610         }
1611     }
1612
1613     if (bTPS)
1614     {
1615         done_top(top);
1616         sfree(top);
1617     }
1618     sfree(xp);
1619     sfree(xmem);
1620     sfree(vmem);
1621     sfree(fmem);
1622     sfree(grpnm);
1623     sfree(index);
1624     sfree(cindex);
1625     done_frame(&fr);
1626
1627     do_view(oenv, out_file, nullptr);
1628
1629     output_env_done(oenv);
1630     return 0;
1631 }