Python detection consolidation.
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_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 <cmath>
40 #include <cstdlib>
41 #include <cstring>
42
43 #include <algorithm>
44 #include <memory>
45
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/commandline/viewit.h"
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/fileio/g96io.h"
50 #include "gromacs/fileio/gmxfio.h"
51 #include "gromacs/fileio/groio.h"
52 #include "gromacs/fileio/pdbio.h"
53 #include "gromacs/fileio/tngio.h"
54 #include "gromacs/fileio/tpxio.h"
55 #include "gromacs/fileio/trrio.h"
56 #include "gromacs/fileio/trxio.h"
57 #include "gromacs/fileio/xtcio.h"
58 #include "gromacs/fileio/xvgr.h"
59 #include "gromacs/gmxana/gmx_ana.h"
60 #include "gromacs/math/do_fit.h"
61 #include "gromacs/math/functions.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdtypes/md_enums.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/pbcutil/rmpbc.h"
66 #include "gromacs/topology/index.h"
67 #include "gromacs/topology/topology.h"
68 #include "gromacs/trajectory/trajectoryframe.h"
69 #include "gromacs/utility/arrayref.h"
70 #include "gromacs/utility/arraysize.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/futil.h"
73 #include "gromacs/utility/smalloc.h"
74
75 enum {
76     euSel, euRect, euTric, euCompact, euNR
77 };
78
79
80 static void calc_pbc_cluster(int ecenter, int nrefat, t_topology *top, int ePBC,
81                              rvec x[], const int index[], matrix box)
82 {
83     int       m, i, j, j0, j1, jj, ai, aj;
84     int       imin, jmin;
85     real      fac, min_dist2;
86     rvec      dx, xtest, box_center;
87     int       nmol, imol_center;
88     int      *molind;
89     gmx_bool *bMol, *bTmp;
90     rvec     *m_com, *m_shift;
91     t_pbc     pbc;
92     int      *cluster;
93     int      *added;
94     int       ncluster, nadded;
95     real      tmp_r2;
96
97     calc_box_center(ecenter, box, box_center);
98
99     /* Initiate the pbc structure */
100     std::memset(&pbc, 0, sizeof(pbc));
101     set_pbc(&pbc, ePBC, box);
102
103     /* Convert atom index to molecular */
104     nmol   = top->mols.nr;
105     molind = top->mols.index;
106     snew(bMol, nmol);
107     snew(m_com, nmol);
108     snew(m_shift, nmol);
109     snew(cluster, nmol);
110     snew(added, nmol);
111     snew(bTmp, top->atoms.nr);
112
113     for (i = 0; (i < nrefat); i++)
114     {
115         /* Mark all molecules in the index */
116         ai       = index[i];
117         bTmp[ai] = TRUE;
118         /* Binary search assuming the molecules are sorted */
119         j0 = 0;
120         j1 = nmol-1;
121         while (j0 < j1)
122         {
123             if (ai < molind[j0+1])
124             {
125                 j1 = j0;
126             }
127             else if (ai >= molind[j1])
128             {
129                 j0 = j1;
130             }
131             else
132             {
133                 jj = (j0+j1)/2;
134                 if (ai < molind[jj+1])
135                 {
136                     j1 = jj;
137                 }
138                 else
139                 {
140                     j0 = jj;
141                 }
142             }
143         }
144         bMol[j0] = TRUE;
145     }
146     /* Double check whether all atoms in all molecules that are marked are part
147      * of the cluster. Simultaneously compute the center of geometry.
148      */
149     min_dist2   = 10*gmx::square(trace(box));
150     imol_center = -1;
151     ncluster    = 0;
152     for (i = 0; i < nmol; i++)
153     {
154         for (j = molind[i]; j < molind[i+1]; j++)
155         {
156             if (bMol[i] && !bTmp[j])
157             {
158                 gmx_fatal(FARGS, "Molecule %d marked for clustering but not atom %d in it - check your index!", i+1, j+1);
159             }
160             else if (!bMol[i] && bTmp[j])
161             {
162                 gmx_fatal(FARGS, "Atom %d marked for clustering but not molecule %d - this is an internal error...", j+1, i+1);
163             }
164             else if (bMol[i])
165             {
166                 /* Make molecule whole, move 2nd and higher atom to same periodicity as 1st atom in molecule */
167                 if (j > molind[i])
168                 {
169                     pbc_dx(&pbc, x[j], x[j-1], dx);
170                     rvec_add(x[j-1], dx, x[j]);
171                 }
172                 /* Compute center of geometry of molecule - m_com[i] was zeroed when we did snew() on it! */
173                 rvec_inc(m_com[i], x[j]);
174             }
175         }
176         if (bMol[i])
177         {
178             /* Normalize center of geometry */
179             fac = 1.0/(molind[i+1]-molind[i]);
180             for (m = 0; (m < DIM); m++)
181             {
182                 m_com[i][m] *= fac;
183             }
184             /* Determine which molecule is closest to the center of the box */
185             pbc_dx(&pbc, box_center, m_com[i], dx);
186             tmp_r2 = iprod(dx, dx);
187
188             if (tmp_r2 < min_dist2)
189             {
190                 min_dist2   = tmp_r2;
191                 imol_center = i;
192             }
193             cluster[ncluster++] = i;
194         }
195     }
196     sfree(bTmp);
197
198     if (ncluster <= 0)
199     {
200         fprintf(stderr, "No molecules selected in the cluster\n");
201         return;
202     }
203     else if (imol_center == -1)
204     {
205         fprintf(stderr, "No central molecules could be found\n");
206         return;
207     }
208
209     nadded            = 0;
210     added[nadded++]   = imol_center;
211     bMol[imol_center] = FALSE;
212
213     while (nadded < ncluster)
214     {
215         /* Find min distance between cluster molecules and those remaining to be added */
216         min_dist2   = 10*gmx::square(trace(box));
217         imin        = -1;
218         jmin        = -1;
219         /* Loop over added mols */
220         for (i = 0; i < nadded; i++)
221         {
222             ai = added[i];
223             /* Loop over all mols */
224             for (j = 0; j < ncluster; j++)
225             {
226                 aj = cluster[j];
227                 /* check those remaining to be added */
228                 if (bMol[aj])
229                 {
230                     pbc_dx(&pbc, m_com[aj], m_com[ai], dx);
231                     tmp_r2 = iprod(dx, dx);
232                     if (tmp_r2 < min_dist2)
233                     {
234                         min_dist2   = tmp_r2;
235                         imin        = ai;
236                         jmin        = aj;
237                     }
238                 }
239             }
240         }
241
242         /* Add the best molecule */
243         added[nadded++]   = jmin;
244         bMol[jmin]        = FALSE;
245         /* Calculate the shift from the ai molecule */
246         pbc_dx(&pbc, m_com[jmin], m_com[imin], dx);
247         rvec_add(m_com[imin], dx, xtest);
248         rvec_sub(xtest, m_com[jmin], m_shift[jmin]);
249         rvec_inc(m_com[jmin], m_shift[jmin]);
250
251         for (j = molind[jmin]; j < molind[jmin+1]; j++)
252         {
253             rvec_inc(x[j], m_shift[jmin]);
254         }
255         fprintf(stdout, "\rClustering iteration %d of %d...", nadded, ncluster);
256         fflush(stdout);
257     }
258
259     sfree(added);
260     sfree(cluster);
261     sfree(bMol);
262     sfree(m_com);
263     sfree(m_shift);
264
265     fprintf(stdout, "\n");
266 }
267
268 static void put_molecule_com_in_box(int unitcell_enum, int ecenter,
269                                     t_block *mols,
270                                     int natoms, t_atom atom[],
271                                     int ePBC, matrix box, rvec x[])
272 {
273     int     i, j;
274     int     d;
275     rvec    com, shift, box_center;
276     real    m;
277     double  mtot;
278     t_pbc   pbc;
279
280     calc_box_center(ecenter, box, box_center);
281     set_pbc(&pbc, ePBC, box);
282     if (mols->nr <= 0)
283     {
284         gmx_fatal(FARGS, "There are no molecule descriptions. I need a .tpr file for this pbc option.");
285     }
286     for (i = 0; (i < mols->nr); i++)
287     {
288         /* calc COM */
289         clear_rvec(com);
290         mtot = 0;
291         for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
292         {
293             m = atom[j].m;
294             for (d = 0; d < DIM; d++)
295             {
296                 com[d] += m*x[j][d];
297             }
298             mtot += m;
299         }
300         /* calculate final COM */
301         svmul(1.0/mtot, com, com);
302
303         /* check if COM is outside box */
304         gmx::RVec newCom;
305         copy_rvec(com, newCom);
306         auto      newComArrayRef = gmx::arrayRefFromArray(&newCom, 1);
307         switch (unitcell_enum)
308         {
309             case euRect:
310                 put_atoms_in_box(ePBC, box, newComArrayRef);
311                 break;
312             case euTric:
313                 put_atoms_in_triclinic_unitcell(ecenter, box, newComArrayRef);
314                 break;
315             case euCompact:
316                 put_atoms_in_compact_unitcell(ePBC, ecenter, box, newComArrayRef);
317                 break;
318         }
319         rvec_sub(newCom, com, shift);
320         if (norm2(shift) > 0)
321         {
322             if (debug)
323             {
324                 fprintf(debug, "\nShifting position of molecule %d "
325                         "by %8.3f  %8.3f  %8.3f\n", i+1,
326                         shift[XX], shift[YY], shift[ZZ]);
327             }
328             for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
329             {
330                 rvec_inc(x[j], shift);
331             }
332         }
333     }
334 }
335
336 static void put_residue_com_in_box(int unitcell_enum, int ecenter,
337                                    int natoms, t_atom atom[],
338                                    int ePBC, matrix box, rvec x[])
339 {
340     int              i, j, res_start, res_end;
341     int              d, presnr;
342     real             m;
343     double           mtot;
344     rvec             box_center, com, shift;
345     static const int NOTSET = -12347;
346     calc_box_center(ecenter, box, box_center);
347
348     presnr    = NOTSET;
349     res_start = 0;
350     clear_rvec(com);
351     mtot = 0;
352     for (i = 0; i < natoms+1; i++)
353     {
354         if (i == natoms || (presnr != atom[i].resind && presnr != NOTSET))
355         {
356             /* calculate final COM */
357             res_end = i;
358             svmul(1.0/mtot, com, com);
359
360             /* check if COM is outside box */
361             gmx::RVec newCom;
362             copy_rvec(com, newCom);
363             auto      newComArrayRef = gmx::arrayRefFromArray(&newCom, 1);
364             switch (unitcell_enum)
365             {
366                 case euRect:
367                     put_atoms_in_box(ePBC, box, newComArrayRef);
368                     break;
369                 case euTric:
370                     put_atoms_in_triclinic_unitcell(ecenter, box, newComArrayRef);
371                     break;
372                 case euCompact:
373                     put_atoms_in_compact_unitcell(ePBC, ecenter, box, newComArrayRef);
374                     break;
375             }
376             rvec_sub(newCom, com, shift);
377             if (norm2(shift) != 0.0f)
378             {
379                 if (debug)
380                 {
381                     fprintf(debug, "\nShifting position of residue %d (atoms %d-%d) "
382                             "by %g,%g,%g\n", atom[res_start].resind+1,
383                             res_start+1, res_end+1, shift[XX], shift[YY], shift[ZZ]);
384                 }
385                 for (j = res_start; j < res_end; j++)
386                 {
387                     rvec_inc(x[j], shift);
388                 }
389             }
390             clear_rvec(com);
391             mtot = 0;
392
393             /* remember start of new residue */
394             res_start = i;
395         }
396         if (i < natoms)
397         {
398             /* calc COM */
399             m = atom[i].m;
400             for (d = 0; d < DIM; d++)
401             {
402                 com[d] += m*x[i][d];
403             }
404             mtot += m;
405
406             presnr = atom[i].resind;
407         }
408     }
409 }
410
411 static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, const int ci[])
412 {
413     int  i, m, ai;
414     rvec cmin, cmax, box_center, dx;
415
416     if (nc > 0)
417     {
418         copy_rvec(x[ci[0]], cmin);
419         copy_rvec(x[ci[0]], cmax);
420         for (i = 0; i < nc; i++)
421         {
422             ai = ci[i];
423             for (m = 0; m < DIM; m++)
424             {
425                 if (x[ai][m] < cmin[m])
426                 {
427                     cmin[m] = x[ai][m];
428                 }
429                 else if (x[ai][m] > cmax[m])
430                 {
431                     cmax[m] = x[ai][m];
432                 }
433             }
434         }
435         calc_box_center(ecenter, box, box_center);
436         for (m = 0; m < DIM; m++)
437         {
438             dx[m] = box_center[m]-(cmin[m]+cmax[m])*0.5;
439         }
440
441         for (i = 0; i < n; i++)
442         {
443             rvec_inc(x[i], dx);
444         }
445     }
446 }
447
448 static void mk_filenm(char *base, const char *ext, int ndigit, int file_nr,
449                       char out_file[])
450 {
451     char nbuf[128];
452     int  nd = 0, fnr;
453
454     std::strcpy(out_file, base);
455     fnr = file_nr;
456     do
457     {
458         fnr /= 10;
459         nd++;
460     }
461     while (fnr > 0);
462
463     if (nd < ndigit)
464     {
465         std::strncat(out_file, "00000000000", ndigit-nd);
466     }
467     sprintf(nbuf, "%d.", file_nr);
468     std::strcat(out_file, nbuf);
469     std::strcat(out_file, ext);
470 }
471
472 static void check_trr(const char *fn)
473 {
474     if (fn2ftp(fn) != efTRR)
475     {
476         gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
477     }
478 }
479
480 static void do_trunc(const char *fn, real t0)
481 {
482     t_fileio        *in;
483     FILE            *fp;
484     gmx_bool         bStop, bOK;
485     gmx_trr_header_t sh;
486     gmx_off_t        fpos;
487     char             yesno[256];
488     int              j;
489     real             t = 0;
490
491     if (t0 == -1)
492     {
493         gmx_fatal(FARGS, "You forgot to set the truncation time");
494     }
495
496     /* Check whether this is a .trr file */
497     check_trr(fn);
498
499     in   = gmx_trr_open(fn, "r");
500     fp   = gmx_fio_getfp(in);
501     if (fp == nullptr)
502     {
503         fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
504         gmx_trr_close(in);
505     }
506     else
507     {
508         j     = 0;
509         fpos  = gmx_fio_ftell(in);
510         bStop = FALSE;
511         while (!bStop && gmx_trr_read_frame_header(in, &sh, &bOK))
512         {
513             gmx_trr_read_frame_data(in, &sh, nullptr, nullptr, nullptr, nullptr);
514             fpos = gmx_ftell(fp);
515             t    = sh.t;
516             if (t >= t0)
517             {
518                 gmx_fseek(fp, fpos, SEEK_SET);
519                 bStop = TRUE;
520             }
521         }
522         if (bStop)
523         {
524             fprintf(stderr, "Do you REALLY want to truncate this trajectory (%s) at:\n"
525                     "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
526                     fn, j, t, static_cast<long int>(fpos));
527             if (1 != scanf("%s", yesno))
528             {
529                 gmx_fatal(FARGS, "Error reading user input");
530             }
531             if (std::strcmp(yesno, "YES") == 0)
532             {
533                 fprintf(stderr, "Once again, I'm gonna DO this...\n");
534                 gmx_trr_close(in);
535                 if (0 != gmx_truncate(fn, fpos))
536                 {
537                     gmx_fatal(FARGS, "Error truncating file %s", fn);
538                 }
539             }
540             else
541             {
542                 fprintf(stderr, "Ok, I'll forget about it\n");
543             }
544         }
545         else
546         {
547             fprintf(stderr, "Already at end of file (t=%g)...\n", t);
548             gmx_trr_close(in);
549         }
550     }
551 }
552
553 /*! \brief Read a full molecular topology if useful and available.
554  *
555  * If the input trajectory file is not in TNG format, and the output
556  * file is in TNG format, then we want to try to read a full topology
557  * (if available), so that we can write molecule information to the
558  * output file. The full topology provides better molecule information
559  * than is available from the normal t_topology data used by GROMACS
560  * tools.
561  *
562  * Also, the t_topology is only read under (different) particular
563  * conditions. If both apply, then a .tpr file might be read
564  * twice. Trying to fix this redundancy while trjconv is still an
565  * all-purpose tool does not seem worthwhile.
566  *
567  * Because of the way gmx_prepare_tng_writing is implemented, the case
568  * where the input TNG file has no molecule information will never
569  * lead to an output TNG file having molecule information. Since
570  * molecule information will generally be present if the input TNG
571  * file was written by a GROMACS tool, this seems like reasonable
572  * behaviour. */
573 static std::unique_ptr<gmx_mtop_t>
574 read_mtop_for_tng(const char *tps_file,
575                   const char *input_file,
576                   const char *output_file)
577 {
578     std::unique_ptr<gmx_mtop_t> mtop;
579
580     if (fn2bTPX(tps_file) &&
581         efTNG != fn2ftp(input_file) &&
582         efTNG == fn2ftp(output_file))
583     {
584         int temp_natoms = -1;
585         mtop = std::make_unique<gmx_mtop_t>();
586         read_tpx(tps_file, nullptr, nullptr, &temp_natoms,
587                  nullptr, nullptr, mtop.get());
588     }
589
590     return mtop;
591 }
592
593 int gmx_trjconv(int argc, char *argv[])
594 {
595     const char *desc[] = {
596         "[THISMODULE] can convert trajectory files in many ways:",
597         "",
598         "* from one format to another",
599         "* select a subset of atoms",
600         "* change the periodicity representation",
601         "* keep multimeric molecules together",
602         "* center atoms in the box",
603         "* fit atoms to reference structure",
604         "* reduce the number of frames",
605         "* change the timestamps of the frames ([TT]-t0[tt] and [TT]-timestep[tt])",
606         "* cut the trajectory in small subtrajectories according",
607         "  to information in an index file. This allows subsequent analysis of",
608         "  the subtrajectories that could, for example, be the result of a",
609         "  cluster analysis. Use option [TT]-sub[tt].",
610         "  This assumes that the entries in the index file are frame numbers and",
611         "  dumps each group in the index file to a separate trajectory file.",
612         "* select frames within a certain range of a quantity given",
613         "  in an [REF].xvg[ref] file.",
614         "",
615         "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
616         "[PAR]",
617
618         "The following formats are supported for input and output:",
619         "[REF].xtc[ref], [REF].trr[ref], [REF].gro[ref], [TT].g96[tt]",
620         "and [REF].pdb[ref].",
621         "The file formats are detected from the file extension.",
622         "The precision of the [REF].xtc[ref] output is taken from the",
623         "input file for [REF].xtc[ref], [REF].gro[ref] and [REF].pdb[ref],",
624         "and from the [TT]-ndec[tt] option for other input formats. The precision",
625         "is always taken from [TT]-ndec[tt], when this option is set.",
626         "All other formats have fixed precision. [REF].trr[ref]",
627         "output can be single or double precision, depending on the precision",
628         "of the [THISMODULE] binary.",
629         "Note that velocities are only supported in",
630         "[REF].trr[ref], [REF].gro[ref] and [TT].g96[tt] files.[PAR]",
631
632         "Option [TT]-sep[tt] can be used to write every frame to a separate",
633         "[TT].gro, .g96[tt] or [REF].pdb[ref] file. By default, all frames all written to one file.",
634         "[REF].pdb[ref] files with all frames concatenated can be viewed with",
635         "[TT]rasmol -nmrpdb[tt].[PAR]",
636
637         "It is possible to select part of your trajectory and write it out",
638         "to a new trajectory file in order to save disk space, e.g. for leaving",
639         "out the water from a trajectory of a protein in water.",
640         "[BB]ALWAYS[bb] put the original trajectory on tape!",
641         "We recommend to use the portable [REF].xtc[ref] format for your analysis",
642         "to save disk space and to have portable files.[PAR]",
643
644         "There are two options for fitting the trajectory to a reference",
645         "either for essential dynamics analysis, etc.",
646         "The first option is just plain fitting to a reference structure",
647         "in the structure file. The second option is a progressive fit",
648         "in which the first timeframe is fitted to the reference structure ",
649         "in the structure file to obtain and each subsequent timeframe is ",
650         "fitted to the previously fitted structure. This way a continuous",
651         "trajectory is generated, which might not be the case when using the",
652         "regular fit method, e.g. when your protein undergoes large",
653         "conformational transitions.[PAR]",
654
655         "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
656         "treatment:",
657         "",
658         " * [TT]mol[tt] puts the center of mass of molecules in the box,",
659         "   and requires a run input file to be supplied with [TT]-s[tt].",
660         " * [TT]res[tt] puts the center of mass of residues in the box.",
661         " * [TT]atom[tt] puts all the atoms in the box.",
662         " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
663         "   them back. This has the effect that all molecules",
664         "   will remain whole (provided they were whole in the initial",
665         "   conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
666         "   molecules may diffuse out of the box. The starting configuration",
667         "   for this procedure is taken from the structure file, if one is",
668         "   supplied, otherwise it is the first frame.",
669         " * [TT]cluster[tt] clusters all the atoms in the selected index",
670         "   such that they are all closest to the center of mass of the cluster,",
671         "   which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
672         "   results if you in fact have a cluster. Luckily that can be checked",
673         "   afterwards using a trajectory viewer. Note also that if your molecules",
674         "   are broken this will not work either.",
675         " * [TT]whole[tt] only makes broken molecules whole.",
676         "",
677
678         "Option [TT]-ur[tt] sets the unit cell representation for options",
679         "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
680         "All three options give different results for triclinic boxes and",
681         "identical results for rectangular boxes.",
682         "[TT]rect[tt] is the ordinary brick shape.",
683         "[TT]tric[tt] is the triclinic unit cell.",
684         "[TT]compact[tt] puts all atoms at the closest distance from the center",
685         "of the box. This can be useful for visualizing e.g. truncated octahedra",
686         "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
687         "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
688         "is set differently.[PAR]",
689
690         "Option [TT]-center[tt] centers the system in the box. The user can",
691         "select the group which is used to determine the geometrical center.",
692         "Option [TT]-boxcenter[tt] sets the location of the center of the box",
693         "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
694         "[TT]tric[tt]: half of the sum of the box vectors,",
695         "[TT]rect[tt]: half of the box diagonal,",
696         "[TT]zero[tt]: zero.",
697         "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
698         "want all molecules in the box after the centering.[PAR]",
699
700         "Option [TT]-box[tt] sets the size of the new box. This option only works",
701         "for leading dimensions and is thus generally only useful for rectangular boxes.",
702         "If you want to modify only some of the dimensions, e.g. when reading from",
703         "a trajectory, you can use -1 for those dimensions that should stay the same",
704
705         "It is not always possible to use combinations of [TT]-pbc[tt],",
706         "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
707         "you want in one call to [THISMODULE]. Consider using multiple",
708         "calls, and check out the GROMACS website for suggestions.[PAR]",
709
710         "With [TT]-dt[tt], it is possible to reduce the number of ",
711         "frames in the output. This option relies on the accuracy of the times",
712         "in your input trajectory, so if these are inaccurate use the",
713         "[TT]-timestep[tt] option to modify the time (this can be done",
714         "simultaneously). For making smooth movies, the program [gmx-filter]",
715         "can reduce the number of frames while using low-pass frequency",
716         "filtering, this reduces aliasing of high frequency motions.[PAR]",
717
718         "Using [TT]-trunc[tt] [THISMODULE] can truncate [REF].trr[ref] in place, i.e.",
719         "without copying the file. This is useful when a run has crashed",
720         "during disk I/O (i.e. full disk), or when two contiguous",
721         "trajectories must be concatenated without having double frames.[PAR]",
722
723         "Option [TT]-dump[tt] can be used to extract a frame at or near",
724         "one specific time from your trajectory, but only works reliably",
725         "if the time interval between frames is uniform.[PAR]",
726
727         "Option [TT]-drop[tt] reads an [REF].xvg[ref] file with times and values.",
728         "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
729         "frames with a value below and above the value of the respective options",
730         "will not be written."
731     };
732
733     int         pbc_enum;
734     enum
735     {
736         epSel,
737         epNone,
738         epComMol,
739         epComRes,
740         epComAtom,
741         epNojump,
742         epCluster,
743         epWhole,
744         epNR
745     };
746     const char *pbc_opt[epNR + 1] =
747     {
748         nullptr, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
749         nullptr
750     };
751
752     int         unitcell_enum;
753     const char *unitcell_opt[euNR+1] =
754     { nullptr, "rect", "tric", "compact", nullptr };
755
756     enum
757     {
758         ecSel, ecTric, ecRect, ecZero, ecNR
759     };
760     const char *center_opt[ecNR+1] =
761     { nullptr, "tric", "rect", "zero", nullptr };
762     int         ecenter;
763
764     int         fit_enum;
765     enum
766     {
767         efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
768     };
769     const char *fit[efNR + 1] =
770     {
771         nullptr, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
772         "progressive", nullptr
773     };
774
775     static gmx_bool  bSeparate     = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
776     static gmx_bool  bCenter       = FALSE;
777     static int       skip_nr       = 1, ndec = 3, nzero = 0;
778     static real      tzero         = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
779     static rvec      newbox        = {0, 0, 0}, shift = {0, 0, 0}, trans = {0, 0, 0};
780     static char     *exec_command  = nullptr;
781     static real      dropunder     = 0, dropover = 0;
782     static gmx_bool  bRound        = FALSE;
783
784     t_pargs
785         pa[] =
786     {
787         { "-skip", FALSE, etINT,
788           { &skip_nr }, "Only write every nr-th frame" },
789         { "-dt", FALSE, etTIME,
790           { &delta_t },
791           "Only write frame when t MOD dt = first time (%t)" },
792         { "-round", FALSE, etBOOL,
793           { &bRound }, "Round measurements to nearest picosecond"},
794         { "-dump", FALSE, etTIME,
795           { &tdump }, "Dump frame nearest specified time (%t)" },
796         { "-t0", FALSE, etTIME,
797           { &tzero },
798           "Starting time (%t) (default: don't change)" },
799         { "-timestep", FALSE, etTIME,
800           { &timestep },
801           "Change time step between input frames (%t)" },
802         { "-pbc", FALSE, etENUM,
803           { pbc_opt },
804           "PBC treatment (see help text for full description)" },
805         { "-ur", FALSE, etENUM,
806           { unitcell_opt }, "Unit-cell representation" },
807         { "-center", FALSE, etBOOL,
808           { &bCenter }, "Center atoms in box" },
809         { "-boxcenter", FALSE, etENUM,
810           { center_opt }, "Center for -pbc and -center" },
811         { "-box", FALSE, etRVEC,
812           { newbox },
813           "Size for new cubic box (default: read from input)" },
814         { "-trans", FALSE, etRVEC,
815           { trans },
816           "All coordinates will be translated by trans. This "
817           "can advantageously be combined with -pbc mol -ur "
818           "compact." },
819         { "-shift", FALSE, etRVEC,
820           { shift },
821           "All coordinates will be shifted by framenr*shift" },
822         { "-fit", FALSE, etENUM,
823           { fit },
824           "Fit molecule to ref structure in the structure file" },
825         { "-ndec", FALSE, etINT,
826           { &ndec },
827           "Number of decimal places to write to .xtc output" },
828         { "-vel", FALSE, etBOOL,
829           { &bVels }, "Read and write velocities if possible" },
830         { "-force", FALSE, etBOOL,
831           { &bForce }, "Read and write forces if possible" },
832         { "-trunc", FALSE, etTIME,
833           { &ttrunc },
834           "Truncate input trajectory file after this time (%t)" },
835         { "-exec", FALSE, etSTR,
836           { &exec_command },
837           "Execute command for every output frame with the "
838           "frame number as argument" },
839         { "-split", FALSE, etTIME,
840           { &split_t },
841           "Start writing new file when t MOD split = first "
842           "time (%t)" },
843         { "-sep", FALSE, etBOOL,
844           { &bSeparate },
845           "Write each frame to a separate .gro, .g96 or .pdb "
846           "file" },
847         { "-nzero", FALSE, etINT,
848           { &nzero },
849           "If the -sep flag is set, use these many digits "
850           "for the file numbers and prepend zeros as needed" },
851         { "-dropunder", FALSE, etREAL,
852           { &dropunder }, "Drop all frames below this value" },
853         { "-dropover", FALSE, etREAL,
854           { &dropover }, "Drop all frames above this value" },
855         { "-conect", FALSE, etBOOL,
856           { &bCONECT },
857           "Add conect records when writing [REF].pdb[ref] files. Useful "
858           "for visualization of non-standard molecules, e.g. "
859           "coarse grained ones" }
860     };
861 #define NPA asize(pa)
862
863     FILE             *out    = nullptr;
864     t_trxstatus      *trxout = nullptr;
865     t_trxstatus      *trxin;
866     int               file_nr;
867     t_trxframe        fr, frout;
868     int               flags;
869     rvec             *xmem  = nullptr, *vmem = nullptr, *fmem = nullptr;
870     rvec             *xp    = nullptr, x_shift, hbox;
871     real             *w_rls = nullptr;
872     int               m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
873 #define SKIP 10
874     t_topology       *top   = nullptr;
875     gmx_conect        gc    = nullptr;
876     int               ePBC  = -1;
877     t_atoms          *atoms = nullptr, useatoms;
878     matrix            top_box;
879     int              *index = nullptr, *cindex = nullptr;
880     char             *grpnm = nullptr;
881     int              *frindex, nrfri;
882     char             *frname;
883     int               ifit, my_clust = -1;
884     int              *ind_fit;
885     char             *gn_fit;
886     t_cluster_ndx    *clust           = nullptr;
887     t_trxstatus     **clust_status    = nullptr;
888     int              *clust_status_id = nullptr;
889     int               ntrxopen        = 0;
890     int              *nfwritten       = nullptr;
891     int               ndrop           = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
892     double          **dropval;
893     real              tshift = 0, dt = -1, prec;
894     gmx_bool          bFit, bPFit, bReset;
895     int               nfitdim;
896     gmx_rmpbc_t       gpbc = nullptr;
897     gmx_bool          bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
898     gmx_bool          bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
899     gmx_bool          bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetXtcPrec, bNeedPrec;
900     gmx_bool          bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
901     gmx_bool          bSubTraj = FALSE, bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
902     gmx_bool          bWriteFrame, bSplitHere;
903     const char       *top_file, *in_file, *out_file = nullptr;
904     char              out_file2[256], *charpt;
905     char             *outf_base = nullptr;
906     const char       *outf_ext  = nullptr;
907     char              top_title[256], timestr[32], stepstr[32], filemode[5];
908     gmx_output_env_t *oenv;
909
910     t_filenm          fnm[] = {
911         { efTRX, "-f",   nullptr,      ffREAD  },
912         { efTRO, "-o",   nullptr,      ffWRITE },
913         { efTPS, nullptr,   nullptr,      ffOPTRD },
914         { efNDX, nullptr,   nullptr,      ffOPTRD },
915         { efNDX, "-fr",  "frames",  ffOPTRD },
916         { efNDX, "-sub", "cluster", ffOPTRD },
917         { efXVG, "-drop", "drop",    ffOPTRD }
918     };
919 #define NFILE asize(fnm)
920
921     if (!parse_common_args(&argc, argv,
922                            PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW |
923                            PCA_TIME_UNIT,
924                            NFILE, fnm, NPA, pa, asize(desc), desc,
925                            0, nullptr, &oenv))
926     {
927         return 0;
928     }
929     fprintf(stdout, "Note that major changes are planned in future for "
930             "trjconv, to improve usability and utility.\n");
931
932     top_file = ftp2fn(efTPS, NFILE, fnm);
933
934     /* Check command line */
935     in_file = opt2fn("-f", NFILE, fnm);
936
937     if (ttrunc != -1)
938     {
939         do_trunc(in_file, ttrunc);
940     }
941     else
942     {
943         /* mark active cmdline options */
944         bSetBox     = opt2parg_bSet("-box", NPA, pa);
945         bSetTime    = opt2parg_bSet("-t0", NPA, pa);
946         bSetXtcPrec = opt2parg_bSet("-ndec", NPA, pa);
947         bSetUR      = opt2parg_bSet("-ur", NPA, pa);
948         bExec       = opt2parg_bSet("-exec", NPA, pa);
949         bTimeStep   = opt2parg_bSet("-timestep", NPA, pa);
950         bTDump      = opt2parg_bSet("-dump", NPA, pa);
951         bDropUnder  = opt2parg_bSet("-dropunder", NPA, pa);
952         bDropOver   = opt2parg_bSet("-dropover", NPA, pa);
953         bTrans      = opt2parg_bSet("-trans", NPA, pa);
954         bSplit      = (split_t != 0);
955
956         /* parse enum options */
957         fit_enum      = nenum(fit);
958         bFit          = (fit_enum == efFit || fit_enum == efFitXY);
959         bReset        = (fit_enum == efReset || fit_enum == efResetXY);
960         bPFit         = fit_enum == efPFit;
961         pbc_enum      = nenum(pbc_opt);
962         bPBCWhole     = pbc_enum == epWhole;
963         bPBCcomRes    = pbc_enum == epComRes;
964         bPBCcomMol    = pbc_enum == epComMol;
965         bPBCcomAtom   = pbc_enum == epComAtom;
966         bNoJump       = pbc_enum == epNojump;
967         bCluster      = pbc_enum == epCluster;
968         bPBC          = pbc_enum != epNone;
969         unitcell_enum = nenum(unitcell_opt);
970         ecenter       = nenum(center_opt) - ecTric;
971
972         /* set and check option dependencies */
973         if (bPFit)
974         {
975             bFit = TRUE;        /* for pfit, fit *must* be set */
976         }
977         if (bFit)
978         {
979             bReset = TRUE;       /* for fit, reset *must* be set */
980         }
981         nfitdim = 0;
982         if (bFit || bReset)
983         {
984             nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
985         }
986         bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
987
988         if (bSetUR)
989         {
990             if (!(bPBCcomRes || bPBCcomMol ||  bPBCcomAtom))
991             {
992                 fprintf(stderr,
993                         "WARNING: Option for unitcell representation (-ur %s)\n"
994                         "         only has effect in combination with -pbc %s, %s or %s.\n"
995                         "         Ingoring unitcell representation.\n\n",
996                         unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
997             }
998         }
999         if (bFit && bPBC)
1000         {
1001             gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
1002                       "Please do the PBC condition treatment first and then run trjconv in a second step\n"
1003                       "for the rotational fit.\n"
1004                       "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1005                       "results!");
1006         }
1007
1008         /* ndec for XTC writing is in nr of decimal places, prec is a multiplication factor: */
1009         prec = 1;
1010         for (i = 0; i < ndec; i++)
1011         {
1012             prec *= 10;
1013         }
1014
1015         bIndex = ftp2bSet(efNDX, NFILE, fnm);
1016
1017
1018         /* Determine output type */
1019         out_file = opt2fn("-o", NFILE, fnm);
1020         int ftp  = fn2ftp(out_file);
1021         fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
1022         bNeedPrec = (ftp == efXTC);
1023         int ftpin = fn2ftp(in_file);
1024         if (bVels)
1025         {
1026             /* check if velocities are possible in input and output files */
1027             bVels = (ftp == efTRR || ftp == efGRO ||
1028                      ftp == efG96 || ftp == efTNG)
1029                 && (ftpin == efTRR || ftpin == efGRO ||
1030                     ftpin == efG96 || ftpin == efTNG || ftpin == efCPT);
1031         }
1032         if (bSeparate || bSplit)
1033         {
1034             outf_ext = std::strrchr(out_file, '.');
1035             if (outf_ext == nullptr)
1036             {
1037                 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1038             }
1039             outf_base = gmx_strdup(out_file);
1040             outf_base[outf_ext - out_file] = '\0';
1041         }
1042
1043         bSubTraj = opt2bSet("-sub", NFILE, fnm);
1044         if (bSubTraj)
1045         {
1046             if ((ftp != efXTC) && (ftp != efTRR))
1047             {
1048                 /* It seems likely that other trajectory file types
1049                  * could work here. */
1050                 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1051                           "xtc and trr");
1052             }
1053             clust = cluster_index(nullptr, opt2fn("-sub", NFILE, fnm));
1054
1055             /* Check for number of files disabled, as FOPEN_MAX is not the correct
1056              * number to check for. In my linux box it is only 16.
1057              */
1058             if (/* DISABLES CODE */ (false))
1059             //if (clust->clust->nr > FOPEN_MAX-4)
1060             {
1061                 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1062                           " trajectories.\ntry splitting the index file in %d parts.\n"
1063                           "FOPEN_MAX = %d",
1064                           clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1065             }
1066             gmx_warning("The -sub option could require as many open output files as there are\n"
1067                         "index groups in the file (%d). If you get I/O errors opening new files,\n"
1068                         "try reducing the number of index groups in the file, and perhaps\n"
1069                         "using trjconv -sub several times on different chunks of your index file.\n",
1070                         clust->clust->nr);
1071
1072             snew(clust_status, clust->clust->nr);
1073             snew(clust_status_id, clust->clust->nr);
1074             snew(nfwritten, clust->clust->nr);
1075             for (i = 0; (i < clust->clust->nr); i++)
1076             {
1077                 clust_status[i]    = nullptr;
1078                 clust_status_id[i] = -1;
1079             }
1080             bSeparate = bSplit = FALSE;
1081         }
1082         /* skipping */
1083         if (skip_nr <= 0)
1084         {
1085             gmx_fatal(FARGS, "Argument for -skip (%d) needs to be greater or equal to 1.", skip_nr);
1086         }
1087
1088         std::unique_ptr<gmx_mtop_t> mtop = read_mtop_for_tng(top_file, in_file, out_file);
1089
1090         /* Determine whether to read a topology */
1091         bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1092                 bRmPBC || bReset || bPBCcomMol || bCluster ||
1093                 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1094
1095         /* Determine if when can read index groups */
1096         bIndex = (bIndex || bTPS);
1097
1098         if (bTPS)
1099         {
1100             snew(top, 1);
1101             read_tps_conf(top_file, top, &ePBC, &xp, nullptr, top_box,
1102                           bReset || bPBCcomRes);
1103             std::strncpy(top_title, *top->name, 255);
1104             top_title[255] = '\0';
1105             atoms          = &top->atoms;
1106
1107             if (0 == top->mols.nr && (bCluster || bPBCcomMol))
1108             {
1109                 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1110             }
1111
1112             /* top_title is only used for gro and pdb,
1113              * the header in such a file is top_title, followed by
1114              * t= ... and/or step= ...
1115              * to prevent double t= or step=, remove it from top_title.
1116              * From GROMACS-2018 we only write t/step when the frame actually
1117              * has a valid time/step, so we need to check for both separately.
1118              */
1119             if ((charpt = std::strstr(top_title, " t= ")))
1120             {
1121                 charpt[0] = '\0';
1122             }
1123             if ((charpt = std::strstr(top_title, " step= ")))
1124             {
1125                 charpt[0] = '\0';
1126             }
1127
1128             if (bCONECT)
1129             {
1130                 gc = gmx_conect_generate(top);
1131             }
1132             if (bRmPBC)
1133             {
1134                 gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
1135             }
1136         }
1137
1138         /* get frame number index */
1139         frindex = nullptr;
1140         if (opt2bSet("-fr", NFILE, fnm))
1141         {
1142             printf("Select groups of frame number indices:\n");
1143             rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, &frindex, &frname);
1144             if (debug)
1145             {
1146                 for (i = 0; i < nrfri; i++)
1147                 {
1148                     fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1149                 }
1150             }
1151         }
1152
1153         /* get index groups etc. */
1154         if (bReset)
1155         {
1156             printf("Select group for %s fit\n",
1157                    bFit ? "least squares" : "translational");
1158             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1159                       1, &ifit, &ind_fit, &gn_fit);
1160
1161             if (bFit)
1162             {
1163                 if (ifit < 2)
1164                 {
1165                     gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1166                 }
1167                 else if (ifit == 3)
1168                 {
1169                     fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1170                 }
1171             }
1172         }
1173         else if (bCluster)
1174         {
1175             printf("Select group for clustering\n");
1176             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1177                       1, &ifit, &ind_fit, &gn_fit);
1178         }
1179
1180         if (bIndex)
1181         {
1182             if (bCenter)
1183             {
1184                 printf("Select group for centering\n");
1185                 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1186                           1, &ncent, &cindex, &grpnm);
1187             }
1188             printf("Select group for output\n");
1189             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1190                       1, &nout, &index, &grpnm);
1191         }
1192         else
1193         {
1194             /* no index file, so read natoms from TRX */
1195             if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
1196             {
1197                 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1198             }
1199             natoms = fr.natoms;
1200             close_trx(trxin);
1201             sfree(fr.x);
1202             snew(index, natoms);
1203             for (i = 0; i < natoms; i++)
1204             {
1205                 index[i] = i;
1206             }
1207             nout = natoms;
1208             if (bCenter)
1209             {
1210                 ncent  = nout;
1211                 cindex = index;
1212             }
1213         }
1214
1215         if (bReset)
1216         {
1217             snew(w_rls, atoms->nr);
1218             for (i = 0; (i < ifit); i++)
1219             {
1220                 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1221             }
1222
1223             /* Restore reference structure and set to origin,
1224                store original location (to put structure back) */
1225             if (bRmPBC)
1226             {
1227                 gmx_rmpbc(gpbc, top->atoms.nr, top_box, xp);
1228             }
1229             copy_rvec(xp[index[0]], x_shift);
1230             reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, nullptr, xp, w_rls);
1231             rvec_dec(x_shift, xp[index[0]]);
1232         }
1233         else
1234         {
1235             clear_rvec(x_shift);
1236         }
1237
1238         if (bDropUnder || bDropOver)
1239         {
1240             /* Read the .xvg file with the drop values */
1241             fprintf(stderr, "\nReading drop file ...");
1242             ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1243             fprintf(stderr, " %d time points\n", ndrop);
1244             if (ndrop == 0 || ncol < 2)
1245             {
1246                 gmx_fatal(FARGS, "Found no data points in %s",
1247                           opt2fn("-drop", NFILE, fnm));
1248             }
1249             drop0 = 0;
1250             drop1 = 0;
1251         }
1252
1253         /* Make atoms struct for output in GRO or PDB files */
1254         if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1255         {
1256             /* get memory for stuff to go in .pdb file, and initialize
1257              * the pdbinfo structure part if the input has it.
1258              */
1259             init_t_atoms(&useatoms, atoms->nr, atoms->havePdbInfo);
1260             sfree(useatoms.resinfo);
1261             useatoms.resinfo = atoms->resinfo;
1262             for (i = 0; (i < nout); i++)
1263             {
1264                 useatoms.atomname[i] = atoms->atomname[index[i]];
1265                 useatoms.atom[i]     = atoms->atom[index[i]];
1266                 if (atoms->havePdbInfo)
1267                 {
1268                     useatoms.pdbinfo[i]  = atoms->pdbinfo[index[i]];
1269                 }
1270                 useatoms.nres        = std::max(useatoms.nres, useatoms.atom[i].resind+1);
1271             }
1272             useatoms.nr = nout;
1273         }
1274         /* select what to read */
1275         if (ftp == efTRR)
1276         {
1277             flags = TRX_READ_X;
1278         }
1279         else
1280         {
1281             flags = TRX_NEED_X;
1282         }
1283         if (bVels)
1284         {
1285             flags = flags | TRX_READ_V;
1286         }
1287         if (bForce)
1288         {
1289             flags = flags | TRX_READ_F;
1290         }
1291
1292         /* open trx file for reading */
1293         bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
1294         if (fr.bPrec)
1295         {
1296             fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1297         }
1298         if (bNeedPrec)
1299         {
1300             if (bSetXtcPrec || !fr.bPrec)
1301             {
1302                 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1303             }
1304             else
1305             {
1306                 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1307             }
1308         }
1309
1310         if (bHaveFirstFrame)
1311         {
1312             if (bTDump)
1313             {
1314                 // Determine timestep (assuming constant spacing for now) if we
1315                 // need to dump frames based on time. This is required so we do not
1316                 // skip the first frame if that was the one that should have been dumped
1317                 double firstFrameTime = fr.time;
1318                 if (read_next_frame(oenv, trxin, &fr))
1319                 {
1320                     dt     = fr.time - firstFrameTime;
1321                     bDTset = TRUE;
1322                     if (dt <= 0)
1323                     {
1324                         fprintf(stderr, "Warning: Frame times are not incrementing - will dump first frame.\n");
1325                     }
1326                 }
1327                 // Now close and reopen so we are at first frame again
1328                 close_trx(trxin);
1329                 done_frame(&fr);
1330                 // Reopen at first frame (We already know it exists if we got here)
1331                 read_first_frame(oenv, &trxin, in_file, &fr, flags);
1332             }
1333
1334             set_trxframe_ePBC(&fr, ePBC);
1335             natoms = fr.natoms;
1336
1337             if (bSetTime)
1338             {
1339                 tshift = tzero-fr.time;
1340             }
1341             else
1342             {
1343                 tzero = fr.time;
1344             }
1345
1346             bCopy = FALSE;
1347             if (bIndex)
1348             {
1349                 /* check if index is meaningful */
1350                 for (i = 0; i < nout; i++)
1351                 {
1352                     if (index[i] >= natoms)
1353                     {
1354                         gmx_fatal(FARGS,
1355                                   "Index[%d] %d is larger than the number of atoms in the\n"
1356                                   "trajectory file (%d). There is a mismatch in the contents\n"
1357                                   "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1358                     }
1359                     bCopy = bCopy || (i != index[i]);
1360                 }
1361             }
1362
1363             /* open output for writing */
1364             std::strcpy(filemode, "w");
1365             switch (ftp)
1366             {
1367                 case efTNG:
1368                     trxout = trjtools_gmx_prepare_tng_writing(out_file,
1369                                                               filemode[0],
1370                                                               trxin,
1371                                                               nullptr,
1372                                                               nout,
1373                                                               mtop.get(),
1374                                                               gmx::arrayRefFromArray(index, nout),
1375                                                               grpnm);
1376                     break;
1377                 case efXTC:
1378                 case efTRR:
1379                     out = nullptr;
1380                     if (!bSplit && !bSubTraj)
1381                     {
1382                         trxout = open_trx(out_file, filemode);
1383                     }
1384                     break;
1385                 case efGRO:
1386                 case efG96:
1387                 case efPDB:
1388                     if (( !bSeparate && !bSplit ) && !bSubTraj)
1389                     {
1390                         out = gmx_ffopen(out_file, filemode);
1391                     }
1392                     break;
1393                 default:
1394                     gmx_incons("Illegal output file format");
1395             }
1396
1397             if (bCopy)
1398             {
1399                 snew(xmem, nout);
1400                 if (bVels)
1401                 {
1402                     snew(vmem, nout);
1403                 }
1404                 if (bForce)
1405                 {
1406                     snew(fmem, nout);
1407                 }
1408             }
1409
1410             /* Start the big loop over frames */
1411             file_nr  =  0;
1412             frame    =  0;
1413             outframe =  0;
1414             model_nr =  0;
1415
1416             /* Main loop over frames */
1417             do
1418             {
1419                 if (!fr.bStep)
1420                 {
1421                     /* set the step */
1422                     fr.step = newstep;
1423                     newstep++;
1424                 }
1425                 if (bSubTraj)
1426                 {
1427                     /*if (frame >= clust->clust->nra)
1428                        gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1429                     if (frame > clust->maxframe)
1430                     {
1431                         my_clust = -1;
1432                     }
1433                     else
1434                     {
1435                         my_clust = clust->inv_clust[frame];
1436                     }
1437                     if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1438                         (my_clust == -1))
1439                     {
1440                         my_clust = -1;
1441                     }
1442                 }
1443
1444                 if (bSetBox)
1445                 {
1446                     /* generate new box */
1447                     if (!fr.bBox)
1448                     {
1449                         clear_mat(fr.box);
1450                     }
1451                     for (m = 0; m < DIM; m++)
1452                     {
1453                         if (newbox[m] >= 0)
1454                         {
1455                             fr.box[m][m] = newbox[m];
1456                         }
1457                         else
1458                         {
1459                             if (!fr.bBox)
1460                             {
1461                                 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1462                             }
1463                         }
1464                     }
1465                 }
1466
1467                 if (bTrans)
1468                 {
1469                     for (i = 0; i < natoms; i++)
1470                     {
1471                         rvec_inc(fr.x[i], trans);
1472                     }
1473                 }
1474
1475                 if (bTDump)
1476                 {
1477                     // If we could not read two frames or times are not incrementing
1478                     // we have almost no idea what to do,
1479                     // but dump the first frame so output is not broken.
1480                     if (dt <= 0 || !bDTset)
1481                     {
1482                         bDumpFrame = true;
1483                     }
1484                     else
1485                     {
1486                         // Dump the frame if we are less than half a frame time
1487                         // below it. This will also ensure we at least dump a
1488                         // somewhat reasonable frame if the spacing is unequal
1489                         // and we have overrun the frame time. Once we dump one
1490                         // frame based on time we quit, so it does not matter
1491                         // that this might be true for all subsequent frames too.
1492                         bDumpFrame = (fr.time > tdump-0.5*dt);
1493                     }
1494                 }
1495                 else
1496                 {
1497                     bDumpFrame = FALSE;
1498                 }
1499
1500                 /* determine if an atom jumped across the box and reset it if so */
1501                 if (bNoJump && (bTPS || frame != 0))
1502                 {
1503                     for (d = 0; d < DIM; d++)
1504                     {
1505                         hbox[d] = 0.5*fr.box[d][d];
1506                     }
1507                     for (i = 0; i < natoms; i++)
1508                     {
1509                         if (bReset)
1510                         {
1511                             rvec_dec(fr.x[i], x_shift);
1512                         }
1513                         for (m = DIM-1; m >= 0; m--)
1514                         {
1515                             if (hbox[m] > 0)
1516                             {
1517                                 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1518                                 {
1519                                     for (d = 0; d <= m; d++)
1520                                     {
1521                                         fr.x[i][d] += fr.box[m][d];
1522                                     }
1523                                 }
1524                                 while (fr.x[i][m]-xp[i][m] > hbox[m])
1525                                 {
1526                                     for (d = 0; d <= m; d++)
1527                                     {
1528                                         fr.x[i][d] -= fr.box[m][d];
1529                                     }
1530                                 }
1531                             }
1532                         }
1533                     }
1534                 }
1535                 else if (bCluster)
1536                 {
1537                     calc_pbc_cluster(ecenter, ifit, top, ePBC, fr.x, ind_fit, fr.box);
1538                 }
1539
1540                 if (bPFit)
1541                 {
1542                     /* Now modify the coords according to the flags,
1543                        for normal fit, this is only done for output frames */
1544                     if (bRmPBC)
1545                     {
1546                         gmx_rmpbc_trxfr(gpbc, &fr);
1547                     }
1548
1549                     reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1550                     do_fit(natoms, w_rls, xp, fr.x);
1551                 }
1552
1553                 /* store this set of coordinates for future use */
1554                 if (bPFit || bNoJump)
1555                 {
1556                     if (xp == nullptr)
1557                     {
1558                         snew(xp, natoms);
1559                     }
1560                     for (i = 0; (i < natoms); i++)
1561                     {
1562                         copy_rvec(fr.x[i], xp[i]);
1563                         rvec_inc(fr.x[i], x_shift);
1564                     }
1565                 }
1566
1567                 if (frindex)
1568                 {
1569                     /* see if we have a frame from the frame index group */
1570                     for (i = 0; i < nrfri && !bDumpFrame; i++)
1571                     {
1572                         bDumpFrame = frame == frindex[i];
1573                     }
1574                 }
1575                 if (debug && bDumpFrame)
1576                 {
1577                     fprintf(debug, "dumping %d\n", frame);
1578                 }
1579
1580                 bWriteFrame =
1581                     ( ( !bTDump && (frindex == nullptr) && frame % skip_nr == 0 ) || bDumpFrame );
1582
1583                 if (bWriteFrame && (bDropUnder || bDropOver))
1584                 {
1585                     while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1586                     {
1587                         drop0 = drop1;
1588                         drop1++;
1589                     }
1590                     if (std::abs(dropval[0][drop0] - fr.time)
1591                         < std::abs(dropval[0][drop1] - fr.time))
1592                     {
1593                         dropuse = drop0;
1594                     }
1595                     else
1596                     {
1597                         dropuse = drop1;
1598                     }
1599                     if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1600                         (bDropOver && dropval[1][dropuse] > dropover))
1601                     {
1602                         bWriteFrame = FALSE;
1603                     }
1604                 }
1605
1606                 if (bWriteFrame)
1607                 {
1608                     /* We should avoid modifying the input frame,
1609                      * but since here we don't have the output frame yet,
1610                      * we introduce a temporary output frame time variable.
1611                      */
1612                     real frout_time;
1613
1614                     frout_time = fr.time;
1615
1616                     /* calc new time */
1617                     if (bTimeStep)
1618                     {
1619                         frout_time = tzero + frame*timestep;
1620                     }
1621                     else if (bSetTime)
1622                     {
1623                         frout_time += tshift;
1624                     }
1625
1626                     if (bTDump)
1627                     {
1628                         fprintf(stderr, "\nDumping frame at t= %g %s\n",
1629                                 output_env_conv_time(oenv, frout_time), output_env_get_time_unit(oenv).c_str());
1630                     }
1631
1632                     /* check for writing at each delta_t */
1633                     bDoIt = (delta_t == 0);
1634                     if (!bDoIt)
1635                     {
1636                         if (!bRound)
1637                         {
1638                             bDoIt = bRmod(frout_time, tzero, delta_t);
1639                         }
1640                         else
1641                         {
1642                             /* round() is not C89 compatible, so we do this:  */
1643                             bDoIt = bRmod(std::floor(frout_time+0.5), std::floor(tzero+0.5),
1644                                           std::floor(delta_t+0.5));
1645                         }
1646                     }
1647
1648                     if (bDoIt || bTDump)
1649                     {
1650                         /* print sometimes */
1651                         if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1652                         {
1653                             fprintf(stderr, " ->  frame %6d time %8.3f      \r",
1654                                     outframe, output_env_conv_time(oenv, frout_time));
1655                             fflush(stderr);
1656                         }
1657
1658                         if (!bPFit)
1659                         {
1660                             /* Now modify the coords according to the flags,
1661                                for PFit we did this already! */
1662
1663                             if (bRmPBC)
1664                             {
1665                                 gmx_rmpbc_trxfr(gpbc, &fr);
1666                             }
1667
1668                             if (bReset)
1669                             {
1670                                 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1671                                 if (bFit)
1672                                 {
1673                                     do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1674                                 }
1675                                 if (!bCenter)
1676                                 {
1677                                     for (i = 0; i < natoms; i++)
1678                                     {
1679                                         rvec_inc(fr.x[i], x_shift);
1680                                     }
1681                                 }
1682                             }
1683
1684                             if (bCenter)
1685                             {
1686                                 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1687                             }
1688                         }
1689
1690                         auto positionsArrayRef = gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec *>(fr.x), natoms);
1691                         if (bPBCcomAtom)
1692                         {
1693                             switch (unitcell_enum)
1694                             {
1695                                 case euRect:
1696                                     put_atoms_in_box(ePBC, fr.box, positionsArrayRef);
1697                                     break;
1698                                 case euTric:
1699                                     put_atoms_in_triclinic_unitcell(ecenter, fr.box, positionsArrayRef);
1700                                     break;
1701                                 case euCompact:
1702                                     put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1703                                                                   positionsArrayRef);
1704                                     break;
1705                             }
1706                         }
1707                         if (bPBCcomRes)
1708                         {
1709                             put_residue_com_in_box(unitcell_enum, ecenter,
1710                                                    natoms, atoms->atom, ePBC, fr.box, fr.x);
1711                         }
1712                         if (bPBCcomMol)
1713                         {
1714                             put_molecule_com_in_box(unitcell_enum, ecenter,
1715                                                     &top->mols,
1716                                                     natoms, atoms->atom, ePBC, fr.box, fr.x);
1717                         }
1718                         /* Copy the input trxframe struct to the output trxframe struct */
1719                         frout        = fr;
1720                         frout.time   = frout_time;
1721                         frout.bV     = (frout.bV && bVels);
1722                         frout.bF     = (frout.bF && bForce);
1723                         frout.natoms = nout;
1724                         if (bNeedPrec && (bSetXtcPrec || !fr.bPrec))
1725                         {
1726                             frout.bPrec = TRUE;
1727                             frout.prec  = prec;
1728                         }
1729                         if (bCopy)
1730                         {
1731                             frout.x = xmem;
1732                             if (frout.bV)
1733                             {
1734                                 frout.v = vmem;
1735                             }
1736                             if (frout.bF)
1737                             {
1738                                 frout.f = fmem;
1739                             }
1740                             for (i = 0; i < nout; i++)
1741                             {
1742                                 copy_rvec(fr.x[index[i]], frout.x[i]);
1743                                 if (frout.bV)
1744                                 {
1745                                     copy_rvec(fr.v[index[i]], frout.v[i]);
1746                                 }
1747                                 if (frout.bF)
1748                                 {
1749                                     copy_rvec(fr.f[index[i]], frout.f[i]);
1750                                 }
1751                             }
1752                         }
1753
1754                         if (opt2parg_bSet("-shift", NPA, pa))
1755                         {
1756                             for (i = 0; i < nout; i++)
1757                             {
1758                                 for (d = 0; d < DIM; d++)
1759                                 {
1760                                     frout.x[i][d] += outframe*shift[d];
1761                                 }
1762                             }
1763                         }
1764
1765                         if (!bRound)
1766                         {
1767                             bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1768                         }
1769                         else
1770                         {
1771                             /* round() is not C89 compatible, so we do this: */
1772                             bSplitHere = bSplit && bRmod(std::floor(frout.time+0.5),
1773                                                          std::floor(tzero+0.5),
1774                                                          std::floor(split_t+0.5));
1775                         }
1776                         if (bSeparate || bSplitHere)
1777                         {
1778                             mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1779                         }
1780
1781                         std::string title;
1782                         switch (ftp)
1783                         {
1784                             case efTNG:
1785                                 write_tng_frame(trxout, &frout);
1786                                 // TODO when trjconv behaves better: work how to read and write lambda
1787                                 break;
1788                             case efTRR:
1789                             case efXTC:
1790                                 if (bSplitHere)
1791                                 {
1792                                     if (trxout)
1793                                     {
1794                                         close_trx(trxout);
1795                                     }
1796                                     trxout = open_trx(out_file2, filemode);
1797                                 }
1798                                 if (bSubTraj)
1799                                 {
1800                                     if (my_clust != -1)
1801                                     {
1802                                         char buf[STRLEN];
1803                                         if (clust_status_id[my_clust] == -1)
1804                                         {
1805                                             sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1806                                             clust_status[my_clust]    = open_trx(buf, "w");
1807                                             clust_status_id[my_clust] = 1;
1808                                             ntrxopen++;
1809                                         }
1810                                         else if (clust_status_id[my_clust] == -2)
1811                                         {
1812                                             gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1813                                                       clust->grpname[my_clust], ntrxopen, frame,
1814                                                       my_clust);
1815                                         }
1816                                         write_trxframe(clust_status[my_clust], &frout, gc);
1817                                         nfwritten[my_clust]++;
1818                                         if (nfwritten[my_clust] ==
1819                                             (clust->clust->index[my_clust+1]-
1820                                              clust->clust->index[my_clust]))
1821                                         {
1822                                             close_trx(clust_status[my_clust]);
1823                                             clust_status[my_clust]    = nullptr;
1824                                             clust_status_id[my_clust] = -2;
1825                                             ntrxopen--;
1826                                             if (ntrxopen < 0)
1827                                             {
1828                                                 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1829                                             }
1830                                         }
1831                                     }
1832                                 }
1833                                 else
1834                                 {
1835                                     write_trxframe(trxout, &frout, gc);
1836                                 }
1837                                 break;
1838                             case efGRO:
1839                             case efG96:
1840                             case efPDB:
1841                                 // Only add a generator statement if title is empty,
1842                                 // to avoid multiple generated-by statements from various programs
1843                                 if (std::strlen(top_title) == 0)
1844                                 {
1845                                     sprintf(top_title, "Generated by trjconv");
1846                                 }
1847                                 if (frout.bTime)
1848                                 {
1849                                     sprintf(timestr, " t= %9.5f", frout.time);
1850                                 }
1851                                 else
1852                                 {
1853                                     std::strcpy(timestr, "");
1854                                 }
1855                                 if (frout.bStep)
1856                                 {
1857                                     sprintf(stepstr, " step= %" PRId64, frout.step);
1858                                 }
1859                                 else
1860                                 {
1861                                     std::strcpy(stepstr, "");
1862                                 }
1863                                 title = gmx::formatString("%s%s%s", top_title, timestr, stepstr);
1864                                 if (bSeparate || bSplitHere)
1865                                 {
1866                                     out = gmx_ffopen(out_file2, "w");
1867                                 }
1868                                 switch (ftp)
1869                                 {
1870                                     case efGRO:
1871                                         write_hconf_p(out, title.c_str(), &useatoms,
1872                                                       frout.x, frout.bV ? frout.v : nullptr, frout.box);
1873                                         break;
1874                                     case efPDB:
1875                                         fprintf(out, "REMARK    GENERATED BY TRJCONV\n");
1876                                         /* if reading from pdb, we want to keep the original
1877                                            model numbering else we write the output frame
1878                                            number plus one, because model 0 is not allowed in pdb */
1879                                         if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1880                                         {
1881                                             model_nr = fr.step;
1882                                         }
1883                                         else
1884                                         {
1885                                             model_nr++;
1886                                         }
1887                                         write_pdbfile(out, title.c_str(), &useatoms, frout.x,
1888                                                       frout.ePBC, frout.box, ' ', model_nr, gc);
1889                                         break;
1890                                     case efG96:
1891                                         const char *outputTitle = "";
1892                                         if (bSeparate || bTDump)
1893                                         {
1894                                             outputTitle = title.c_str();
1895                                             if (bTPS)
1896                                             {
1897                                                 frout.bAtoms = TRUE;
1898                                             }
1899                                             frout.atoms  = &useatoms;
1900                                             frout.bStep  = FALSE;
1901                                             frout.bTime  = FALSE;
1902                                         }
1903                                         else
1904                                         {
1905                                             if (outframe == 0)
1906                                             {
1907                                                 outputTitle = title.c_str();
1908                                             }
1909                                             frout.bAtoms = FALSE;
1910                                             frout.bStep  = TRUE;
1911                                             frout.bTime  = TRUE;
1912                                         }
1913                                         write_g96_conf(out, outputTitle, &frout, -1, nullptr);
1914                                 }
1915                                 if (bSeparate || bSplitHere)
1916                                 {
1917                                     gmx_ffclose(out);
1918                                     out = nullptr;
1919                                 }
1920                                 break;
1921                             default:
1922                                 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1923                         }
1924                         if (bSeparate || bSplitHere)
1925                         {
1926                             file_nr++;
1927                         }
1928
1929                         /* execute command */
1930                         if (bExec)
1931                         {
1932                             char c[255];
1933                             sprintf(c, "%s  %d", exec_command, file_nr-1);
1934                             /*fprintf(stderr,"Executing '%s'\n",c);*/
1935                             if (0 != system(c))
1936                             {
1937                                 gmx_fatal(FARGS, "Error executing command: %s", c);
1938                             }
1939                         }
1940                         outframe++;
1941                     }
1942                 }
1943                 frame++;
1944                 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1945             }
1946             while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1947         }
1948
1949         if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1950         {
1951             fprintf(stderr, "\nWARNING no output, "
1952                     "last frame read at t=%g\n", fr.time);
1953         }
1954         fprintf(stderr, "\n");
1955
1956         close_trx(trxin);
1957         sfree(outf_base);
1958
1959         if (bRmPBC)
1960         {
1961             gmx_rmpbc_done(gpbc);
1962         }
1963
1964         if (trxout)
1965         {
1966             close_trx(trxout);
1967         }
1968         else if (out != nullptr)
1969         {
1970             gmx_ffclose(out);
1971         }
1972         if (bSubTraj)
1973         {
1974             for (i = 0; (i < clust->clust->nr); i++)
1975             {
1976                 if (clust_status_id[i] >= 0)
1977                 {
1978                     close_trx(clust_status[i]);
1979                 }
1980             }
1981         }
1982     }
1983
1984     if (bTPS)
1985     {
1986         done_top(top);
1987         sfree(top);
1988     }
1989     sfree(xp);
1990     sfree(xmem);
1991     sfree(vmem);
1992     sfree(fmem);
1993     sfree(grpnm);
1994     sfree(index);
1995     sfree(cindex);
1996     done_frame(&fr);
1997
1998     do_view(oenv, out_file, nullptr);
1999
2000     output_env_done(oenv);
2001     return 0;
2002 }