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