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