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