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