73100255d9ed5e4d3e30e5f866cb792f411e9ebf
[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 #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:[BR]",
588         "* from one format to another[BR]",
589         "* select a subset of atoms[BR]",
590         "* change the periodicity representation[BR]",
591         "* keep multimeric molecules together[BR]",
592         "* center atoms in the box[BR]",
593         "* fit atoms to reference structure[BR]",
594         "* reduce the number of frames[BR]",
595         "* change the timestamps of the frames ",
596         "([TT]-t0[tt] and [TT]-timestep[tt])[BR]",
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.[BR]",
603         "* select frames within a certain range of a quantity given",
604         "in an [TT].xvg[tt] file.[PAR]",
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         "[TT].xtc[tt], [TT].trr[tt], [TT].gro[tt], [TT].g96[tt]",
611         "and [TT].pdb[tt].",
612         "The file formats are detected from the file extension.",
613         "The precision of [TT].xtc[tt] and [TT].gro[tt] output is taken from the",
614         "input file for [TT].xtc[tt], [TT].gro[tt] and [TT].pdb[tt],",
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. [TT].trr[tt]",
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         "[TT].trr[tt], [TT].gro[tt] 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 [TT].pdb[tt] file. By default, all frames all written to one file.",
625         "[TT].pdb[tt] 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 [TT].xtc[tt] 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:[BR]",
648         "[TT]* mol[tt] puts the center of mass of molecules in the box,",
649         "and requires a run input file to be supplied with [TT]-s[tt].[BR]",
650         "[TT]* res[tt] puts the center of mass of residues in the box.[BR]",
651         "[TT]* atom[tt] puts all the atoms in the box.[BR]",
652         "[TT]* nojump[tt] checks if atoms jump across the box and then puts",
653         "them back. This has the effect that all molecules",
654         "will remain whole (provided they were whole in the initial",
655         "conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
656         "molecules may diffuse out of the box. The starting configuration",
657         "for this procedure is taken from the structure file, if one is",
658         "supplied, otherwise it is the first frame.[BR]",
659         "[TT]* cluster[tt] clusters all the atoms in the selected index",
660         "such that they are all closest to the center of mass of the cluster,",
661         "which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
662         "results if you in fact have a cluster. Luckily that can be checked",
663         "afterwards using a trajectory viewer. Note also that if your molecules",
664         "are broken this will not work either.[BR]",
665         "The separate option [TT]-clustercenter[tt] can be used to specify an",
666         "approximate center for the cluster. This is useful e.g. if you have",
667         "two big vesicles, and you want to maintain their relative positions.[BR]",
668         "[TT]* whole[tt] only makes broken molecules whole.[PAR]",
669
670         "Option [TT]-ur[tt] sets the unit cell representation for options",
671         "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
672         "All three options give different results for triclinic boxes and",
673         "identical results for rectangular boxes.",
674         "[TT]rect[tt] is the ordinary brick shape.",
675         "[TT]tric[tt] is the triclinic unit cell.",
676         "[TT]compact[tt] puts all atoms at the closest distance from the center",
677         "of the box. This can be useful for visualizing e.g. truncated octahedra",
678         "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
679         "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
680         "is set differently.[PAR]",
681
682         "Option [TT]-center[tt] centers the system in the box. The user can",
683         "select the group which is used to determine the geometrical center.",
684         "Option [TT]-boxcenter[tt] sets the location of the center of the box",
685         "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
686         "[TT]tric[tt]: half of the sum of the box vectors,",
687         "[TT]rect[tt]: half of the box diagonal,",
688         "[TT]zero[tt]: zero.",
689         "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
690         "want all molecules in the box after the centering.[PAR]",
691
692         "Option [TT]-box[tt] sets the size of the new box. If you want to"
693         "modify only some of the dimensions, e.g. when reading from a trajectory,"
694         "you can use -1 for those dimensions that should stay the same"
695
696         "It is not always possible to use combinations of [TT]-pbc[tt],",
697         "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
698         "you want in one call to [THISMODULE]. Consider using multiple",
699         "calls, and check out the GROMACS website for suggestions.[PAR]",
700
701         "With [TT]-dt[tt], it is possible to reduce the number of ",
702         "frames in the output. This option relies on the accuracy of the times",
703         "in your input trajectory, so if these are inaccurate use the",
704         "[TT]-timestep[tt] option to modify the time (this can be done",
705         "simultaneously). For making smooth movies, the program [gmx-filter]",
706         "can reduce the number of frames while using low-pass frequency",
707         "filtering, this reduces aliasing of high frequency motions.[PAR]",
708
709         "Using [TT]-trunc[tt] [THISMODULE] can truncate [TT].trr[tt] in place, i.e.",
710         "without copying the file. This is useful when a run has crashed",
711         "during disk I/O (i.e. full disk), or when two contiguous",
712         "trajectories must be concatenated without having double frames.[PAR]",
713
714         "Option [TT]-dump[tt] can be used to extract a frame at or near",
715         "one specific time from your trajectory.[PAR]",
716
717         "Option [TT]-drop[tt] reads an [TT].xvg[tt] file with times and values.",
718         "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
719         "frames with a value below and above the value of the respective options",
720         "will not be written."
721     };
722
723     int         pbc_enum;
724     enum
725     {
726         epSel,
727         epNone,
728         epComMol,
729         epComRes,
730         epComAtom,
731         epNojump,
732         epCluster,
733         epWhole,
734         epNR
735     };
736     const char *pbc_opt[epNR + 1] =
737     {
738         NULL, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
739         NULL
740     };
741
742     int         unitcell_enum;
743     const char *unitcell_opt[euNR+1] =
744     { NULL, "rect", "tric", "compact", NULL };
745
746     enum
747     {
748         ecSel, ecTric, ecRect, ecZero, ecNR
749     };
750     const char *center_opt[ecNR+1] =
751     { NULL, "tric", "rect", "zero", NULL };
752     int         ecenter;
753
754     int         fit_enum;
755     enum
756     {
757         efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
758     };
759     const char *fit[efNR + 1] =
760     {
761         NULL, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
762         "progressive", NULL
763     };
764
765     static gmx_bool  bSeparate     = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
766     static gmx_bool  bCenter       = FALSE;
767     static int       skip_nr       = 1, ndec = 3, nzero = 0;
768     static real      tzero         = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
769     static rvec      newbox        = {0, 0, 0}, shift = {0, 0, 0}, trans = {0, 0, 0};
770     static char     *exec_command  = NULL;
771     static real      dropunder     = 0, dropover = 0;
772     static gmx_bool  bRound        = FALSE;
773
774     t_pargs
775         pa[] =
776     {
777         { "-skip", FALSE, etINT,
778           { &skip_nr }, "Only write every nr-th frame" },
779         { "-dt", FALSE, etTIME,
780           { &delta_t },
781           "Only write frame when t MOD dt = first time (%t)" },
782         { "-round", FALSE, etBOOL,
783           { &bRound }, "Round measurements to nearest picosecond"},
784         { "-dump", FALSE, etTIME,
785           { &tdump }, "Dump frame nearest specified time (%t)" },
786         { "-t0", FALSE, etTIME,
787           { &tzero },
788           "Starting time (%t) (default: don't change)" },
789         { "-timestep", FALSE, etTIME,
790           { &timestep },
791           "Change time step between input frames (%t)" },
792         { "-pbc", FALSE, etENUM,
793           { pbc_opt },
794           "PBC treatment (see help text for full description)" },
795         { "-ur", FALSE, etENUM,
796           { unitcell_opt }, "Unit-cell representation" },
797         { "-center", FALSE, etBOOL,
798           { &bCenter }, "Center atoms in box" },
799         { "-boxcenter", FALSE, etENUM,
800           { center_opt }, "Center for -pbc and -center" },
801         { "-box", FALSE, etRVEC,
802           { newbox },
803           "Size for new cubic box (default: read from input)" },
804         { "-trans", FALSE, etRVEC,
805           { trans },
806           "All coordinates will be translated by trans. This "
807           "can advantageously be combined with -pbc mol -ur "
808           "compact." },
809         { "-shift", FALSE, etRVEC,
810           { shift },
811           "All coordinates will be shifted by framenr*shift" },
812         { "-fit", FALSE, etENUM,
813           { fit },
814           "Fit molecule to ref structure in the structure file" },
815         { "-ndec", FALSE, etINT,
816           { &ndec },
817           "Precision for .xtc and .gro writing in number of "
818           "decimal places" },
819         { "-vel", FALSE, etBOOL,
820           { &bVels }, "Read and write velocities if possible" },
821         { "-force", FALSE, etBOOL,
822           { &bForce }, "Read and write forces if possible" },
823         { "-trunc", FALSE, etTIME,
824           { &ttrunc },
825           "Truncate input trajectory file after this time (%t)" },
826         { "-exec", FALSE, etSTR,
827           { &exec_command },
828           "Execute command for every output frame with the "
829           "frame number as argument" },
830         { "-split", FALSE, etTIME,
831           { &split_t },
832           "Start writing new file when t MOD split = first "
833           "time (%t)" },
834         { "-sep", FALSE, etBOOL,
835           { &bSeparate },
836           "Write each frame to a separate .gro, .g96 or .pdb "
837           "file" },
838         { "-nzero", FALSE, etINT,
839           { &nzero },
840           "If the -sep flag is set, use these many digits "
841           "for the file numbers and prepend zeros as needed" },
842         { "-dropunder", FALSE, etREAL,
843           { &dropunder }, "Drop all frames below this value" },
844         { "-dropover", FALSE, etREAL,
845           { &dropover }, "Drop all frames above this value" },
846         { "-conect", FALSE, etBOOL,
847           { &bCONECT },
848           "Add conect records when writing [TT].pdb[tt] files. Useful "
849           "for visualization of non-standard molecules, e.g. "
850           "coarse grained ones" }
851     };
852 #define NPA asize(pa)
853
854     FILE            *out    = NULL;
855     t_trxstatus     *trxout = NULL;
856     t_trxstatus     *trxin;
857     int              ftp, ftpin = 0, file_nr;
858     t_trxframe       fr, frout;
859     int              flags;
860     rvec            *xmem = NULL, *vmem = NULL, *fmem = NULL;
861     rvec            *xp   = NULL, x_shift, hbox, box_center, dx;
862     real             xtcpr, lambda, *w_rls = NULL;
863     int              m, i, d, frame, outframe, natoms, nout, ncent, nre, newstep = 0, model_nr;
864 #define SKIP 10
865     t_topology       top;
866     gmx_mtop_t      *mtop  = NULL;
867     gmx_conect       gc    = NULL;
868     int              ePBC  = -1;
869     t_atoms         *atoms = NULL, useatoms;
870     matrix           top_box;
871     atom_id         *index, *cindex;
872     char            *grpnm;
873     int             *frindex, nrfri;
874     char            *frname;
875     int              ifit, irms, my_clust = -1;
876     atom_id         *ind_fit, *ind_rms;
877     char            *gn_fit, *gn_rms;
878     t_cluster_ndx   *clust           = NULL;
879     t_trxstatus    **clust_status    = NULL;
880     int             *clust_status_id = NULL;
881     int              ntrxopen        = 0;
882     int             *nfwritten       = NULL;
883     int              ndrop           = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
884     double         **dropval;
885     real             tshift = 0, t0 = -1, dt = 0.001, prec;
886     gmx_bool         bFit, bFitXY, bPFit, bReset;
887     int              nfitdim;
888     gmx_rmpbc_t      gpbc = NULL;
889     gmx_bool         bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
890     gmx_bool         bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
891     gmx_bool         bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetPrec, bNeedPrec;
892     gmx_bool         bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
893     gmx_bool         bSubTraj = FALSE, bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
894     gmx_bool         bWriteFrame, bSplitHere;
895     const char      *top_file, *in_file, *out_file = NULL;
896     char             out_file2[256], *charpt;
897     char            *outf_base = NULL;
898     const char      *outf_ext  = NULL;
899     char             top_title[256], title[256], command[256], filemode[5];
900     int              xdr          = 0;
901     gmx_bool         bWarnCompact = FALSE;
902     const char      *warn;
903     output_env_t     oenv;
904
905     t_filenm         fnm[] = {
906         { efTRX, "-f",   NULL,      ffREAD  },
907         { efTRO, "-o",   NULL,      ffWRITE },
908         { efTPS, NULL,   NULL,      ffOPTRD },
909         { efNDX, NULL,   NULL,      ffOPTRD },
910         { efNDX, "-fr",  "frames",  ffOPTRD },
911         { efNDX, "-sub", "cluster", ffOPTRD },
912         { efXVG, "-drop", "drop",    ffOPTRD }
913     };
914 #define NFILE asize(fnm)
915
916     if (!parse_common_args(&argc, argv,
917                            PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW |
918                            PCA_TIME_UNIT,
919                            NFILE, fnm, NPA, pa, asize(desc), desc,
920                            0, NULL, &oenv))
921     {
922         return 0;
923     }
924
925     top_file = ftp2fn(efTPS, NFILE, fnm);
926     init_top(&top);
927
928     /* Check command line */
929     in_file = opt2fn("-f", NFILE, fnm);
930
931     if (ttrunc != -1)
932     {
933         do_trunc(in_file, ttrunc);
934     }
935     else
936     {
937         /* mark active cmdline options */
938         bSetBox    = opt2parg_bSet("-box", NPA, pa);
939         bSetTime   = opt2parg_bSet("-t0", NPA, pa);
940         bSetPrec   = opt2parg_bSet("-ndec", NPA, pa);
941         bSetUR     = opt2parg_bSet("-ur", NPA, pa);
942         bExec      = opt2parg_bSet("-exec", NPA, pa);
943         bTimeStep  = opt2parg_bSet("-timestep", NPA, pa);
944         bTDump     = opt2parg_bSet("-dump", NPA, pa);
945         bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
946         bDropOver  = opt2parg_bSet("-dropover", NPA, pa);
947         bTrans     = opt2parg_bSet("-trans", NPA, pa);
948         bSplit     = (split_t != 0);
949
950         /* parse enum options */
951         fit_enum      = nenum(fit);
952         bFit          = (fit_enum == efFit || fit_enum == efFitXY);
953         bFitXY        = fit_enum == efFitXY;
954         bReset        = (fit_enum == efReset || fit_enum == efResetXY);
955         bPFit         = fit_enum == efPFit;
956         pbc_enum      = nenum(pbc_opt);
957         bPBCWhole     = pbc_enum == epWhole;
958         bPBCcomRes    = pbc_enum == epComRes;
959         bPBCcomMol    = pbc_enum == epComMol;
960         bPBCcomAtom   = pbc_enum == epComAtom;
961         bNoJump       = pbc_enum == epNojump;
962         bCluster      = pbc_enum == epCluster;
963         bPBC          = pbc_enum != epNone;
964         unitcell_enum = nenum(unitcell_opt);
965         ecenter       = nenum(center_opt) - ecTric;
966
967         /* set and check option dependencies */
968         if (bPFit)
969         {
970             bFit = TRUE;        /* for pfit, fit *must* be set */
971         }
972         if (bFit)
973         {
974             bReset = TRUE;       /* for fit, reset *must* be set */
975         }
976         nfitdim = 0;
977         if (bFit || bReset)
978         {
979             nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
980         }
981         bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
982
983         if (bSetUR)
984         {
985             if (!(bPBCcomRes || bPBCcomMol ||  bPBCcomAtom))
986             {
987                 fprintf(stderr,
988                         "WARNING: Option for unitcell representation (-ur %s)\n"
989                         "         only has effect in combination with -pbc %s, %s or %s.\n"
990                         "         Ingoring unitcell representation.\n\n",
991                         unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
992                 bSetUR = FALSE;
993             }
994         }
995         if (bFit && bPBC)
996         {
997             gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
998                       "Please do the PBC condition treatment first and then run trjconv in a second step\n"
999                       "for the rotational fit.\n"
1000                       "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1001                       "results!");
1002         }
1003
1004         /* ndec is in nr of decimal places, prec is a multiplication factor: */
1005         prec = 1;
1006         for (i = 0; i < ndec; i++)
1007         {
1008             prec *= 10;
1009         }
1010
1011         bIndex = ftp2bSet(efNDX, NFILE, fnm);
1012
1013
1014         /* Determine output type */
1015         out_file = opt2fn("-o", NFILE, fnm);
1016         ftp      = fn2ftp(out_file);
1017         fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
1018         bNeedPrec = (ftp == efXTC || ftp == efGRO);
1019         if (bVels)
1020         {
1021             /* check if velocities are possible in input and output files */
1022             ftpin = fn2ftp(in_file);
1023             bVels = (ftp == efTRR || ftp == efGRO ||
1024                      ftp == efG96 || ftp == efTNG)
1025                 && (ftpin == efTRR || ftpin == efGRO ||
1026                     ftpin == efG96 || ftpin == efTNG || ftpin == efCPT);
1027         }
1028         if (bSeparate || bSplit)
1029         {
1030             outf_ext = strrchr(out_file, '.');
1031             if (outf_ext == NULL)
1032             {
1033                 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1034             }
1035             outf_base = gmx_strdup(out_file);
1036             outf_base[outf_ext - out_file] = '\0';
1037         }
1038
1039         bSubTraj = opt2bSet("-sub", NFILE, fnm);
1040         if (bSubTraj)
1041         {
1042             if ((ftp != efXTC) && (ftp != efTRR))
1043             {
1044                 /* It seems likely that other trajectory file types
1045                  * could work here. */
1046                 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1047                           "xtc and trr");
1048             }
1049             clust = cluster_index(NULL, opt2fn("-sub", NFILE, fnm));
1050
1051             /* Check for number of files disabled, as FOPEN_MAX is not the correct
1052              * number to check for. In my linux box it is only 16.
1053              */
1054             if (0 && (clust->clust->nr > FOPEN_MAX-4))
1055             {
1056                 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1057                           " trajectories.\ntry splitting the index file in %d parts.\n"
1058                           "FOPEN_MAX = %d",
1059                           clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1060             }
1061             gmx_warning("The -sub option could require as many open output files as there are\n"
1062                         "index groups in the file (%d). If you get I/O errors opening new files,\n"
1063                         "try reducing the number of index groups in the file, and perhaps\n"
1064                         "using trjconv -sub several times on different chunks of your index file.\n",
1065                         clust->clust->nr);
1066
1067             snew(clust_status, clust->clust->nr);
1068             snew(clust_status_id, clust->clust->nr);
1069             snew(nfwritten, clust->clust->nr);
1070             for (i = 0; (i < clust->clust->nr); i++)
1071             {
1072                 clust_status[i]    = NULL;
1073                 clust_status_id[i] = -1;
1074             }
1075             bSeparate = bSplit = FALSE;
1076         }
1077         /* skipping */
1078         if (skip_nr <= 0)
1079         {
1080         }
1081
1082         mtop = read_mtop_for_tng(top_file, in_file, out_file);
1083
1084         /* Determine whether to read a topology */
1085         bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1086                 bRmPBC || bReset || bPBCcomMol || bCluster ||
1087                 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1088
1089         /* Determine if when can read index groups */
1090         bIndex = (bIndex || bTPS);
1091
1092         if (bTPS)
1093         {
1094             read_tps_conf(top_file, top_title, &top, &ePBC, &xp, NULL, top_box,
1095                           bReset || bPBCcomRes);
1096             atoms = &top.atoms;
1097
1098             if (0 == top.mols.nr && (bCluster || bPBCcomMol))
1099             {
1100                 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1101             }
1102
1103             /* top_title is only used for gro and pdb,
1104              * the header in such a file is top_title t= ...
1105              * to prevent a double t=, remove it from top_title
1106              */
1107             if ((charpt = strstr(top_title, " t= ")))
1108             {
1109                 charpt[0] = '\0';
1110             }
1111
1112             if (bCONECT)
1113             {
1114                 gc = gmx_conect_generate(&top);
1115             }
1116             if (bRmPBC)
1117             {
1118                 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
1119             }
1120         }
1121
1122         /* get frame number index */
1123         frindex = NULL;
1124         if (opt2bSet("-fr", NFILE, fnm))
1125         {
1126             printf("Select groups of frame number indices:\n");
1127             rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, (atom_id **)&frindex, &frname);
1128             if (debug)
1129             {
1130                 for (i = 0; i < nrfri; i++)
1131                 {
1132                     fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1133                 }
1134             }
1135         }
1136
1137         /* get index groups etc. */
1138         if (bReset)
1139         {
1140             printf("Select group for %s fit\n",
1141                    bFit ? "least squares" : "translational");
1142             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1143                       1, &ifit, &ind_fit, &gn_fit);
1144
1145             if (bFit)
1146             {
1147                 if (ifit < 2)
1148                 {
1149                     gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1150                 }
1151                 else if (ifit == 3)
1152                 {
1153                     fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1154                 }
1155             }
1156         }
1157         else if (bCluster)
1158         {
1159             printf("Select group for clustering\n");
1160             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1161                       1, &ifit, &ind_fit, &gn_fit);
1162         }
1163
1164         if (bIndex)
1165         {
1166             if (bCenter)
1167             {
1168                 printf("Select group for centering\n");
1169                 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1170                           1, &ncent, &cindex, &grpnm);
1171             }
1172             printf("Select group for output\n");
1173             get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1174                       1, &nout, &index, &grpnm);
1175         }
1176         else
1177         {
1178             /* no index file, so read natoms from TRX */
1179             if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
1180             {
1181                 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1182             }
1183             natoms = fr.natoms;
1184             close_trj(trxin);
1185             sfree(fr.x);
1186             snew(index, natoms);
1187             for (i = 0; i < natoms; i++)
1188             {
1189                 index[i] = i;
1190             }
1191             nout = natoms;
1192             if (bCenter)
1193             {
1194                 ncent  = nout;
1195                 cindex = index;
1196             }
1197         }
1198
1199         if (bReset)
1200         {
1201             snew(w_rls, atoms->nr);
1202             for (i = 0; (i < ifit); i++)
1203             {
1204                 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1205             }
1206
1207             /* Restore reference structure and set to origin,
1208                store original location (to put structure back) */
1209             if (bRmPBC)
1210             {
1211                 gmx_rmpbc(gpbc, top.atoms.nr, top_box, xp);
1212             }
1213             copy_rvec(xp[index[0]], x_shift);
1214             reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, NULL, xp, w_rls);
1215             rvec_dec(x_shift, xp[index[0]]);
1216         }
1217         else
1218         {
1219             clear_rvec(x_shift);
1220         }
1221
1222         if (bDropUnder || bDropOver)
1223         {
1224             /* Read the .xvg file with the drop values */
1225             fprintf(stderr, "\nReading drop file ...");
1226             ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1227             fprintf(stderr, " %d time points\n", ndrop);
1228             if (ndrop == 0 || ncol < 2)
1229             {
1230                 gmx_fatal(FARGS, "Found no data points in %s",
1231                           opt2fn("-drop", NFILE, fnm));
1232             }
1233             drop0 = 0;
1234             drop1 = 0;
1235         }
1236
1237         /* Make atoms struct for output in GRO or PDB files */
1238         if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1239         {
1240             /* get memory for stuff to go in .pdb file, and initialize
1241              * the pdbinfo structure part if the input has it.
1242              */
1243             init_t_atoms(&useatoms, atoms->nr, (atoms->pdbinfo != NULL));
1244             sfree(useatoms.resinfo);
1245             useatoms.resinfo = atoms->resinfo;
1246             for (i = 0; (i < nout); i++)
1247             {
1248                 useatoms.atomname[i] = atoms->atomname[index[i]];
1249                 useatoms.atom[i]     = atoms->atom[index[i]];
1250                 if (atoms->pdbinfo != NULL)
1251                 {
1252                     useatoms.pdbinfo[i]  = atoms->pdbinfo[index[i]];
1253                 }
1254                 useatoms.nres        = max(useatoms.nres, useatoms.atom[i].resind+1);
1255             }
1256             useatoms.nr = nout;
1257         }
1258         /* select what to read */
1259         if (ftp == efTRR)
1260         {
1261             flags = TRX_READ_X;
1262         }
1263         else
1264         {
1265             flags = TRX_NEED_X;
1266         }
1267         if (bVels)
1268         {
1269             flags = flags | TRX_READ_V;
1270         }
1271         if (bForce)
1272         {
1273             flags = flags | TRX_READ_F;
1274         }
1275
1276         /* open trx file for reading */
1277         bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
1278         if (fr.bPrec)
1279         {
1280             fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1281         }
1282         if (bNeedPrec)
1283         {
1284             if (bSetPrec || !fr.bPrec)
1285             {
1286                 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1287             }
1288             else
1289             {
1290                 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1291             }
1292         }
1293
1294         if (bHaveFirstFrame)
1295         {
1296             set_trxframe_ePBC(&fr, ePBC);
1297
1298             natoms = fr.natoms;
1299
1300             if (bSetTime)
1301             {
1302                 tshift = tzero-fr.time;
1303             }
1304             else
1305             {
1306                 tzero = fr.time;
1307             }
1308
1309             bCopy = FALSE;
1310             if (bIndex)
1311             {
1312                 /* check if index is meaningful */
1313                 for (i = 0; i < nout; i++)
1314                 {
1315                     if (index[i] >= natoms)
1316                     {
1317                         gmx_fatal(FARGS,
1318                                   "Index[%d] %d is larger than the number of atoms in the\n"
1319                                   "trajectory file (%d). There is a mismatch in the contents\n"
1320                                   "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1321                     }
1322                     bCopy = bCopy || (i != index[i]);
1323                 }
1324             }
1325
1326             /* open output for writing */
1327             strcpy(filemode, "w");
1328             switch (ftp)
1329             {
1330                 case efTNG:
1331                     trjtools_gmx_prepare_tng_writing(out_file,
1332                                                      filemode[0],
1333                                                      trxin,
1334                                                      &trxout,
1335                                                      NULL,
1336                                                      nout,
1337                                                      mtop,
1338                                                      index,
1339                                                      grpnm);
1340                     break;
1341                 case efXTC:
1342                 case efTRR:
1343                     out = NULL;
1344                     if (!bSplit && !bSubTraj)
1345                     {
1346                         trxout = open_trx(out_file, filemode);
1347                     }
1348                     break;
1349                 case efGRO:
1350                 case efG96:
1351                 case efPDB:
1352                     if (( !bSeparate && !bSplit ) && !bSubTraj)
1353                     {
1354                         out = gmx_ffopen(out_file, filemode);
1355                     }
1356                     break;
1357                 default:
1358                     gmx_incons("Illegal output file format");
1359             }
1360
1361             if (bCopy)
1362             {
1363                 snew(xmem, nout);
1364                 if (bVels)
1365                 {
1366                     snew(vmem, nout);
1367                 }
1368                 if (bForce)
1369                 {
1370                     snew(fmem, nout);
1371                 }
1372             }
1373
1374             /* Start the big loop over frames */
1375             file_nr  =  0;
1376             frame    =  0;
1377             outframe =  0;
1378             model_nr =  0;
1379             bDTset   = FALSE;
1380
1381             /* Main loop over frames */
1382             do
1383             {
1384                 if (!fr.bStep)
1385                 {
1386                     /* set the step */
1387                     fr.step = newstep;
1388                     newstep++;
1389                 }
1390                 if (bSubTraj)
1391                 {
1392                     /*if (frame >= clust->clust->nra)
1393                        gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1394                     if (frame > clust->maxframe)
1395                     {
1396                         my_clust = -1;
1397                     }
1398                     else
1399                     {
1400                         my_clust = clust->inv_clust[frame];
1401                     }
1402                     if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1403                         (my_clust == NO_ATID))
1404                     {
1405                         my_clust = -1;
1406                     }
1407                 }
1408
1409                 if (bSetBox)
1410                 {
1411                     /* generate new box */
1412                     clear_mat(fr.box);
1413                     for (m = 0; m < DIM; m++)
1414                     {
1415                         if (newbox[m] >= 0)
1416                         {
1417                             fr.box[m][m] = newbox[m];
1418                         }
1419                     }
1420                 }
1421
1422                 if (bTrans)
1423                 {
1424                     for (i = 0; i < natoms; i++)
1425                     {
1426                         rvec_inc(fr.x[i], trans);
1427                     }
1428                 }
1429
1430                 if (bTDump)
1431                 {
1432                     /* determine timestep */
1433                     if (t0 == -1)
1434                     {
1435                         t0 = fr.time;
1436                     }
1437                     else
1438                     {
1439                         if (!bDTset)
1440                         {
1441                             dt     = fr.time-t0;
1442                             bDTset = TRUE;
1443                         }
1444                     }
1445                     /* This is not very elegant, as one can not dump a frame after
1446                      * a timestep with is more than twice as small as the first one. */
1447                     bDumpFrame = (fr.time > tdump-0.5*dt) && (fr.time <= tdump+0.5*dt);
1448                 }
1449                 else
1450                 {
1451                     bDumpFrame = FALSE;
1452                 }
1453
1454                 /* determine if an atom jumped across the box and reset it if so */
1455                 if (bNoJump && (bTPS || frame != 0))
1456                 {
1457                     for (d = 0; d < DIM; d++)
1458                     {
1459                         hbox[d] = 0.5*fr.box[d][d];
1460                     }
1461                     for (i = 0; i < natoms; i++)
1462                     {
1463                         if (bReset)
1464                         {
1465                             rvec_dec(fr.x[i], x_shift);
1466                         }
1467                         for (m = DIM-1; m >= 0; m--)
1468                         {
1469                             if (hbox[m] > 0)
1470                             {
1471                                 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1472                                 {
1473                                     for (d = 0; d <= m; d++)
1474                                     {
1475                                         fr.x[i][d] += fr.box[m][d];
1476                                     }
1477                                 }
1478                                 while (fr.x[i][m]-xp[i][m] > hbox[m])
1479                                 {
1480                                     for (d = 0; d <= m; d++)
1481                                     {
1482                                         fr.x[i][d] -= fr.box[m][d];
1483                                     }
1484                                 }
1485                             }
1486                         }
1487                     }
1488                 }
1489                 else if (bCluster)
1490                 {
1491                     calc_pbc_cluster(ecenter, ifit, &top, ePBC, fr.x, ind_fit, fr.box);
1492                 }
1493
1494                 if (bPFit)
1495                 {
1496                     /* Now modify the coords according to the flags,
1497                        for normal fit, this is only done for output frames */
1498                     if (bRmPBC)
1499                     {
1500                         gmx_rmpbc_trxfr(gpbc, &fr);
1501                     }
1502
1503                     reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1504                     do_fit(natoms, w_rls, xp, fr.x);
1505                 }
1506
1507                 /* store this set of coordinates for future use */
1508                 if (bPFit || bNoJump)
1509                 {
1510                     if (xp == NULL)
1511                     {
1512                         snew(xp, natoms);
1513                     }
1514                     for (i = 0; (i < natoms); i++)
1515                     {
1516                         copy_rvec(fr.x[i], xp[i]);
1517                         rvec_inc(fr.x[i], x_shift);
1518                     }
1519                 }
1520
1521                 if (frindex)
1522                 {
1523                     /* see if we have a frame from the frame index group */
1524                     for (i = 0; i < nrfri && !bDumpFrame; i++)
1525                     {
1526                         bDumpFrame = frame == frindex[i];
1527                     }
1528                 }
1529                 if (debug && bDumpFrame)
1530                 {
1531                     fprintf(debug, "dumping %d\n", frame);
1532                 }
1533
1534                 bWriteFrame =
1535                     ( ( !bTDump && !frindex && frame % skip_nr == 0 ) || bDumpFrame );
1536
1537                 if (bWriteFrame && (bDropUnder || bDropOver))
1538                 {
1539                     while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1540                     {
1541                         drop0 = drop1;
1542                         drop1++;
1543                     }
1544                     if (fabs(dropval[0][drop0] - fr.time)
1545                         < fabs(dropval[0][drop1] - fr.time))
1546                     {
1547                         dropuse = drop0;
1548                     }
1549                     else
1550                     {
1551                         dropuse = drop1;
1552                     }
1553                     if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1554                         (bDropOver && dropval[1][dropuse] > dropover))
1555                     {
1556                         bWriteFrame = FALSE;
1557                     }
1558                 }
1559
1560                 if (bWriteFrame)
1561                 {
1562                     /* We should avoid modifying the input frame,
1563                      * but since here we don't have the output frame yet,
1564                      * we introduce a temporary output frame time variable.
1565                      */
1566                     real frout_time;
1567
1568                     frout_time = fr.time;
1569
1570                     /* calc new time */
1571                     if (bTimeStep)
1572                     {
1573                         frout_time = tzero + frame*timestep;
1574                     }
1575                     else
1576                     if (bSetTime)
1577                     {
1578                         frout_time += tshift;
1579                     }
1580
1581                     if (bTDump)
1582                     {
1583                         fprintf(stderr, "\nDumping frame at t= %g %s\n",
1584                                 output_env_conv_time(oenv, frout_time), output_env_get_time_unit(oenv));
1585                     }
1586
1587                     /* check for writing at each delta_t */
1588                     bDoIt = (delta_t == 0);
1589                     if (!bDoIt)
1590                     {
1591                         if (!bRound)
1592                         {
1593                             bDoIt = bRmod(frout_time, tzero, delta_t);
1594                         }
1595                         else
1596                         {
1597                             /* round() is not C89 compatible, so we do this:  */
1598                             bDoIt = bRmod(floor(frout_time+0.5), floor(tzero+0.5),
1599                                           floor(delta_t+0.5));
1600                         }
1601                     }
1602
1603                     if (bDoIt || bTDump)
1604                     {
1605                         /* print sometimes */
1606                         if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1607                         {
1608                             fprintf(stderr, " ->  frame %6d time %8.3f      \r",
1609                                     outframe, output_env_conv_time(oenv, frout_time));
1610                         }
1611
1612                         if (!bPFit)
1613                         {
1614                             /* Now modify the coords according to the flags,
1615                                for PFit we did this already! */
1616
1617                             if (bRmPBC)
1618                             {
1619                                 gmx_rmpbc_trxfr(gpbc, &fr);
1620                             }
1621
1622                             if (bReset)
1623                             {
1624                                 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1625                                 if (bFit)
1626                                 {
1627                                     do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1628                                 }
1629                                 if (!bCenter)
1630                                 {
1631                                     for (i = 0; i < natoms; i++)
1632                                     {
1633                                         rvec_inc(fr.x[i], x_shift);
1634                                     }
1635                                 }
1636                             }
1637
1638                             if (bCenter)
1639                             {
1640                                 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1641                             }
1642                         }
1643
1644                         if (bPBCcomAtom)
1645                         {
1646                             switch (unitcell_enum)
1647                             {
1648                                 case euRect:
1649                                     put_atoms_in_box(ePBC, fr.box, natoms, fr.x);
1650                                     break;
1651                                 case euTric:
1652                                     put_atoms_in_triclinic_unitcell(ecenter, fr.box, natoms, fr.x);
1653                                     break;
1654                                 case euCompact:
1655                                     warn = put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1656                                                                          natoms, fr.x);
1657                                     if (warn && !bWarnCompact)
1658                                     {
1659                                         fprintf(stderr, "\n%s\n", warn);
1660                                         bWarnCompact = TRUE;
1661                                     }
1662                                     break;
1663                             }
1664                         }
1665                         if (bPBCcomRes)
1666                         {
1667                             put_residue_com_in_box(unitcell_enum, ecenter,
1668                                                    natoms, atoms->atom, ePBC, fr.box, fr.x);
1669                         }
1670                         if (bPBCcomMol)
1671                         {
1672                             put_molecule_com_in_box(unitcell_enum, ecenter,
1673                                                     &top.mols,
1674                                                     natoms, atoms->atom, ePBC, fr.box, fr.x);
1675                         }
1676                         /* Copy the input trxframe struct to the output trxframe struct */
1677                         frout        = fr;
1678                         frout.time   = frout_time;
1679                         frout.bV     = (frout.bV && bVels);
1680                         frout.bF     = (frout.bF && bForce);
1681                         frout.natoms = nout;
1682                         if (bNeedPrec && (bSetPrec || !fr.bPrec))
1683                         {
1684                             frout.bPrec = TRUE;
1685                             frout.prec  = prec;
1686                         }
1687                         if (bCopy)
1688                         {
1689                             frout.x = xmem;
1690                             if (frout.bV)
1691                             {
1692                                 frout.v = vmem;
1693                             }
1694                             if (frout.bF)
1695                             {
1696                                 frout.f = fmem;
1697                             }
1698                             for (i = 0; i < nout; i++)
1699                             {
1700                                 copy_rvec(fr.x[index[i]], frout.x[i]);
1701                                 if (frout.bV)
1702                                 {
1703                                     copy_rvec(fr.v[index[i]], frout.v[i]);
1704                                 }
1705                                 if (frout.bF)
1706                                 {
1707                                     copy_rvec(fr.f[index[i]], frout.f[i]);
1708                                 }
1709                             }
1710                         }
1711
1712                         if (opt2parg_bSet("-shift", NPA, pa))
1713                         {
1714                             for (i = 0; i < nout; i++)
1715                             {
1716                                 for (d = 0; d < DIM; d++)
1717                                 {
1718                                     frout.x[i][d] += outframe*shift[d];
1719                                 }
1720                             }
1721                         }
1722
1723                         if (!bRound)
1724                         {
1725                             bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1726                         }
1727                         else
1728                         {
1729                             /* round() is not C89 compatible, so we do this: */
1730                             bSplitHere = bSplit && bRmod(floor(frout.time+0.5),
1731                                                          floor(tzero+0.5),
1732                                                          floor(split_t+0.5));
1733                         }
1734                         if (bSeparate || bSplitHere)
1735                         {
1736                             mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1737                         }
1738
1739                         switch (ftp)
1740                         {
1741                             case efTNG:
1742                                 write_tng_frame(trxout, &frout);
1743                                 // TODO when trjconv behaves better: work how to read and write lambda
1744                                 break;
1745                             case efTRR:
1746                             case efXTC:
1747                                 if (bSplitHere)
1748                                 {
1749                                     if (trxout)
1750                                     {
1751                                         close_trx(trxout);
1752                                     }
1753                                     trxout = open_trx(out_file2, filemode);
1754                                 }
1755                                 if (bSubTraj)
1756                                 {
1757                                     if (my_clust != -1)
1758                                     {
1759                                         char buf[STRLEN];
1760                                         if (clust_status_id[my_clust] == -1)
1761                                         {
1762                                             sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1763                                             clust_status[my_clust]    = open_trx(buf, "w");
1764                                             clust_status_id[my_clust] = 1;
1765                                             ntrxopen++;
1766                                         }
1767                                         else if (clust_status_id[my_clust] == -2)
1768                                         {
1769                                             gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1770                                                       clust->grpname[my_clust], ntrxopen, frame,
1771                                                       my_clust);
1772                                         }
1773                                         write_trxframe(clust_status[my_clust], &frout, gc);
1774                                         nfwritten[my_clust]++;
1775                                         if (nfwritten[my_clust] ==
1776                                             (clust->clust->index[my_clust+1]-
1777                                              clust->clust->index[my_clust]))
1778                                         {
1779                                             close_trx(clust_status[my_clust]);
1780                                             clust_status[my_clust]    = NULL;
1781                                             clust_status_id[my_clust] = -2;
1782                                             ntrxopen--;
1783                                             if (ntrxopen < 0)
1784                                             {
1785                                                 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1786                                             }
1787                                         }
1788                                     }
1789                                 }
1790                                 else
1791                                 {
1792                                     write_trxframe(trxout, &frout, gc);
1793                                 }
1794                                 break;
1795                             case efGRO:
1796                             case efG96:
1797                             case efPDB:
1798                                 sprintf(title, "Generated by trjconv : %s t= %9.5f",
1799                                         top_title, frout.time);
1800                                 if (bSeparate || bSplitHere)
1801                                 {
1802                                     out = gmx_ffopen(out_file2, "w");
1803                                 }
1804                                 switch (ftp)
1805                                 {
1806                                     case efGRO:
1807                                         write_hconf_p(out, title, &useatoms, prec2ndec(frout.prec),
1808                                                       frout.x, frout.bV ? frout.v : NULL, frout.box);
1809                                         break;
1810                                     case efPDB:
1811                                         fprintf(out, "REMARK    GENERATED BY TRJCONV\n");
1812                                         sprintf(title, "%s t= %9.5f", top_title, frout.time);
1813                                         /* if reading from pdb, we want to keep the original
1814                                            model numbering else we write the output frame
1815                                            number plus one, because model 0 is not allowed in pdb */
1816                                         if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1817                                         {
1818                                             model_nr = fr.step;
1819                                         }
1820                                         else
1821                                         {
1822                                             model_nr++;
1823                                         }
1824                                         write_pdbfile(out, title, &useatoms, frout.x,
1825                                                       frout.ePBC, frout.box, ' ', model_nr, gc, TRUE);
1826                                         break;
1827                                     case efG96:
1828                                         frout.title = title;
1829                                         if (bSeparate || bTDump)
1830                                         {
1831                                             frout.bTitle = TRUE;
1832                                             if (bTPS)
1833                                             {
1834                                                 frout.bAtoms = TRUE;
1835                                             }
1836                                             frout.atoms  = &useatoms;
1837                                             frout.bStep  = FALSE;
1838                                             frout.bTime  = FALSE;
1839                                         }
1840                                         else
1841                                         {
1842                                             frout.bTitle = (outframe == 0);
1843                                             frout.bAtoms = FALSE;
1844                                             frout.bStep  = TRUE;
1845                                             frout.bTime  = TRUE;
1846                                         }
1847                                         write_g96_conf(out, &frout, -1, NULL);
1848                                 }
1849                                 if (bSeparate || bSplitHere)
1850                                 {
1851                                     gmx_ffclose(out);
1852                                     out = NULL;
1853                                 }
1854                                 break;
1855                             default:
1856                                 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1857                         }
1858                         if (bSeparate || bSplitHere)
1859                         {
1860                             file_nr++;
1861                         }
1862
1863                         /* execute command */
1864                         if (bExec)
1865                         {
1866                             char c[255];
1867                             sprintf(c, "%s  %d", exec_command, file_nr-1);
1868                             /*fprintf(stderr,"Executing '%s'\n",c);*/
1869                             if (0 != system(c))
1870                             {
1871                                 gmx_fatal(FARGS, "Error executing command: %s", c);
1872                             }
1873                         }
1874                         outframe++;
1875                     }
1876                 }
1877                 frame++;
1878                 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1879             }
1880             while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1881         }
1882
1883         if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1884         {
1885             fprintf(stderr, "\nWARNING no output, "
1886                     "last frame read at t=%g\n", fr.time);
1887         }
1888         fprintf(stderr, "\n");
1889
1890         close_trj(trxin);
1891         sfree(outf_base);
1892
1893         if (bRmPBC)
1894         {
1895             gmx_rmpbc_done(gpbc);
1896         }
1897
1898         if (trxout)
1899         {
1900             close_trx(trxout);
1901         }
1902         else if (out != NULL)
1903         {
1904             gmx_ffclose(out);
1905         }
1906         if (bSubTraj)
1907         {
1908             for (i = 0; (i < clust->clust->nr); i++)
1909             {
1910                 if (clust_status_id[i] >= 0)
1911                 {
1912                     close_trx(clust_status[i]);
1913                 }
1914             }
1915         }
1916     }
1917
1918     sfree(mtop);
1919
1920     do_view(oenv, out_file, NULL);
1921
1922     return 0;
1923 }