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