Moved second half of gmxana tools to C++
[gromacs.git] / src / gromacs / gmxana / gmx_trjconv.cpp
bloba5e9913ee663bffb0412fd101389c9cff780e5f6
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/fileio/confio.h"
47 #include "gromacs/fileio/g96io.h"
48 #include "gromacs/fileio/gmxfio.h"
49 #include "gromacs/fileio/groio.h"
50 #include "gromacs/fileio/pdbio.h"
51 #include "gromacs/fileio/tngio_for_tools.h"
52 #include "gromacs/fileio/tpxio.h"
53 #include "gromacs/fileio/trrio.h"
54 #include "gromacs/fileio/trxio.h"
55 #include "gromacs/fileio/xtcio.h"
56 #include "gromacs/fileio/xvgr.h"
57 #include "gromacs/gmxana/gmx_ana.h"
58 #include "gromacs/legacyheaders/copyrite.h"
59 #include "gromacs/legacyheaders/macros.h"
60 #include "gromacs/legacyheaders/names.h"
61 #include "gromacs/legacyheaders/typedefs.h"
62 #include "gromacs/legacyheaders/viewit.h"
63 #include "gromacs/math/do_fit.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/pbcutil/rmpbc.h"
67 #include "gromacs/topology/index.h"
68 #include "gromacs/topology/topology.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/futil.h"
71 #include "gromacs/utility/smalloc.h"
73 enum {
74 euSel, euRect, euTric, euCompact, euNR
78 static void calc_pbc_cluster(int ecenter, int nrefat, t_topology *top, int ePBC,
79 rvec x[], atom_id index[], matrix box)
81 int m, i, j, j0, j1, jj, ai, aj;
82 int imin, jmin;
83 real fac, min_dist2;
84 rvec dx, xtest, box_center;
85 int nmol, imol_center;
86 atom_id *molind;
87 gmx_bool *bMol, *bTmp;
88 rvec *m_com, *m_shift;
89 t_pbc pbc;
90 int *cluster;
91 int *added;
92 int ncluster, nadded;
93 real tmp_r2;
95 calc_box_center(ecenter, box, box_center);
97 /* Initiate the pbc structure */
98 std::memset(&pbc, 0, sizeof(pbc));
99 set_pbc(&pbc, ePBC, box);
101 /* Convert atom index to molecular */
102 nmol = top->mols.nr;
103 molind = top->mols.index;
104 snew(bMol, nmol);
105 snew(m_com, nmol);
106 snew(m_shift, nmol);
107 snew(cluster, nmol);
108 snew(added, nmol);
109 snew(bTmp, top->atoms.nr);
111 for (i = 0; (i < nrefat); i++)
113 /* Mark all molecules in the index */
114 ai = index[i];
115 bTmp[ai] = TRUE;
116 /* Binary search assuming the molecules are sorted */
117 j0 = 0;
118 j1 = nmol-1;
119 while (j0 < j1)
121 if (ai < molind[j0+1])
123 j1 = j0;
125 else if (ai >= molind[j1])
127 j0 = j1;
129 else
131 jj = (j0+j1)/2;
132 if (ai < molind[jj+1])
134 j1 = jj;
136 else
138 j0 = jj;
142 bMol[j0] = TRUE;
144 /* Double check whether all atoms in all molecules that are marked are part
145 * of the cluster. Simultaneously compute the center of geometry.
147 min_dist2 = 10*sqr(trace(box));
148 imol_center = -1;
149 ncluster = 0;
150 for (i = 0; i < nmol; i++)
152 for (j = molind[i]; j < molind[i+1]; j++)
154 if (bMol[i] && !bTmp[j])
156 gmx_fatal(FARGS, "Molecule %d marked for clustering but not atom %d in it - check your index!", i+1, j+1);
158 else if (!bMol[i] && bTmp[j])
160 gmx_fatal(FARGS, "Atom %d marked for clustering but not molecule %d - this is an internal error...", j+1, i+1);
162 else if (bMol[i])
164 /* Make molecule whole, move 2nd and higher atom to same periodicity as 1st atom in molecule */
165 if (j > molind[i])
167 pbc_dx(&pbc, x[j], x[j-1], dx);
168 rvec_add(x[j-1], dx, x[j]);
170 /* Compute center of geometry of molecule - m_com[i] was zeroed when we did snew() on it! */
171 rvec_inc(m_com[i], x[j]);
174 if (bMol[i])
176 /* Normalize center of geometry */
177 fac = 1.0/(molind[i+1]-molind[i]);
178 for (m = 0; (m < DIM); m++)
180 m_com[i][m] *= fac;
182 /* Determine which molecule is closest to the center of the box */
183 pbc_dx(&pbc, box_center, m_com[i], dx);
184 tmp_r2 = iprod(dx, dx);
186 if (tmp_r2 < min_dist2)
188 min_dist2 = tmp_r2;
189 imol_center = i;
191 cluster[ncluster++] = i;
194 sfree(bTmp);
196 if (ncluster <= 0)
198 fprintf(stderr, "No molecules selected in the cluster\n");
199 return;
201 else if (imol_center == -1)
203 fprintf(stderr, "No central molecules could be found\n");
204 return;
207 nadded = 0;
208 added[nadded++] = imol_center;
209 bMol[imol_center] = FALSE;
211 while (nadded < ncluster)
213 /* Find min distance between cluster molecules and those remaining to be added */
214 min_dist2 = 10*sqr(trace(box));
215 imin = -1;
216 jmin = -1;
217 /* Loop over added mols */
218 for (i = 0; i < nadded; i++)
220 ai = added[i];
221 /* Loop over all mols */
222 for (j = 0; j < ncluster; j++)
224 aj = cluster[j];
225 /* check those remaining to be added */
226 if (bMol[aj])
228 pbc_dx(&pbc, m_com[aj], m_com[ai], dx);
229 tmp_r2 = iprod(dx, dx);
230 if (tmp_r2 < min_dist2)
232 min_dist2 = tmp_r2;
233 imin = ai;
234 jmin = aj;
240 /* Add the best molecule */
241 added[nadded++] = jmin;
242 bMol[jmin] = FALSE;
243 /* Calculate the shift from the ai molecule */
244 pbc_dx(&pbc, m_com[jmin], m_com[imin], dx);
245 rvec_add(m_com[imin], dx, xtest);
246 rvec_sub(xtest, m_com[jmin], m_shift[jmin]);
247 rvec_inc(m_com[jmin], m_shift[jmin]);
249 for (j = molind[jmin]; j < molind[jmin+1]; j++)
251 rvec_inc(x[j], m_shift[jmin]);
253 fprintf(stdout, "\rClustering iteration %d of %d...", nadded, ncluster);
254 fflush(stdout);
257 sfree(added);
258 sfree(cluster);
259 sfree(bMol);
260 sfree(m_com);
261 sfree(m_shift);
263 fprintf(stdout, "\n");
266 static void put_molecule_com_in_box(int unitcell_enum, int ecenter,
267 t_block *mols,
268 int natoms, t_atom atom[],
269 int ePBC, matrix box, rvec x[])
271 atom_id i, j;
272 int d;
273 rvec com, new_com, shift, box_center;
274 real m;
275 double mtot;
276 t_pbc pbc;
278 calc_box_center(ecenter, box, box_center);
279 set_pbc(&pbc, ePBC, box);
280 if (mols->nr <= 0)
282 gmx_fatal(FARGS, "There are no molecule descriptions. I need a .tpr file for this pbc option.");
284 for (i = 0; (i < mols->nr); i++)
286 /* calc COM */
287 clear_rvec(com);
288 mtot = 0;
289 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
291 m = atom[j].m;
292 for (d = 0; d < DIM; d++)
294 com[d] += m*x[j][d];
296 mtot += m;
298 /* calculate final COM */
299 svmul(1.0/mtot, com, com);
301 /* check if COM is outside box */
302 copy_rvec(com, new_com);
303 switch (unitcell_enum)
305 case euRect:
306 put_atoms_in_box(ePBC, box, 1, &new_com);
307 break;
308 case euTric:
309 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
310 break;
311 case euCompact:
312 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
313 break;
315 rvec_sub(new_com, com, shift);
316 if (norm2(shift) > 0)
318 if (debug)
320 fprintf(debug, "\nShifting position of molecule %d "
321 "by %8.3f %8.3f %8.3f\n", i+1,
322 shift[XX], shift[YY], shift[ZZ]);
324 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
326 rvec_inc(x[j], shift);
332 static void put_residue_com_in_box(int unitcell_enum, int ecenter,
333 int natoms, t_atom atom[],
334 int ePBC, matrix box, rvec x[])
336 atom_id i, j, res_start, res_end;
337 int d, presnr;
338 real m;
339 double mtot;
340 rvec box_center, com, new_com, shift;
342 calc_box_center(ecenter, box, box_center);
344 presnr = NOTSET;
345 res_start = 0;
346 clear_rvec(com);
347 mtot = 0;
348 for (i = 0; i < natoms+1; i++)
350 if (i == natoms || (presnr != atom[i].resind && presnr != NOTSET))
352 /* calculate final COM */
353 res_end = i;
354 svmul(1.0/mtot, com, com);
356 /* check if COM is outside box */
357 copy_rvec(com, new_com);
358 switch (unitcell_enum)
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;
370 rvec_sub(new_com, com, shift);
371 if (norm2(shift))
373 if (debug)
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]);
379 for (j = res_start; j < res_end; j++)
381 rvec_inc(x[j], shift);
384 clear_rvec(com);
385 mtot = 0;
387 /* remember start of new residue */
388 res_start = i;
390 if (i < natoms)
392 /* calc COM */
393 m = atom[i].m;
394 for (d = 0; d < DIM; d++)
396 com[d] += m*x[i][d];
398 mtot += m;
400 presnr = atom[i].resind;
405 static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, atom_id ci[])
407 int i, m, ai;
408 rvec cmin, cmax, box_center, dx;
410 if (nc > 0)
412 copy_rvec(x[ci[0]], cmin);
413 copy_rvec(x[ci[0]], cmax);
414 for (i = 0; i < nc; i++)
416 ai = ci[i];
417 for (m = 0; m < DIM; m++)
419 if (x[ai][m] < cmin[m])
421 cmin[m] = x[ai][m];
423 else if (x[ai][m] > cmax[m])
425 cmax[m] = x[ai][m];
429 calc_box_center(ecenter, box, box_center);
430 for (m = 0; m < DIM; m++)
432 dx[m] = box_center[m]-(cmin[m]+cmax[m])*0.5;
435 for (i = 0; i < n; i++)
437 rvec_inc(x[i], dx);
442 static void mk_filenm(char *base, const char *ext, int ndigit, int file_nr,
443 char out_file[])
445 char nbuf[128];
446 int nd = 0, fnr;
448 std::strcpy(out_file, base);
449 fnr = file_nr;
452 fnr /= 10;
453 nd++;
455 while (fnr > 0);
457 if (nd < ndigit)
459 std::strncat(out_file, "00000000000", ndigit-nd);
461 sprintf(nbuf, "%d.", file_nr);
462 std::strcat(out_file, nbuf);
463 std::strcat(out_file, ext);
466 void check_trr(const char *fn)
468 if (fn2ftp(fn) != efTRR)
470 gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
474 void do_trunc(const char *fn, real t0)
476 t_fileio *in;
477 FILE *fp;
478 gmx_bool bStop, bOK;
479 gmx_trr_header_t sh;
480 gmx_off_t fpos;
481 char yesno[256];
482 int j;
483 real t = 0;
485 if (t0 == -1)
487 gmx_fatal(FARGS, "You forgot to set the truncation time");
490 /* Check whether this is a .trr file */
491 check_trr(fn);
493 in = gmx_trr_open(fn, "r");
494 fp = gmx_fio_getfp(in);
495 if (fp == NULL)
497 fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
498 gmx_trr_close(in);
500 else
502 j = 0;
503 fpos = gmx_fio_ftell(in);
504 bStop = FALSE;
505 while (!bStop && gmx_trr_read_frame_header(in, &sh, &bOK))
507 gmx_trr_read_frame_data(in, &sh, NULL, NULL, NULL, NULL);
508 fpos = gmx_ftell(fp);
509 t = sh.t;
510 if (t >= t0)
512 gmx_fseek(fp, fpos, SEEK_SET);
513 bStop = TRUE;
516 if (bStop)
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))
523 gmx_fatal(FARGS, "Error reading user input");
525 if (std::strcmp(yesno, "YES") == 0)
527 fprintf(stderr, "Once again, I'm gonna DO this...\n");
528 gmx_trr_close(in);
529 if (0 != gmx_truncate(fn, fpos))
531 gmx_fatal(FARGS, "Error truncating file %s", fn);
534 else
536 fprintf(stderr, "Ok, I'll forget about it\n");
539 else
541 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
542 gmx_trr_close(in);
547 /*! \brief Read a full molecular topology if useful and available.
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.
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.
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)
571 gmx_mtop_t *mtop = NULL;
573 if (fn2bTPX(tps_file) &&
574 efTNG != fn2ftp(input_file) &&
575 efTNG == fn2ftp(output_file))
577 int temp_natoms = -1;
578 snew(mtop, 1);
579 read_tpx(tps_file, NULL, NULL, &temp_natoms,
580 NULL, NULL, NULL, mtop);
583 return mtop;
586 int gmx_trjconv(int argc, char *argv[])
588 const char *desc[] = {
589 "[THISMODULE] can convert trajectory files in many ways:",
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.",
608 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
609 "[PAR]",
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]",
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]",
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]",
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]",
648 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
649 "treatment:",
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.",
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.",
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]",
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]",
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",
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]",
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]",
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]",
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]",
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."
730 int pbc_enum;
731 enum
733 epSel,
734 epNone,
735 epComMol,
736 epComRes,
737 epComAtom,
738 epNojump,
739 epCluster,
740 epWhole,
741 epNR
743 const char *pbc_opt[epNR + 1] =
745 NULL, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
746 NULL
749 int unitcell_enum;
750 const char *unitcell_opt[euNR+1] =
751 { NULL, "rect", "tric", "compact", NULL };
753 enum
755 ecSel, ecTric, ecRect, ecZero, ecNR
757 const char *center_opt[ecNR+1] =
758 { NULL, "tric", "rect", "zero", NULL };
759 int ecenter;
761 int fit_enum;
762 enum
764 efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
766 const char *fit[efNR + 1] =
768 NULL, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
769 "progressive", NULL
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;
781 t_pargs
782 pa[] =
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" }
859 #define NPA asize(pa)
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;
869 real *w_rls = NULL;
870 int m, i, d, frame, outframe, natoms, nout, ncent, 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, my_clust = -1;
883 atom_id *ind_fit;
884 char *gn_fit;
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, 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], filemode[5];
907 gmx_bool bWarnCompact = FALSE;
908 const char *warn;
909 output_env_t oenv;
911 t_filenm fnm[] = {
912 { efTRX, "-f", NULL, ffREAD },
913 { efTRO, "-o", NULL, ffWRITE },
914 { efTPS, NULL, NULL, ffOPTRD },
915 { efNDX, NULL, NULL, ffOPTRD },
916 { efNDX, "-fr", "frames", ffOPTRD },
917 { efNDX, "-sub", "cluster", ffOPTRD },
918 { efXVG, "-drop", "drop", ffOPTRD }
920 #define NFILE asize(fnm)
922 if (!parse_common_args(&argc, argv,
923 PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW |
924 PCA_TIME_UNIT,
925 NFILE, fnm, NPA, pa, asize(desc), desc,
926 0, NULL, &oenv))
928 return 0;
931 top_file = ftp2fn(efTPS, NFILE, fnm);
932 init_top(&top);
934 /* Check command line */
935 in_file = opt2fn("-f", NFILE, fnm);
937 if (ttrunc != -1)
939 do_trunc(in_file, ttrunc);
941 else
943 /* mark active cmdline options */
944 bSetBox = opt2parg_bSet("-box", NPA, pa);
945 bSetTime = opt2parg_bSet("-t0", NPA, pa);
946 bSetPrec = opt2parg_bSet("-ndec", NPA, pa);
947 bSetUR = opt2parg_bSet("-ur", NPA, pa);
948 bExec = opt2parg_bSet("-exec", NPA, pa);
949 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
950 bTDump = opt2parg_bSet("-dump", NPA, pa);
951 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
952 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
953 bTrans = opt2parg_bSet("-trans", NPA, pa);
954 bSplit = (split_t != 0);
956 /* parse enum options */
957 fit_enum = nenum(fit);
958 bFit = (fit_enum == efFit || fit_enum == efFitXY);
959 bReset = (fit_enum == efReset || fit_enum == efResetXY);
960 bPFit = fit_enum == efPFit;
961 pbc_enum = nenum(pbc_opt);
962 bPBCWhole = pbc_enum == epWhole;
963 bPBCcomRes = pbc_enum == epComRes;
964 bPBCcomMol = pbc_enum == epComMol;
965 bPBCcomAtom = pbc_enum == epComAtom;
966 bNoJump = pbc_enum == epNojump;
967 bCluster = pbc_enum == epCluster;
968 bPBC = pbc_enum != epNone;
969 unitcell_enum = nenum(unitcell_opt);
970 ecenter = nenum(center_opt) - ecTric;
972 /* set and check option dependencies */
973 if (bPFit)
975 bFit = TRUE; /* for pfit, fit *must* be set */
977 if (bFit)
979 bReset = TRUE; /* for fit, reset *must* be set */
981 nfitdim = 0;
982 if (bFit || bReset)
984 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
986 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
988 if (bSetUR)
990 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
992 fprintf(stderr,
993 "WARNING: Option for unitcell representation (-ur %s)\n"
994 " only has effect in combination with -pbc %s, %s or %s.\n"
995 " Ingoring unitcell representation.\n\n",
996 unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
999 if (bFit && bPBC)
1001 gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
1002 "Please do the PBC condition treatment first and then run trjconv in a second step\n"
1003 "for the rotational fit.\n"
1004 "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1005 "results!");
1008 /* ndec is in nr of decimal places, prec is a multiplication factor: */
1009 prec = 1;
1010 for (i = 0; i < ndec; i++)
1012 prec *= 10;
1015 bIndex = ftp2bSet(efNDX, NFILE, fnm);
1018 /* Determine output type */
1019 out_file = opt2fn("-o", NFILE, fnm);
1020 ftp = fn2ftp(out_file);
1021 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
1022 bNeedPrec = (ftp == efXTC || ftp == efGRO);
1023 if (bVels)
1025 /* check if velocities are possible in input and output files */
1026 ftpin = fn2ftp(in_file);
1027 bVels = (ftp == efTRR || ftp == efGRO ||
1028 ftp == efG96 || ftp == efTNG)
1029 && (ftpin == efTRR || ftpin == efGRO ||
1030 ftpin == efG96 || ftpin == efTNG || ftpin == efCPT);
1032 if (bSeparate || bSplit)
1034 outf_ext = std::strrchr(out_file, '.');
1035 if (outf_ext == NULL)
1037 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1039 outf_base = gmx_strdup(out_file);
1040 outf_base[outf_ext - out_file] = '\0';
1043 bSubTraj = opt2bSet("-sub", NFILE, fnm);
1044 if (bSubTraj)
1046 if ((ftp != efXTC) && (ftp != efTRR))
1048 /* It seems likely that other trajectory file types
1049 * could work here. */
1050 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1051 "xtc and trr");
1053 clust = cluster_index(NULL, opt2fn("-sub", NFILE, fnm));
1055 /* Check for number of files disabled, as FOPEN_MAX is not the correct
1056 * number to check for. In my linux box it is only 16.
1058 if (0 && (clust->clust->nr > FOPEN_MAX-4))
1060 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1061 " trajectories.\ntry splitting the index file in %d parts.\n"
1062 "FOPEN_MAX = %d",
1063 clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1065 gmx_warning("The -sub option could require as many open output files as there are\n"
1066 "index groups in the file (%d). If you get I/O errors opening new files,\n"
1067 "try reducing the number of index groups in the file, and perhaps\n"
1068 "using trjconv -sub several times on different chunks of your index file.\n",
1069 clust->clust->nr);
1071 snew(clust_status, clust->clust->nr);
1072 snew(clust_status_id, clust->clust->nr);
1073 snew(nfwritten, clust->clust->nr);
1074 for (i = 0; (i < clust->clust->nr); i++)
1076 clust_status[i] = NULL;
1077 clust_status_id[i] = -1;
1079 bSeparate = bSplit = FALSE;
1081 /* skipping */
1082 if (skip_nr <= 0)
1086 mtop = read_mtop_for_tng(top_file, in_file, out_file);
1088 /* Determine whether to read a topology */
1089 bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1090 bRmPBC || bReset || bPBCcomMol || bCluster ||
1091 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1093 /* Determine if when can read index groups */
1094 bIndex = (bIndex || bTPS);
1096 if (bTPS)
1098 read_tps_conf(top_file, top_title, &top, &ePBC, &xp, NULL, top_box,
1099 bReset || bPBCcomRes);
1100 atoms = &top.atoms;
1102 if (0 == top.mols.nr && (bCluster || bPBCcomMol))
1104 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1107 /* top_title is only used for gro and pdb,
1108 * the header in such a file is top_title t= ...
1109 * to prevent a double t=, remove it from top_title
1111 if ((charpt = std::strstr(top_title, " t= ")))
1113 charpt[0] = '\0';
1116 if (bCONECT)
1118 gc = gmx_conect_generate(&top);
1120 if (bRmPBC)
1122 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
1126 /* get frame number index */
1127 frindex = NULL;
1128 if (opt2bSet("-fr", NFILE, fnm))
1130 printf("Select groups of frame number indices:\n");
1131 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, (atom_id **)&frindex, &frname);
1132 if (debug)
1134 for (i = 0; i < nrfri; i++)
1136 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1141 /* get index groups etc. */
1142 if (bReset)
1144 printf("Select group for %s fit\n",
1145 bFit ? "least squares" : "translational");
1146 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1147 1, &ifit, &ind_fit, &gn_fit);
1149 if (bFit)
1151 if (ifit < 2)
1153 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1155 else if (ifit == 3)
1157 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1161 else if (bCluster)
1163 printf("Select group for clustering\n");
1164 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1165 1, &ifit, &ind_fit, &gn_fit);
1168 if (bIndex)
1170 if (bCenter)
1172 printf("Select group for centering\n");
1173 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1174 1, &ncent, &cindex, &grpnm);
1176 printf("Select group for output\n");
1177 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1178 1, &nout, &index, &grpnm);
1180 else
1182 /* no index file, so read natoms from TRX */
1183 if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
1185 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1187 natoms = fr.natoms;
1188 close_trj(trxin);
1189 sfree(fr.x);
1190 snew(index, natoms);
1191 for (i = 0; i < natoms; i++)
1193 index[i] = i;
1195 nout = natoms;
1196 if (bCenter)
1198 ncent = nout;
1199 cindex = index;
1203 if (bReset)
1205 snew(w_rls, atoms->nr);
1206 for (i = 0; (i < ifit); i++)
1208 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1211 /* Restore reference structure and set to origin,
1212 store original location (to put structure back) */
1213 if (bRmPBC)
1215 gmx_rmpbc(gpbc, top.atoms.nr, top_box, xp);
1217 copy_rvec(xp[index[0]], x_shift);
1218 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, NULL, xp, w_rls);
1219 rvec_dec(x_shift, xp[index[0]]);
1221 else
1223 clear_rvec(x_shift);
1226 if (bDropUnder || bDropOver)
1228 /* Read the .xvg file with the drop values */
1229 fprintf(stderr, "\nReading drop file ...");
1230 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1231 fprintf(stderr, " %d time points\n", ndrop);
1232 if (ndrop == 0 || ncol < 2)
1234 gmx_fatal(FARGS, "Found no data points in %s",
1235 opt2fn("-drop", NFILE, fnm));
1237 drop0 = 0;
1238 drop1 = 0;
1241 /* Make atoms struct for output in GRO or PDB files */
1242 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1244 /* get memory for stuff to go in .pdb file, and initialize
1245 * the pdbinfo structure part if the input has it.
1247 init_t_atoms(&useatoms, atoms->nr, (atoms->pdbinfo != NULL));
1248 sfree(useatoms.resinfo);
1249 useatoms.resinfo = atoms->resinfo;
1250 for (i = 0; (i < nout); i++)
1252 useatoms.atomname[i] = atoms->atomname[index[i]];
1253 useatoms.atom[i] = atoms->atom[index[i]];
1254 if (atoms->pdbinfo != NULL)
1256 useatoms.pdbinfo[i] = atoms->pdbinfo[index[i]];
1258 useatoms.nres = std::max(useatoms.nres, useatoms.atom[i].resind+1);
1260 useatoms.nr = nout;
1262 /* select what to read */
1263 if (ftp == efTRR)
1265 flags = TRX_READ_X;
1267 else
1269 flags = TRX_NEED_X;
1271 if (bVels)
1273 flags = flags | TRX_READ_V;
1275 if (bForce)
1277 flags = flags | TRX_READ_F;
1280 /* open trx file for reading */
1281 bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
1282 if (fr.bPrec)
1284 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1286 if (bNeedPrec)
1288 if (bSetPrec || !fr.bPrec)
1290 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1292 else
1294 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1298 if (bHaveFirstFrame)
1300 set_trxframe_ePBC(&fr, ePBC);
1302 natoms = fr.natoms;
1304 if (bSetTime)
1306 tshift = tzero-fr.time;
1308 else
1310 tzero = fr.time;
1313 bCopy = FALSE;
1314 if (bIndex)
1316 /* check if index is meaningful */
1317 for (i = 0; i < nout; i++)
1319 if (index[i] >= natoms)
1321 gmx_fatal(FARGS,
1322 "Index[%d] %d is larger than the number of atoms in the\n"
1323 "trajectory file (%d). There is a mismatch in the contents\n"
1324 "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1326 bCopy = bCopy || (i != index[i]);
1330 /* open output for writing */
1331 std::strcpy(filemode, "w");
1332 switch (ftp)
1334 case efTNG:
1335 trjtools_gmx_prepare_tng_writing(out_file,
1336 filemode[0],
1337 trxin,
1338 &trxout,
1339 NULL,
1340 nout,
1341 mtop,
1342 index,
1343 grpnm);
1344 break;
1345 case efXTC:
1346 case efTRR:
1347 out = NULL;
1348 if (!bSplit && !bSubTraj)
1350 trxout = open_trx(out_file, filemode);
1352 break;
1353 case efGRO:
1354 case efG96:
1355 case efPDB:
1356 if (( !bSeparate && !bSplit ) && !bSubTraj)
1358 out = gmx_ffopen(out_file, filemode);
1360 break;
1361 default:
1362 gmx_incons("Illegal output file format");
1365 if (bCopy)
1367 snew(xmem, nout);
1368 if (bVels)
1370 snew(vmem, nout);
1372 if (bForce)
1374 snew(fmem, nout);
1378 /* Start the big loop over frames */
1379 file_nr = 0;
1380 frame = 0;
1381 outframe = 0;
1382 model_nr = 0;
1383 bDTset = FALSE;
1385 /* Main loop over frames */
1388 if (!fr.bStep)
1390 /* set the step */
1391 fr.step = newstep;
1392 newstep++;
1394 if (bSubTraj)
1396 /*if (frame >= clust->clust->nra)
1397 gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1398 if (frame > clust->maxframe)
1400 my_clust = -1;
1402 else
1404 my_clust = clust->inv_clust[frame];
1406 if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1407 (my_clust == NO_ATID))
1409 my_clust = -1;
1413 if (bSetBox)
1415 /* generate new box */
1416 if (fr.bBox == FALSE)
1418 clear_mat(fr.box);
1420 for (m = 0; m < DIM; m++)
1422 if (newbox[m] >= 0)
1424 fr.box[m][m] = newbox[m];
1426 else
1428 if (fr.bBox == FALSE)
1430 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1436 if (bTrans)
1438 for (i = 0; i < natoms; i++)
1440 rvec_inc(fr.x[i], trans);
1444 if (bTDump)
1446 /* determine timestep */
1447 if (t0 == -1)
1449 t0 = fr.time;
1451 else
1453 if (!bDTset)
1455 dt = fr.time-t0;
1456 bDTset = TRUE;
1459 /* This is not very elegant, as one can not dump a frame after
1460 * a timestep with is more than twice as small as the first one. */
1461 bDumpFrame = (fr.time > tdump-0.5*dt) && (fr.time <= tdump+0.5*dt);
1463 else
1465 bDumpFrame = FALSE;
1468 /* determine if an atom jumped across the box and reset it if so */
1469 if (bNoJump && (bTPS || frame != 0))
1471 for (d = 0; d < DIM; d++)
1473 hbox[d] = 0.5*fr.box[d][d];
1475 for (i = 0; i < natoms; i++)
1477 if (bReset)
1479 rvec_dec(fr.x[i], x_shift);
1481 for (m = DIM-1; m >= 0; m--)
1483 if (hbox[m] > 0)
1485 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1487 for (d = 0; d <= m; d++)
1489 fr.x[i][d] += fr.box[m][d];
1492 while (fr.x[i][m]-xp[i][m] > hbox[m])
1494 for (d = 0; d <= m; d++)
1496 fr.x[i][d] -= fr.box[m][d];
1503 else if (bCluster)
1505 calc_pbc_cluster(ecenter, ifit, &top, ePBC, fr.x, ind_fit, fr.box);
1508 if (bPFit)
1510 /* Now modify the coords according to the flags,
1511 for normal fit, this is only done for output frames */
1512 if (bRmPBC)
1514 gmx_rmpbc_trxfr(gpbc, &fr);
1517 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1518 do_fit(natoms, w_rls, xp, fr.x);
1521 /* store this set of coordinates for future use */
1522 if (bPFit || bNoJump)
1524 if (xp == NULL)
1526 snew(xp, natoms);
1528 for (i = 0; (i < natoms); i++)
1530 copy_rvec(fr.x[i], xp[i]);
1531 rvec_inc(fr.x[i], x_shift);
1535 if (frindex)
1537 /* see if we have a frame from the frame index group */
1538 for (i = 0; i < nrfri && !bDumpFrame; i++)
1540 bDumpFrame = frame == frindex[i];
1543 if (debug && bDumpFrame)
1545 fprintf(debug, "dumping %d\n", frame);
1548 bWriteFrame =
1549 ( ( !bTDump && !frindex && frame % skip_nr == 0 ) || bDumpFrame );
1551 if (bWriteFrame && (bDropUnder || bDropOver))
1553 while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1555 drop0 = drop1;
1556 drop1++;
1558 if (std::abs(dropval[0][drop0] - fr.time)
1559 < std::abs(dropval[0][drop1] - fr.time))
1561 dropuse = drop0;
1563 else
1565 dropuse = drop1;
1567 if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1568 (bDropOver && dropval[1][dropuse] > dropover))
1570 bWriteFrame = FALSE;
1574 if (bWriteFrame)
1576 /* We should avoid modifying the input frame,
1577 * but since here we don't have the output frame yet,
1578 * we introduce a temporary output frame time variable.
1580 real frout_time;
1582 frout_time = fr.time;
1584 /* calc new time */
1585 if (bTimeStep)
1587 frout_time = tzero + frame*timestep;
1589 else
1590 if (bSetTime)
1592 frout_time += tshift;
1595 if (bTDump)
1597 fprintf(stderr, "\nDumping frame at t= %g %s\n",
1598 output_env_conv_time(oenv, frout_time), output_env_get_time_unit(oenv));
1601 /* check for writing at each delta_t */
1602 bDoIt = (delta_t == 0);
1603 if (!bDoIt)
1605 if (!bRound)
1607 bDoIt = bRmod(frout_time, tzero, delta_t);
1609 else
1611 /* round() is not C89 compatible, so we do this: */
1612 bDoIt = bRmod(std::floor(frout_time+0.5), std::floor(tzero+0.5),
1613 std::floor(delta_t+0.5));
1617 if (bDoIt || bTDump)
1619 /* print sometimes */
1620 if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1622 fprintf(stderr, " -> frame %6d time %8.3f \r",
1623 outframe, output_env_conv_time(oenv, frout_time));
1626 if (!bPFit)
1628 /* Now modify the coords according to the flags,
1629 for PFit we did this already! */
1631 if (bRmPBC)
1633 gmx_rmpbc_trxfr(gpbc, &fr);
1636 if (bReset)
1638 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1639 if (bFit)
1641 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1643 if (!bCenter)
1645 for (i = 0; i < natoms; i++)
1647 rvec_inc(fr.x[i], x_shift);
1652 if (bCenter)
1654 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1658 if (bPBCcomAtom)
1660 switch (unitcell_enum)
1662 case euRect:
1663 put_atoms_in_box(ePBC, fr.box, natoms, fr.x);
1664 break;
1665 case euTric:
1666 put_atoms_in_triclinic_unitcell(ecenter, fr.box, natoms, fr.x);
1667 break;
1668 case euCompact:
1669 warn = put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1670 natoms, fr.x);
1671 if (warn && !bWarnCompact)
1673 fprintf(stderr, "\n%s\n", warn);
1674 bWarnCompact = TRUE;
1676 break;
1679 if (bPBCcomRes)
1681 put_residue_com_in_box(unitcell_enum, ecenter,
1682 natoms, atoms->atom, ePBC, fr.box, fr.x);
1684 if (bPBCcomMol)
1686 put_molecule_com_in_box(unitcell_enum, ecenter,
1687 &top.mols,
1688 natoms, atoms->atom, ePBC, fr.box, fr.x);
1690 /* Copy the input trxframe struct to the output trxframe struct */
1691 frout = fr;
1692 frout.time = frout_time;
1693 frout.bV = (frout.bV && bVels);
1694 frout.bF = (frout.bF && bForce);
1695 frout.natoms = nout;
1696 if (bNeedPrec && (bSetPrec || !fr.bPrec))
1698 frout.bPrec = TRUE;
1699 frout.prec = prec;
1701 if (bCopy)
1703 frout.x = xmem;
1704 if (frout.bV)
1706 frout.v = vmem;
1708 if (frout.bF)
1710 frout.f = fmem;
1712 for (i = 0; i < nout; i++)
1714 copy_rvec(fr.x[index[i]], frout.x[i]);
1715 if (frout.bV)
1717 copy_rvec(fr.v[index[i]], frout.v[i]);
1719 if (frout.bF)
1721 copy_rvec(fr.f[index[i]], frout.f[i]);
1726 if (opt2parg_bSet("-shift", NPA, pa))
1728 for (i = 0; i < nout; i++)
1730 for (d = 0; d < DIM; d++)
1732 frout.x[i][d] += outframe*shift[d];
1737 if (!bRound)
1739 bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1741 else
1743 /* round() is not C89 compatible, so we do this: */
1744 bSplitHere = bSplit && bRmod(std::floor(frout.time+0.5),
1745 std::floor(tzero+0.5),
1746 std::floor(split_t+0.5));
1748 if (bSeparate || bSplitHere)
1750 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1753 switch (ftp)
1755 case efTNG:
1756 write_tng_frame(trxout, &frout);
1757 // TODO when trjconv behaves better: work how to read and write lambda
1758 break;
1759 case efTRR:
1760 case efXTC:
1761 if (bSplitHere)
1763 if (trxout)
1765 close_trx(trxout);
1767 trxout = open_trx(out_file2, filemode);
1769 if (bSubTraj)
1771 if (my_clust != -1)
1773 char buf[STRLEN];
1774 if (clust_status_id[my_clust] == -1)
1776 sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1777 clust_status[my_clust] = open_trx(buf, "w");
1778 clust_status_id[my_clust] = 1;
1779 ntrxopen++;
1781 else if (clust_status_id[my_clust] == -2)
1783 gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1784 clust->grpname[my_clust], ntrxopen, frame,
1785 my_clust);
1787 write_trxframe(clust_status[my_clust], &frout, gc);
1788 nfwritten[my_clust]++;
1789 if (nfwritten[my_clust] ==
1790 (clust->clust->index[my_clust+1]-
1791 clust->clust->index[my_clust]))
1793 close_trx(clust_status[my_clust]);
1794 clust_status[my_clust] = NULL;
1795 clust_status_id[my_clust] = -2;
1796 ntrxopen--;
1797 if (ntrxopen < 0)
1799 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1804 else
1806 write_trxframe(trxout, &frout, gc);
1808 break;
1809 case efGRO:
1810 case efG96:
1811 case efPDB:
1812 sprintf(title, "Generated by trjconv : %s t= %9.5f",
1813 top_title, frout.time);
1814 if (bSeparate || bSplitHere)
1816 out = gmx_ffopen(out_file2, "w");
1818 switch (ftp)
1820 case efGRO:
1821 write_hconf_p(out, title, &useatoms, prec2ndec(frout.prec),
1822 frout.x, frout.bV ? frout.v : NULL, frout.box);
1823 break;
1824 case efPDB:
1825 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1826 sprintf(title, "%s t= %9.5f", top_title, frout.time);
1827 /* if reading from pdb, we want to keep the original
1828 model numbering else we write the output frame
1829 number plus one, because model 0 is not allowed in pdb */
1830 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1832 model_nr = fr.step;
1834 else
1836 model_nr++;
1838 write_pdbfile(out, title, &useatoms, frout.x,
1839 frout.ePBC, frout.box, ' ', model_nr, gc, TRUE);
1840 break;
1841 case efG96:
1842 frout.title = title;
1843 if (bSeparate || bTDump)
1845 frout.bTitle = TRUE;
1846 if (bTPS)
1848 frout.bAtoms = TRUE;
1850 frout.atoms = &useatoms;
1851 frout.bStep = FALSE;
1852 frout.bTime = FALSE;
1854 else
1856 frout.bTitle = (outframe == 0);
1857 frout.bAtoms = FALSE;
1858 frout.bStep = TRUE;
1859 frout.bTime = TRUE;
1861 write_g96_conf(out, &frout, -1, NULL);
1863 if (bSeparate || bSplitHere)
1865 gmx_ffclose(out);
1866 out = NULL;
1868 break;
1869 default:
1870 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1872 if (bSeparate || bSplitHere)
1874 file_nr++;
1877 /* execute command */
1878 if (bExec)
1880 char c[255];
1881 sprintf(c, "%s %d", exec_command, file_nr-1);
1882 /*fprintf(stderr,"Executing '%s'\n",c);*/
1883 if (0 != system(c))
1885 gmx_fatal(FARGS, "Error executing command: %s", c);
1888 outframe++;
1891 frame++;
1892 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1894 while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1897 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1899 fprintf(stderr, "\nWARNING no output, "
1900 "last frame read at t=%g\n", fr.time);
1902 fprintf(stderr, "\n");
1904 close_trj(trxin);
1905 sfree(outf_base);
1907 if (bRmPBC)
1909 gmx_rmpbc_done(gpbc);
1912 if (trxout)
1914 close_trx(trxout);
1916 else if (out != NULL)
1918 gmx_ffclose(out);
1920 if (bSubTraj)
1922 for (i = 0; (i < clust->clust->nr); i++)
1924 if (clust_status_id[i] >= 0)
1926 close_trx(clust_status[i]);
1932 sfree(mtop);
1934 do_view(oenv, out_file, NULL);
1936 return 0;