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