55771c8e17bccb0378c6208a2d69dbaaeedff215
[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, 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
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/commandline/viewit.h"
47 #include "gromacs/compat/make_unique.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))
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 = gmx::compat::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           "Precision for .xtc and .gro writing in number of "
828           "decimal places" },
829         { "-vel", FALSE, etBOOL,
830           { &bVels }, "Read and write velocities if possible" },
831         { "-force", FALSE, etBOOL,
832           { &bForce }, "Read and write forces if possible" },
833         { "-trunc", FALSE, etTIME,
834           { &ttrunc },
835           "Truncate input trajectory file after this time (%t)" },
836         { "-exec", FALSE, etSTR,
837           { &exec_command },
838           "Execute command for every output frame with the "
839           "frame number as argument" },
840         { "-split", FALSE, etTIME,
841           { &split_t },
842           "Start writing new file when t MOD split = first "
843           "time (%t)" },
844         { "-sep", FALSE, etBOOL,
845           { &bSeparate },
846           "Write each frame to a separate .gro, .g96 or .pdb "
847           "file" },
848         { "-nzero", FALSE, etINT,
849           { &nzero },
850           "If the -sep flag is set, use these many digits "
851           "for the file numbers and prepend zeros as needed" },
852         { "-dropunder", FALSE, etREAL,
853           { &dropunder }, "Drop all frames below this value" },
854         { "-dropover", FALSE, etREAL,
855           { &dropover }, "Drop all frames above this value" },
856         { "-conect", FALSE, etBOOL,
857           { &bCONECT },
858           "Add conect records when writing [REF].pdb[ref] files. Useful "
859           "for visualization of non-standard molecules, e.g. "
860           "coarse grained ones" }
861     };
862 #define NPA asize(pa)
863
864     FILE             *out    = nullptr;
865     t_trxstatus      *trxout = nullptr;
866     t_trxstatus      *trxin;
867     int               file_nr;
868     t_trxframe        fr, frout;
869     int               flags;
870     rvec             *xmem  = nullptr, *vmem = nullptr, *fmem = nullptr;
871     rvec             *xp    = nullptr, x_shift, hbox;
872     real             *w_rls = nullptr;
873     int               m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
874 #define SKIP 10
875     t_topology        top;
876     gmx_conect        gc    = nullptr;
877     int               ePBC  = -1;
878     t_atoms          *atoms = nullptr, useatoms;
879     matrix            top_box;
880     int              *index = nullptr, *cindex = nullptr;
881     char             *grpnm = nullptr;
882     int              *frindex, nrfri;
883     char             *frname;
884     int               ifit, my_clust = -1;
885     int              *ind_fit;
886     char             *gn_fit;
887     t_cluster_ndx    *clust           = nullptr;
888     t_trxstatus     **clust_status    = nullptr;
889     int              *clust_status_id = nullptr;
890     int               ntrxopen        = 0;
891     int              *nfwritten       = nullptr;
892     int               ndrop           = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
893     double          **dropval;
894     real              tshift = 0, dt = -1, prec;
895     gmx_bool          bFit, bPFit, bReset;
896     int               nfitdim;
897     gmx_rmpbc_t       gpbc = nullptr;
898     gmx_bool          bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
899     gmx_bool          bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
900     gmx_bool          bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetPrec, bNeedPrec;
901     gmx_bool          bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
902     gmx_bool          bSubTraj = FALSE, bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
903     gmx_bool          bWriteFrame, bSplitHere;
904     const char       *top_file, *in_file, *out_file = nullptr;
905     char              out_file2[256], *charpt;
906     char             *outf_base = nullptr;
907     const char       *outf_ext  = nullptr;
908     char              top_title[256], title[256], timestr[32], stepstr[32], filemode[5];
909     gmx_output_env_t *oenv;
910
911     t_filenm          fnm[] = {
912         { efTRX, "-f",   nullptr,      ffREAD  },
913         { efTRO, "-o",   nullptr,      ffWRITE },
914         { efTPS, nullptr,   nullptr,      ffOPTRD },
915         { efNDX, nullptr,   nullptr,      ffOPTRD },
916         { efNDX, "-fr",  "frames",  ffOPTRD },
917         { efNDX, "-sub", "cluster", ffOPTRD },
918         { efXVG, "-drop", "drop",    ffOPTRD }
919     };
920 #define NFILE asize(fnm)
921
922     if (!parse_common_args(&argc, argv,
923                            PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW |
924                            PCA_TIME_UNIT,
925                            NFILE, fnm, NPA, pa, asize(desc), desc,
926                            0, nullptr, &oenv))
927     {
928         return 0;
929     }
930
931     top_file = ftp2fn(efTPS, NFILE, fnm);
932
933     /* Check command line */
934     in_file = opt2fn("-f", NFILE, fnm);
935
936     if (ttrunc != -1)
937     {
938         do_trunc(in_file, ttrunc);
939     }
940     else
941     {
942         /* mark active cmdline options */
943         bSetBox    = opt2parg_bSet("-box", NPA, pa);
944         bSetTime   = opt2parg_bSet("-t0", NPA, pa);
945         bSetPrec   = opt2parg_bSet("-ndec", NPA, pa);
946         bSetUR     = opt2parg_bSet("-ur", NPA, pa);
947         bExec      = opt2parg_bSet("-exec", NPA, pa);
948         bTimeStep  = opt2parg_bSet("-timestep", NPA, pa);
949         bTDump     = opt2parg_bSet("-dump", NPA, pa);
950         bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
951         bDropOver  = opt2parg_bSet("-dropover", NPA, pa);
952         bTrans     = opt2parg_bSet("-trans", NPA, pa);
953         bSplit     = (split_t != 0);
954
955         /* parse enum options */
956         fit_enum      = nenum(fit);
957         bFit          = (fit_enum == efFit || fit_enum == efFitXY);
958         bReset        = (fit_enum == efReset || fit_enum == efResetXY);
959         bPFit         = fit_enum == efPFit;
960         pbc_enum      = nenum(pbc_opt);
961         bPBCWhole     = pbc_enum == epWhole;
962         bPBCcomRes    = pbc_enum == epComRes;
963         bPBCcomMol    = pbc_enum == epComMol;
964         bPBCcomAtom   = pbc_enum == epComAtom;
965         bNoJump       = pbc_enum == epNojump;
966         bCluster      = pbc_enum == epCluster;
967         bPBC          = pbc_enum != epNone;
968         unitcell_enum = nenum(unitcell_opt);
969         ecenter       = nenum(center_opt) - ecTric;
970
971         /* set and check option dependencies */
972         if (bPFit)
973         {
974             bFit = TRUE;        /* for pfit, fit *must* be set */
975         }
976         if (bFit)
977         {
978             bReset = TRUE;       /* for fit, reset *must* be set */
979         }
980         nfitdim = 0;
981         if (bFit || bReset)
982         {
983             nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
984         }
985         bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
986
987         if (bSetUR)
988         {
989             if (!(bPBCcomRes || bPBCcomMol ||  bPBCcomAtom))
990             {
991                 fprintf(stderr,
992                         "WARNING: Option for unitcell representation (-ur %s)\n"
993                         "         only has effect in combination with -pbc %s, %s or %s.\n"
994                         "         Ingoring unitcell representation.\n\n",
995                         unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
996             }
997         }
998         if (bFit && bPBC)
999         {
1000             gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
1001                       "Please do the PBC condition treatment first and then run trjconv in a second step\n"
1002                       "for the rotational fit.\n"
1003                       "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1004                       "results!");
1005         }
1006
1007         /* ndec is in nr of decimal places, prec is a multiplication factor: */
1008         prec = 1;
1009         for (i = 0; i < ndec; i++)
1010         {
1011             prec *= 10;
1012         }
1013
1014         bIndex = ftp2bSet(efNDX, NFILE, fnm);
1015
1016
1017         /* Determine output type */
1018         out_file = opt2fn("-o", NFILE, fnm);
1019         int ftp  = fn2ftp(out_file);
1020         fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
1021         bNeedPrec = (ftp == efXTC || ftp == efGRO);
1022         int ftpin = fn2ftp(in_file);
1023         if (bVels)
1024         {
1025             /* check if velocities are possible in input and output files */
1026             bVels = (ftp == efTRR || ftp == efGRO ||
1027                      ftp == efG96 || ftp == efTNG)
1028                 && (ftpin == efTRR || ftpin == efGRO ||
1029                     ftpin == efG96 || ftpin == efTNG || ftpin == efCPT);
1030         }
1031         if (bSeparate || bSplit)
1032         {
1033             outf_ext = std::strrchr(out_file, '.');
1034             if (outf_ext == nullptr)
1035             {
1036                 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1037             }
1038             outf_base = gmx_strdup(out_file);
1039             outf_base[outf_ext - out_file] = '\0';
1040         }
1041
1042         bSubTraj = opt2bSet("-sub", NFILE, fnm);
1043         if (bSubTraj)
1044         {
1045             if ((ftp != efXTC) && (ftp != efTRR))
1046             {
1047                 /* It seems likely that other trajectory file types
1048                  * could work here. */
1049                 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1050                           "xtc and trr");
1051             }
1052             clust = cluster_index(nullptr, opt2fn("-sub", NFILE, fnm));
1053
1054             /* Check for number of files disabled, as FOPEN_MAX is not the correct
1055              * number to check for. In my linux box it is only 16.
1056              */
1057             if (/* DISABLES CODE */ (false))
1058             //if (clust->clust->nr > FOPEN_MAX-4)
1059             {
1060                 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1061                           " trajectories.\ntry splitting the index file in %d parts.\n"
1062                           "FOPEN_MAX = %d",
1063                           clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1064             }
1065             gmx_warning("The -sub option could require as many open output files as there are\n"
1066                         "index groups in the file (%d). If you get I/O errors opening new files,\n"
1067                         "try reducing the number of index groups in the file, and perhaps\n"
1068                         "using trjconv -sub several times on different chunks of your index file.\n",
1069                         clust->clust->nr);
1070
1071             snew(clust_status, clust->clust->nr);
1072             snew(clust_status_id, clust->clust->nr);
1073             snew(nfwritten, clust->clust->nr);
1074             for (i = 0; (i < clust->clust->nr); i++)
1075             {
1076                 clust_status[i]    = nullptr;
1077                 clust_status_id[i] = -1;
1078             }
1079             bSeparate = bSplit = FALSE;
1080         }
1081         /* skipping */
1082         if (skip_nr <= 0)
1083         {
1084             gmx_fatal(FARGS, "Argument for -skip (%d) needs to be greater or equal to 1.", skip_nr);
1085         }
1086
1087         std::unique_ptr<gmx_mtop_t> mtop = read_mtop_for_tng(top_file, in_file, out_file);
1088
1089         /* Determine whether to read a topology */
1090         bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1091                 bRmPBC || bReset || bPBCcomMol || bCluster ||
1092                 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1093
1094         /* Determine if when can read index groups */
1095         bIndex = (bIndex || bTPS);
1096
1097         if (bTPS)
1098         {
1099             read_tps_conf(top_file, &top, &ePBC, &xp, nullptr, top_box,
1100                           bReset || bPBCcomRes);
1101             std::strncpy(top_title, *top.name, 255);
1102             top_title[255] = '\0';
1103             atoms          = &top.atoms;
1104
1105             if (0 == top.mols.nr && (bCluster || bPBCcomMol))
1106             {
1107                 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1108             }
1109
1110             /* top_title is only used for gro and pdb,
1111              * the header in such a file is top_title, followed by
1112              * t= ... and/or step= ...
1113              * to prevent double t= or step=, remove it from top_title.
1114              * From GROMACS-2018 we only write t/step when the frame actually
1115              * has a valid time/step, so we need to check for both separately.
1116              */
1117             if ((charpt = std::strstr(top_title, " t= ")))
1118             {
1119                 charpt[0] = '\0';
1120             }
1121             if ((charpt = std::strstr(top_title, " step= ")))
1122             {
1123                 charpt[0] = '\0';
1124             }
1125
1126             if (bCONECT)
1127             {
1128                 gc = gmx_conect_generate(&top);
1129             }
1130             if (bRmPBC)
1131             {
1132                 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
1133             }
1134         }
1135
1136         /* get frame number index */
1137         frindex = nullptr;
1138         if (opt2bSet("-fr", NFILE, fnm))
1139         {
1140             printf("Select groups of frame number indices:\n");
1141             rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, &frindex, &frname);
1142             if (debug)
1143             {
1144                 for (i = 0; i < nrfri; i++)
1145                 {
1146                     fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1147                 }
1148             }
1149         }
1150
1151         /* get index groups etc. */
1152         if (bReset)
1153         {
1154             printf("Select group for %s fit\n",
1155                    bFit ? "least squares" : "translational");
1156             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1157                       1, &ifit, &ind_fit, &gn_fit);
1158
1159             if (bFit)
1160             {
1161                 if (ifit < 2)
1162                 {
1163                     gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1164                 }
1165                 else if (ifit == 3)
1166                 {
1167                     fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1168                 }
1169             }
1170         }
1171         else if (bCluster)
1172         {
1173             printf("Select group for clustering\n");
1174             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1175                       1, &ifit, &ind_fit, &gn_fit);
1176         }
1177
1178         if (bIndex)
1179         {
1180             if (bCenter)
1181             {
1182                 printf("Select group for centering\n");
1183                 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1184                           1, &ncent, &cindex, &grpnm);
1185             }
1186             printf("Select group for output\n");
1187             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1188                       1, &nout, &index, &grpnm);
1189         }
1190         else
1191         {
1192             /* no index file, so read natoms from TRX */
1193             if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
1194             {
1195                 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1196             }
1197             natoms = fr.natoms;
1198             close_trx(trxin);
1199             sfree(fr.x);
1200             snew(index, natoms);
1201             for (i = 0; i < natoms; i++)
1202             {
1203                 index[i] = i;
1204             }
1205             nout = natoms;
1206             if (bCenter)
1207             {
1208                 ncent  = nout;
1209                 cindex = index;
1210             }
1211         }
1212
1213         if (bReset)
1214         {
1215             snew(w_rls, atoms->nr);
1216             for (i = 0; (i < ifit); i++)
1217             {
1218                 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1219             }
1220
1221             /* Restore reference structure and set to origin,
1222                store original location (to put structure back) */
1223             if (bRmPBC)
1224             {
1225                 gmx_rmpbc(gpbc, top.atoms.nr, top_box, xp);
1226             }
1227             copy_rvec(xp[index[0]], x_shift);
1228             reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, nullptr, xp, w_rls);
1229             rvec_dec(x_shift, xp[index[0]]);
1230         }
1231         else
1232         {
1233             clear_rvec(x_shift);
1234         }
1235
1236         if (bDropUnder || bDropOver)
1237         {
1238             /* Read the .xvg file with the drop values */
1239             fprintf(stderr, "\nReading drop file ...");
1240             ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1241             fprintf(stderr, " %d time points\n", ndrop);
1242             if (ndrop == 0 || ncol < 2)
1243             {
1244                 gmx_fatal(FARGS, "Found no data points in %s",
1245                           opt2fn("-drop", NFILE, fnm));
1246             }
1247             drop0 = 0;
1248             drop1 = 0;
1249         }
1250
1251         /* Make atoms struct for output in GRO or PDB files */
1252         if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1253         {
1254             /* get memory for stuff to go in .pdb file, and initialize
1255              * the pdbinfo structure part if the input has it.
1256              */
1257             init_t_atoms(&useatoms, atoms->nr, atoms->havePdbInfo);
1258             sfree(useatoms.resinfo);
1259             useatoms.resinfo = atoms->resinfo;
1260             for (i = 0; (i < nout); i++)
1261             {
1262                 useatoms.atomname[i] = atoms->atomname[index[i]];
1263                 useatoms.atom[i]     = atoms->atom[index[i]];
1264                 if (atoms->havePdbInfo)
1265                 {
1266                     useatoms.pdbinfo[i]  = atoms->pdbinfo[index[i]];
1267                 }
1268                 useatoms.nres        = std::max(useatoms.nres, useatoms.atom[i].resind+1);
1269             }
1270             useatoms.nr = nout;
1271         }
1272         /* select what to read */
1273         if (ftp == efTRR)
1274         {
1275             flags = TRX_READ_X;
1276         }
1277         else
1278         {
1279             flags = TRX_NEED_X;
1280         }
1281         if (bVels)
1282         {
1283             flags = flags | TRX_READ_V;
1284         }
1285         if (bForce)
1286         {
1287             flags = flags | TRX_READ_F;
1288         }
1289
1290         /* open trx file for reading */
1291         bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
1292         if (fr.bPrec)
1293         {
1294             fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1295         }
1296         if (bNeedPrec)
1297         {
1298             if (bSetPrec || !fr.bPrec)
1299             {
1300                 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1301             }
1302             else
1303             {
1304                 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1305             }
1306         }
1307
1308         if (bHaveFirstFrame)
1309         {
1310             if (bTDump)
1311             {
1312                 // Determine timestep (assuming constant spacing for now) if we
1313                 // need to dump frames based on time. This is required so we do not
1314                 // skip the first frame if that was the one that should have been dumped
1315                 double firstFrameTime = fr.time;
1316                 if (read_next_frame(oenv, trxin, &fr))
1317                 {
1318                     dt     = fr.time - firstFrameTime;
1319                     bDTset = TRUE;
1320                     if (dt <= 0)
1321                     {
1322                         fprintf(stderr, "Warning: Frame times are not incrementing - will dump first frame.\n");
1323                     }
1324                 }
1325                 // Now close and reopen so we are at first frame again
1326                 close_trx(trxin);
1327                 done_frame(&fr);
1328                 // Reopen at first frame (We already know it exists if we got here)
1329                 read_first_frame(oenv, &trxin, in_file, &fr, flags);
1330             }
1331
1332             set_trxframe_ePBC(&fr, ePBC);
1333             natoms = fr.natoms;
1334
1335             if (bSetTime)
1336             {
1337                 tshift = tzero-fr.time;
1338             }
1339             else
1340             {
1341                 tzero = fr.time;
1342             }
1343
1344             bCopy = FALSE;
1345             if (bIndex)
1346             {
1347                 /* check if index is meaningful */
1348                 for (i = 0; i < nout; i++)
1349                 {
1350                     if (index[i] >= natoms)
1351                     {
1352                         gmx_fatal(FARGS,
1353                                   "Index[%d] %d is larger than the number of atoms in the\n"
1354                                   "trajectory file (%d). There is a mismatch in the contents\n"
1355                                   "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1356                     }
1357                     bCopy = bCopy || (i != index[i]);
1358                 }
1359             }
1360
1361             /* open output for writing */
1362             std::strcpy(filemode, "w");
1363             switch (ftp)
1364             {
1365                 case efTNG:
1366                     trxout = trjtools_gmx_prepare_tng_writing(out_file,
1367                                                               filemode[0],
1368                                                               trxin,
1369                                                               nullptr,
1370                                                               nout,
1371                                                               mtop.get(),
1372                                                               index,
1373                                                               grpnm);
1374                     break;
1375                 case efXTC:
1376                 case efTRR:
1377                     out = nullptr;
1378                     if (!bSplit && !bSubTraj)
1379                     {
1380                         trxout = open_trx(out_file, filemode);
1381                     }
1382                     break;
1383                 case efGRO:
1384                 case efG96:
1385                 case efPDB:
1386                     if (( !bSeparate && !bSplit ) && !bSubTraj)
1387                     {
1388                         out = gmx_ffopen(out_file, filemode);
1389                     }
1390                     break;
1391                 default:
1392                     gmx_incons("Illegal output file format");
1393             }
1394
1395             if (bCopy)
1396             {
1397                 snew(xmem, nout);
1398                 if (bVels)
1399                 {
1400                     snew(vmem, nout);
1401                 }
1402                 if (bForce)
1403                 {
1404                     snew(fmem, nout);
1405                 }
1406             }
1407
1408             /* Start the big loop over frames */
1409             file_nr  =  0;
1410             frame    =  0;
1411             outframe =  0;
1412             model_nr =  0;
1413
1414             /* Main loop over frames */
1415             do
1416             {
1417                 if (!fr.bStep)
1418                 {
1419                     /* set the step */
1420                     fr.step = newstep;
1421                     newstep++;
1422                 }
1423                 if (bSubTraj)
1424                 {
1425                     /*if (frame >= clust->clust->nra)
1426                        gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1427                     if (frame > clust->maxframe)
1428                     {
1429                         my_clust = -1;
1430                     }
1431                     else
1432                     {
1433                         my_clust = clust->inv_clust[frame];
1434                     }
1435                     if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1436                         (my_clust == -1))
1437                     {
1438                         my_clust = -1;
1439                     }
1440                 }
1441
1442                 if (bSetBox)
1443                 {
1444                     /* generate new box */
1445                     if (fr.bBox == FALSE)
1446                     {
1447                         clear_mat(fr.box);
1448                     }
1449                     for (m = 0; m < DIM; m++)
1450                     {
1451                         if (newbox[m] >= 0)
1452                         {
1453                             fr.box[m][m] = newbox[m];
1454                         }
1455                         else
1456                         {
1457                             if (fr.bBox == FALSE)
1458                             {
1459                                 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1460                             }
1461                         }
1462                     }
1463                 }
1464
1465                 if (bTrans)
1466                 {
1467                     for (i = 0; i < natoms; i++)
1468                     {
1469                         rvec_inc(fr.x[i], trans);
1470                     }
1471                 }
1472
1473                 if (bTDump)
1474                 {
1475                     // If we could not read two frames or times are not incrementing
1476                     // we have almost no idea what to do,
1477                     // but dump the first frame so output is not broken.
1478                     if (dt <= 0 || !bDTset)
1479                     {
1480                         bDumpFrame = true;
1481                     }
1482                     else
1483                     {
1484                         // Dump the frame if we are less than half a frame time
1485                         // below it. This will also ensure we at least dump a
1486                         // somewhat reasonable frame if the spacing is unequal
1487                         // and we have overrun the frame time. Once we dump one
1488                         // frame based on time we quit, so it does not matter
1489                         // that this might be true for all subsequent frames too.
1490                         bDumpFrame = (fr.time > tdump-0.5*dt);
1491                     }
1492                 }
1493                 else
1494                 {
1495                     bDumpFrame = FALSE;
1496                 }
1497
1498                 /* determine if an atom jumped across the box and reset it if so */
1499                 if (bNoJump && (bTPS || frame != 0))
1500                 {
1501                     for (d = 0; d < DIM; d++)
1502                     {
1503                         hbox[d] = 0.5*fr.box[d][d];
1504                     }
1505                     for (i = 0; i < natoms; i++)
1506                     {
1507                         if (bReset)
1508                         {
1509                             rvec_dec(fr.x[i], x_shift);
1510                         }
1511                         for (m = DIM-1; m >= 0; m--)
1512                         {
1513                             if (hbox[m] > 0)
1514                             {
1515                                 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1516                                 {
1517                                     for (d = 0; d <= m; d++)
1518                                     {
1519                                         fr.x[i][d] += fr.box[m][d];
1520                                     }
1521                                 }
1522                                 while (fr.x[i][m]-xp[i][m] > hbox[m])
1523                                 {
1524                                     for (d = 0; d <= m; d++)
1525                                     {
1526                                         fr.x[i][d] -= fr.box[m][d];
1527                                     }
1528                                 }
1529                             }
1530                         }
1531                     }
1532                 }
1533                 else if (bCluster)
1534                 {
1535                     calc_pbc_cluster(ecenter, ifit, &top, ePBC, fr.x, ind_fit, fr.box);
1536                 }
1537
1538                 if (bPFit)
1539                 {
1540                     /* Now modify the coords according to the flags,
1541                        for normal fit, this is only done for output frames */
1542                     if (bRmPBC)
1543                     {
1544                         gmx_rmpbc_trxfr(gpbc, &fr);
1545                     }
1546
1547                     reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1548                     do_fit(natoms, w_rls, xp, fr.x);
1549                 }
1550
1551                 /* store this set of coordinates for future use */
1552                 if (bPFit || bNoJump)
1553                 {
1554                     if (xp == nullptr)
1555                     {
1556                         snew(xp, natoms);
1557                     }
1558                     for (i = 0; (i < natoms); i++)
1559                     {
1560                         copy_rvec(fr.x[i], xp[i]);
1561                         rvec_inc(fr.x[i], x_shift);
1562                     }
1563                 }
1564
1565                 if (frindex)
1566                 {
1567                     /* see if we have a frame from the frame index group */
1568                     for (i = 0; i < nrfri && !bDumpFrame; i++)
1569                     {
1570                         bDumpFrame = frame == frindex[i];
1571                     }
1572                 }
1573                 if (debug && bDumpFrame)
1574                 {
1575                     fprintf(debug, "dumping %d\n", frame);
1576                 }
1577
1578                 bWriteFrame =
1579                     ( ( !bTDump && !frindex && frame % skip_nr == 0 ) || bDumpFrame );
1580
1581                 if (bWriteFrame && (bDropUnder || bDropOver))
1582                 {
1583                     while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1584                     {
1585                         drop0 = drop1;
1586                         drop1++;
1587                     }
1588                     if (std::abs(dropval[0][drop0] - fr.time)
1589                         < std::abs(dropval[0][drop1] - fr.time))
1590                     {
1591                         dropuse = drop0;
1592                     }
1593                     else
1594                     {
1595                         dropuse = drop1;
1596                     }
1597                     if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1598                         (bDropOver && dropval[1][dropuse] > dropover))
1599                     {
1600                         bWriteFrame = FALSE;
1601                     }
1602                 }
1603
1604                 if (bWriteFrame)
1605                 {
1606                     /* We should avoid modifying the input frame,
1607                      * but since here we don't have the output frame yet,
1608                      * we introduce a temporary output frame time variable.
1609                      */
1610                     real frout_time;
1611
1612                     frout_time = fr.time;
1613
1614                     /* calc new time */
1615                     if (bTimeStep)
1616                     {
1617                         frout_time = tzero + frame*timestep;
1618                     }
1619                     else if (bSetTime)
1620                     {
1621                         frout_time += tshift;
1622                     }
1623
1624                     if (bTDump)
1625                     {
1626                         fprintf(stderr, "\nDumping frame at t= %g %s\n",
1627                                 output_env_conv_time(oenv, frout_time), output_env_get_time_unit(oenv).c_str());
1628                     }
1629
1630                     /* check for writing at each delta_t */
1631                     bDoIt = (delta_t == 0);
1632                     if (!bDoIt)
1633                     {
1634                         if (!bRound)
1635                         {
1636                             bDoIt = bRmod(frout_time, tzero, delta_t);
1637                         }
1638                         else
1639                         {
1640                             /* round() is not C89 compatible, so we do this:  */
1641                             bDoIt = bRmod(std::floor(frout_time+0.5), std::floor(tzero+0.5),
1642                                           std::floor(delta_t+0.5));
1643                         }
1644                     }
1645
1646                     if (bDoIt || bTDump)
1647                     {
1648                         /* print sometimes */
1649                         if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1650                         {
1651                             fprintf(stderr, " ->  frame %6d time %8.3f      \r",
1652                                     outframe, output_env_conv_time(oenv, frout_time));
1653                             fflush(stderr);
1654                         }
1655
1656                         if (!bPFit)
1657                         {
1658                             /* Now modify the coords according to the flags,
1659                                for PFit we did this already! */
1660
1661                             if (bRmPBC)
1662                             {
1663                                 gmx_rmpbc_trxfr(gpbc, &fr);
1664                             }
1665
1666                             if (bReset)
1667                             {
1668                                 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1669                                 if (bFit)
1670                                 {
1671                                     do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1672                                 }
1673                                 if (!bCenter)
1674                                 {
1675                                     for (i = 0; i < natoms; i++)
1676                                     {
1677                                         rvec_inc(fr.x[i], x_shift);
1678                                     }
1679                                 }
1680                             }
1681
1682                             if (bCenter)
1683                             {
1684                                 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1685                             }
1686                         }
1687
1688                         auto positionsArrayRef = gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec *>(fr.x), natoms);
1689                         if (bPBCcomAtom)
1690                         {
1691                             switch (unitcell_enum)
1692                             {
1693                                 case euRect:
1694                                     put_atoms_in_box(ePBC, fr.box, positionsArrayRef);
1695                                     break;
1696                                 case euTric:
1697                                     put_atoms_in_triclinic_unitcell(ecenter, fr.box, positionsArrayRef);
1698                                     break;
1699                                 case euCompact:
1700                                     put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1701                                                                   positionsArrayRef);
1702                                     break;
1703                             }
1704                         }
1705                         if (bPBCcomRes)
1706                         {
1707                             put_residue_com_in_box(unitcell_enum, ecenter,
1708                                                    natoms, atoms->atom, ePBC, fr.box, fr.x);
1709                         }
1710                         if (bPBCcomMol)
1711                         {
1712                             put_molecule_com_in_box(unitcell_enum, ecenter,
1713                                                     &top.mols,
1714                                                     natoms, atoms->atom, ePBC, fr.box, fr.x);
1715                         }
1716                         /* Copy the input trxframe struct to the output trxframe struct */
1717                         frout        = fr;
1718                         frout.time   = frout_time;
1719                         frout.bV     = (frout.bV && bVels);
1720                         frout.bF     = (frout.bF && bForce);
1721                         frout.natoms = nout;
1722                         if (bNeedPrec && (bSetPrec || !fr.bPrec))
1723                         {
1724                             frout.bPrec = TRUE;
1725                             frout.prec  = prec;
1726                         }
1727                         if (bCopy)
1728                         {
1729                             frout.x = xmem;
1730                             if (frout.bV)
1731                             {
1732                                 frout.v = vmem;
1733                             }
1734                             if (frout.bF)
1735                             {
1736                                 frout.f = fmem;
1737                             }
1738                             for (i = 0; i < nout; i++)
1739                             {
1740                                 copy_rvec(fr.x[index[i]], frout.x[i]);
1741                                 if (frout.bV)
1742                                 {
1743                                     copy_rvec(fr.v[index[i]], frout.v[i]);
1744                                 }
1745                                 if (frout.bF)
1746                                 {
1747                                     copy_rvec(fr.f[index[i]], frout.f[i]);
1748                                 }
1749                             }
1750                         }
1751
1752                         if (opt2parg_bSet("-shift", NPA, pa))
1753                         {
1754                             for (i = 0; i < nout; i++)
1755                             {
1756                                 for (d = 0; d < DIM; d++)
1757                                 {
1758                                     frout.x[i][d] += outframe*shift[d];
1759                                 }
1760                             }
1761                         }
1762
1763                         if (!bRound)
1764                         {
1765                             bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1766                         }
1767                         else
1768                         {
1769                             /* round() is not C89 compatible, so we do this: */
1770                             bSplitHere = bSplit && bRmod(std::floor(frout.time+0.5),
1771                                                          std::floor(tzero+0.5),
1772                                                          std::floor(split_t+0.5));
1773                         }
1774                         if (bSeparate || bSplitHere)
1775                         {
1776                             mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1777                         }
1778
1779                         switch (ftp)
1780                         {
1781                             case efTNG:
1782                                 write_tng_frame(trxout, &frout);
1783                                 // TODO when trjconv behaves better: work how to read and write lambda
1784                                 break;
1785                             case efTRR:
1786                             case efXTC:
1787                                 if (bSplitHere)
1788                                 {
1789                                     if (trxout)
1790                                     {
1791                                         close_trx(trxout);
1792                                     }
1793                                     trxout = open_trx(out_file2, filemode);
1794                                 }
1795                                 if (bSubTraj)
1796                                 {
1797                                     if (my_clust != -1)
1798                                     {
1799                                         char buf[STRLEN];
1800                                         if (clust_status_id[my_clust] == -1)
1801                                         {
1802                                             sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1803                                             clust_status[my_clust]    = open_trx(buf, "w");
1804                                             clust_status_id[my_clust] = 1;
1805                                             ntrxopen++;
1806                                         }
1807                                         else if (clust_status_id[my_clust] == -2)
1808                                         {
1809                                             gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1810                                                       clust->grpname[my_clust], ntrxopen, frame,
1811                                                       my_clust);
1812                                         }
1813                                         write_trxframe(clust_status[my_clust], &frout, gc);
1814                                         nfwritten[my_clust]++;
1815                                         if (nfwritten[my_clust] ==
1816                                             (clust->clust->index[my_clust+1]-
1817                                              clust->clust->index[my_clust]))
1818                                         {
1819                                             close_trx(clust_status[my_clust]);
1820                                             clust_status[my_clust]    = nullptr;
1821                                             clust_status_id[my_clust] = -2;
1822                                             ntrxopen--;
1823                                             if (ntrxopen < 0)
1824                                             {
1825                                                 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1826                                             }
1827                                         }
1828                                     }
1829                                 }
1830                                 else
1831                                 {
1832                                     write_trxframe(trxout, &frout, gc);
1833                                 }
1834                                 break;
1835                             case efGRO:
1836                             case efG96:
1837                             case efPDB:
1838                                 // Only add a generator statement if title is empty,
1839                                 // to avoid multiple generated-by statements from various programs
1840                                 if (std::strlen(top_title) == 0)
1841                                 {
1842                                     sprintf(top_title, "Generated by trjconv");
1843                                 }
1844                                 if (frout.bTime)
1845                                 {
1846                                     sprintf(timestr, " t= %9.5f", frout.time);
1847                                 }
1848                                 else
1849                                 {
1850                                     std::strcpy(timestr, "");
1851                                 }
1852                                 if (frout.bStep)
1853                                 {
1854                                     sprintf(stepstr, " step= %" PRId64, frout.step);
1855                                 }
1856                                 else
1857                                 {
1858                                     std::strcpy(stepstr, "");
1859                                 }
1860                                 snprintf(title, 256, "%s%s%s", top_title, timestr, stepstr);
1861                                 if (bSeparate || bSplitHere)
1862                                 {
1863                                     out = gmx_ffopen(out_file2, "w");
1864                                 }
1865                                 switch (ftp)
1866                                 {
1867                                     case efGRO:
1868                                         write_hconf_p(out, title, &useatoms,
1869                                                       frout.x, frout.bV ? frout.v : nullptr, frout.box);
1870                                         break;
1871                                     case efPDB:
1872                                         fprintf(out, "REMARK    GENERATED BY TRJCONV\n");
1873                                         /* if reading from pdb, we want to keep the original
1874                                            model numbering else we write the output frame
1875                                            number plus one, because model 0 is not allowed in pdb */
1876                                         if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1877                                         {
1878                                             model_nr = fr.step;
1879                                         }
1880                                         else
1881                                         {
1882                                             model_nr++;
1883                                         }
1884                                         write_pdbfile(out, title, &useatoms, frout.x,
1885                                                       frout.ePBC, frout.box, ' ', model_nr, gc, TRUE);
1886                                         break;
1887                                     case efG96:
1888                                         const char *outputTitle = "";
1889                                         if (bSeparate || bTDump)
1890                                         {
1891                                             outputTitle = title;
1892                                             if (bTPS)
1893                                             {
1894                                                 frout.bAtoms = TRUE;
1895                                             }
1896                                             frout.atoms  = &useatoms;
1897                                             frout.bStep  = FALSE;
1898                                             frout.bTime  = FALSE;
1899                                         }
1900                                         else
1901                                         {
1902                                             if (outframe == 0)
1903                                             {
1904                                                 outputTitle = title;
1905                                             }
1906                                             frout.bAtoms = FALSE;
1907                                             frout.bStep  = TRUE;
1908                                             frout.bTime  = TRUE;
1909                                         }
1910                                         write_g96_conf(out, outputTitle, &frout, -1, nullptr);
1911                                 }
1912                                 if (bSeparate || bSplitHere)
1913                                 {
1914                                     gmx_ffclose(out);
1915                                     out = nullptr;
1916                                 }
1917                                 break;
1918                             default:
1919                                 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1920                         }
1921                         if (bSeparate || bSplitHere)
1922                         {
1923                             file_nr++;
1924                         }
1925
1926                         /* execute command */
1927                         if (bExec)
1928                         {
1929                             char c[255];
1930                             sprintf(c, "%s  %d", exec_command, file_nr-1);
1931                             /*fprintf(stderr,"Executing '%s'\n",c);*/
1932                             if (0 != system(c))
1933                             {
1934                                 gmx_fatal(FARGS, "Error executing command: %s", c);
1935                             }
1936                         }
1937                         outframe++;
1938                     }
1939                 }
1940                 frame++;
1941                 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1942             }
1943             while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1944         }
1945
1946         if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1947         {
1948             fprintf(stderr, "\nWARNING no output, "
1949                     "last frame read at t=%g\n", fr.time);
1950         }
1951         fprintf(stderr, "\n");
1952
1953         close_trx(trxin);
1954         sfree(outf_base);
1955
1956         if (bRmPBC)
1957         {
1958             gmx_rmpbc_done(gpbc);
1959         }
1960
1961         if (trxout)
1962         {
1963             close_trx(trxout);
1964         }
1965         else if (out != nullptr)
1966         {
1967             gmx_ffclose(out);
1968         }
1969         if (bSubTraj)
1970         {
1971             for (i = 0; (i < clust->clust->nr); i++)
1972             {
1973                 if (clust_status_id[i] >= 0)
1974                 {
1975                     close_trx(clust_status[i]);
1976                 }
1977             }
1978         }
1979     }
1980
1981     if (bTPS)
1982     {
1983         done_top(&top);
1984     }
1985     sfree(xp);
1986     sfree(xmem);
1987     sfree(vmem);
1988     sfree(fmem);
1989     sfree(grpnm);
1990     sfree(index);
1991     sfree(cindex);
1992     done_frame(&fr);
1993
1994     do_view(oenv, out_file, nullptr);
1995
1996     output_env_done(oenv);
1997     return 0;
1998 }