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