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