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