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