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