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