Create separate module for trajectory data
[gromacs.git] / src / gromacs / gmxana / gmx_trjconv.cpp
blob59fa86a8d91ec3b5d3b76dfac026fe58299035f6
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdlib>
41 #include <cstring>
43 #include <algorithm>
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/commandline/viewit.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/fileio/g96io.h"
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/fileio/groio.h"
51 #include "gromacs/fileio/pdbio.h"
52 #include "gromacs/fileio/tngio_for_tools.h"
53 #include "gromacs/fileio/tpxio.h"
54 #include "gromacs/fileio/trrio.h"
55 #include "gromacs/fileio/trxio.h"
56 #include "gromacs/fileio/xtcio.h"
57 #include "gromacs/fileio/xvgr.h"
58 #include "gromacs/gmxana/gmx_ana.h"
59 #include "gromacs/math/do_fit.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/mdtypes/md_enums.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/pbcutil/rmpbc.h"
64 #include "gromacs/topology/index.h"
65 #include "gromacs/topology/topology.h"
66 #include "gromacs/trajectory/trajectoryframe.h"
67 #include "gromacs/utility/arraysize.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/futil.h"
70 #include "gromacs/utility/smalloc.h"
72 enum {
73 euSel, euRect, euTric, euCompact, euNR
77 static void calc_pbc_cluster(int ecenter, int nrefat, t_topology *top, int ePBC,
78 rvec x[], int index[], matrix box)
80 int m, i, j, j0, j1, jj, ai, aj;
81 int imin, jmin;
82 real fac, min_dist2;
83 rvec dx, xtest, box_center;
84 int nmol, imol_center;
85 int *molind;
86 gmx_bool *bMol, *bTmp;
87 rvec *m_com, *m_shift;
88 t_pbc pbc;
89 int *cluster;
90 int *added;
91 int ncluster, nadded;
92 real tmp_r2;
94 calc_box_center(ecenter, box, box_center);
96 /* Initiate the pbc structure */
97 std::memset(&pbc, 0, sizeof(pbc));
98 set_pbc(&pbc, ePBC, box);
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);
110 for (i = 0; (i < nrefat); i++)
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)
120 if (ai < molind[j0+1])
122 j1 = j0;
124 else if (ai >= molind[j1])
126 j0 = j1;
128 else
130 jj = (j0+j1)/2;
131 if (ai < molind[jj+1])
133 j1 = jj;
135 else
137 j0 = jj;
141 bMol[j0] = TRUE;
143 /* Double check whether all atoms in all molecules that are marked are part
144 * of the cluster. Simultaneously compute the center of geometry.
146 min_dist2 = 10*sqr(trace(box));
147 imol_center = -1;
148 ncluster = 0;
149 for (i = 0; i < nmol; i++)
151 for (j = molind[i]; j < molind[i+1]; j++)
153 if (bMol[i] && !bTmp[j])
155 gmx_fatal(FARGS, "Molecule %d marked for clustering but not atom %d in it - check your index!", i+1, j+1);
157 else if (!bMol[i] && bTmp[j])
159 gmx_fatal(FARGS, "Atom %d marked for clustering but not molecule %d - this is an internal error...", j+1, i+1);
161 else if (bMol[i])
163 /* Make molecule whole, move 2nd and higher atom to same periodicity as 1st atom in molecule */
164 if (j > molind[i])
166 pbc_dx(&pbc, x[j], x[j-1], dx);
167 rvec_add(x[j-1], dx, x[j]);
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]);
173 if (bMol[i])
175 /* Normalize center of geometry */
176 fac = 1.0/(molind[i+1]-molind[i]);
177 for (m = 0; (m < DIM); m++)
179 m_com[i][m] *= fac;
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);
185 if (tmp_r2 < min_dist2)
187 min_dist2 = tmp_r2;
188 imol_center = i;
190 cluster[ncluster++] = i;
193 sfree(bTmp);
195 if (ncluster <= 0)
197 fprintf(stderr, "No molecules selected in the cluster\n");
198 return;
200 else if (imol_center == -1)
202 fprintf(stderr, "No central molecules could be found\n");
203 return;
206 nadded = 0;
207 added[nadded++] = imol_center;
208 bMol[imol_center] = FALSE;
210 while (nadded < ncluster)
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++)
219 ai = added[i];
220 /* Loop over all mols */
221 for (j = 0; j < ncluster; j++)
223 aj = cluster[j];
224 /* check those remaining to be added */
225 if (bMol[aj])
227 pbc_dx(&pbc, m_com[aj], m_com[ai], dx);
228 tmp_r2 = iprod(dx, dx);
229 if (tmp_r2 < min_dist2)
231 min_dist2 = tmp_r2;
232 imin = ai;
233 jmin = aj;
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]);
248 for (j = molind[jmin]; j < molind[jmin+1]; j++)
250 rvec_inc(x[j], m_shift[jmin]);
252 fprintf(stdout, "\rClustering iteration %d of %d...", nadded, ncluster);
253 fflush(stdout);
256 sfree(added);
257 sfree(cluster);
258 sfree(bMol);
259 sfree(m_com);
260 sfree(m_shift);
262 fprintf(stdout, "\n");
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[])
270 int i, j;
271 int d;
272 rvec com, new_com, shift, box_center;
273 real m;
274 double mtot;
275 t_pbc pbc;
277 calc_box_center(ecenter, box, box_center);
278 set_pbc(&pbc, ePBC, box);
279 if (mols->nr <= 0)
281 gmx_fatal(FARGS, "There are no molecule descriptions. I need a .tpr file for this pbc option.");
283 for (i = 0; (i < mols->nr); i++)
285 /* calc COM */
286 clear_rvec(com);
287 mtot = 0;
288 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
290 m = atom[j].m;
291 for (d = 0; d < DIM; d++)
293 com[d] += m*x[j][d];
295 mtot += m;
297 /* calculate final COM */
298 svmul(1.0/mtot, com, com);
300 /* check if COM is outside box */
301 copy_rvec(com, new_com);
302 switch (unitcell_enum)
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;
314 rvec_sub(new_com, com, shift);
315 if (norm2(shift) > 0)
317 if (debug)
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]);
323 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
325 rvec_inc(x[j], shift);
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[])
335 int i, j, res_start, res_end;
336 int d, presnr;
337 real m;
338 double mtot;
339 rvec box_center, com, new_com, shift;
340 static const int NOTSET = -12347;
341 calc_box_center(ecenter, box, box_center);
343 presnr = NOTSET;
344 res_start = 0;
345 clear_rvec(com);
346 mtot = 0;
347 for (i = 0; i < natoms+1; i++)
349 if (i == natoms || (presnr != atom[i].resind && presnr != NOTSET))
351 /* calculate final COM */
352 res_end = i;
353 svmul(1.0/mtot, com, com);
355 /* check if COM is outside box */
356 copy_rvec(com, new_com);
357 switch (unitcell_enum)
359 case euRect:
360 put_atoms_in_box(ePBC, box, 1, &new_com);
361 break;
362 case euTric:
363 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
364 break;
365 case euCompact:
366 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
367 break;
369 rvec_sub(new_com, com, shift);
370 if (norm2(shift))
372 if (debug)
374 fprintf(debug, "\nShifting position of residue %d (atoms %d-%d) "
375 "by %g,%g,%g\n", atom[res_start].resind+1,
376 res_start+1, res_end+1, shift[XX], shift[YY], shift[ZZ]);
378 for (j = res_start; j < res_end; j++)
380 rvec_inc(x[j], shift);
383 clear_rvec(com);
384 mtot = 0;
386 /* remember start of new residue */
387 res_start = i;
389 if (i < natoms)
391 /* calc COM */
392 m = atom[i].m;
393 for (d = 0; d < DIM; d++)
395 com[d] += m*x[i][d];
397 mtot += m;
399 presnr = atom[i].resind;
404 static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, int ci[])
406 int i, m, ai;
407 rvec cmin, cmax, box_center, dx;
409 if (nc > 0)
411 copy_rvec(x[ci[0]], cmin);
412 copy_rvec(x[ci[0]], cmax);
413 for (i = 0; i < nc; i++)
415 ai = ci[i];
416 for (m = 0; m < DIM; m++)
418 if (x[ai][m] < cmin[m])
420 cmin[m] = x[ai][m];
422 else if (x[ai][m] > cmax[m])
424 cmax[m] = x[ai][m];
428 calc_box_center(ecenter, box, box_center);
429 for (m = 0; m < DIM; m++)
431 dx[m] = box_center[m]-(cmin[m]+cmax[m])*0.5;
434 for (i = 0; i < n; i++)
436 rvec_inc(x[i], dx);
441 static void mk_filenm(char *base, const char *ext, int ndigit, int file_nr,
442 char out_file[])
444 char nbuf[128];
445 int nd = 0, fnr;
447 std::strcpy(out_file, base);
448 fnr = file_nr;
451 fnr /= 10;
452 nd++;
454 while (fnr > 0);
456 if (nd < ndigit)
458 std::strncat(out_file, "00000000000", ndigit-nd);
460 sprintf(nbuf, "%d.", file_nr);
461 std::strcat(out_file, nbuf);
462 std::strcat(out_file, ext);
465 void check_trr(const char *fn)
467 if (fn2ftp(fn) != efTRR)
469 gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
473 void do_trunc(const char *fn, real t0)
475 t_fileio *in;
476 FILE *fp;
477 gmx_bool bStop, bOK;
478 gmx_trr_header_t sh;
479 gmx_off_t fpos;
480 char yesno[256];
481 int j;
482 real t = 0;
484 if (t0 == -1)
486 gmx_fatal(FARGS, "You forgot to set the truncation time");
489 /* Check whether this is a .trr file */
490 check_trr(fn);
492 in = gmx_trr_open(fn, "r");
493 fp = gmx_fio_getfp(in);
494 if (fp == NULL)
496 fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
497 gmx_trr_close(in);
499 else
501 j = 0;
502 fpos = gmx_fio_ftell(in);
503 bStop = FALSE;
504 while (!bStop && gmx_trr_read_frame_header(in, &sh, &bOK))
506 gmx_trr_read_frame_data(in, &sh, NULL, NULL, NULL, NULL);
507 fpos = gmx_ftell(fp);
508 t = sh.t;
509 if (t >= t0)
511 gmx_fseek(fp, fpos, SEEK_SET);
512 bStop = TRUE;
515 if (bStop)
517 fprintf(stderr, "Do you REALLY want to truncate this trajectory (%s) at:\n"
518 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
519 fn, j, t, (long int)fpos);
520 if (1 != scanf("%s", yesno))
522 gmx_fatal(FARGS, "Error reading user input");
524 if (std::strcmp(yesno, "YES") == 0)
526 fprintf(stderr, "Once again, I'm gonna DO this...\n");
527 gmx_trr_close(in);
528 if (0 != gmx_truncate(fn, fpos))
530 gmx_fatal(FARGS, "Error truncating file %s", fn);
533 else
535 fprintf(stderr, "Ok, I'll forget about it\n");
538 else
540 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
541 gmx_trr_close(in);
546 /*! \brief Read a full molecular topology if useful and available.
548 * If the input trajectory file is not in TNG format, and the output
549 * file is in TNG format, then we want to try to read a full topology
550 * (if available), so that we can write molecule information to the
551 * output file. The full topology provides better molecule information
552 * than is available from the normal t_topology data used by GROMACS
553 * tools.
555 * Also, the t_topology is only read under (different) particular
556 * conditions. If both apply, then a .tpr file might be read
557 * twice. Trying to fix this redundancy while trjconv is still an
558 * all-purpose tool does not seem worthwhile.
560 * Because of the way gmx_prepare_tng_writing is implemented, the case
561 * where the input TNG file has no molecule information will never
562 * lead to an output TNG file having molecule information. Since
563 * molecule information will generally be present if the input TNG
564 * file was written by a GROMACS tool, this seems like reasonable
565 * behaviour. */
566 static gmx_mtop_t *read_mtop_for_tng(const char *tps_file,
567 const char *input_file,
568 const char *output_file)
570 gmx_mtop_t *mtop = NULL;
572 if (fn2bTPX(tps_file) &&
573 efTNG != fn2ftp(input_file) &&
574 efTNG == fn2ftp(output_file))
576 int temp_natoms = -1;
577 snew(mtop, 1);
578 read_tpx(tps_file, NULL, NULL, &temp_natoms,
579 NULL, NULL, mtop);
582 return mtop;
585 int gmx_trjconv(int argc, char *argv[])
587 const char *desc[] = {
588 "[THISMODULE] can convert trajectory files in many ways:",
590 "* from one format to another",
591 "* select a subset of atoms",
592 "* change the periodicity representation",
593 "* keep multimeric molecules together",
594 "* center atoms in the box",
595 "* fit atoms to reference structure",
596 "* reduce the number of frames",
597 "* change the timestamps of the frames ([TT]-t0[tt] and [TT]-timestep[tt])",
598 "* cut the trajectory in small subtrajectories according",
599 " to information in an index file. This allows subsequent analysis of",
600 " the subtrajectories that could, for example, be the result of a",
601 " cluster analysis. Use option [TT]-sub[tt].",
602 " This assumes that the entries in the index file are frame numbers and",
603 " dumps each group in the index file to a separate trajectory file.",
604 "* select frames within a certain range of a quantity given",
605 " in an [REF].xvg[ref] file.",
607 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
608 "[PAR]",
610 "The following formats are supported for input and output:",
611 "[REF].xtc[ref], [REF].trr[ref], [REF].gro[ref], [TT].g96[tt]",
612 "and [REF].pdb[ref].",
613 "The file formats are detected from the file extension.",
614 "The precision of [REF].xtc[ref] and [REF].gro[ref] output is taken from the",
615 "input file for [REF].xtc[ref], [REF].gro[ref] and [REF].pdb[ref],",
616 "and from the [TT]-ndec[tt] option for other input formats. The precision",
617 "is always taken from [TT]-ndec[tt], when this option is set.",
618 "All other formats have fixed precision. [REF].trr[ref]",
619 "output can be single or double precision, depending on the precision",
620 "of the [THISMODULE] binary.",
621 "Note that velocities are only supported in",
622 "[REF].trr[ref], [REF].gro[ref] and [TT].g96[tt] files.[PAR]",
624 "Option [TT]-sep[tt] can be used to write every frame to a separate",
625 "[TT].gro, .g96[tt] or [REF].pdb[ref] file. By default, all frames all written to one file.",
626 "[REF].pdb[ref] files with all frames concatenated can be viewed with",
627 "[TT]rasmol -nmrpdb[tt].[PAR]",
629 "It is possible to select part of your trajectory and write it out",
630 "to a new trajectory file in order to save disk space, e.g. for leaving",
631 "out the water from a trajectory of a protein in water.",
632 "[BB]ALWAYS[bb] put the original trajectory on tape!",
633 "We recommend to use the portable [REF].xtc[ref] format for your analysis",
634 "to save disk space and to have portable files.[PAR]",
636 "There are two options for fitting the trajectory to a reference",
637 "either for essential dynamics analysis, etc.",
638 "The first option is just plain fitting to a reference structure",
639 "in the structure file. The second option is a progressive fit",
640 "in which the first timeframe is fitted to the reference structure ",
641 "in the structure file to obtain and each subsequent timeframe is ",
642 "fitted to the previously fitted structure. This way a continuous",
643 "trajectory is generated, which might not be the case when using the",
644 "regular fit method, e.g. when your protein undergoes large",
645 "conformational transitions.[PAR]",
647 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
648 "treatment:",
650 " * [TT]mol[tt] puts the center of mass of molecules in the box,",
651 " and requires a run input file to be supplied with [TT]-s[tt].",
652 " * [TT]res[tt] puts the center of mass of residues in the box.",
653 " * [TT]atom[tt] puts all the atoms in the box.",
654 " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
655 " them back. This has the effect that all molecules",
656 " will remain whole (provided they were whole in the initial",
657 " conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
658 " molecules may diffuse out of the box. The starting configuration",
659 " for this procedure is taken from the structure file, if one is",
660 " supplied, otherwise it is the first frame.",
661 " * [TT]cluster[tt] clusters all the atoms in the selected index",
662 " such that they are all closest to the center of mass of the cluster,",
663 " which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
664 " results if you in fact have a cluster. Luckily that can be checked",
665 " afterwards using a trajectory viewer. Note also that if your molecules",
666 " are broken this will not work either.",
668 " The separate option [TT]-clustercenter[tt] can be used to specify an",
669 " approximate center for the cluster. This is useful e.g. if you have",
670 " two big vesicles, and you want to maintain their relative positions.",
671 " * [TT]whole[tt] only makes broken molecules whole.",
674 "Option [TT]-ur[tt] sets the unit cell representation for options",
675 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
676 "All three options give different results for triclinic boxes and",
677 "identical results for rectangular boxes.",
678 "[TT]rect[tt] is the ordinary brick shape.",
679 "[TT]tric[tt] is the triclinic unit cell.",
680 "[TT]compact[tt] puts all atoms at the closest distance from the center",
681 "of the box. This can be useful for visualizing e.g. truncated octahedra",
682 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
683 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
684 "is set differently.[PAR]",
686 "Option [TT]-center[tt] centers the system in the box. The user can",
687 "select the group which is used to determine the geometrical center.",
688 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
689 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
690 "[TT]tric[tt]: half of the sum of the box vectors,",
691 "[TT]rect[tt]: half of the box diagonal,",
692 "[TT]zero[tt]: zero.",
693 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
694 "want all molecules in the box after the centering.[PAR]",
696 "Option [TT]-box[tt] sets the size of the new box. This option only works",
697 "for leading dimensions and is thus generally only useful for rectangular boxes.",
698 "If you want to modify only some of the dimensions, e.g. when reading from",
699 "a trajectory, you can use -1 for those dimensions that should stay the same",
701 "It is not always possible to use combinations of [TT]-pbc[tt],",
702 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
703 "you want in one call to [THISMODULE]. Consider using multiple",
704 "calls, and check out the GROMACS website for suggestions.[PAR]",
706 "With [TT]-dt[tt], it is possible to reduce the number of ",
707 "frames in the output. This option relies on the accuracy of the times",
708 "in your input trajectory, so if these are inaccurate use the",
709 "[TT]-timestep[tt] option to modify the time (this can be done",
710 "simultaneously). For making smooth movies, the program [gmx-filter]",
711 "can reduce the number of frames while using low-pass frequency",
712 "filtering, this reduces aliasing of high frequency motions.[PAR]",
714 "Using [TT]-trunc[tt] [THISMODULE] can truncate [REF].trr[ref] in place, i.e.",
715 "without copying the file. This is useful when a run has crashed",
716 "during disk I/O (i.e. full disk), or when two contiguous",
717 "trajectories must be concatenated without having double frames.[PAR]",
719 "Option [TT]-dump[tt] can be used to extract a frame at or near",
720 "one specific time from your trajectory, but only works reliably",
721 "if the time interval between frames is uniform.[PAR]",
723 "Option [TT]-drop[tt] reads an [REF].xvg[ref] file with times and values.",
724 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
725 "frames with a value below and above the value of the respective options",
726 "will not be written."
729 int pbc_enum;
730 enum
732 epSel,
733 epNone,
734 epComMol,
735 epComRes,
736 epComAtom,
737 epNojump,
738 epCluster,
739 epWhole,
740 epNR
742 const char *pbc_opt[epNR + 1] =
744 NULL, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
745 NULL
748 int unitcell_enum;
749 const char *unitcell_opt[euNR+1] =
750 { NULL, "rect", "tric", "compact", NULL };
752 enum
754 ecSel, ecTric, ecRect, ecZero, ecNR
756 const char *center_opt[ecNR+1] =
757 { NULL, "tric", "rect", "zero", NULL };
758 int ecenter;
760 int fit_enum;
761 enum
763 efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
765 const char *fit[efNR + 1] =
767 NULL, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
768 "progressive", NULL
771 static gmx_bool bSeparate = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
772 static gmx_bool bCenter = FALSE;
773 static int skip_nr = 1, ndec = 3, nzero = 0;
774 static real tzero = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
775 static rvec newbox = {0, 0, 0}, shift = {0, 0, 0}, trans = {0, 0, 0};
776 static char *exec_command = NULL;
777 static real dropunder = 0, dropover = 0;
778 static gmx_bool bRound = FALSE;
780 t_pargs
781 pa[] =
783 { "-skip", FALSE, etINT,
784 { &skip_nr }, "Only write every nr-th frame" },
785 { "-dt", FALSE, etTIME,
786 { &delta_t },
787 "Only write frame when t MOD dt = first time (%t)" },
788 { "-round", FALSE, etBOOL,
789 { &bRound }, "Round measurements to nearest picosecond"},
790 { "-dump", FALSE, etTIME,
791 { &tdump }, "Dump frame nearest specified time (%t)" },
792 { "-t0", FALSE, etTIME,
793 { &tzero },
794 "Starting time (%t) (default: don't change)" },
795 { "-timestep", FALSE, etTIME,
796 { &timestep },
797 "Change time step between input frames (%t)" },
798 { "-pbc", FALSE, etENUM,
799 { pbc_opt },
800 "PBC treatment (see help text for full description)" },
801 { "-ur", FALSE, etENUM,
802 { unitcell_opt }, "Unit-cell representation" },
803 { "-center", FALSE, etBOOL,
804 { &bCenter }, "Center atoms in box" },
805 { "-boxcenter", FALSE, etENUM,
806 { center_opt }, "Center for -pbc and -center" },
807 { "-box", FALSE, etRVEC,
808 { newbox },
809 "Size for new cubic box (default: read from input)" },
810 { "-trans", FALSE, etRVEC,
811 { trans },
812 "All coordinates will be translated by trans. This "
813 "can advantageously be combined with -pbc mol -ur "
814 "compact." },
815 { "-shift", FALSE, etRVEC,
816 { shift },
817 "All coordinates will be shifted by framenr*shift" },
818 { "-fit", FALSE, etENUM,
819 { fit },
820 "Fit molecule to ref structure in the structure file" },
821 { "-ndec", FALSE, etINT,
822 { &ndec },
823 "Precision for .xtc and .gro writing in number of "
824 "decimal places" },
825 { "-vel", FALSE, etBOOL,
826 { &bVels }, "Read and write velocities if possible" },
827 { "-force", FALSE, etBOOL,
828 { &bForce }, "Read and write forces if possible" },
829 { "-trunc", FALSE, etTIME,
830 { &ttrunc },
831 "Truncate input trajectory file after this time (%t)" },
832 { "-exec", FALSE, etSTR,
833 { &exec_command },
834 "Execute command for every output frame with the "
835 "frame number as argument" },
836 { "-split", FALSE, etTIME,
837 { &split_t },
838 "Start writing new file when t MOD split = first "
839 "time (%t)" },
840 { "-sep", FALSE, etBOOL,
841 { &bSeparate },
842 "Write each frame to a separate .gro, .g96 or .pdb "
843 "file" },
844 { "-nzero", FALSE, etINT,
845 { &nzero },
846 "If the -sep flag is set, use these many digits "
847 "for the file numbers and prepend zeros as needed" },
848 { "-dropunder", FALSE, etREAL,
849 { &dropunder }, "Drop all frames below this value" },
850 { "-dropover", FALSE, etREAL,
851 { &dropover }, "Drop all frames above this value" },
852 { "-conect", FALSE, etBOOL,
853 { &bCONECT },
854 "Add conect records when writing [REF].pdb[ref] files. Useful "
855 "for visualization of non-standard molecules, e.g. "
856 "coarse grained ones" }
858 #define NPA asize(pa)
860 FILE *out = NULL;
861 t_trxstatus *trxout = NULL;
862 t_trxstatus *trxin;
863 int ftp, ftpin = 0, file_nr;
864 t_trxframe fr, frout;
865 int flags;
866 rvec *xmem = NULL, *vmem = NULL, *fmem = NULL;
867 rvec *xp = NULL, x_shift, hbox;
868 real *w_rls = NULL;
869 int m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
870 #define SKIP 10
871 t_topology top;
872 gmx_mtop_t *mtop = NULL;
873 gmx_conect gc = NULL;
874 int ePBC = -1;
875 t_atoms *atoms = NULL, useatoms;
876 matrix top_box;
877 int *index, *cindex;
878 char *grpnm;
879 int *frindex, nrfri;
880 char *frname;
881 int ifit, my_clust = -1;
882 int *ind_fit;
883 char *gn_fit;
884 t_cluster_ndx *clust = NULL;
885 t_trxstatus **clust_status = NULL;
886 int *clust_status_id = NULL;
887 int ntrxopen = 0;
888 int *nfwritten = NULL;
889 int ndrop = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
890 double **dropval;
891 real tshift = 0, t0 = -1, dt = 0.001, prec;
892 gmx_bool bFit, bPFit, bReset;
893 int nfitdim;
894 gmx_rmpbc_t gpbc = NULL;
895 gmx_bool bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
896 gmx_bool bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
897 gmx_bool bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetPrec, bNeedPrec;
898 gmx_bool bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
899 gmx_bool bSubTraj = FALSE, bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
900 gmx_bool bWriteFrame, bSplitHere;
901 const char *top_file, *in_file, *out_file = NULL;
902 char out_file2[256], *charpt;
903 char *outf_base = NULL;
904 const char *outf_ext = NULL;
905 char top_title[256], title[256], filemode[5];
906 gmx_output_env_t *oenv;
908 t_filenm fnm[] = {
909 { efTRX, "-f", NULL, ffREAD },
910 { efTRO, "-o", NULL, ffWRITE },
911 { efTPS, NULL, NULL, ffOPTRD },
912 { efNDX, NULL, NULL, ffOPTRD },
913 { efNDX, "-fr", "frames", ffOPTRD },
914 { efNDX, "-sub", "cluster", ffOPTRD },
915 { efXVG, "-drop", "drop", ffOPTRD }
917 #define NFILE asize(fnm)
919 if (!parse_common_args(&argc, argv,
920 PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW |
921 PCA_TIME_UNIT,
922 NFILE, fnm, NPA, pa, asize(desc), desc,
923 0, NULL, &oenv))
925 return 0;
928 top_file = ftp2fn(efTPS, NFILE, fnm);
929 init_top(&top);
931 /* Check command line */
932 in_file = opt2fn("-f", NFILE, fnm);
934 if (ttrunc != -1)
936 do_trunc(in_file, ttrunc);
938 else
940 /* mark active cmdline options */
941 bSetBox = opt2parg_bSet("-box", NPA, pa);
942 bSetTime = opt2parg_bSet("-t0", NPA, pa);
943 bSetPrec = opt2parg_bSet("-ndec", NPA, pa);
944 bSetUR = opt2parg_bSet("-ur", NPA, pa);
945 bExec = opt2parg_bSet("-exec", NPA, pa);
946 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
947 bTDump = opt2parg_bSet("-dump", NPA, pa);
948 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
949 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
950 bTrans = opt2parg_bSet("-trans", NPA, pa);
951 bSplit = (split_t != 0);
953 /* parse enum options */
954 fit_enum = nenum(fit);
955 bFit = (fit_enum == efFit || fit_enum == efFitXY);
956 bReset = (fit_enum == efReset || fit_enum == efResetXY);
957 bPFit = fit_enum == efPFit;
958 pbc_enum = nenum(pbc_opt);
959 bPBCWhole = pbc_enum == epWhole;
960 bPBCcomRes = pbc_enum == epComRes;
961 bPBCcomMol = pbc_enum == epComMol;
962 bPBCcomAtom = pbc_enum == epComAtom;
963 bNoJump = pbc_enum == epNojump;
964 bCluster = pbc_enum == epCluster;
965 bPBC = pbc_enum != epNone;
966 unitcell_enum = nenum(unitcell_opt);
967 ecenter = nenum(center_opt) - ecTric;
969 /* set and check option dependencies */
970 if (bPFit)
972 bFit = TRUE; /* for pfit, fit *must* be set */
974 if (bFit)
976 bReset = TRUE; /* for fit, reset *must* be set */
978 nfitdim = 0;
979 if (bFit || bReset)
981 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
983 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
985 if (bSetUR)
987 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
989 fprintf(stderr,
990 "WARNING: Option for unitcell representation (-ur %s)\n"
991 " only has effect in combination with -pbc %s, %s or %s.\n"
992 " Ingoring unitcell representation.\n\n",
993 unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
996 if (bFit && bPBC)
998 gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
999 "Please do the PBC condition treatment first and then run trjconv in a second step\n"
1000 "for the rotational fit.\n"
1001 "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1002 "results!");
1005 /* ndec is in nr of decimal places, prec is a multiplication factor: */
1006 prec = 1;
1007 for (i = 0; i < ndec; i++)
1009 prec *= 10;
1012 bIndex = ftp2bSet(efNDX, NFILE, fnm);
1015 /* Determine output type */
1016 out_file = opt2fn("-o", NFILE, fnm);
1017 ftp = fn2ftp(out_file);
1018 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
1019 bNeedPrec = (ftp == efXTC || ftp == efGRO);
1020 if (bVels)
1022 /* check if velocities are possible in input and output files */
1023 ftpin = fn2ftp(in_file);
1024 bVels = (ftp == efTRR || ftp == efGRO ||
1025 ftp == efG96 || ftp == efTNG)
1026 && (ftpin == efTRR || ftpin == efGRO ||
1027 ftpin == efG96 || ftpin == efTNG || ftpin == efCPT);
1029 if (bSeparate || bSplit)
1031 outf_ext = std::strrchr(out_file, '.');
1032 if (outf_ext == NULL)
1034 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1036 outf_base = gmx_strdup(out_file);
1037 outf_base[outf_ext - out_file] = '\0';
1040 bSubTraj = opt2bSet("-sub", NFILE, fnm);
1041 if (bSubTraj)
1043 if ((ftp != efXTC) && (ftp != efTRR))
1045 /* It seems likely that other trajectory file types
1046 * could work here. */
1047 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1048 "xtc and trr");
1050 clust = cluster_index(NULL, opt2fn("-sub", NFILE, fnm));
1052 /* Check for number of files disabled, as FOPEN_MAX is not the correct
1053 * number to check for. In my linux box it is only 16.
1055 if (0 && (clust->clust->nr > FOPEN_MAX-4))
1057 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1058 " trajectories.\ntry splitting the index file in %d parts.\n"
1059 "FOPEN_MAX = %d",
1060 clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1062 gmx_warning("The -sub option could require as many open output files as there are\n"
1063 "index groups in the file (%d). If you get I/O errors opening new files,\n"
1064 "try reducing the number of index groups in the file, and perhaps\n"
1065 "using trjconv -sub several times on different chunks of your index file.\n",
1066 clust->clust->nr);
1068 snew(clust_status, clust->clust->nr);
1069 snew(clust_status_id, clust->clust->nr);
1070 snew(nfwritten, clust->clust->nr);
1071 for (i = 0; (i < clust->clust->nr); i++)
1073 clust_status[i] = NULL;
1074 clust_status_id[i] = -1;
1076 bSeparate = bSplit = FALSE;
1078 /* skipping */
1079 if (skip_nr <= 0)
1083 mtop = read_mtop_for_tng(top_file, in_file, out_file);
1085 /* Determine whether to read a topology */
1086 bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1087 bRmPBC || bReset || bPBCcomMol || bCluster ||
1088 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1090 /* Determine if when can read index groups */
1091 bIndex = (bIndex || bTPS);
1093 if (bTPS)
1095 read_tps_conf(top_file, &top, &ePBC, &xp, NULL, top_box,
1096 bReset || bPBCcomRes);
1097 std::strncpy(top_title, *top.name, 255);
1098 top_title[255] = '\0';
1099 atoms = &top.atoms;
1101 if (0 == top.mols.nr && (bCluster || bPBCcomMol))
1103 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1106 /* top_title is only used for gro and pdb,
1107 * the header in such a file is top_title t= ...
1108 * to prevent a double t=, remove it from top_title
1110 if ((charpt = std::strstr(top_title, " t= ")))
1112 charpt[0] = '\0';
1115 if (bCONECT)
1117 gc = gmx_conect_generate(&top);
1119 if (bRmPBC)
1121 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
1125 /* get frame number index */
1126 frindex = NULL;
1127 if (opt2bSet("-fr", NFILE, fnm))
1129 printf("Select groups of frame number indices:\n");
1130 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, (int **)&frindex, &frname);
1131 if (debug)
1133 for (i = 0; i < nrfri; i++)
1135 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1140 /* get index groups etc. */
1141 if (bReset)
1143 printf("Select group for %s fit\n",
1144 bFit ? "least squares" : "translational");
1145 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1146 1, &ifit, &ind_fit, &gn_fit);
1148 if (bFit)
1150 if (ifit < 2)
1152 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1154 else if (ifit == 3)
1156 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1160 else if (bCluster)
1162 printf("Select group for clustering\n");
1163 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1164 1, &ifit, &ind_fit, &gn_fit);
1167 if (bIndex)
1169 if (bCenter)
1171 printf("Select group for centering\n");
1172 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1173 1, &ncent, &cindex, &grpnm);
1175 printf("Select group for output\n");
1176 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1177 1, &nout, &index, &grpnm);
1179 else
1181 /* no index file, so read natoms from TRX */
1182 if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
1184 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1186 natoms = fr.natoms;
1187 close_trj(trxin);
1188 sfree(fr.x);
1189 snew(index, natoms);
1190 for (i = 0; i < natoms; i++)
1192 index[i] = i;
1194 nout = natoms;
1195 if (bCenter)
1197 ncent = nout;
1198 cindex = index;
1202 if (bReset)
1204 snew(w_rls, atoms->nr);
1205 for (i = 0; (i < ifit); i++)
1207 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1210 /* Restore reference structure and set to origin,
1211 store original location (to put structure back) */
1212 if (bRmPBC)
1214 gmx_rmpbc(gpbc, top.atoms.nr, top_box, xp);
1216 copy_rvec(xp[index[0]], x_shift);
1217 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, NULL, xp, w_rls);
1218 rvec_dec(x_shift, xp[index[0]]);
1220 else
1222 clear_rvec(x_shift);
1225 if (bDropUnder || bDropOver)
1227 /* Read the .xvg file with the drop values */
1228 fprintf(stderr, "\nReading drop file ...");
1229 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1230 fprintf(stderr, " %d time points\n", ndrop);
1231 if (ndrop == 0 || ncol < 2)
1233 gmx_fatal(FARGS, "Found no data points in %s",
1234 opt2fn("-drop", NFILE, fnm));
1236 drop0 = 0;
1237 drop1 = 0;
1240 /* Make atoms struct for output in GRO or PDB files */
1241 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1243 /* get memory for stuff to go in .pdb file, and initialize
1244 * the pdbinfo structure part if the input has it.
1246 init_t_atoms(&useatoms, atoms->nr, (atoms->pdbinfo != NULL));
1247 sfree(useatoms.resinfo);
1248 useatoms.resinfo = atoms->resinfo;
1249 for (i = 0; (i < nout); i++)
1251 useatoms.atomname[i] = atoms->atomname[index[i]];
1252 useatoms.atom[i] = atoms->atom[index[i]];
1253 if (atoms->pdbinfo != NULL)
1255 useatoms.pdbinfo[i] = atoms->pdbinfo[index[i]];
1257 useatoms.nres = std::max(useatoms.nres, useatoms.atom[i].resind+1);
1259 useatoms.nr = nout;
1261 /* select what to read */
1262 if (ftp == efTRR)
1264 flags = TRX_READ_X;
1266 else
1268 flags = TRX_NEED_X;
1270 if (bVels)
1272 flags = flags | TRX_READ_V;
1274 if (bForce)
1276 flags = flags | TRX_READ_F;
1279 /* open trx file for reading */
1280 bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
1281 if (fr.bPrec)
1283 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1285 if (bNeedPrec)
1287 if (bSetPrec || !fr.bPrec)
1289 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1291 else
1293 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1297 if (bHaveFirstFrame)
1299 set_trxframe_ePBC(&fr, ePBC);
1301 natoms = fr.natoms;
1303 if (bSetTime)
1305 tshift = tzero-fr.time;
1307 else
1309 tzero = fr.time;
1312 bCopy = FALSE;
1313 if (bIndex)
1315 /* check if index is meaningful */
1316 for (i = 0; i < nout; i++)
1318 if (index[i] >= natoms)
1320 gmx_fatal(FARGS,
1321 "Index[%d] %d is larger than the number of atoms in the\n"
1322 "trajectory file (%d). There is a mismatch in the contents\n"
1323 "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1325 bCopy = bCopy || (i != index[i]);
1329 /* open output for writing */
1330 std::strcpy(filemode, "w");
1331 switch (ftp)
1333 case efTNG:
1334 trjtools_gmx_prepare_tng_writing(out_file,
1335 filemode[0],
1336 trxin,
1337 &trxout,
1338 NULL,
1339 nout,
1340 mtop,
1341 index,
1342 grpnm);
1343 break;
1344 case efXTC:
1345 case efTRR:
1346 out = NULL;
1347 if (!bSplit && !bSubTraj)
1349 trxout = open_trx(out_file, filemode);
1351 break;
1352 case efGRO:
1353 case efG96:
1354 case efPDB:
1355 if (( !bSeparate && !bSplit ) && !bSubTraj)
1357 out = gmx_ffopen(out_file, filemode);
1359 break;
1360 default:
1361 gmx_incons("Illegal output file format");
1364 if (bCopy)
1366 snew(xmem, nout);
1367 if (bVels)
1369 snew(vmem, nout);
1371 if (bForce)
1373 snew(fmem, nout);
1377 /* Start the big loop over frames */
1378 file_nr = 0;
1379 frame = 0;
1380 outframe = 0;
1381 model_nr = 0;
1382 bDTset = FALSE;
1384 /* Main loop over frames */
1387 if (!fr.bStep)
1389 /* set the step */
1390 fr.step = newstep;
1391 newstep++;
1393 if (bSubTraj)
1395 /*if (frame >= clust->clust->nra)
1396 gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1397 if (frame > clust->maxframe)
1399 my_clust = -1;
1401 else
1403 my_clust = clust->inv_clust[frame];
1405 if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1406 (my_clust == -1))
1408 my_clust = -1;
1412 if (bSetBox)
1414 /* generate new box */
1415 if (fr.bBox == FALSE)
1417 clear_mat(fr.box);
1419 for (m = 0; m < DIM; m++)
1421 if (newbox[m] >= 0)
1423 fr.box[m][m] = newbox[m];
1425 else
1427 if (fr.bBox == FALSE)
1429 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1435 if (bTrans)
1437 for (i = 0; i < natoms; i++)
1439 rvec_inc(fr.x[i], trans);
1443 if (bTDump)
1445 /* determine timestep */
1446 if (t0 == -1)
1448 t0 = fr.time;
1450 else
1452 if (!bDTset)
1454 dt = fr.time-t0;
1455 bDTset = TRUE;
1458 /* This is not very elegant, as one can not dump a frame after
1459 * a timestep with is more than twice as small as the first one. */
1460 bDumpFrame = (fr.time > tdump-0.5*dt) && (fr.time <= tdump+0.5*dt);
1462 else
1464 bDumpFrame = FALSE;
1467 /* determine if an atom jumped across the box and reset it if so */
1468 if (bNoJump && (bTPS || frame != 0))
1470 for (d = 0; d < DIM; d++)
1472 hbox[d] = 0.5*fr.box[d][d];
1474 for (i = 0; i < natoms; i++)
1476 if (bReset)
1478 rvec_dec(fr.x[i], x_shift);
1480 for (m = DIM-1; m >= 0; m--)
1482 if (hbox[m] > 0)
1484 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1486 for (d = 0; d <= m; d++)
1488 fr.x[i][d] += fr.box[m][d];
1491 while (fr.x[i][m]-xp[i][m] > hbox[m])
1493 for (d = 0; d <= m; d++)
1495 fr.x[i][d] -= fr.box[m][d];
1502 else if (bCluster)
1504 calc_pbc_cluster(ecenter, ifit, &top, ePBC, fr.x, ind_fit, fr.box);
1507 if (bPFit)
1509 /* Now modify the coords according to the flags,
1510 for normal fit, this is only done for output frames */
1511 if (bRmPBC)
1513 gmx_rmpbc_trxfr(gpbc, &fr);
1516 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1517 do_fit(natoms, w_rls, xp, fr.x);
1520 /* store this set of coordinates for future use */
1521 if (bPFit || bNoJump)
1523 if (xp == NULL)
1525 snew(xp, natoms);
1527 for (i = 0; (i < natoms); i++)
1529 copy_rvec(fr.x[i], xp[i]);
1530 rvec_inc(fr.x[i], x_shift);
1534 if (frindex)
1536 /* see if we have a frame from the frame index group */
1537 for (i = 0; i < nrfri && !bDumpFrame; i++)
1539 bDumpFrame = frame == frindex[i];
1542 if (debug && bDumpFrame)
1544 fprintf(debug, "dumping %d\n", frame);
1547 bWriteFrame =
1548 ( ( !bTDump && !frindex && frame % skip_nr == 0 ) || bDumpFrame );
1550 if (bWriteFrame && (bDropUnder || bDropOver))
1552 while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1554 drop0 = drop1;
1555 drop1++;
1557 if (std::abs(dropval[0][drop0] - fr.time)
1558 < std::abs(dropval[0][drop1] - fr.time))
1560 dropuse = drop0;
1562 else
1564 dropuse = drop1;
1566 if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1567 (bDropOver && dropval[1][dropuse] > dropover))
1569 bWriteFrame = FALSE;
1573 if (bWriteFrame)
1575 /* We should avoid modifying the input frame,
1576 * but since here we don't have the output frame yet,
1577 * we introduce a temporary output frame time variable.
1579 real frout_time;
1581 frout_time = fr.time;
1583 /* calc new time */
1584 if (bTimeStep)
1586 frout_time = tzero + frame*timestep;
1588 else
1589 if (bSetTime)
1591 frout_time += tshift;
1594 if (bTDump)
1596 fprintf(stderr, "\nDumping frame at t= %g %s\n",
1597 output_env_conv_time(oenv, frout_time), output_env_get_time_unit(oenv));
1600 /* check for writing at each delta_t */
1601 bDoIt = (delta_t == 0);
1602 if (!bDoIt)
1604 if (!bRound)
1606 bDoIt = bRmod(frout_time, tzero, delta_t);
1608 else
1610 /* round() is not C89 compatible, so we do this: */
1611 bDoIt = bRmod(std::floor(frout_time+0.5), std::floor(tzero+0.5),
1612 std::floor(delta_t+0.5));
1616 if (bDoIt || bTDump)
1618 /* print sometimes */
1619 if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1621 fprintf(stderr, " -> frame %6d time %8.3f \r",
1622 outframe, output_env_conv_time(oenv, frout_time));
1625 if (!bPFit)
1627 /* Now modify the coords according to the flags,
1628 for PFit we did this already! */
1630 if (bRmPBC)
1632 gmx_rmpbc_trxfr(gpbc, &fr);
1635 if (bReset)
1637 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1638 if (bFit)
1640 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1642 if (!bCenter)
1644 for (i = 0; i < natoms; i++)
1646 rvec_inc(fr.x[i], x_shift);
1651 if (bCenter)
1653 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1657 if (bPBCcomAtom)
1659 switch (unitcell_enum)
1661 case euRect:
1662 put_atoms_in_box(ePBC, fr.box, natoms, fr.x);
1663 break;
1664 case euTric:
1665 put_atoms_in_triclinic_unitcell(ecenter, fr.box, natoms, fr.x);
1666 break;
1667 case euCompact:
1668 put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1669 natoms, fr.x);
1670 break;
1673 if (bPBCcomRes)
1675 put_residue_com_in_box(unitcell_enum, ecenter,
1676 natoms, atoms->atom, ePBC, fr.box, fr.x);
1678 if (bPBCcomMol)
1680 put_molecule_com_in_box(unitcell_enum, ecenter,
1681 &top.mols,
1682 natoms, atoms->atom, ePBC, fr.box, fr.x);
1684 /* Copy the input trxframe struct to the output trxframe struct */
1685 frout = fr;
1686 frout.time = frout_time;
1687 frout.bV = (frout.bV && bVels);
1688 frout.bF = (frout.bF && bForce);
1689 frout.natoms = nout;
1690 if (bNeedPrec && (bSetPrec || !fr.bPrec))
1692 frout.bPrec = TRUE;
1693 frout.prec = prec;
1695 if (bCopy)
1697 frout.x = xmem;
1698 if (frout.bV)
1700 frout.v = vmem;
1702 if (frout.bF)
1704 frout.f = fmem;
1706 for (i = 0; i < nout; i++)
1708 copy_rvec(fr.x[index[i]], frout.x[i]);
1709 if (frout.bV)
1711 copy_rvec(fr.v[index[i]], frout.v[i]);
1713 if (frout.bF)
1715 copy_rvec(fr.f[index[i]], frout.f[i]);
1720 if (opt2parg_bSet("-shift", NPA, pa))
1722 for (i = 0; i < nout; i++)
1724 for (d = 0; d < DIM; d++)
1726 frout.x[i][d] += outframe*shift[d];
1731 if (!bRound)
1733 bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1735 else
1737 /* round() is not C89 compatible, so we do this: */
1738 bSplitHere = bSplit && bRmod(std::floor(frout.time+0.5),
1739 std::floor(tzero+0.5),
1740 std::floor(split_t+0.5));
1742 if (bSeparate || bSplitHere)
1744 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1747 switch (ftp)
1749 case efTNG:
1750 write_tng_frame(trxout, &frout);
1751 // TODO when trjconv behaves better: work how to read and write lambda
1752 break;
1753 case efTRR:
1754 case efXTC:
1755 if (bSplitHere)
1757 if (trxout)
1759 close_trx(trxout);
1761 trxout = open_trx(out_file2, filemode);
1763 if (bSubTraj)
1765 if (my_clust != -1)
1767 char buf[STRLEN];
1768 if (clust_status_id[my_clust] == -1)
1770 sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1771 clust_status[my_clust] = open_trx(buf, "w");
1772 clust_status_id[my_clust] = 1;
1773 ntrxopen++;
1775 else if (clust_status_id[my_clust] == -2)
1777 gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1778 clust->grpname[my_clust], ntrxopen, frame,
1779 my_clust);
1781 write_trxframe(clust_status[my_clust], &frout, gc);
1782 nfwritten[my_clust]++;
1783 if (nfwritten[my_clust] ==
1784 (clust->clust->index[my_clust+1]-
1785 clust->clust->index[my_clust]))
1787 close_trx(clust_status[my_clust]);
1788 clust_status[my_clust] = NULL;
1789 clust_status_id[my_clust] = -2;
1790 ntrxopen--;
1791 if (ntrxopen < 0)
1793 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1798 else
1800 write_trxframe(trxout, &frout, gc);
1802 break;
1803 case efGRO:
1804 case efG96:
1805 case efPDB:
1806 sprintf(title, "Generated by trjconv : %s t= %9.5f",
1807 top_title, frout.time);
1808 if (bSeparate || bSplitHere)
1810 out = gmx_ffopen(out_file2, "w");
1812 switch (ftp)
1814 case efGRO:
1815 write_hconf_p(out, title, &useatoms, prec2ndec(frout.prec),
1816 frout.x, frout.bV ? frout.v : NULL, frout.box);
1817 break;
1818 case efPDB:
1819 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1820 sprintf(title, "%s t= %9.5f", top_title, frout.time);
1821 /* if reading from pdb, we want to keep the original
1822 model numbering else we write the output frame
1823 number plus one, because model 0 is not allowed in pdb */
1824 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1826 model_nr = fr.step;
1828 else
1830 model_nr++;
1832 write_pdbfile(out, title, &useatoms, frout.x,
1833 frout.ePBC, frout.box, ' ', model_nr, gc, TRUE);
1834 break;
1835 case efG96:
1836 frout.title = title;
1837 if (bSeparate || bTDump)
1839 frout.bTitle = TRUE;
1840 if (bTPS)
1842 frout.bAtoms = TRUE;
1844 frout.atoms = &useatoms;
1845 frout.bStep = FALSE;
1846 frout.bTime = FALSE;
1848 else
1850 frout.bTitle = (outframe == 0);
1851 frout.bAtoms = FALSE;
1852 frout.bStep = TRUE;
1853 frout.bTime = TRUE;
1855 write_g96_conf(out, &frout, -1, NULL);
1857 if (bSeparate || bSplitHere)
1859 gmx_ffclose(out);
1860 out = NULL;
1862 break;
1863 default:
1864 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1866 if (bSeparate || bSplitHere)
1868 file_nr++;
1871 /* execute command */
1872 if (bExec)
1874 char c[255];
1875 sprintf(c, "%s %d", exec_command, file_nr-1);
1876 /*fprintf(stderr,"Executing '%s'\n",c);*/
1877 if (0 != system(c))
1879 gmx_fatal(FARGS, "Error executing command: %s", c);
1882 outframe++;
1885 frame++;
1886 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1888 while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1891 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1893 fprintf(stderr, "\nWARNING no output, "
1894 "last frame read at t=%g\n", fr.time);
1896 fprintf(stderr, "\n");
1898 close_trj(trxin);
1899 sfree(outf_base);
1901 if (bRmPBC)
1903 gmx_rmpbc_done(gpbc);
1906 if (trxout)
1908 close_trx(trxout);
1910 else if (out != NULL)
1912 gmx_ffclose(out);
1914 if (bSubTraj)
1916 for (i = 0; (i < clust->clust->nr); i++)
1918 if (clust_status_id[i] >= 0)
1920 close_trx(clust_status[i]);
1926 sfree(mtop);
1928 do_view(oenv, out_file, NULL);
1930 return 0;