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