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